Go to the documentation of this file.00001
00002
00003
00004
00005
00006
00007 #include "Adapt_NodalDataBlock.hpp"
00008 #include "Epetra_Import.h"
00009
00010 Adapt::NodalDataBlock::NodalDataBlock() :
00011 nodeContainer(Teuchos::rcp(new Albany::NodeFieldContainer)),
00012 blocksize(0),
00013 mapsHaveChanged(false)
00014 {
00015 }
00016
00017 void
00018 Adapt::NodalDataBlock::resizeOverlapMap(const std::vector<int>& overlap_nodeGIDs, const Epetra_Comm& comm){
00019
00020
00021 overlap_node_map = Teuchos::rcp(new Epetra_BlockMap(-1,
00022 overlap_nodeGIDs.size(),
00023 &overlap_nodeGIDs[0],
00024 blocksize,
00025 0,
00026 comm));
00027
00028
00029 overlap_node_vec = Teuchos::rcp(new Epetra_Vector(*overlap_node_map, false));
00030
00031 mapsHaveChanged = true;
00032
00033 }
00034
00035 void
00036 Adapt::NodalDataBlock::resizeLocalMap(const std::vector<int>& local_nodeGIDs, const Epetra_Comm& comm){
00037
00038
00039
00040
00041 local_node_map = Teuchos::rcp(new Epetra_BlockMap(-1,
00042 local_nodeGIDs.size(),
00043 &local_nodeGIDs[0],
00044 blocksize,
00045 0,
00046 comm));
00047
00048
00049 local_node_vec = Teuchos::rcp(new Epetra_Vector(*local_node_map, false));
00050
00051 mapsHaveChanged = true;
00052
00053 }
00054
00055 void
00056 Adapt::NodalDataBlock::initializeExport(){
00057
00058 if(mapsHaveChanged){
00059
00060 importer = Teuchos::rcp(new Epetra_Import(*overlap_node_map, *local_node_map));
00061 mapsHaveChanged = false;
00062
00063 }
00064
00065 }
00066
00067 void
00068 Adapt::NodalDataBlock::exportAddNodalDataBlock(){
00069
00070 overlap_node_vec->Import(*local_node_vec, *importer, Add);
00071
00072 }
00073
00074 void
00075 Adapt::NodalDataBlock::registerState(const std::string &stateName, int ndofs){
00076
00077
00078
00079 NodeFieldSizeMap::const_iterator it;
00080 it = nodeBlockMap.find(stateName);
00081
00082 TEUCHOS_TEST_FOR_EXCEPTION((it != nodeBlockMap.end()), std::logic_error,
00083 std::endl << "Error: found duplicate entry " << stateName << " in NodalDataBlock" << std::endl);
00084
00085 NodeFieldSize size;
00086 size.name = stateName;
00087 size.offset = blocksize;
00088 size.ndofs = ndofs;
00089
00090 nodeBlockMap[stateName] = nodeBlockLayout.size();
00091 nodeBlockLayout.push_back(size);
00092
00093 blocksize += ndofs;
00094
00095 }
00096
00097 void
00098 Adapt::NodalDataBlock::getNDofsAndOffset(const std::string &stateName, int& offset, int& ndofs) const {
00099
00100 NodeFieldSizeMap::const_iterator it;
00101 it = nodeBlockMap.find(stateName);
00102
00103 TEUCHOS_TEST_FOR_EXCEPTION((it == nodeBlockMap.end()), std::logic_error,
00104 std::endl << "Error: cannot find state " << stateName << " in NodalDataBlock" << std::endl);
00105
00106 std::size_t value = it->second;
00107
00108 offset = nodeBlockLayout[value].offset;
00109 ndofs = nodeBlockLayout[value].ndofs;
00110
00111 }
00112
00113 void
00114 Adapt::NodalDataBlock::saveNodalDataState() const {
00115
00116
00117 for(NodeFieldSizeVector::const_iterator i = nodeBlockLayout.begin(); i != nodeBlockLayout.end(); ++i){
00118
00119
00120
00121 (*nodeContainer)[i->name]->saveField(overlap_node_vec, i->offset);
00122
00123 }
00124
00125 }
00126
00127 void
00128 Adapt::NodalDataBlock::saveEpetraNodalDataVector(const std::string& name,
00129 const Teuchos::RCP<const Epetra_Vector>& overlap_node_vec,
00130 int offset, int blocksize) const {
00131
00132 Albany::NodeFieldContainer::const_iterator it;
00133 it = nodeContainer->find(name);
00134
00135 TEUCHOS_TEST_FOR_EXCEPTION((it == nodeContainer->end()), std::logic_error,
00136 std::endl << "Error: Cannot locate nodal field " << name << " in NodalDataBlock" << std::endl);
00137
00138
00139
00140 (*nodeContainer)[name]->saveField(overlap_node_vec, offset, blocksize);
00141
00142 }
00143
00144