11#ifndef PANZER_UNIQUEGLOBALINDEXER_EPETRAUTILITIES_IMPL_HPP
12#define PANZER_UNIQUEGLOBALINDEXER_EPETRAUTILITIES_IMPL_HPP
17#include "Teuchos_FancyOStream.hpp"
18#include "Teuchos_ArrayView.hpp"
19#include "Teuchos_CommHelpers.hpp"
21#include "Epetra_Map.h"
22#include "Epetra_IntVector.h"
23#include "Epetra_MultiVector.h"
24#include "Epetra_Import.h"
25#include "Epetra_MpiComm.h"
35template <
typename ScalarT,
typename ArrayT>
42 TEUCHOS_TEST_FOR_EXCEPTION(!ugi.
fieldInBlock(fieldName,blockId),std::runtime_error,
43 "panzer::updateGhostedDataReducedVector: field name = \""+fieldName+
"\" is not in element block = \"" +blockId +
"\"!");
45 Teuchos::RCP<const Map> dataMap = Teuchos::rcpFromRef(dataVector.
Map());
48 const std::vector<panzer::LocalOrdinal> & elements = ugi.
getElementBlock(blockId);
51 TEUCHOS_TEST_FOR_EXCEPTION(data.extent(0)!=elements.size(),std::runtime_error,
52 "panzer::updateGhostedDataReducedVector: data cell dimension does not match up with block cell count");
54 int rank = data.rank();
58 std::vector<panzer::GlobalOrdinal> gids;
59 for(std::size_t e=0;e<elements.size();e++) {
62 for(std::size_t f=0;f<fieldOffsets.size();f++) {
63 std::size_t localIndex = dataMap->LID(gids[fieldOffsets[f]]);
69 std::size_t entries = data.extent(2);
71 TEUCHOS_TEST_FOR_EXCEPTION(dataVector.
NumVectors()!=
static_cast<int>(entries),std::runtime_error,
72 "panzer::updateGhostedDataReducedVector: number of columns in data vector inconsistent with data array");
75 std::vector<panzer::GlobalOrdinal> gids;
76 for(std::size_t e=0;e<elements.size();e++) {
79 for(std::size_t f=0;f<fieldOffsets.size();f++) {
80 std::size_t localIndex = dataMap->LID(gids[fieldOffsets[f]]);
81 for(std::size_t v=0;v<entries;v++)
87 TEUCHOS_TEST_FOR_EXCEPTION(
true,std::runtime_error,
88 "panzer::updateGhostedDataReducedVector: data array rank must be 2 or 3");
91template <
typename ScalarT,
typename ArrayT>
92Teuchos::RCP<Epetra_MultiVector>
94 const std::map<std::string,ArrayT> & data)
const
96 TEUCHOS_ASSERT(data.size()>0);
98 int fieldNum =
ugi_->getFieldNum(fieldName);
99 std::vector<std::string> blockIds;
100 ugi_->getElementBlockIds(blockIds);
103 int rank = data.begin()->second.rank();
108 numCols = data.begin()->second.extent(2);
110 TEUCHOS_TEST_FOR_EXCEPTION(
true,std::runtime_error,
111 "ArrayToFieldVectorEpetra::getGhostedDataVector: data array must have rank 2 or 3. This array has rank " << rank <<
".");
125 Teuchos::RCP<MultiVector> finalReducedVec
126 = Teuchos::rcp(
new MultiVector(*reducedMap,numCols));
127 for(std::size_t b=0;b<blockIds.size();b++) {
128 std::string block = blockIds[b];
131 if(!
ugi_->fieldInBlock(fieldName,block))
135 typename std::map<std::string,ArrayT>::const_iterator blockItr = data.find(block);
136 TEUCHOS_TEST_FOR_EXCEPTION(blockItr==data.end(),std::runtime_error,
137 "ArrayToFieldVectorEpetra::getDataVector: can not find block \""+block+
"\".");
139 const ArrayT & d = blockItr->second;
140 updateGhostedDataReducedVectorEpetra<ScalarT,ArrayT>(fieldName,block,*
ugi_,d,*finalReducedVec);
152 Teuchos::RCP<MultiVector> finalVec
157 finalVec->Import(*finalReducedVec,importer,
Insert);
162template <
typename ScalarT,
typename ArrayT>
163Teuchos::RCP<Epetra_MultiVector>
165 const std::map<std::string,ArrayT> & data)
const
171 Teuchos::RCP<const MultiVector> sourceVec
172 = getGhostedDataVector<ScalarT,ArrayT>(fieldName,data);
175 int fieldNum =
ugi_->getFieldNum(fieldName);
176 Teuchos::RCP<const Map> destMap =
fieldMaps_[fieldNum];
182 Teuchos::RCP<MultiVector> destVec
183 = Teuchos::rcp(
new MultiVector(*destMap,sourceVec->NumVectors()));
187 destVec->Import(*sourceVec,importer,
Insert);
const Epetra_BlockMap & Map() const
int ReplaceMyValue(int MyRow, int VectorIndex, double ScalarValue)
Teuchos::RCP< const GlobalIndexer > ugi_
DOF mapping.
std::map< int, Teuchos::RCP< const Map > > gh_reducedFieldMaps_
ghosted field vector
Teuchos::RCP< const IntVector > gh_reducedFieldVector_
std::map< int, Teuchos::RCP< const Map > > gh_fieldMaps_
Maps for each field (as needed)
void buildFieldVector(const Epetra_IntVector &source) const
build unghosted field vector from ghosted field vector
Teuchos::RCP< const IntVector > gh_fieldVector_
ghosted reduced field vector
Teuchos::RCP< Epetra_MultiVector > getGhostedDataVector(const std::string &fieldName, const std::map< std::string, ArrayT > &data) const
Epetra_MultiVector MultiVector
Teuchos::RCP< const IntVector > fieldVector_
Maps for each field (as needed)
std::map< int, Teuchos::RCP< const Map > > fieldMaps_
(unghosted) field vector (as needed)
Teuchos::RCP< Epetra_MultiVector > getDataVector(const std::string &fieldName, const std::map< std::string, ArrayT > &data) const
virtual const std::vector< panzer::LocalOrdinal > & getElementBlock(const std::string &blockId) const =0
virtual const std::vector< int > & getGIDFieldOffsets(const std::string &blockId, int fieldNum) const =0
Use the field pattern so that you can find a particular field in the GIDs array.
virtual void getElementGIDs(panzer::LocalOrdinal localElmtId, std::vector< panzer::GlobalOrdinal > &gids, const std::string &blockIdHint="") const =0
Get the global IDs for a particular element. This function overwrites the gids variable.
virtual int getFieldNum(const std::string &str) const =0
Get the number used for access to this field.
virtual bool fieldInBlock(const std::string &field, const std::string &block) const =0
void updateGhostedDataReducedVectorEpetra(const std::string &fieldName, const std::string blockId, const GlobalIndexer &ugi, const ArrayT &data, Epetra_MultiVector &dataVector)
Teuchos::RCP< const Epetra_BlockMap > getFieldMapEpetra(int fieldNum, const Epetra_IntVector &fieldVector)