Panzer  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Panzer_STK_Interface.hpp
Go to the documentation of this file.
1 // @HEADER
2 // ***********************************************************************
3 //
4 // Panzer: A partial differential equation assembly
5 // engine for strongly coupled complex multiphysics systems
6 // Copyright (2011) Sandia Corporation
7 //
8 // Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
9 // the U.S. Government retains certain rights in this software.
10 //
11 // Redistribution and use in source and binary forms, with or without
12 // modification, are permitted provided that the following conditions are
13 // met:
14 //
15 // 1. Redistributions of source code must retain the above copyright
16 // notice, this list of conditions and the following disclaimer.
17 //
18 // 2. Redistributions in binary form must reproduce the above copyright
19 // notice, this list of conditions and the following disclaimer in the
20 // documentation and/or other materials provided with the distribution.
21 //
22 // 3. Neither the name of the Corporation nor the names of the
23 // contributors may be used to endorse or promote products derived from
24 // this software without specific prior written permission.
25 //
26 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
27 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
28 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
29 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
30 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
31 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
32 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
33 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
34 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
35 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
36 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
37 //
38 // Questions? Contact Roger P. Pawlowski (rppawlo@sandia.gov) and
39 // Eric C. Cyr (eccyr@sandia.gov)
40 // ***********************************************************************
41 // @HEADER
42 
43 #ifndef __Panzer_STK_Interface_hpp__
44 #define __Panzer_STK_Interface_hpp__
45 
46 #include <Teuchos_RCP.hpp>
48 
49 #include <stk_mesh/base/Types.hpp>
50 #include <stk_mesh/base/MetaData.hpp>
51 #include <stk_mesh/base/BulkData.hpp>
52 #include <stk_mesh/base/Field.hpp>
53 #include <stk_mesh/base/FieldBase.hpp>
54 #include <stk_mesh/base/MetaData.hpp>
55 #include <stk_mesh/base/CoordinateSystems.hpp>
56 
57 #include "Kokkos_Core.hpp"
58 
59 #include <Shards_CellTopology.hpp>
60 #include <Shards_CellTopologyData.h>
61 
62 #include <PanzerAdaptersSTK_config.hpp>
63 #include <Kokkos_ViewFactory.hpp>
64 
65 #include <unordered_map>
66 
67 #ifdef PANZER_HAVE_IOSS
68 #include <stk_io/StkMeshIoBroker.hpp>
69 #endif
70 
71 #ifdef PANZER_HAVE_PERCEPT
72 namespace percept {
73  class PerceptMesh;
74  class URP_Heterogeneous_3D;
75 }
76 #endif
77 
78 namespace panzer_stk {
79 
80 class PeriodicBC_MatcherBase;
81 
86 public:
87  ElementDescriptor(stk::mesh::EntityId gid,const std::vector<stk::mesh::EntityId> & nodes);
88  virtual ~ElementDescriptor();
89 
90  stk::mesh::EntityId getGID() const { return gid_; }
91  const std::vector<stk::mesh::EntityId> & getNodes() const { return nodes_; }
92 protected:
93  stk::mesh::EntityId gid_;
94  std::vector<stk::mesh::EntityId> nodes_;
95 
97 };
98 
102 buildElementDescriptor(stk::mesh::EntityId elmtId,std::vector<stk::mesh::EntityId> & nodes);
103 
105 public:
106  typedef double ProcIdData; // ECC: Not sure why?
107  typedef stk::mesh::Field<double> SolutionFieldType;
108  typedef stk::mesh::Field<double,stk::mesh::Cartesian> VectorFieldType;
109  typedef stk::mesh::Field<ProcIdData> ProcIdFieldType;
110 
111  // some simple exception classes
112  struct ElementBlockException : public std::logic_error
113  { ElementBlockException(const std::string & what) : std::logic_error(what) {} };
114 
115  struct SidesetException : public std::logic_error
116  { SidesetException(const std::string & what) : std::logic_error(what) {} };
117 
118  STK_Interface();
119 
122  STK_Interface(unsigned dim);
123 
125 
126  // functions called before initialize
128 
131  void addElementBlock(const std::string & name,const CellTopologyData * ctData);
132 
135  void addSideset(const std::string & name,const CellTopologyData * ctData);
136 
139  void addNodeset(const std::string & name);
140 
143  void addSolutionField(const std::string & fieldName,const std::string & blockId);
144 
147  void addCellField(const std::string & fieldName,const std::string & blockId);
148 
157  void addMeshCoordFields(const std::string & blockId,
158  const std::vector<std::string> & coordField,
159  const std::string & dispPrefix);
160 
162 
174  void initialize(stk::ParallelMachine parallelMach,bool setupIO=true,
175  const bool buildRefinementSupport = false);
176 
182  void instantiateBulkData(stk::ParallelMachine parallelMach);
183 
184  // functions to manage and manipulate bulk data
186 
189  void beginModification();
190 
193  void endModification();
194 
200  void addNode(stk::mesh::EntityId gid, const std::vector<double> & coord);
201 
202  void addElement(const Teuchos::RCP<ElementDescriptor> & ed,stk::mesh::Part * block);
203 
204  void addEdges();
205 
206  void addFaces();
207 
210  void addEntityToSideset(stk::mesh::Entity entity,stk::mesh::Part * sideset);
211 
214  void addEntityToNodeset(stk::mesh::Entity entity,stk::mesh::Part * nodeset);
215 
216  // Methods to interrogate the mesh topology and structure
218 
222  { return *coordinatesField_; }
223 
227  { return *edgesField_; }
228 
230  { return *facesField_; }
231 
234  const double * getNodeCoordinates(stk::mesh::EntityId nodeId) const;
235 
238  const double * getNodeCoordinates(stk::mesh::Entity node) const;
239 
242  void getSubcellIndices(unsigned entityRank,stk::mesh::EntityId elementId,
243  std::vector<stk::mesh::EntityId> & subcellIds) const;
244 
247  void getMyElements(std::vector<stk::mesh::Entity> & elements) const;
248 
251  void getMyElements(const std::string & blockID,std::vector<stk::mesh::Entity> & elements) const;
252 
256  void getNeighborElements(std::vector<stk::mesh::Entity> & elements) const;
257 
260  void getNeighborElements(const std::string & blockID,std::vector<stk::mesh::Entity> & elements) const;
261 
269  void getMySides(const std::string & sideName,std::vector<stk::mesh::Entity> & sides) const;
270 
279  void getMySides(const std::string & sideName,const std::string & blockName,std::vector<stk::mesh::Entity> & sides) const;
280 
288  void getAllSides(const std::string & sideName,std::vector<stk::mesh::Entity> & sides) const;
289 
299  void getAllSides(const std::string & sideName,const std::string & blockName,std::vector<stk::mesh::Entity> & sides) const;
300 
310  void getMyNodes(const std::string & sideName,const std::string & blockName,std::vector<stk::mesh::Entity> & nodes) const;
311 
320  stk::mesh::Entity findConnectivityById(stk::mesh::Entity src, stk::mesh::EntityRank tgt_rank, unsigned rel_id) const;
321 
322  // Utility functions
324 
332  void
334  const std::string& filename);
335 
350  void
352  const std::string& filename);
353 
366  void
368  double timestep);
369 
389  void
391  const std::string& key,
392  const int& value);
393 
413  void
415  const std::string& key,
416  const double& value);
417 
437  void
439  const std::string& key,
440  const std::vector<int>& value);
441 
461  void
463  const std::string& key,
464  const std::vector<double>& value);
465 
466  // Accessor functions
468 
471 
474 
475  bool isWritable() const;
476 
477  bool isModifiable() const
478  { if(bulkData_==Teuchos::null) return false;
479  return bulkData_->in_modifiable_state(); }
480 
482  unsigned getDimension() const
483  { return dimension_; }
484 
486  std::size_t getNumElementBlocks() const
487  { return elementBlocks_.size(); }
488 
496  void getElementBlockNames(std::vector<std::string> & names) const;
497 
505  void getSidesetNames(std::vector<std::string> & name) const;
506 
514  void getNodesetNames(std::vector<std::string> & name) const;
515 
517  stk::mesh::Part * getOwnedPart() const
518  { return &getMetaData()->locally_owned_part(); } // I don't like the pointer access here, but it will do for now!
519 
521  stk::mesh::Part * getElementBlockPart(const std::string & name) const
522  {
523  std::map<std::string, stk::mesh::Part*>::const_iterator itr = elementBlocks_.find(name); // Element blocks
524  if(itr==elementBlocks_.end()) return 0;
525  return itr->second;
526  }
527 
529  std::size_t getNumSidesets() const
530  { return sidesets_.size(); }
531 
532  stk::mesh::Part * getSideset(const std::string & name) const
533  {
534  auto itr = sidesets_.find(name);
535  return (itr != sidesets_.end()) ? itr->second : nullptr;
536  }
537 
539  std::size_t getNumNodesets() const
540  { return nodesets_.size(); }
541 
542  stk::mesh::Part * getNodeset(const std::string & name) const
543  {
544  auto itr = nodesets_.find(name);
545  return (itr != nodesets_.end()) ? itr->second : nullptr;
546  }
547 
549  std::size_t getEntityCounts(unsigned entityRank) const;
550 
552  stk::mesh::EntityId getMaxEntityId(unsigned entityRank) const;
553 
554  // Utilities
556 
558  void getElementsSharingNode(stk::mesh::EntityId nodeId,std::vector<stk::mesh::Entity> & elements) const;
559 
561  void getNodeIdsForElement(stk::mesh::Entity element,std::vector<stk::mesh::EntityId> & nodeIds) const;
562 
565  void getOwnedElementsSharingNode(stk::mesh::Entity node,std::vector<stk::mesh::Entity> & elements,
566  std::vector<int> & relIds) const;
567 
570  void getOwnedElementsSharingNode(stk::mesh::EntityId nodeId,std::vector<stk::mesh::Entity> & elements,
571  std::vector<int> & relIds, unsigned int matchType) const;
572 
573 
575  void getElementsSharingNodes(const std::vector<stk::mesh::EntityId> nodeId,std::vector<stk::mesh::Entity> & elements) const;
576 
578  void buildSubcells();
579 
582  std::size_t elementLocalId(stk::mesh::Entity elmt) const;
583 
586  std::size_t elementLocalId(stk::mesh::EntityId gid) const;
587 
590  inline stk::mesh::EntityId elementGlobalId(std::size_t lid) const
591  { return bulkData_->identifier((*orderedElementVector_)[lid]); }
592 
595  inline stk::mesh::EntityId elementGlobalId(stk::mesh::Entity elmt) const
596  { return bulkData_->identifier(elmt); }
597 
600  inline unsigned entityOwnerRank(stk::mesh::Entity entity) const
601  { return bulkData_->parallel_owner_rank(entity); }
602 
605  inline bool isValid(stk::mesh::Entity entity) const
606  { return bulkData_->is_valid(entity); }
607 
610  std::string containingBlockId(stk::mesh::Entity elmt) const;
611 
616  stk::mesh::Field<double> * getSolutionField(const std::string & fieldName,
617  const std::string & blockId) const;
618 
623  stk::mesh::Field<double> * getCellField(const std::string & fieldName,
624  const std::string & blockId) const;
625 
627 
629  bool isInitialized() const { return initialized_; }
630 
635 
650  template <typename ArrayT>
651  void setSolutionFieldData(const std::string & fieldName,const std::string & blockId,
652  const std::vector<std::size_t> & localElementIds,const ArrayT & solutionValues,double scaleValue=1.0);
653 
668  template <typename ArrayT>
669  void getSolutionFieldData(const std::string & fieldName,const std::string & blockId,
670  const std::vector<std::size_t> & localElementIds,ArrayT & solutionValues) const;
671 
686  template <typename ArrayT>
687  void setCellFieldData(const std::string & fieldName,const std::string & blockId,
688  const std::vector<std::size_t> & localElementIds,const ArrayT & solutionValues,double scaleValue=1.0);
689 
698  template <typename ArrayT>
699  void getElementVertices(const std::vector<std::size_t> & localIds, ArrayT & vertices) const;
700 
709  template <typename ArrayT>
710  void getElementVertices(const std::vector<stk::mesh::Entity> & elements, ArrayT & vertices) const;
711 
721  template <typename ArrayT>
722  void getElementVertices(const std::vector<std::size_t> & localIds,const std::string & eBlock, ArrayT & vertices) const;
723 
733  template <typename ArrayT>
734  void getElementVertices(const std::vector<stk::mesh::Entity> & elements,const std::string & eBlock, ArrayT & vertices) const;
735 
744  template <typename ArrayT>
745  void getElementVerticesNoResize(const std::vector<std::size_t> & localIds, ArrayT & vertices) const;
746 
755  template <typename ArrayT>
756  void getElementVerticesNoResize(const std::vector<stk::mesh::Entity> & elements, ArrayT & vertices) const;
757 
767  template <typename ArrayT>
768  void getElementVerticesNoResize(const std::vector<std::size_t> & localIds,const std::string & eBlock, ArrayT & vertices) const;
769 
779  template <typename ArrayT>
780  void getElementVerticesNoResize(const std::vector<stk::mesh::Entity> & elements,const std::string & eBlock, ArrayT & vertices) const;
781 
782  // const stk::mesh::FEMInterface & getFEMInterface() const
783  // { return *femPtr_; }
784 
785  stk::mesh::EntityRank getElementRank() const { return stk::topology::ELEMENT_RANK; }
786  stk::mesh::EntityRank getSideRank() const { return metaData_->side_rank(); }
787  stk::mesh::EntityRank getFaceRank() const { return stk::topology::FACE_RANK; }
788  stk::mesh::EntityRank getEdgeRank() const { return stk::topology::EDGE_RANK; }
789  stk::mesh::EntityRank getNodeRank() const { return stk::topology::NODE_RANK; }
790 
793  void initializeFromMetaData();
794 
797  void buildLocalElementIDs();
798 
801  const std::vector<Teuchos::RCP<const PeriodicBC_MatcherBase> > &
803  { return periodicBCs_; }
804 
807  std::vector<Teuchos::RCP<const PeriodicBC_MatcherBase> > &
809  { return periodicBCs_; }
810 
817  { periodicBCs_.push_back(bc); }
818 
825  { periodicBCs_.insert(periodicBCs_.end(),bc_vec.begin(),bc_vec.end()); }
826 
827  std::pair<Teuchos::RCP<std::vector<std::pair<std::size_t,std::size_t> > >, Teuchos::RCP<std::vector<unsigned int> > >
828  getPeriodicNodePairing() const;
829 
832  bool validBlockId(const std::string & blockId) const;
833 
836  void print(std::ostream & os) const;
837 
840  void printMetaData(std::ostream & os) const;
841 
844  Teuchos::RCP<const shards::CellTopology> getCellTopology(const std::string & eBlock) const;
845 
850  double getCurrentStateTime() const { return currentStateTime_; }
851 
857  double getInitialStateTime() const { return initialStateTime_; }
858 
863  void setInitialStateTime(double value) { initialStateTime_ = value; }
864 
867  void rebalance(const Teuchos::ParameterList & params);
868 
872  void setBlockWeight(const std::string & blockId,double weight)
873  { blockWeights_[blockId] = weight; }
874 
881  void setUseFieldCoordinates(bool useFieldCoordinates)
882  { useFieldCoordinates_ = useFieldCoordinates; }
883 
886  { return useFieldCoordinates_; }
887 
889  void setUseLowerCaseForIO(bool useLowerCase)
890  { useLowerCase_ = useLowerCase; }
891 
893  bool getUseLowerCaseForIO() const
894  { return useLowerCase_; }
895 
906  template <typename ArrayT>
907  void getElementVertices_FromField(const std::vector<stk::mesh::Entity> & elements,const std::string & eBlock, ArrayT & vertices) const;
908 
909  template <typename ArrayT>
910  void getElementVertices_FromFieldNoResize(const std::vector<stk::mesh::Entity> & elements,
911  const std::string & eBlock, ArrayT & vertices) const;
912 
922  template <typename ArrayT>
923  void getElementVertices_FromCoords(const std::vector<stk::mesh::Entity> & elements, ArrayT & vertices) const;
924 
925  template <typename ArrayT>
926  void getElementVertices_FromCoordsNoResize(const std::vector<stk::mesh::Entity> & elements, ArrayT & vertices) const;
927 
933  void refineMesh(const int numberOfLevels, const bool deleteParentElements);
934 
935 public: // static operations
936  static const std::string coordsString;
937  static const std::string nodesString;
938  static const std::string edgesString;
939  static const std::string facesString;
940 
941 protected:
942 
945  void buildEntityCounts();
946 
949  void buildMaxEntityIds();
950 
954  void initializeFieldsInSTK(const std::map<std::pair<std::string,std::string>,SolutionFieldType*> & nameToField,
955  bool setupIO);
956 
961  Teuchos::RCP<Teuchos::MpiComm<int> > getSafeCommunicator(stk::ParallelMachine parallelMach) const;
962 
969 
973  bool isMeshCoordField(const std::string & eBlock,const std::string & fieldName,int & axis) const;
974 
990  template <typename ArrayT>
991  void setDispFieldData(const std::string & fieldName,const std::string & blockId,int axis,
992  const std::vector<std::size_t> & localElementIds,const ArrayT & solutionValues);
993 
994  std::vector<Teuchos::RCP<const PeriodicBC_MatcherBase> > periodicBCs_;
995 
998 #ifdef PANZER_HAVE_PERCEPT
1001 #endif
1002 
1003  std::map<std::string, stk::mesh::Part*> elementBlocks_; // Element blocks
1004  std::map<std::string, stk::mesh::Part*> sidesets_; // Side sets
1005  std::map<std::string, stk::mesh::Part*> nodesets_; // Node sets
1006  std::map<std::string, Teuchos::RCP<shards::CellTopology> > elementBlockCT_;
1007 
1008  // for storing/accessing nodes
1009  stk::mesh::Part * nodesPart_;
1010  std::vector<stk::mesh::Part*> nodesPartVec_;
1011  stk::mesh::Part * edgesPart_;
1012  std::vector<stk::mesh::Part*> edgesPartVec_;
1013  stk::mesh::Part * facesPart_;
1014  std::vector<stk::mesh::Part*> facesPartVec_;
1015 
1021 
1022  // maps field names to solution field stk mesh handles
1023  std::map<std::pair<std::string,std::string>,SolutionFieldType*> fieldNameToSolution_;
1024  std::map<std::pair<std::string,std::string>,SolutionFieldType*> fieldNameToCellField_;
1025 
1026  unsigned dimension_;
1027 
1029 
1030  // how many elements, faces, edges, and nodes are there globally
1031  std::vector<std::size_t> entityCounts_;
1032 
1033  // what is maximum entity ID
1034  std::vector<stk::mesh::EntityId> maxEntityId_;
1035 
1036  unsigned procRank_;
1037  std::size_t currentLocalId_;
1038 
1040 
1041  double initialStateTime_; // the time stamp at the time this object was constructed (default 0.0)
1042  double currentStateTime_; // the time stamp set by the user when writeToExodus is called (default 0.0)
1043 
1044 #ifdef PANZER_HAVE_IOSS
1045  // I/O support
1047  int meshIndex_;
1048 
1053  enum class GlobalVariable
1054  {
1055  ADD,
1056  WRITE
1057  }; // end of enum class GlobalVariable
1058 
1075  void
1076  globalToExodus(
1077  const GlobalVariable& flag);
1078 
1082  Teuchos::ParameterList globalData_;
1083 #endif
1084 
1085  // uses lazy evaluation
1087 
1088  // for element block weights
1089  std::map<std::string,double> blockWeights_;
1090 
1091  std::unordered_map<stk::mesh::EntityId,std::size_t> localIDHash_;
1092 
1093  // Store mesh displacement fields by element block. This map
1094  // goes like this meshCoordFields_[eBlock][axis_index] => coordinate FieldName
1095  // goes like this meshDispFields_[eBlock][axis_index] => displacement FieldName
1096  std::map<std::string,std::vector<std::string> > meshCoordFields_; // coordinate fields written by user
1097  std::map<std::string,std::vector<std::string> > meshDispFields_; // displacement fields, output to exodus
1098 
1100 
1102 
1103  // Object describing how to sort a vector of elements using
1104  // local ID as the key, very short lived object
1106  public:
1107  LocalIdCompare(const STK_Interface * mesh) : mesh_(mesh) {}
1108 
1109  // Compares two stk mesh entities based on local ID
1110  bool operator() (stk::mesh::Entity a,stk::mesh::Entity b)
1111  { return mesh_->elementLocalId(a) < mesh_->elementLocalId(b);}
1112 
1113  private:
1115  };
1116 };
1117 
1118 template <typename ArrayT>
1119 void STK_Interface::setSolutionFieldData(const std::string & fieldName,const std::string & blockId,
1120  const std::vector<std::size_t> & localElementIds,const ArrayT & solutionValues,double scaleValue)
1121 {
1122  const std::vector<stk::mesh::Entity> & elements = *(this->getElementsOrderedByLID());
1123 
1124  int field_axis = -1;
1125  if(isMeshCoordField(blockId,fieldName,field_axis)) {
1126  setDispFieldData(fieldName,blockId,field_axis,localElementIds,solutionValues);
1127  return;
1128  }
1129 
1130  SolutionFieldType * field = this->getSolutionField(fieldName,blockId);
1131 
1132  for(std::size_t cell=0;cell<localElementIds.size();cell++) {
1133  std::size_t localId = localElementIds[cell];
1134  stk::mesh::Entity element = elements[localId];
1135 
1136  // loop over nodes set solution values
1137  const size_t num_nodes = bulkData_->num_nodes(element);
1138  stk::mesh::Entity const* nodes = bulkData_->begin_nodes(element);
1139  for(std::size_t i=0; i<num_nodes; ++i) {
1140  stk::mesh::Entity node = nodes[i];
1141 
1142  double * solnData = stk::mesh::field_data(*field,node);
1143  // TEUCHOS_ASSERT(solnData!=0); // only needed if blockId is not specified
1144  solnData[0] = scaleValue*solutionValues(cell,i);
1145  }
1146  }
1147 }
1148 
1149 template <typename ArrayT>
1150 void STK_Interface::setDispFieldData(const std::string & fieldName,const std::string & blockId,int axis,
1151  const std::vector<std::size_t> & localElementIds,const ArrayT & dispValues)
1152 {
1153  TEUCHOS_ASSERT(axis>=0); // sanity check
1154 
1155  const std::vector<stk::mesh::Entity> & elements = *(this->getElementsOrderedByLID());
1156 
1157  SolutionFieldType * field = this->getSolutionField(fieldName,blockId);
1158  const VectorFieldType & coord_field = this->getCoordinatesField();
1159 
1160  for(std::size_t cell=0;cell<localElementIds.size();cell++) {
1161  std::size_t localId = localElementIds[cell];
1162  stk::mesh::Entity element = elements[localId];
1163 
1164  // loop over nodes set solution values
1165  const size_t num_nodes = bulkData_->num_nodes(element);
1166  stk::mesh::Entity const* nodes = bulkData_->begin_nodes(element);
1167  for(std::size_t i=0; i<num_nodes; ++i) {
1168  stk::mesh::Entity node = nodes[i];
1169 
1170  double * solnData = stk::mesh::field_data(*field,node);
1171  double * coordData = stk::mesh::field_data(coord_field,node);
1172  // TEUCHOS_ASSERT(solnData!=0); // only needed if blockId is not specified
1173  solnData[0] = dispValues(cell,i)-coordData[axis];
1174  }
1175  }
1176 }
1177 
1178 template <typename ArrayT>
1179 void STK_Interface::getSolutionFieldData(const std::string & fieldName,const std::string & blockId,
1180  const std::vector<std::size_t> & localElementIds,ArrayT & solutionValues) const
1181 {
1182  const std::vector<stk::mesh::Entity> & elements = *(this->getElementsOrderedByLID());
1183 
1184  solutionValues = Kokkos::createDynRankView(solutionValues,
1185  "solutionValues",
1186  localElementIds.size(),
1187  bulkData_->num_nodes(elements[localElementIds[0]]));
1188 
1189  SolutionFieldType * field = this->getSolutionField(fieldName,blockId);
1190 
1191  for(std::size_t cell=0;cell<localElementIds.size();cell++) {
1192  std::size_t localId = localElementIds[cell];
1193  stk::mesh::Entity element = elements[localId];
1194 
1195  // loop over nodes set solution values
1196  const size_t num_nodes = bulkData_->num_nodes(element);
1197  stk::mesh::Entity const* nodes = bulkData_->begin_nodes(element);
1198  for(std::size_t i=0; i<num_nodes; ++i) {
1199  stk::mesh::Entity node = nodes[i];
1200 
1201  double * solnData = stk::mesh::field_data(*field,node);
1202  // TEUCHOS_ASSERT(solnData!=0); // only needed if blockId is not specified
1203  solutionValues(cell,i) = solnData[0];
1204  }
1205  }
1206 }
1207 
1208 template <typename ArrayT>
1209 void STK_Interface::setCellFieldData(const std::string & fieldName,const std::string & blockId,
1210  const std::vector<std::size_t> & localElementIds,const ArrayT & solutionValues,double scaleValue)
1211 {
1212  const std::vector<stk::mesh::Entity> & elements = *(this->getElementsOrderedByLID());
1213 
1214  SolutionFieldType * field = this->getCellField(fieldName,blockId);
1215 
1216  for(std::size_t cell=0;cell<localElementIds.size();cell++) {
1217  std::size_t localId = localElementIds[cell];
1218  stk::mesh::Entity element = elements[localId];
1219 
1220  double * solnData = stk::mesh::field_data(*field,element);
1221  TEUCHOS_ASSERT(solnData!=0); // only needed if blockId is not specified
1222  solnData[0] = scaleValue*solutionValues.access(cell,0);
1223  }
1224 }
1225 
1226 template <typename ArrayT>
1227 void STK_Interface::getElementVertices(const std::vector<std::size_t> & localElementIds, ArrayT & vertices) const
1228 {
1229  if(!useFieldCoordinates_) {
1230  //
1231  // gather from the intrinsic mesh coordinates (non-lagrangian)
1232  //
1233 
1234  const std::vector<stk::mesh::Entity> & elements = *(this->getElementsOrderedByLID());
1235 
1236  // convert to a vector of entity objects
1237  std::vector<stk::mesh::Entity> selected_elements;
1238  for(std::size_t cell=0;cell<localElementIds.size();cell++)
1239  selected_elements.push_back(elements[localElementIds[cell]]);
1240 
1241  getElementVertices_FromCoords(selected_elements,vertices);
1242  }
1243  else {
1244  TEUCHOS_TEST_FOR_EXCEPTION(true,std::invalid_argument,
1245  "STK_Interface::getElementVertices: Cannot call this method when field coordinates are used "
1246  "without specifying an element block.");
1247  }
1248 }
1249 
1250 template <typename ArrayT>
1251 void STK_Interface::getElementVertices(const std::vector<stk::mesh::Entity> & elements, ArrayT & vertices) const
1252 {
1253  if(!useFieldCoordinates_) {
1254  getElementVertices_FromCoords(elements,vertices);
1255  }
1256  else {
1257  TEUCHOS_TEST_FOR_EXCEPTION(true,std::invalid_argument,
1258  "STK_Interface::getElementVertices: Cannot call this method when field coordinates are used "
1259  "without specifying an element block.");
1260  }
1261 }
1262 
1263 template <typename ArrayT>
1264 void STK_Interface::getElementVertices(const std::vector<stk::mesh::Entity> & elements,const std::string & eBlock, ArrayT & vertices) const
1265 {
1266  if(!useFieldCoordinates_) {
1267  getElementVertices_FromCoords(elements,vertices);
1268  }
1269  else {
1270  getElementVertices_FromField(elements,eBlock,vertices);
1271  }
1272 }
1273 
1274 template <typename ArrayT>
1275 void STK_Interface::getElementVertices(const std::vector<std::size_t> & localElementIds,const std::string & eBlock, ArrayT & vertices) const
1276 {
1277  const std::vector<stk::mesh::Entity> & elements = *(this->getElementsOrderedByLID());
1278 
1279  // convert to a vector of entity objects
1280  std::vector<stk::mesh::Entity> selected_elements;
1281  for(std::size_t cell=0;cell<localElementIds.size();cell++)
1282  selected_elements.push_back(elements[localElementIds[cell]]);
1283 
1284  if(!useFieldCoordinates_) {
1285  getElementVertices_FromCoords(selected_elements,vertices);
1286  }
1287  else {
1288  getElementVertices_FromField(selected_elements,eBlock,vertices);
1289  }
1290 }
1291 
1292 template <typename ArrayT>
1293 void STK_Interface::getElementVerticesNoResize(const std::vector<std::size_t> & localElementIds, ArrayT & vertices) const
1294 {
1295  if(!useFieldCoordinates_) {
1296  //
1297  // gather from the intrinsic mesh coordinates (non-lagrangian)
1298  //
1299 
1300  const std::vector<stk::mesh::Entity> & elements = *(this->getElementsOrderedByLID());
1301 
1302  // convert to a vector of entity objects
1303  std::vector<stk::mesh::Entity> selected_elements;
1304  for(std::size_t cell=0;cell<localElementIds.size();cell++)
1305  selected_elements.push_back(elements[localElementIds[cell]]);
1306 
1307  getElementVertices_FromCoordsNoResize(selected_elements,vertices);
1308  }
1309  else {
1310  TEUCHOS_TEST_FOR_EXCEPTION(true,std::invalid_argument,
1311  "STK_Interface::getElementVerticesNoResize: Cannot call this method when field coordinates are used "
1312  "without specifying an element block.");
1313  }
1314 }
1315 
1316 template <typename ArrayT>
1317 void STK_Interface::getElementVerticesNoResize(const std::vector<stk::mesh::Entity> & elements, ArrayT & vertices) const
1318 {
1319  if(!useFieldCoordinates_) {
1320  getElementVertices_FromCoordsNoResize(elements,vertices);
1321  }
1322  else {
1323  TEUCHOS_TEST_FOR_EXCEPTION(true,std::invalid_argument,
1324  "STK_Interface::getElementVerticesNoResize: Cannot call this method when field coordinates are used "
1325  "without specifying an element block.");
1326  }
1327 }
1328 
1329 template <typename ArrayT>
1330 void STK_Interface::getElementVerticesNoResize(const std::vector<stk::mesh::Entity> & elements,const std::string & eBlock, ArrayT & vertices) const
1331 {
1332  if(!useFieldCoordinates_) {
1333  getElementVertices_FromCoordsNoResize(elements,vertices);
1334  }
1335  else {
1336  getElementVertices_FromFieldNoResize(elements,eBlock,vertices);
1337  }
1338 }
1339 
1340 template <typename ArrayT>
1341 void STK_Interface::getElementVerticesNoResize(const std::vector<std::size_t> & localElementIds,const std::string & eBlock, ArrayT & vertices) const
1342 {
1343  const std::vector<stk::mesh::Entity> & elements = *(this->getElementsOrderedByLID());
1344 
1345  // convert to a vector of entity objects
1346  std::vector<stk::mesh::Entity> selected_elements;
1347  for(std::size_t cell=0;cell<localElementIds.size();cell++)
1348  selected_elements.push_back(elements[localElementIds[cell]]);
1349 
1350  if(!useFieldCoordinates_) {
1351  getElementVertices_FromCoordsNoResize(selected_elements,vertices);
1352  }
1353  else {
1354  getElementVertices_FromFieldNoResize(selected_elements,eBlock,vertices);
1355  }
1356 }
1357 
1358 template <typename ArrayT>
1359 void STK_Interface::getElementVertices_FromCoords(const std::vector<stk::mesh::Entity> & elements, ArrayT & vertices) const
1360 {
1361  // nothing to do! silently return
1362  if(elements.size()==0) {
1363  vertices = Kokkos::createDynRankView(vertices,"vertices",0,0,0);
1364  return;
1365  }
1366 
1367  //
1368  // gather from the intrinsic mesh coordinates (non-lagrangian)
1369  //
1370 
1371  // get *master* cell toplogy...(belongs to first element)
1372  unsigned masterVertexCount
1373  = stk::mesh::get_cell_topology(bulkData_->bucket(elements[0]).topology()).getCellTopologyData()->vertex_count;
1374 
1375  // allocate space
1376  vertices = Kokkos::createDynRankView(vertices,"vertices",elements.size(),masterVertexCount,getDimension());
1377 
1378  // loop over each requested element
1379  unsigned dim = getDimension();
1380  for(std::size_t cell=0;cell<elements.size();cell++) {
1381  stk::mesh::Entity element = elements[cell];
1382  TEUCHOS_ASSERT(element!=0);
1383 
1384  unsigned vertexCount
1385  = stk::mesh::get_cell_topology(bulkData_->bucket(element).topology()).getCellTopologyData()->vertex_count;
1386  TEUCHOS_TEST_FOR_EXCEPTION(vertexCount!=masterVertexCount,std::runtime_error,
1387  "In call to STK_Interface::getElementVertices all elements "
1388  "must have the same vertex count!");
1389 
1390  // loop over all element nodes
1391  const size_t num_nodes = bulkData_->num_nodes(element);
1392  stk::mesh::Entity const* nodes = bulkData_->begin_nodes(element);
1393  TEUCHOS_TEST_FOR_EXCEPTION(num_nodes!=masterVertexCount,std::runtime_error,
1394  "In call to STK_Interface::getElementVertices cardinality of "
1395  "element node relations must be the vertex count!");
1396  for(std::size_t node=0; node<num_nodes; ++node) {
1397  const double * coord = getNodeCoordinates(nodes[node]);
1398 
1399  // set each dimension of the coordinate
1400  for(unsigned d=0;d<dim;d++)
1401  vertices(cell,node,d) = coord[d];
1402  }
1403  }
1404 }
1405 
1406 template <typename ArrayT>
1407 void STK_Interface::getElementVertices_FromCoordsNoResize(const std::vector<stk::mesh::Entity> & elements, ArrayT & vertices) const
1408 {
1409  // nothing to do! silently return
1410  if(elements.size()==0) {
1411  return;
1412  }
1413 
1414  //
1415  // gather from the intrinsic mesh coordinates (non-lagrangian)
1416  //
1417 
1418  // get *master* cell toplogy...(belongs to first element)
1419  unsigned masterVertexCount
1420  = stk::mesh::get_cell_topology(bulkData_->bucket(elements[0]).topology()).getCellTopologyData()->vertex_count;
1421 
1422  // loop over each requested element
1423  unsigned dim = getDimension();
1424  for(std::size_t cell=0;cell<elements.size();cell++) {
1425  stk::mesh::Entity element = elements[cell];
1426  TEUCHOS_ASSERT(element!=0);
1427 
1428  unsigned vertexCount
1429  = stk::mesh::get_cell_topology(bulkData_->bucket(element).topology()).getCellTopologyData()->vertex_count;
1430  TEUCHOS_TEST_FOR_EXCEPTION(vertexCount!=masterVertexCount,std::runtime_error,
1431  "In call to STK_Interface::getElementVertices all elements "
1432  "must have the same vertex count!");
1433 
1434  // loop over all element nodes
1435  const size_t num_nodes = bulkData_->num_nodes(element);
1436  stk::mesh::Entity const* nodes = bulkData_->begin_nodes(element);
1437  TEUCHOS_TEST_FOR_EXCEPTION(num_nodes!=masterVertexCount,std::runtime_error,
1438  "In call to STK_Interface::getElementVertices cardinality of "
1439  "element node relations must be the vertex count!");
1440  for(std::size_t node=0; node<num_nodes; ++node) {
1441  const double * coord = getNodeCoordinates(nodes[node]);
1442 
1443  // set each dimension of the coordinate
1444  for(unsigned d=0;d<dim;d++)
1445  vertices(cell,node,d) = coord[d];
1446  }
1447  }
1448 }
1449 
1450 template <typename ArrayT>
1451 void STK_Interface::getElementVertices_FromField(const std::vector<stk::mesh::Entity> & elements,const std::string & eBlock, ArrayT & vertices) const
1452 {
1454 
1455  // nothing to do! silently return
1456  if(elements.size()==0) {
1457  vertices = Kokkos::createDynRankView(vertices,"vertices",0,0,0);
1458  return;
1459  }
1460 
1461  // get *master* cell toplogy...(belongs to first element)
1462  unsigned masterVertexCount
1463  = stk::mesh::get_cell_topology(bulkData_->bucket(elements[0]).topology()).getCellTopologyData()->vertex_count;
1464 
1465  // allocate space
1466  vertices = Kokkos::createDynRankView(vertices,"vertices",elements.size(),masterVertexCount,getDimension());
1467 
1468  std::map<std::string,std::vector<std::string> >::const_iterator itr = meshCoordFields_.find(eBlock);
1469  if(itr==meshCoordFields_.end()) {
1470  // no coordinate field set for this element block
1471  TEUCHOS_ASSERT(false);
1472  }
1473 
1474  const std::vector<std::string> & coordField = itr->second;
1475  std::vector<SolutionFieldType*> fields(getDimension());
1476  for(std::size_t d=0;d<fields.size();d++) {
1477  fields[d] = this->getSolutionField(coordField[d],eBlock);
1478  }
1479 
1480  for(std::size_t cell=0;cell<elements.size();cell++) {
1481  stk::mesh::Entity element = elements[cell];
1482 
1483  // loop over nodes set solution values
1484  const size_t num_nodes = bulkData_->num_nodes(element);
1485  stk::mesh::Entity const* nodes = bulkData_->begin_nodes(element);
1486  for(std::size_t i=0; i<num_nodes; ++i) {
1487  stk::mesh::Entity node = nodes[i];
1488 
1489  const double * coord = getNodeCoordinates(node);
1490 
1491  for(unsigned d=0;d<getDimension();d++) {
1492  double * solnData = stk::mesh::field_data(*fields[d],node);
1493 
1494  // recall mesh field coordinates are stored as displacements
1495  // from the mesh coordinates, make sure to add them together
1496  vertices(cell,i,d) = solnData[0]+coord[d];
1497  }
1498  }
1499  }
1500 }
1501 
1502 template <typename ArrayT>
1503 void STK_Interface::getElementVertices_FromFieldNoResize(const std::vector<stk::mesh::Entity> & elements,
1504  const std::string & eBlock, ArrayT & vertices) const
1505 {
1507 
1508  // nothing to do! silently return
1509  if(elements.size()==0) {
1510  return;
1511  }
1512 
1513  std::map<std::string,std::vector<std::string> >::const_iterator itr = meshCoordFields_.find(eBlock);
1514  if(itr==meshCoordFields_.end()) {
1515  // no coordinate field set for this element block
1516  TEUCHOS_ASSERT(false);
1517  }
1518 
1519  const std::vector<std::string> & coordField = itr->second;
1520  std::vector<SolutionFieldType*> fields(getDimension());
1521  for(std::size_t d=0;d<fields.size();d++) {
1522  fields[d] = this->getSolutionField(coordField[d],eBlock);
1523  }
1524 
1525  for(std::size_t cell=0;cell<elements.size();cell++) {
1526  stk::mesh::Entity element = elements[cell];
1527 
1528  // loop over nodes set solution values
1529  const size_t num_nodes = bulkData_->num_nodes(element);
1530  stk::mesh::Entity const* nodes = bulkData_->begin_nodes(element);
1531  for(std::size_t i=0; i<num_nodes; ++i) {
1532  stk::mesh::Entity node = nodes[i];
1533 
1534  const double * coord = getNodeCoordinates(node);
1535 
1536  for(unsigned d=0;d<getDimension();d++) {
1537  double * solnData = stk::mesh::field_data(*fields[d],node);
1538 
1539  // recall mesh field coordinates are stored as displacements
1540  // from the mesh coordinates, make sure to add them together
1541  vertices(cell,i,d) = solnData[0]+coord[d];
1542  }
1543  }
1544  }
1545 }
1546 
1547 }
1548 
1549 #endif
Teuchos::RCP< ElementDescriptor > buildElementDescriptor(stk::mesh::EntityId elmtId, std::vector< stk::mesh::EntityId > &nodes)
Teuchos::RCP< Teuchos::MpiComm< int > > getSafeCommunicator(stk::ParallelMachine parallelMach) const
void getSidesetNames(std::vector< std::string > &name) const
std::pair< Teuchos::RCP< std::vector< std::pair< std::size_t, std::size_t > > >, Teuchos::RCP< std::vector< unsigned int > > > getPeriodicNodePairing() const
Kokkos::DynRankView< typename InputArray::value_type, PHX::Device > createDynRankView(const InputArray &a, const std::string &name, const DimensionPack...dims)
Wrapper to simplify Panzer use of Sacado ViewFactory.
void addNodeset(const std::string &name)
Teuchos::RCP< stk::mesh::BulkData > bulkData_
void print(std::ostream &os) const
unsigned entityOwnerRank(stk::mesh::Entity entity) const
Teuchos::RCP< const std::vector< stk::mesh::Entity > > getElementsOrderedByLID() const
void setSolutionFieldData(const std::string &fieldName, const std::string &blockId, const std::vector< std::size_t > &localElementIds, const ArrayT &solutionValues, double scaleValue=1.0)
std::vector< stk::mesh::Part * > facesPartVec_
std::vector< stk::mesh::EntityId > nodes_
std::map< std::pair< std::string, std::string >, SolutionFieldType * > fieldNameToSolution_
void getElementVerticesNoResize(const std::vector< std::size_t > &localIds, ArrayT &vertices) const
void addEntityToNodeset(stk::mesh::Entity entity, stk::mesh::Part *nodeset)
void addGlobalToExodus(const std::string &key, const int &value)
Add an int global variable to the information to be written to the Exodus output file.
void setDispFieldData(const std::string &fieldName, const std::string &blockId, int axis, const std::vector< std::size_t > &localElementIds, const ArrayT &solutionValues)
void getElementBlockNames(std::vector< std::string > &names) const
void addPeriodicBC(const Teuchos::RCP< const PeriodicBC_MatcherBase > &bc)
void getSubcellIndices(unsigned entityRank, stk::mesh::EntityId elementId, std::vector< stk::mesh::EntityId > &subcellIds) const
std::vector< Teuchos::RCP< const PeriodicBC_MatcherBase > > & getPeriodicBCVector()
void getElementsSharingNodes(const std::vector< stk::mesh::EntityId > nodeId, std::vector< stk::mesh::Entity > &elements) const
get a set of elements sharing multiple nodes
stk::mesh::Field< double > SolutionFieldType
stk::mesh::Part * getElementBlockPart(const std::string &name) const
get the block count
std::vector< stk::mesh::Part * > edgesPartVec_
std::map< std::string, double > blockWeights_
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
std::size_t getEntityCounts(unsigned entityRank) const
get the global counts for the entity of specified rank
void printMetaData(std::ostream &os) const
void getElementVertices(const std::vector< std::size_t > &localIds, ArrayT &vertices) const
void addEntityToSideset(stk::mesh::Entity entity, stk::mesh::Part *sideset)
void getElementVertices_FromCoords(const std::vector< stk::mesh::Entity > &elements, ArrayT &vertices) const
void addSolutionField(const std::string &fieldName, const std::string &blockId)
std::map< std::string, Teuchos::RCP< shards::CellTopology > > elementBlockCT_
stk::mesh::EntityId getMaxEntityId(unsigned entityRank) const
get max entity ID of type entityRank
std::map< std::string, stk::mesh::Part * > sidesets_
std::vector< std::size_t > entityCounts_
bool validBlockId(const std::string &blockId) const
stk::mesh::Part * getNodeset(const std::string &name) const
stk::mesh::Field< double > * getSolutionField(const std::string &fieldName, const std::string &blockId) const
const double * getNodeCoordinates(stk::mesh::EntityId nodeId) const
stk::mesh::Part * getSideset(const std::string &name) const
static const std::string nodesString
void instantiateBulkData(stk::ParallelMachine parallelMach)
void getElementVertices_FromFieldNoResize(const std::vector< stk::mesh::Entity > &elements, const std::string &eBlock, ArrayT &vertices) const
void getMySides(const std::string &sideName, std::vector< stk::mesh::Entity > &sides) const
std::size_t getNumSidesets() const
get the side set count
stk::mesh::EntityId elementGlobalId(stk::mesh::Entity elmt) const
std::size_t elementLocalId(stk::mesh::Entity elmt) const
void getElementVertices_FromCoordsNoResize(const std::vector< stk::mesh::Entity > &elements, ArrayT &vertices) const
unsigned getDimension() const
get the dimension
void addElement(const Teuchos::RCP< ElementDescriptor > &ed, stk::mesh::Part *block)
std::map< std::string, std::vector< std::string > > meshDispFields_
void setUseFieldCoordinates(bool useFieldCoordinates)
stk::mesh::EntityRank getSideRank() const
void initialize(stk::ParallelMachine parallelMach, bool setupIO=true, const bool buildRefinementSupport=false)
static const std::string coordsString
stk::mesh::Field< double, stk::mesh::Cartesian > VectorFieldType
std::size_t getNumElementBlocks() const
get the block count
Teuchos::RCP< stk::mesh::BulkData > getBulkData() const
void getElementVertices_FromField(const std::vector< stk::mesh::Entity > &elements, const std::string &eBlock, ArrayT &vertices) const
void writeToExodus(const std::string &filename)
Write this mesh and associated fields to the given output file.
stk::mesh::EntityId elementGlobalId(std::size_t lid) const
std::map< std::string, stk::mesh::Part * > elementBlocks_
std::vector< stk::mesh::Part * > nodesPartVec_
static const std::string edgesString
std::string containingBlockId(stk::mesh::Entity elmt) const
stk::mesh::EntityRank getFaceRank() const
std::vector< Teuchos::RCP< const PeriodicBC_MatcherBase > > periodicBCs_
void getOwnedElementsSharingNode(stk::mesh::Entity node, std::vector< stk::mesh::Entity > &elements, std::vector< int > &relIds) const
void addNode(stk::mesh::EntityId gid, const std::vector< double > &coord)
void setBlockWeight(const std::string &blockId, double weight)
void buildSubcells()
force the mesh to build subcells: edges and faces
Teuchos::RCP< Teuchos::MpiComm< int > > mpiComm_
void setInitialStateTime(double value)
Teuchos::RCP< stk::mesh::MetaData > metaData_
Teuchos::RCP< stk::mesh::MetaData > getMetaData() const
Teuchos::RCP< std::vector< stk::mesh::Entity > > orderedElementVector_
std::size_t getNumNodesets() const
get the side set count
bool isInitialized() const
Has initialize been called on this mesh object?
void getMyNodes(const std::string &sideName, const std::string &blockName, std::vector< stk::mesh::Entity > &nodes) const
void getSolutionFieldData(const std::string &fieldName, const std::string &blockId, const std::vector< std::size_t > &localElementIds, ArrayT &solutionValues) const
void addMeshCoordFields(const std::string &blockId, const std::vector< std::string > &coordField, const std::string &dispPrefix)
stk::mesh::Entity findConnectivityById(stk::mesh::Entity src, stk::mesh::EntityRank tgt_rank, unsigned rel_id) const
void setUseLowerCaseForIO(bool useLowerCase)
bool isMeshCoordField(const std::string &eBlock, const std::string &fieldName, int &axis) const
void getAllSides(const std::string &sideName, std::vector< stk::mesh::Entity > &sides) const
void addSideset(const std::string &name, const CellTopologyData *ctData)
void setCellFieldData(const std::string &fieldName, const std::string &blockId, const std::vector< std::size_t > &localElementIds, const ArrayT &solutionValues, double scaleValue=1.0)
const VectorFieldType & getEdgesField() const
ProcIdFieldType * getProcessorIdField()
PHX::MDField< ScalarT, panzer::Cell, panzer::BASIS > field
A field to which we&#39;ll contribute, or in which we&#39;ll store, the result of computing this integral...
void rebalance(const Teuchos::ParameterList &params)
std::map< std::string, std::vector< std::string > > meshCoordFields_
stk::mesh::EntityRank getEdgeRank() const
static const std::string facesString
void initializeFieldsInSTK(const std::map< std::pair< std::string, std::string >, SolutionFieldType * > &nameToField, bool setupIO)
Teuchos::RCP< const Teuchos::Comm< int > > getComm() const
get the comm associated with this mesh
const std::vector< Teuchos::RCP< const PeriodicBC_MatcherBase > > & getPeriodicBCVector() const
const VectorFieldType & getFacesField() const
void refineMesh(const int numberOfLevels, const bool deleteParentElements)
void setupExodusFile(const std::string &filename)
Set up an output Exodus file for writing results.
std::unordered_map< stk::mesh::EntityId, std::size_t > localIDHash_
bool isValid(stk::mesh::Entity entity) const
std::vector< stk::mesh::EntityId > maxEntityId_
stk::mesh::Field< ProcIdData > ProcIdFieldType
void addPeriodicBCs(const std::vector< Teuchos::RCP< const PeriodicBC_MatcherBase > > &bc_vec)
void getElementsSharingNode(stk::mesh::EntityId nodeId, std::vector< stk::mesh::Entity > &elements) const
get a set of elements sharing a single node
#define TEUCHOS_ASSERT(assertion_test)
void getMyElements(std::vector< stk::mesh::Entity > &elements) const
stk::mesh::EntityRank getNodeRank() const
void addCellField(const std::string &fieldName, const std::string &blockId)
void getNodeIdsForElement(stk::mesh::Entity element, std::vector< stk::mesh::EntityId > &nodeIds) const
get a list of node ids for nodes connected to an element
bool operator()(stk::mesh::Entity a, stk::mesh::Entity b)
stk::mesh::EntityRank getElementRank() const
void addElementBlock(const std::string &name, const CellTopologyData *ctData)
const std::vector< stk::mesh::EntityId > & getNodes() const
stk::mesh::Field< double > * getCellField(const std::string &fieldName, const std::string &blockId) const
const VectorFieldType & getCoordinatesField() const
stk::mesh::Part * getOwnedPart() const
Get a pointer to the locally owned part.
stk::mesh::EntityId getGID() const
void getNeighborElements(std::vector< stk::mesh::Entity > &elements) const
void getNodesetNames(std::vector< std::string > &name) const
Teuchos::RCP< const shards::CellTopology > getCellTopology(const std::string &eBlock) const
std::map< std::pair< std::string, std::string >, SolutionFieldType * > fieldNameToCellField_
std::map< std::string, stk::mesh::Part * > nodesets_