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"
35 template <
typename ScalarT,
typename ArrayT>
43 "panzer::updateGhostedDataReducedVector: field name = \""+fieldName+
"\" is not in element block = \"" +blockId +
"\"!");
48 const std::vector<panzer::LocalOrdinal> & elements = ugi.
getElementBlock(blockId);
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]]);
64 dataVector.ReplaceMyValue(localIndex,0,data(e,f));
69 std::size_t entries = data.extent(2);
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++)
82 dataVector.ReplaceMyValue(localIndex,v,data(e,f,v));
88 "panzer::updateGhostedDataReducedVector: data array rank must be 2 or 3");
91 template <
typename ScalarT,
typename ArrayT>
94 const std::map<std::string,ArrayT> & data)
const
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);
111 "ArrayToFieldVectorEpetra::getGhostedDataVector: data array must have rank 2 or 3. This array has rank " << rank <<
".");
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);
137 "ArrayToFieldVectorEpetra::getDataVector: can not find block \""+block+
"\".");
139 const ArrayT & d = blockItr->second;
140 updateGhostedDataReducedVectorEpetra<ScalarT,ArrayT>(fieldName,block,*
ugi_,d,*finalReducedVec);
157 finalVec->Import(*finalReducedVec,importer,
Insert);
162 template <
typename ScalarT,
typename ArrayT>
165 const std::map<std::string,ArrayT> & data)
const
172 = getGhostedDataVector<ScalarT,ArrayT>(fieldName,data);
175 int fieldNum =
ugi_->getFieldNum(fieldName);
187 destVec->Import(*sourceVec,importer,
Insert);
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
Teuchos::RCP< const IntVector > fieldVector_
Maps for each field (as needed)
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
virtual const std::vector< panzer::LocalOrdinal > & getElementBlock(const std::string &blockId) const =0
void buildFieldVector(const Epetra_IntVector &source) const
build unghosted field vector from ghosted field vector
Epetra_MultiVector MultiVector
std::map< int, Teuchos::RCP< const Map > > gh_fieldMaps_
Maps for each field (as needed)
Teuchos::RCP< const IntVector > gh_fieldVector_
ghosted reduced field vector
virtual bool fieldInBlock(const std::string &field, const std::string &block) const =0
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
Teuchos::RCP< const IntVector > gh_reducedFieldVector_
virtual int getFieldNum(const std::string &str) const =0
Get the number used for access to this field.
Teuchos::RCP< const GlobalIndexer > ugi_
DOF mapping.
Teuchos::RCP< Epetra_MultiVector > getGhostedDataVector(const std::string &fieldName, const std::map< std::string, ArrayT > &data) const
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.
std::map< int, Teuchos::RCP< const Map > > gh_reducedFieldMaps_
ghosted field vector
#define TEUCHOS_ASSERT(assertion_test)
Teuchos::RCP< const Epetra_BlockMap > getFieldMapEpetra(int fieldNum, const Epetra_IntVector &fieldVector)
void updateGhostedDataReducedVectorEpetra(const std::string &fieldName, const std::string blockId, const GlobalIndexer &ugi, const ArrayT &data, Epetra_MultiVector &dataVector)