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 // Panzer: A partial differential equation assembly
4 // engine for strongly coupled complex multiphysics systems
5 //
6 // Copyright 2011 NTESS and the Panzer contributors.
7 // SPDX-License-Identifier: BSD-3-Clause
8 // *****************************************************************************
9 // @HEADER
10 
11 #ifndef __Panzer_STK_Interface_hpp__
12 #define __Panzer_STK_Interface_hpp__
13 
14 #include <Teuchos_RCP.hpp>
16 
17 #include <stk_mesh/base/Types.hpp>
18 #include <stk_mesh/base/MetaData.hpp>
19 #include <stk_mesh/base/BulkData.hpp>
20 #include <stk_mesh/base/Field.hpp>
21 #include <stk_mesh/base/FieldBase.hpp>
22 
23 #include "Kokkos_Core.hpp"
24 
25 #include <Shards_CellTopology.hpp>
26 #include <Shards_CellTopologyData.h>
27 
28 #include <PanzerAdaptersSTK_config.hpp>
29 #include <Kokkos_ViewFactory.hpp>
30 
31 #include <unordered_map>
32 
33 #ifdef PANZER_HAVE_IOSS
34 #include <stk_io/StkMeshIoBroker.hpp>
35 #endif
36 
37 #ifdef PANZER_HAVE_PERCEPT
38 namespace percept {
39  class PerceptMesh;
40  class URP_Heterogeneous_3D;
41 }
42 #endif
43 
44 namespace panzer_stk {
45 
46 class PeriodicBC_MatcherBase;
47 
52 public:
53  ElementDescriptor(stk::mesh::EntityId gid,const std::vector<stk::mesh::EntityId> & nodes);
54  virtual ~ElementDescriptor();
55 
56  stk::mesh::EntityId getGID() const { return gid_; }
57  const std::vector<stk::mesh::EntityId> & getNodes() const { return nodes_; }
58 protected:
59  stk::mesh::EntityId gid_;
60  std::vector<stk::mesh::EntityId> nodes_;
61 
63 };
64 
68 buildElementDescriptor(stk::mesh::EntityId elmtId,std::vector<stk::mesh::EntityId> & nodes);
69 
71 public:
72  typedef double ProcIdData; // ECC: Not sure why?
73  typedef stk::mesh::Field<double> SolutionFieldType;
74  typedef stk::mesh::Field<double> VectorFieldType;
75  typedef stk::mesh::Field<ProcIdData> ProcIdFieldType;
76 
77  // some simple exception classes
78  struct ElementBlockException : public std::logic_error
79  { ElementBlockException(const std::string & what) : std::logic_error(what) {} };
80 
81  struct SidesetException : public std::logic_error
82  { SidesetException(const std::string & what) : std::logic_error(what) {} };
83 
84  struct EdgeBlockException : public std::logic_error
85  { EdgeBlockException(const std::string & what) : std::logic_error(what) {} };
86 
87  struct FaceBlockException : public std::logic_error
88  { FaceBlockException(const std::string & what) : std::logic_error(what) {} };
89 
90  STK_Interface();
91 
94  STK_Interface(unsigned dim);
95 
97 
98  // functions called before initialize
100 
103  void addElementBlock(const std::string & name,const CellTopologyData * ctData);
104 
107 // void addEdgeBlock(const std::string & name,const CellTopologyData * ctData);
108  void addEdgeBlock(const std::string & elemBlockName,
109  const std::string & edgeBlockName,
110  const stk::topology & topology);
111  void addEdgeBlock(const std::string & elemBlockName,
112  const std::string & edgeBlockName,
113  const CellTopologyData * ctData);
114 
117  void addFaceBlock(const std::string & elemBlockName,
118  const std::string & faceBlockName,
119  const stk::topology & topology);
120  void addFaceBlock(const std::string & elemBlockName,
121  const std::string & faceBlockName,
122  const CellTopologyData * ctData);
123 
126  void addSideset(const std::string & name,const CellTopologyData * ctData);
127 
130  void addNodeset(const std::string & name);
131 
134  void addSolutionField(const std::string & fieldName,const std::string & blockId);
135 
138  void addCellField(const std::string & fieldName,const std::string & blockId);
139 
142  void addEdgeField(const std::string & fieldName,const std::string & blockId);
143 
146  void addFaceField(const std::string & fieldName,const std::string & blockId);
147 
156  void addMeshCoordFields(const std::string & blockId,
157  const std::vector<std::string> & coordField,
158  const std::string & dispPrefix);
159 
164  void addInformationRecords(const std::vector<std::string> & info_records);
165 
167 
179  void initialize(stk::ParallelMachine parallelMach,bool setupIO=true,
180  const bool buildRefinementSupport = false);
181 
187  void instantiateBulkData(stk::ParallelMachine parallelMach);
188 
189  // functions to manage and manipulate bulk data
191 
194  void beginModification();
195 
198  void endModification();
199 
205  void addNode(stk::mesh::EntityId gid, const std::vector<double> & coord);
206 
207  void addElement(const Teuchos::RCP<ElementDescriptor> & ed,stk::mesh::Part * block);
208 
209  void addEdges();
210 
211  void addFaces();
212 
215  void addEntityToSideset(stk::mesh::Entity entity,stk::mesh::Part * sideset);
216 
219  void addEntityToNodeset(stk::mesh::Entity entity,stk::mesh::Part * nodeset);
220 
223  void addEntityToEdgeBlock(stk::mesh::Entity entity,stk::mesh::Part * edgeblock);
226  void addEntitiesToEdgeBlock(std::vector<stk::mesh::Entity> entities,stk::mesh::Part * edgeblock);
227 
230  void addEntityToFaceBlock(stk::mesh::Entity entity,stk::mesh::Part * faceblock);
233  void addEntitiesToFaceBlock(std::vector<stk::mesh::Entity> entities,stk::mesh::Part * faceblock);
234 
235  // Methods to interrogate the mesh topology and structure
237 
241  { return *coordinatesField_; }
242 
246  { return *edgesField_; }
247 
249  { return *facesField_; }
250 
253  const double * getNodeCoordinates(stk::mesh::EntityId nodeId) const;
254 
257  const double * getNodeCoordinates(stk::mesh::Entity node) const;
258 
261  void getSubcellIndices(unsigned entityRank,stk::mesh::EntityId elementId,
262  std::vector<stk::mesh::EntityId> & subcellIds) const;
263 
266  void getMyElements(std::vector<stk::mesh::Entity> & elements) const;
267 
270  void getMyElements(const std::string & blockID,std::vector<stk::mesh::Entity> & elements) const;
271 
275  void getNeighborElements(std::vector<stk::mesh::Entity> & elements) const;
276 
279  void getNeighborElements(const std::string & blockID,std::vector<stk::mesh::Entity> & elements) const;
280 
283  void getMyEdges(std::vector<stk::mesh::Entity> & edges) const;
284 
292  void getMyEdges(const std::string & edgeBlockName,std::vector<stk::mesh::Entity> & edges) const;
293 
302  void getMyEdges(const std::string & edgeBlockName,const std::string & blockName,std::vector<stk::mesh::Entity> & edges) const;
303 
311  void getAllEdges(const std::string & edgeBlockName,std::vector<stk::mesh::Entity> & edges) const;
312 
321  void getAllEdges(const std::string & edgeBlockName,const std::string & blockName,std::vector<stk::mesh::Entity> & edges) const;
322 
325  void getMyFaces(std::vector<stk::mesh::Entity> & faces) const;
326 
334  void getMyFaces(const std::string & faceBlockName,std::vector<stk::mesh::Entity> & faces) const;
335 
344  void getMyFaces(const std::string & faceBlockName,const std::string & blockName,std::vector<stk::mesh::Entity> & faces) const;
345 
353  void getAllFaces(const std::string & faceBlockName,std::vector<stk::mesh::Entity> & faces) const;
354 
363  void getAllFaces(const std::string & faceBlockName,const std::string & blockName,std::vector<stk::mesh::Entity> & faces) const;
364 
372  void getMySides(const std::string & sideName,std::vector<stk::mesh::Entity> & sides) const;
373 
382  void getMySides(const std::string & sideName,const std::string & blockName,std::vector<stk::mesh::Entity> & sides) const;
383 
391  void getAllSides(const std::string & sideName,std::vector<stk::mesh::Entity> & sides) const;
392 
401  void getAllSides(const std::string & sideName,const std::string & blockName,std::vector<stk::mesh::Entity> & sides) const;
402 
412  void getMyNodes(const std::string & sideName,const std::string & blockName,std::vector<stk::mesh::Entity> & nodes) const;
413 
422  stk::mesh::Entity findConnectivityById(stk::mesh::Entity src, stk::mesh::EntityRank tgt_rank, unsigned rel_id) const;
423 
424  // Utility functions
426 
435  void
436  writeToExodus(const std::string& filename,
437  const bool append = false);
438 
456  void
457  setupExodusFile(const std::string& filename,
458  const bool append = false,
459  const bool append_after_restart_time = false,
460  const double restart_time = 0.0);
461 
474  void
476  double timestep);
477 
497  void
499  const std::string& key,
500  const int& value);
501 
521  void
523  const std::string& key,
524  const double& value);
525 
545  void
547  const std::string& key,
548  const std::vector<int>& value);
549 
569  void
571  const std::string& key,
572  const std::vector<double>& value);
573 
574  // Accessor functions
576 
579 
582 
583 #ifdef PANZER_HAVE_PERCEPT
584  Teuchos::RCP<percept::PerceptMesh> getRefinedMesh() const
586  { TEUCHOS_ASSERT(Teuchos::nonnull(refinedMesh_)); return refinedMesh_; }
587 #endif
588 
589  bool isWritable() const;
590 
591  bool isModifiable() const
592  { if(bulkData_==Teuchos::null) return false;
593  return bulkData_->in_modifiable_state(); }
594 
596  unsigned getDimension() const
597  { return dimension_; }
598 
600  std::size_t getNumElementBlocks() const
601  { return elementBlocks_.size(); }
602 
610  void getElementBlockNames(std::vector<std::string> & names) const;
611 
619  void getSidesetNames(std::vector<std::string> & name) const;
620 
628  void getNodesetNames(std::vector<std::string> & name) const;
629 
637  void getEdgeBlockNames(std::vector<std::string> & names) const;
638 
646  void getFaceBlockNames(std::vector<std::string> & names) const;
647 
649  stk::mesh::Part * getOwnedPart() const
650  { return &getMetaData()->locally_owned_part(); } // I don't like the pointer access here, but it will do for now!
651 
653  stk::mesh::Part * getElementBlockPart(const std::string & name) const
654  {
655  std::map<std::string, stk::mesh::Part*>::const_iterator itr = elementBlocks_.find(name); // Element blocks
656  if(itr==elementBlocks_.end()) return 0;
657  return itr->second;
658  }
659 
661  stk::mesh::Part * getEdgeBlock(const std::string & name) const
662  {
663  std::map<std::string, stk::mesh::Part*>::const_iterator itr = edgeBlocks_.find(name); // edge blocks
664  if(itr==edgeBlocks_.end()) return 0;
665  return itr->second;
666  }
667 
669  stk::mesh::Part * getFaceBlock(const std::string & name) const
670  {
671  std::map<std::string, stk::mesh::Part*>::const_iterator itr = faceBlocks_.find(name); // face blocks
672  if(itr==faceBlocks_.end()) return 0;
673  return itr->second;
674  }
675 
677  std::size_t getNumSidesets() const
678  { return sidesets_.size(); }
679 
680  stk::mesh::Part * getSideset(const std::string & name) const
681  {
682  auto itr = sidesets_.find(name);
683  return (itr != sidesets_.end()) ? itr->second : nullptr;
684  }
685 
687  std::size_t getNumNodesets() const
688  { return nodesets_.size(); }
689 
690  stk::mesh::Part * getNodeset(const std::string & name) const
691  {
692  auto itr = nodesets_.find(name);
693  return (itr != nodesets_.end()) ? itr->second : nullptr;
694  }
695 
697  std::size_t getEntityCounts(unsigned entityRank) const;
698 
700  stk::mesh::EntityId getMaxEntityId(unsigned entityRank) const;
701 
702  // Utilities
704 
706  void getElementsSharingNode(stk::mesh::EntityId nodeId,std::vector<stk::mesh::Entity> & elements) const;
707 
709  void getNodeIdsForElement(stk::mesh::Entity element,std::vector<stk::mesh::EntityId> & nodeIds) const;
710 
713  void getOwnedElementsSharingNode(stk::mesh::Entity node,std::vector<stk::mesh::Entity> & elements,
714  std::vector<int> & relIds) const;
715 
718  void getOwnedElementsSharingNode(stk::mesh::EntityId nodeId,std::vector<stk::mesh::Entity> & elements,
719  std::vector<int> & relIds, unsigned int matchType) const;
720 
721 
723  void getElementsSharingNodes(const std::vector<stk::mesh::EntityId> nodeId,std::vector<stk::mesh::Entity> & elements) const;
724 
726  void buildSubcells();
727 
730  std::size_t elementLocalId(stk::mesh::Entity elmt) const;
731 
734  std::size_t elementLocalId(stk::mesh::EntityId gid) const;
735 
738  inline stk::mesh::EntityId elementGlobalId(std::size_t lid) const
739  { return bulkData_->identifier((*orderedElementVector_)[lid]); }
740 
743  inline stk::mesh::EntityId elementGlobalId(stk::mesh::Entity elmt) const
744  { return bulkData_->identifier(elmt); }
745 
748  bool isEdgeLocal(stk::mesh::Entity edge) const;
749 
752  bool isEdgeLocal(stk::mesh::EntityId gid) const;
753 
756  std::size_t edgeLocalId(stk::mesh::Entity elmt) const;
757 
760  std::size_t edgeLocalId(stk::mesh::EntityId gid) const;
761 
764  inline stk::mesh::EntityId edgeGlobalId(std::size_t lid) const
765  { return bulkData_->identifier((*orderedEdgeVector_)[lid]); }
766 
769  inline stk::mesh::EntityId edgeGlobalId(stk::mesh::Entity edge) const
770  { return bulkData_->identifier(edge); }
771 
774  bool isFaceLocal(stk::mesh::Entity face) const;
775 
778  bool isFaceLocal(stk::mesh::EntityId gid) const;
779 
782  std::size_t faceLocalId(stk::mesh::Entity elmt) const;
783 
786  std::size_t faceLocalId(stk::mesh::EntityId gid) const;
787 
790  inline stk::mesh::EntityId faceGlobalId(std::size_t lid) const
791  { return bulkData_->identifier((*orderedFaceVector_)[lid]); }
792 
795  inline stk::mesh::EntityId faceGlobalId(stk::mesh::Entity face) const
796  { return bulkData_->identifier(face); }
797 
800  inline unsigned entityOwnerRank(stk::mesh::Entity entity) const
801  { return bulkData_->parallel_owner_rank(entity); }
802 
805  inline bool isValid(stk::mesh::Entity entity) const
806  { return bulkData_->is_valid(entity); }
807 
810  std::string containingBlockId(stk::mesh::Entity elmt) const;
811 
816  stk::mesh::Field<double> * getSolutionField(const std::string & fieldName,
817  const std::string & blockId) const;
818 
823  stk::mesh::Field<double> * getCellField(const std::string & fieldName,
824  const std::string & blockId) const;
825 
830  stk::mesh::Field<double> * getEdgeField(const std::string & fieldName,
831  const std::string & blockId) const;
832 
837  stk::mesh::Field<double> * getFaceField(const std::string & fieldName,
838  const std::string & blockId) const;
839 
841 
843  bool isInitialized() const { return initialized_; }
844 
849 
864  template <typename ArrayT>
865  void setSolutionFieldData(const std::string & fieldName,const std::string & blockId,
866  const std::vector<std::size_t> & localElementIds,const ArrayT & solutionValues,double scaleValue=1.0);
867 
882  template <typename ArrayT>
883  void getSolutionFieldData(const std::string & fieldName,const std::string & blockId,
884  const std::vector<std::size_t> & localElementIds,ArrayT & solutionValues) const;
885 
900  template <typename ArrayT>
901  void setCellFieldData(const std::string & fieldName,const std::string & blockId,
902  const std::vector<std::size_t> & localElementIds,const ArrayT & solutionValues,double scaleValue=1.0);
903 
908 
913 
928  template <typename ArrayT>
929  void setEdgeFieldData(const std::string & fieldName,const std::string & blockId,
930  const std::vector<std::size_t> & localEdgeIds,const ArrayT & edgeValues,double scaleValue=1.0);
931 
946  template <typename ArrayT>
947  void setFaceFieldData(const std::string & fieldName,const std::string & blockId,
948  const std::vector<std::size_t> & localFaceIds,const ArrayT & faceValues,double scaleValue=1.0);
949 
951 
960  template <typename ArrayT>
961  void getElementVertices(const std::vector<std::size_t> & localIds, ArrayT & vertices) const;
962 
971  template <typename ArrayT>
972  void getElementVertices(const std::vector<stk::mesh::Entity> & elements, ArrayT & vertices) const;
973 
983  template <typename ArrayT>
984  void getElementVertices(const std::vector<std::size_t> & localIds,const std::string & eBlock, ArrayT & vertices) const;
985 
995  template <typename ArrayT>
996  void getElementVertices(const std::vector<stk::mesh::Entity> & elements,const std::string & eBlock, ArrayT & vertices) const;
997 
1006  template <typename ArrayT>
1007  void getElementVerticesNoResize(const std::vector<std::size_t> & localIds, ArrayT & vertices) const;
1008 
1017  template <typename ArrayT>
1018  void getElementVerticesNoResize(const std::vector<stk::mesh::Entity> & elements, ArrayT & vertices) const;
1019 
1029  template <typename ArrayT>
1030  void getElementVerticesNoResize(const std::vector<std::size_t> & localIds,const std::string & eBlock, ArrayT & vertices) const;
1031 
1041  template <typename ArrayT>
1042  void getElementVerticesNoResize(const std::vector<stk::mesh::Entity> & elements,const std::string & eBlock, ArrayT & vertices) const;
1043 
1045 
1055  template <typename ArrayT>
1056  void getElementNodes(const std::vector<std::size_t> & localIds, ArrayT & nodes) const;
1057 
1066  template <typename ArrayT>
1067  void getElementNodes(const std::vector<stk::mesh::Entity> & elements, ArrayT & nodes) const;
1068 
1078  template <typename ArrayT>
1079  void getElementNodes(const std::vector<std::size_t> & localIds,const std::string & eBlock, ArrayT & nodes) const;
1080 
1090  template <typename ArrayT>
1091  void getElementNodes(const std::vector<stk::mesh::Entity> & elements,const std::string & eBlock, ArrayT & nodes) const;
1092 
1101  template <typename ArrayT>
1102  void getElementNodesNoResize(const std::vector<std::size_t> & localIds, ArrayT & nodes) const;
1103 
1112  template <typename ArrayT>
1113  void getElementNodesNoResize(const std::vector<stk::mesh::Entity> & elements, ArrayT & nodes) const;
1114 
1124  template <typename ArrayT>
1125  void getElementNodesNoResize(const std::vector<std::size_t> & localIds,const std::string & eBlock, ArrayT & nodes) const;
1126 
1136  template <typename ArrayT>
1137  void getElementNodesNoResize(const std::vector<stk::mesh::Entity> & elements,const std::string & eBlock, ArrayT & nodes) const;
1138 
1139  // const stk::mesh::FEMInterface & getFEMInterface() const
1140  // { return *femPtr_; }
1141 
1142  stk::mesh::EntityRank getElementRank() const { return stk::topology::ELEMENT_RANK; }
1143  stk::mesh::EntityRank getSideRank() const { return metaData_->side_rank(); }
1144  stk::mesh::EntityRank getFaceRank() const { return stk::topology::FACE_RANK; }
1145  stk::mesh::EntityRank getEdgeRank() const { return stk::topology::EDGE_RANK; }
1146  stk::mesh::EntityRank getNodeRank() const { return stk::topology::NODE_RANK; }
1147 
1150  void initializeFromMetaData();
1151 
1154  void buildLocalElementIDs();
1155 
1158  void buildLocalEdgeIDs();
1159 
1162  void buildLocalFaceIDs();
1163 
1166  const std::vector<Teuchos::RCP<const PeriodicBC_MatcherBase> > &
1168  { return periodicBCs_; }
1169 
1172  std::vector<Teuchos::RCP<const PeriodicBC_MatcherBase> > &
1174  { return periodicBCs_; }
1175 
1178  const bool & useBoundingBoxSearch() const
1179  { return useBBoxSearch_; }
1180 
1182  void setBoundingBoxSearchFlag(const bool & searchFlag)
1183  { useBBoxSearch_ = searchFlag; return; }
1184 
1191  { periodicBCs_.push_back(bc); }
1192 
1199  { periodicBCs_.insert(periodicBCs_.end(),bc_vec.begin(),bc_vec.end()); }
1200 
1203  std::pair<Teuchos::RCP<std::vector<std::pair<std::size_t,std::size_t> > >, Teuchos::RCP<std::vector<unsigned int> > >
1204  getPeriodicNodePairing() const;
1205 
1208  bool validBlockId(const std::string & blockId) const;
1209 
1212  void print(std::ostream & os) const;
1213 
1216  void printMetaData(std::ostream & os) const;
1217 
1220  Teuchos::RCP<const shards::CellTopology> getCellTopology(const std::string & eBlock) const;
1221 
1226  double getCurrentStateTime() const { return currentStateTime_; }
1227 
1233  double getInitialStateTime() const { return initialStateTime_; }
1234 
1239  void setInitialStateTime(double value) { initialStateTime_ = value; }
1240 
1243  void rebalance(const Teuchos::ParameterList & params);
1244 
1248  void setBlockWeight(const std::string & blockId,double weight)
1249  { blockWeights_[blockId] = weight; }
1250 
1257  void setUseFieldCoordinates(bool useFieldCoordinates)
1258  { useFieldCoordinates_ = useFieldCoordinates; }
1259 
1262  { return useFieldCoordinates_; }
1263 
1265  void setUseLowerCaseForIO(bool useLowerCase)
1266  { useLowerCase_ = useLowerCase; }
1267 
1270  { return useLowerCase_; }
1271 
1273 
1284  template <typename ArrayT>
1285  void getElementVertices_FromField(const std::vector<stk::mesh::Entity> & elements,const std::string & eBlock, ArrayT & vertices) const;
1286 
1287  template <typename ArrayT>
1288  void getElementVertices_FromFieldNoResize(const std::vector<stk::mesh::Entity> & elements,
1289  const std::string & eBlock, ArrayT & vertices) const;
1290 
1300  template <typename ArrayT>
1301  void getElementVertices_FromCoords(const std::vector<stk::mesh::Entity> & elements, ArrayT & vertices) const;
1302 
1303  template <typename ArrayT>
1304  void getElementVertices_FromCoordsNoResize(const std::vector<stk::mesh::Entity> & elements, ArrayT & vertices) const;
1305 
1307 
1318  template <typename ArrayT>
1319  void getElementNodes_FromField(const std::vector<stk::mesh::Entity> & elements,const std::string & eBlock, ArrayT & nodes) const;
1320 
1321  template <typename ArrayT>
1322  void getElementNodes_FromFieldNoResize(const std::vector<stk::mesh::Entity> & elements,
1323  const std::string & eBlock, ArrayT & nodes) const;
1324 
1334  template <typename ArrayT>
1335  void getElementNodes_FromCoords(const std::vector<stk::mesh::Entity> & elements, ArrayT & nodes) const;
1336 
1337  template <typename ArrayT>
1338  void getElementNodes_FromCoordsNoResize(const std::vector<stk::mesh::Entity> & elements, ArrayT & nodes) const;
1339 
1345  void refineMesh(const int numberOfLevels, const bool deleteParentElements);
1346 
1347 public: // static operations
1348  static const std::string coordsString;
1349  static const std::string nodesString;
1350  static const std::string edgesString;
1351  static const std::string edgeBlockString;
1352  static const std::string faceBlockString;
1353  static const std::string facesString;
1354 
1355 protected:
1356 
1359  void buildEntityCounts();
1360 
1363  void buildMaxEntityIds();
1364 
1368  void initializeFieldsInSTK(const std::map<std::pair<std::string,std::string>,SolutionFieldType*> & nameToField,
1369  bool setupIO);
1370 
1375  Teuchos::RCP<Teuchos::MpiComm<int> > getSafeCommunicator(stk::ParallelMachine parallelMach) const;
1376 
1383 
1387  bool isMeshCoordField(const std::string & eBlock,const std::string & fieldName,int & axis) const;
1388 
1404  template <typename ArrayT>
1405  void setDispFieldData(const std::string & fieldName,const std::string & blockId,int axis,
1406  const std::vector<std::size_t> & localElementIds,const ArrayT & solutionValues);
1407 
1408  std::vector<Teuchos::RCP<const PeriodicBC_MatcherBase> > periodicBCs_;
1409  bool useBBoxSearch_ = false; // TODO swap this to change default periodic BC search (see also PeriodicBC_Parser.cpp)
1410 
1413 #ifdef PANZER_HAVE_PERCEPT
1416 #endif
1417 
1418  std::map<std::string, stk::mesh::Part*> elementBlocks_; // Element blocks
1419  std::map<std::string, stk::mesh::Part*> sidesets_; // Side sets
1420  std::map<std::string, stk::mesh::Part*> nodesets_; // Node sets
1421  std::map<std::string, stk::mesh::Part*> edgeBlocks_; // Edge blocks
1422  std::map<std::string, stk::mesh::Part*> faceBlocks_; // Face blocks
1423 
1424  std::map<std::string, Teuchos::RCP<shards::CellTopology> > elementBlockCT_;
1425 
1426  // for storing/accessing nodes
1427  stk::mesh::Part * nodesPart_;
1428  std::vector<stk::mesh::Part*> nodesPartVec_;
1429  stk::mesh::Part * edgesPart_;
1430  std::vector<stk::mesh::Part*> edgesPartVec_;
1431  stk::mesh::Part * facesPart_;
1432  std::vector<stk::mesh::Part*> facesPartVec_;
1433 
1439 
1440  // maps field names to solution field stk mesh handles
1441  std::map<std::pair<std::string,std::string>,SolutionFieldType*> fieldNameToSolution_;
1442  std::map<std::pair<std::string,std::string>,SolutionFieldType*> fieldNameToCellField_;
1443  std::map<std::pair<std::string,std::string>,SolutionFieldType*> fieldNameToEdgeField_;
1444  std::map<std::pair<std::string,std::string>,SolutionFieldType*> fieldNameToFaceField_;
1445 
1446  // use a set to maintain a list of unique information records
1447  std::set<std::string> informationRecords_;
1448 
1449  unsigned dimension_;
1450 
1452 
1453  // how many elements, faces, edges, and nodes are there globally
1454  std::vector<std::size_t> entityCounts_;
1455 
1456  // what is maximum entity ID
1457  std::vector<stk::mesh::EntityId> maxEntityId_;
1458 
1459  unsigned procRank_;
1460  std::size_t currentLocalId_;
1461 
1463 
1464  double initialStateTime_; // the time stamp at the time this object was constructed (default 0.0)
1465  double currentStateTime_; // the time stamp set by the user when writeToExodus is called (default 0.0)
1466 
1467 #ifdef PANZER_HAVE_IOSS
1468  // I/O support
1470  int meshIndex_;
1471 
1476  enum class GlobalVariable
1477  {
1478  ADD,
1479  WRITE
1480  }; // end of enum class GlobalVariable
1481 
1498  void
1499  globalToExodus(
1500  const GlobalVariable& flag);
1501 
1505  Teuchos::ParameterList globalData_;
1506 #endif
1507 
1508  // uses lazy evaluation
1510 
1511  // uses lazy evaluation
1513 
1514  // uses lazy evaluation
1516 
1517  // for element block weights
1518  std::map<std::string,double> blockWeights_;
1519 
1520  std::unordered_map<stk::mesh::EntityId,std::size_t> localIDHash_;
1521  std::unordered_map<stk::mesh::EntityId,std::size_t> localEdgeIDHash_;
1522  std::unordered_map<stk::mesh::EntityId,std::size_t> localFaceIDHash_;
1523 
1524  // Store mesh displacement fields by element block. This map
1525  // goes like this meshCoordFields_[eBlock][axis_index] => coordinate FieldName
1526  // goes like this meshDispFields_[eBlock][axis_index] => displacement FieldName
1527  std::map<std::string,std::vector<std::string> > meshCoordFields_; // coordinate fields written by user
1528  std::map<std::string,std::vector<std::string> > meshDispFields_; // displacement fields, output to exodus
1529 
1531 
1533 
1534  // Object describing how to sort a vector of elements using
1535  // local ID as the key, very short lived object
1537  public:
1538  LocalIdCompare(const STK_Interface * mesh) : mesh_(mesh) {}
1539 
1540  // Compares two stk mesh entities based on local ID
1541  bool operator() (stk::mesh::Entity a,stk::mesh::Entity b)
1542  { return mesh_->elementLocalId(a) < mesh_->elementLocalId(b);}
1543 
1544  private:
1546  };
1547 };
1548 
1549 template <typename ArrayT>
1550 void STK_Interface::setSolutionFieldData(const std::string & fieldName,const std::string & blockId,
1551  const std::vector<std::size_t> & localElementIds,const ArrayT & solutionValues,double scaleValue)
1552 {
1553  const std::vector<stk::mesh::Entity> & elements = *(this->getElementsOrderedByLID());
1554  auto solutionValues_h = Kokkos::create_mirror_view(solutionValues);
1555  Kokkos::deep_copy(solutionValues_h, solutionValues);
1556 
1557  int field_axis = -1;
1558  if(isMeshCoordField(blockId,fieldName,field_axis)) {
1559  setDispFieldData(fieldName,blockId,field_axis,localElementIds,solutionValues_h);
1560  return;
1561  }
1562 
1563  SolutionFieldType * field = this->getSolutionField(fieldName,blockId);
1564 
1565  for(std::size_t cell=0;cell<localElementIds.size();cell++) {
1566  std::size_t localId = localElementIds[cell];
1567  stk::mesh::Entity element = elements[localId];
1568 
1569  // loop over nodes set solution values
1570  const size_t num_nodes = bulkData_->num_nodes(element);
1571  stk::mesh::Entity const* nodes = bulkData_->begin_nodes(element);
1572  for(std::size_t i=0; i<num_nodes; ++i) {
1573  stk::mesh::Entity node = nodes[i];
1574 
1575  double * solnData = stk::mesh::field_data(*field,node);
1576  // TEUCHOS_ASSERT(solnData!=0); // only needed if blockId is not specified
1577  solnData[0] = scaleValue*solutionValues_h(cell,i);
1578  }
1579  }
1580 }
1581 
1582 template <typename ArrayT>
1583 void STK_Interface::setDispFieldData(const std::string & fieldName,const std::string & blockId,int axis,
1584  const std::vector<std::size_t> & localElementIds,const ArrayT & dispValues)
1585 {
1586  TEUCHOS_ASSERT(axis>=0); // sanity check
1587 
1588  const std::vector<stk::mesh::Entity> & elements = *(this->getElementsOrderedByLID());
1589 
1590  SolutionFieldType * field = this->getSolutionField(fieldName,blockId);
1591  const VectorFieldType & coord_field = this->getCoordinatesField();
1592 
1593  for(std::size_t cell=0;cell<localElementIds.size();cell++) {
1594  std::size_t localId = localElementIds[cell];
1595  stk::mesh::Entity element = elements[localId];
1596 
1597  // loop over nodes set solution values
1598  const size_t num_nodes = bulkData_->num_nodes(element);
1599  stk::mesh::Entity const* nodes = bulkData_->begin_nodes(element);
1600  for(std::size_t i=0; i<num_nodes; ++i) {
1601  stk::mesh::Entity node = nodes[i];
1602 
1603  double * solnData = stk::mesh::field_data(*field,node);
1604  double * coordData = stk::mesh::field_data(coord_field,node);
1605  // TEUCHOS_ASSERT(solnData!=0); // only needed if blockId is not specified
1606  solnData[0] = dispValues(cell,i)-coordData[axis];
1607  }
1608  }
1609 }
1610 
1611 template <typename ArrayT>
1612 void STK_Interface::getSolutionFieldData(const std::string & fieldName,const std::string & blockId,
1613  const std::vector<std::size_t> & localElementIds,ArrayT & solutionValues) const
1614 {
1615  const std::vector<stk::mesh::Entity> & elements = *(this->getElementsOrderedByLID());
1616 
1617  solutionValues = Kokkos::createDynRankView(solutionValues,
1618  "solutionValues",
1619  localElementIds.size(),
1620  bulkData_->num_nodes(elements[localElementIds[0]]));
1621 
1622  SolutionFieldType * field = this->getSolutionField(fieldName,blockId);
1623 
1624  for(std::size_t cell=0;cell<localElementIds.size();cell++) {
1625  std::size_t localId = localElementIds[cell];
1626  stk::mesh::Entity element = elements[localId];
1627 
1628  // loop over nodes set solution values
1629  const size_t num_nodes = bulkData_->num_nodes(element);
1630  stk::mesh::Entity const* nodes = bulkData_->begin_nodes(element);
1631  for(std::size_t i=0; i<num_nodes; ++i) {
1632  stk::mesh::Entity node = nodes[i];
1633 
1634  double * solnData = stk::mesh::field_data(*field,node);
1635  // TEUCHOS_ASSERT(solnData!=0); // only needed if blockId is not specified
1636  solutionValues(cell,i) = solnData[0];
1637  }
1638  }
1639 }
1640 
1641 template <typename ArrayT>
1642 void STK_Interface::setCellFieldData(const std::string & fieldName,const std::string & blockId,
1643  const std::vector<std::size_t> & localElementIds,const ArrayT & solutionValues,double scaleValue)
1644 {
1645  const std::vector<stk::mesh::Entity> & elements = *(this->getElementsOrderedByLID());
1646 
1647  SolutionFieldType * field = this->getCellField(fieldName,blockId);
1648 
1649  auto solutionValues_h = Kokkos::create_mirror_view(solutionValues);
1650  Kokkos::deep_copy(solutionValues_h, solutionValues);
1651 
1652  for(std::size_t cell=0;cell<localElementIds.size();cell++) {
1653  std::size_t localId = localElementIds[cell];
1654  stk::mesh::Entity element = elements[localId];
1655 
1656  double * solnData = stk::mesh::field_data(*field,element);
1657  TEUCHOS_ASSERT(solnData!=0); // only needed if blockId is not specified
1658  solnData[0] = scaleValue*solutionValues_h.access(cell,0);
1659  }
1660 }
1661 
1662 template <typename ArrayT>
1663 void STK_Interface::setEdgeFieldData(const std::string & fieldName,const std::string & blockId,
1664  const std::vector<std::size_t> & localEdgeIds,const ArrayT & edgeValues,double scaleValue)
1665 {
1666  const std::vector<stk::mesh::Entity> & edges = *(this->getEdgesOrderedByLID());
1667 
1668  SolutionFieldType * field = this->getEdgeField(fieldName,blockId);
1669 
1670  auto edgeValues_h = Kokkos::create_mirror_view(edgeValues);
1671  Kokkos::deep_copy(edgeValues_h, edgeValues);
1672 
1673  for(std::size_t idx=0;idx<localEdgeIds.size();idx++) {
1674  std::size_t localId = localEdgeIds[idx];
1675  stk::mesh::Entity edge = edges[localId];
1676 
1677  double * solnData = stk::mesh::field_data(*field,edge);
1678  TEUCHOS_ASSERT(solnData!=0); // only needed if blockId is not specified
1679  solnData[0] = scaleValue*edgeValues_h.access(idx,0);
1680  }
1681 }
1682 
1683 template <typename ArrayT>
1684 void STK_Interface::setFaceFieldData(const std::string & fieldName,const std::string & blockId,
1685  const std::vector<std::size_t> & localFaceIds,const ArrayT & faceValues,double scaleValue)
1686 {
1687  const std::vector<stk::mesh::Entity> & faces = *(this->getFacesOrderedByLID());
1688 
1689  SolutionFieldType * field = this->getFaceField(fieldName,blockId);
1690 
1691  auto faceValues_h = Kokkos::create_mirror_view(faceValues);
1692  Kokkos::deep_copy(faceValues_h, faceValues);
1693 
1694  for(std::size_t idx=0;idx<localFaceIds.size();idx++) {
1695  std::size_t localId = localFaceIds[idx];
1696  stk::mesh::Entity face = faces[localId];
1697 
1698  double * solnData = stk::mesh::field_data(*field,face);
1699  TEUCHOS_ASSERT(solnData!=0); // only needed if blockId is not specified
1700  solnData[0] = scaleValue*faceValues_h.access(idx,0);
1701  }
1702 }
1703 
1705 
1706 template <typename ArrayT>
1707 void STK_Interface::getElementVertices(const std::vector<std::size_t> & localElementIds, ArrayT & vertices) const
1708 {
1709  if(!useFieldCoordinates_) {
1710  //
1711  // gather from the intrinsic mesh coordinates (non-lagrangian)
1712  //
1713 
1714  const std::vector<stk::mesh::Entity> & elements = *(this->getElementsOrderedByLID());
1715 
1716  // convert to a vector of entity objects
1717  std::vector<stk::mesh::Entity> selected_elements;
1718  for(std::size_t cell=0;cell<localElementIds.size();cell++)
1719  selected_elements.push_back(elements[localElementIds[cell]]);
1720 
1721  getElementVertices_FromCoords(selected_elements,vertices);
1722  }
1723  else {
1724  TEUCHOS_TEST_FOR_EXCEPTION(true,std::invalid_argument,
1725  "STK_Interface::getElementVertices: Cannot call this method when field coordinates are used "
1726  "without specifying an element block.");
1727  }
1728 }
1729 
1730 template <typename ArrayT>
1731 void STK_Interface::getElementVertices(const std::vector<stk::mesh::Entity> & elements, ArrayT & vertices) const
1732 {
1733  if(!useFieldCoordinates_) {
1734  getElementVertices_FromCoords(elements,vertices);
1735  }
1736  else {
1737  TEUCHOS_TEST_FOR_EXCEPTION(true,std::invalid_argument,
1738  "STK_Interface::getElementVertices: Cannot call this method when field coordinates are used "
1739  "without specifying an element block.");
1740  }
1741 }
1742 
1743 template <typename ArrayT>
1744 void STK_Interface::getElementVertices(const std::vector<stk::mesh::Entity> & elements,const std::string & eBlock, ArrayT & vertices) const
1745 {
1746  if(!useFieldCoordinates_) {
1747  getElementVertices_FromCoords(elements,vertices);
1748  }
1749  else {
1750  getElementVertices_FromField(elements,eBlock,vertices);
1751  }
1752 }
1753 
1754 template <typename ArrayT>
1755 void STK_Interface::getElementVertices(const std::vector<std::size_t> & localElementIds,const std::string & eBlock, ArrayT & vertices) const
1756 {
1757  const std::vector<stk::mesh::Entity> & elements = *(this->getElementsOrderedByLID());
1758 
1759  // convert to a vector of entity objects
1760  std::vector<stk::mesh::Entity> selected_elements;
1761  for(std::size_t cell=0;cell<localElementIds.size();cell++)
1762  selected_elements.push_back(elements[localElementIds[cell]]);
1763 
1764  if(!useFieldCoordinates_) {
1765  getElementVertices_FromCoords(selected_elements,vertices);
1766  }
1767  else {
1768  getElementVertices_FromField(selected_elements,eBlock,vertices);
1769  }
1770 }
1771 
1772 template <typename ArrayT>
1773 void STK_Interface::getElementVerticesNoResize(const std::vector<std::size_t> & localElementIds, ArrayT & vertices) const
1774 {
1775  if(!useFieldCoordinates_) {
1776  //
1777  // gather from the intrinsic mesh coordinates (non-lagrangian)
1778  //
1779 
1780  const std::vector<stk::mesh::Entity> & elements = *(this->getElementsOrderedByLID());
1781 
1782  // convert to a vector of entity objects
1783  std::vector<stk::mesh::Entity> selected_elements;
1784  for(std::size_t cell=0;cell<localElementIds.size();cell++)
1785  selected_elements.push_back(elements[localElementIds[cell]]);
1786 
1787  getElementVertices_FromCoordsNoResize(selected_elements,vertices);
1788  }
1789  else {
1790  TEUCHOS_TEST_FOR_EXCEPTION(true,std::invalid_argument,
1791  "STK_Interface::getElementVerticesNoResize: Cannot call this method when field coordinates are used "
1792  "without specifying an element block.");
1793  }
1794 }
1795 
1796 template <typename ArrayT>
1797 void STK_Interface::getElementVerticesNoResize(const std::vector<stk::mesh::Entity> & elements, ArrayT & vertices) const
1798 {
1799  if(!useFieldCoordinates_) {
1800  getElementVertices_FromCoordsNoResize(elements,vertices);
1801  }
1802  else {
1803  TEUCHOS_TEST_FOR_EXCEPTION(true,std::invalid_argument,
1804  "STK_Interface::getElementVerticesNoResize: Cannot call this method when field coordinates are used "
1805  "without specifying an element block.");
1806  }
1807 }
1808 
1809 template <typename ArrayT>
1810 void STK_Interface::getElementVerticesNoResize(const std::vector<stk::mesh::Entity> & elements,const std::string & eBlock, ArrayT & vertices) const
1811 {
1812  if(!useFieldCoordinates_) {
1813  getElementVertices_FromCoordsNoResize(elements,vertices);
1814  }
1815  else {
1816  getElementVertices_FromFieldNoResize(elements,eBlock,vertices);
1817  }
1818 }
1819 
1820 template <typename ArrayT>
1821 void STK_Interface::getElementVerticesNoResize(const std::vector<std::size_t> & localElementIds,const std::string & eBlock, ArrayT & vertices) const
1822 {
1823  const std::vector<stk::mesh::Entity> & elements = *(this->getElementsOrderedByLID());
1824 
1825  // convert to a vector of entity objects
1826  std::vector<stk::mesh::Entity> selected_elements;
1827  for(std::size_t cell=0;cell<localElementIds.size();cell++)
1828  selected_elements.push_back(elements[localElementIds[cell]]);
1829 
1830  if(!useFieldCoordinates_) {
1831  getElementVertices_FromCoordsNoResize(selected_elements,vertices);
1832  }
1833  else {
1834  getElementVertices_FromFieldNoResize(selected_elements,eBlock,vertices);
1835  }
1836 }
1837 
1838 template <typename ArrayT>
1839 void STK_Interface::getElementVertices_FromCoords(const std::vector<stk::mesh::Entity> & elements, ArrayT & vertices) const
1840 {
1841  // nothing to do! silently return
1842  if(elements.size() == 0) {
1843  vertices = Kokkos::createDynRankView(vertices, "vertices", 0, 0, 0);
1844  return;
1845  }
1846 
1847  //
1848  // gather from the intrinsic mesh coordinates (non-lagrangian)
1849  //
1850 
1851  // get *master* cell toplogy...(belongs to first element)
1852  const auto masterVertexCount
1853  = stk::mesh::get_cell_topology(bulkData_->bucket(elements[0]).topology()).getCellTopologyData()->vertex_count;
1854 
1855  // allocate space
1856  vertices = Kokkos::createDynRankView(vertices, "vertices", elements.size(), masterVertexCount,getDimension());
1857  auto vertices_h = Kokkos::create_mirror_view(vertices);
1858  Kokkos::deep_copy(vertices_h, vertices);
1859 
1860  // loop over each requested element
1861  const auto dim = getDimension();
1862  for(std::size_t cell = 0; cell < elements.size(); cell++) {
1863  const auto element = elements[cell];
1864  TEUCHOS_ASSERT(element != 0);
1865 
1866  const auto vertexCount
1867  = stk::mesh::get_cell_topology(bulkData_->bucket(element).topology()).getCellTopologyData()->vertex_count;
1868  TEUCHOS_TEST_FOR_EXCEPTION(vertexCount != masterVertexCount, std::runtime_error,
1869  "In call to STK_Interface::getElementVertices all elements "
1870  "must have the same vertex count!");
1871 
1872  // loop over all element nodes
1873  const size_t num_nodes = bulkData_->num_nodes(element);
1874  auto const* nodes = bulkData_->begin_nodes(element);
1875  TEUCHOS_TEST_FOR_EXCEPTION(num_nodes!=masterVertexCount,std::runtime_error,
1876  "In call to STK_Interface::getElementVertices cardinality of "
1877  "element node relations must be the vertex count!");
1878  for(std::size_t node = 0; node < num_nodes; ++node) {
1879  const double * coord = getNodeCoordinates(nodes[node]);
1880 
1881  // set each dimension of the coordinate
1882  for(unsigned d=0;d<dim;d++)
1883  vertices_h(cell,node,d) = coord[d];
1884  }
1885  }
1886  Kokkos::deep_copy(vertices, vertices_h);
1887 }
1888 
1889 template <typename ArrayT>
1890 void STK_Interface::getElementVertices_FromCoordsNoResize(const std::vector<stk::mesh::Entity> & elements, ArrayT & vertices) const
1891 {
1892  // nothing to do! silently return
1893  if(elements.size()==0) {
1894  return;
1895  }
1896 
1897  //
1898  // gather from the intrinsic mesh coordinates (non-lagrangian)
1899  //
1900 
1901  // get *master* cell toplogy...(belongs to first element)
1902  unsigned masterVertexCount
1903  = stk::mesh::get_cell_topology(bulkData_->bucket(elements[0]).topology()).getCellTopologyData()->vertex_count;
1904 
1905  // loop over each requested element
1906  unsigned dim = getDimension();
1907  auto vertices_h = Kokkos::create_mirror_view(vertices);
1908  for(std::size_t cell=0;cell<elements.size();cell++) {
1909  stk::mesh::Entity element = elements[cell];
1910  TEUCHOS_ASSERT(element!=0);
1911 
1912  unsigned vertexCount
1913  = stk::mesh::get_cell_topology(bulkData_->bucket(element).topology()).getCellTopologyData()->vertex_count;
1914  TEUCHOS_TEST_FOR_EXCEPTION(vertexCount!=masterVertexCount,std::runtime_error,
1915  "In call to STK_Interface::getElementVertices all elements "
1916  "must have the same vertex count!");
1917 
1918  // loop over all element nodes
1919  const size_t num_nodes = bulkData_->num_nodes(element);
1920  stk::mesh::Entity const* nodes = bulkData_->begin_nodes(element);
1921  TEUCHOS_TEST_FOR_EXCEPTION(num_nodes!=masterVertexCount,std::runtime_error,
1922  "In call to STK_Interface::getElementVertices cardinality of "
1923  "element node relations must be the vertex count!");
1924  for(std::size_t node=0; node<num_nodes; ++node) {
1925  const double * coord = getNodeCoordinates(nodes[node]);
1926 
1927  // set each dimension of the coordinate
1928  for(unsigned d=0;d<dim;d++)
1929  vertices_h(cell,node,d) = coord[d];
1930  }
1931  }
1932  Kokkos::deep_copy(vertices, vertices_h);
1933 }
1934 
1935 template <typename ArrayT>
1936 void STK_Interface::getElementVertices_FromField(const std::vector<stk::mesh::Entity> & elements,const std::string & eBlock, ArrayT & vertices) const
1937 {
1939 
1940  // nothing to do! silently return
1941  if(elements.size()==0) {
1942  vertices = Kokkos::createDynRankView(vertices,"vertices",0,0,0);
1943  return;
1944  }
1945 
1946  // get *master* cell toplogy...(belongs to first element)
1947  unsigned masterVertexCount
1948  = stk::mesh::get_cell_topology(bulkData_->bucket(elements[0]).topology()).getCellTopologyData()->vertex_count;
1949 
1950  // allocate space
1951  vertices = Kokkos::createDynRankView(vertices,"vertices",elements.size(),masterVertexCount,getDimension());
1952  auto vertices_h = Kokkos::create_mirror_view(vertices);
1953  std::map<std::string,std::vector<std::string> >::const_iterator itr = meshCoordFields_.find(eBlock);
1954  if(itr==meshCoordFields_.end()) {
1955  // no coordinate field set for this element block
1956  TEUCHOS_ASSERT(false);
1957  }
1958 
1959  const std::vector<std::string> & coordField = itr->second;
1960  std::vector<SolutionFieldType*> fields(getDimension());
1961  for(std::size_t d=0;d<fields.size();d++) {
1962  fields[d] = this->getSolutionField(coordField[d],eBlock);
1963  }
1964 
1965  for(std::size_t cell=0;cell<elements.size();cell++) {
1966  stk::mesh::Entity element = elements[cell];
1967 
1968  // loop over nodes set solution values
1969  const size_t num_nodes = bulkData_->num_nodes(element);
1970  stk::mesh::Entity const* nodes = bulkData_->begin_nodes(element);
1971  for(std::size_t i=0; i<num_nodes; ++i) {
1972  stk::mesh::Entity node = nodes[i];
1973 
1974  const double * coord = getNodeCoordinates(node);
1975 
1976  for(unsigned d=0;d<getDimension();d++) {
1977  double * solnData = stk::mesh::field_data(*fields[d],node);
1978 
1979  // recall mesh field coordinates are stored as displacements
1980  // from the mesh coordinates, make sure to add them together
1981  vertices_h(cell,i,d) = solnData[0]+coord[d];
1982  }
1983  }
1984  }
1985  Kokkos::deep_copy(vertices, vertices_h);
1986 }
1987 
1988 template <typename ArrayT>
1989 void STK_Interface::getElementVertices_FromFieldNoResize(const std::vector<stk::mesh::Entity> & elements,
1990  const std::string & eBlock, ArrayT & vertices) const
1991 {
1993 
1994  // nothing to do! silently return
1995  if(elements.size()==0) {
1996  return;
1997  }
1998 
1999  std::map<std::string,std::vector<std::string> >::const_iterator itr = meshCoordFields_.find(eBlock);
2000  if(itr==meshCoordFields_.end()) {
2001  // no coordinate field set for this element block
2002  TEUCHOS_ASSERT(false);
2003  }
2004 
2005  const std::vector<std::string> & coordField = itr->second;
2006  std::vector<SolutionFieldType*> fields(getDimension());
2007  for(std::size_t d=0;d<fields.size();d++) {
2008  fields[d] = this->getSolutionField(coordField[d],eBlock);
2009  }
2010 
2011  for(std::size_t cell=0;cell<elements.size();cell++) {
2012  stk::mesh::Entity element = elements[cell];
2013 
2014  // loop over nodes set solution values
2015  const size_t num_nodes = bulkData_->num_nodes(element);
2016  stk::mesh::Entity const* nodes = bulkData_->begin_nodes(element);
2017  for(std::size_t i=0; i<num_nodes; ++i) {
2018  stk::mesh::Entity node = nodes[i];
2019 
2020  const double * coord = getNodeCoordinates(node);
2021 
2022  for(unsigned d=0;d<getDimension();d++) {
2023  double * solnData = stk::mesh::field_data(*fields[d],node);
2024 
2025  // recall mesh field coordinates are stored as displacements
2026  // from the mesh coordinates, make sure to add them together
2027  vertices(cell,i,d) = solnData[0]+coord[d];
2028  }
2029  }
2030  }
2031 }
2032 
2034 
2035 template <typename ArrayT>
2036 void STK_Interface::getElementNodes(const std::vector<std::size_t> & localElementIds, ArrayT & nodes) const
2037 {
2038  if(!useFieldCoordinates_) {
2039  //
2040  // gather from the intrinsic mesh coordinates (non-lagrangian)
2041  //
2042 
2043  const std::vector<stk::mesh::Entity> & elements = *(this->getElementsOrderedByLID());
2044 
2045  // convert to a vector of entity objects
2046  std::vector<stk::mesh::Entity> selected_elements;
2047  for(std::size_t cell=0;cell<localElementIds.size();cell++)
2048  selected_elements.push_back(elements[localElementIds[cell]]);
2049 
2050  getElementNodes_FromCoords(selected_elements,nodes);
2051  }
2052  else {
2053  TEUCHOS_TEST_FOR_EXCEPTION(true,std::invalid_argument,
2054  "STK_Interface::getElementNodes: Cannot call this method when field coordinates are used "
2055  "without specifying an element block.");
2056  }
2057 }
2058 
2059 template <typename ArrayT>
2060 void STK_Interface::getElementNodes(const std::vector<stk::mesh::Entity> & elements, ArrayT & nodes) const
2061 {
2062  if(!useFieldCoordinates_) {
2063  getElementNodes_FromCoords(elements,nodes);
2064  }
2065  else {
2066  TEUCHOS_TEST_FOR_EXCEPTION(true,std::invalid_argument,
2067  "STK_Interface::getElementNodes: Cannot call this method when field coordinates are used "
2068  "without specifying an element block.");
2069  }
2070 }
2071 
2072 template <typename ArrayT>
2073 void STK_Interface::getElementNodes(const std::vector<stk::mesh::Entity> & elements,const std::string & eBlock, ArrayT & nodes) const
2074 {
2075  if(!useFieldCoordinates_) {
2076  getElementNodes_FromCoords(elements,nodes);
2077  }
2078  else {
2079  getElementNodes_FromField(elements,eBlock,nodes);
2080  }
2081 }
2082 
2083 template <typename ArrayT>
2084 void STK_Interface::getElementNodes(const std::vector<std::size_t> & localElementIds,const std::string & eBlock, ArrayT & nodes) const
2085 {
2086  const std::vector<stk::mesh::Entity> & elements = *(this->getElementsOrderedByLID());
2087 
2088  // convert to a vector of entity objects
2089  std::vector<stk::mesh::Entity> selected_elements;
2090  for(std::size_t cell=0;cell<localElementIds.size();cell++)
2091  selected_elements.push_back(elements[localElementIds[cell]]);
2092 
2093  if(!useFieldCoordinates_) {
2094  getElementNodes_FromCoords(selected_elements,nodes);
2095  }
2096  else {
2097  getElementNodes_FromField(selected_elements,eBlock,nodes);
2098  }
2099 }
2100 
2101 template <typename ArrayT>
2102 void STK_Interface::getElementNodesNoResize(const std::vector<std::size_t> & localElementIds, ArrayT & nodes) const
2103 {
2104  if(!useFieldCoordinates_) {
2105  //
2106  // gather from the intrinsic mesh coordinates (non-lagrangian)
2107  //
2108 
2109  const std::vector<stk::mesh::Entity> & elements = *(this->getElementsOrderedByLID());
2110 
2111  // convert to a vector of entity objects
2112  std::vector<stk::mesh::Entity> selected_elements;
2113  for(std::size_t cell=0;cell<localElementIds.size();cell++)
2114  selected_elements.push_back(elements[localElementIds[cell]]);
2115 
2116  getElementNodes_FromCoordsNoResize(selected_elements,nodes);
2117  }
2118  else {
2119  TEUCHOS_TEST_FOR_EXCEPTION(true,std::invalid_argument,
2120  "STK_Interface::getElementNodesNoResize: Cannot call this method when field coordinates are used "
2121  "without specifying an element block.");
2122  }
2123 }
2124 
2125 template <typename ArrayT>
2126 void STK_Interface::getElementNodesNoResize(const std::vector<stk::mesh::Entity> & elements, ArrayT & nodes) const
2127 {
2128  if(!useFieldCoordinates_) {
2129  getElementNodes_FromCoordsNoResize(elements,nodes);
2130  }
2131  else {
2132  TEUCHOS_TEST_FOR_EXCEPTION(true,std::invalid_argument,
2133  "STK_Interface::getElementNodesNoResize: Cannot call this method when field coordinates are used "
2134  "without specifying an element block.");
2135  }
2136 }
2137 
2138 template <typename ArrayT>
2139 void STK_Interface::getElementNodesNoResize(const std::vector<stk::mesh::Entity> & elements,const std::string & eBlock, ArrayT & nodes) const
2140 {
2141  if(!useFieldCoordinates_) {
2142  getElementNodes_FromCoordsNoResize(elements,nodes);
2143  }
2144  else {
2145  getElementNodes_FromFieldNoResize(elements,eBlock,nodes);
2146  }
2147 }
2148 
2149 template <typename ArrayT>
2150 void STK_Interface::getElementNodesNoResize(const std::vector<std::size_t> & localElementIds,const std::string & eBlock, ArrayT & nodes) const
2151 {
2152  const std::vector<stk::mesh::Entity> & elements = *(this->getElementsOrderedByLID());
2153 
2154  // convert to a vector of entity objects
2155  std::vector<stk::mesh::Entity> selected_elements;
2156  for(std::size_t cell=0;cell<localElementIds.size();cell++)
2157  selected_elements.push_back(elements[localElementIds[cell]]);
2158 
2159  if(!useFieldCoordinates_) {
2160  getElementNodes_FromCoordsNoResize(selected_elements,nodes);
2161  }
2162  else {
2163  getElementNodes_FromFieldNoResize(selected_elements,eBlock,nodes);
2164  }
2165 }
2166 
2167 template <typename ArrayT>
2168 void STK_Interface::getElementNodes_FromCoords(const std::vector<stk::mesh::Entity> & elements, ArrayT & nodes) const
2169 {
2170  // nothing to do! silently return
2171  if(elements.size() == 0) {
2172  nodes = Kokkos::createDynRankView(nodes, "nodes", 0, 0, 0);
2173  return;
2174  }
2175 
2176  //
2177  // gather from the intrinsic mesh coordinates (non-lagrangian)
2178  //
2179 
2180  // get *master* cell toplogy...(belongs to first element)
2181  const auto masterNodeCount
2182  = stk::mesh::get_cell_topology(bulkData_->bucket(elements[0]).topology()).getCellTopologyData()->node_count;
2183 
2184  // allocate space
2185  nodes = Kokkos::createDynRankView(nodes, "nodes", elements.size(), masterNodeCount,getDimension());
2186  auto nodes_h = Kokkos::create_mirror_view(nodes);
2187  Kokkos::deep_copy(nodes_h, nodes);
2188 
2189  // loop over each requested element
2190  const auto dim = getDimension();
2191  for(std::size_t cell = 0; cell < elements.size(); cell++) {
2192  const auto element = elements[cell];
2193  TEUCHOS_ASSERT(element != 0);
2194 
2195  const auto nodeCount
2196  = stk::mesh::get_cell_topology(bulkData_->bucket(element).topology()).getCellTopologyData()->node_count;
2197  TEUCHOS_TEST_FOR_EXCEPTION(nodeCount != masterNodeCount, std::runtime_error,
2198  "In call to STK_Interface::getElementNodes all elements "
2199  "must have the same node count!");
2200 
2201  // loop over all element nodes
2202  const size_t num_nodes = bulkData_->num_nodes(element);
2203  auto const* elem_nodes = bulkData_->begin_nodes(element);
2204  TEUCHOS_TEST_FOR_EXCEPTION(num_nodes!=masterNodeCount,std::runtime_error,
2205  "In call to STK_Interface::getElementNodes cardinality of "
2206  "element node relations must be the node count!");
2207  for(std::size_t node = 0; node < num_nodes; ++node) {
2208  const double * coord = getNodeCoordinates(elem_nodes[node]);
2209 
2210  // set each dimension of the coordinate
2211  for(unsigned d=0;d<dim;d++)
2212  nodes_h(cell,node,d) = coord[d];
2213  }
2214  }
2215  Kokkos::deep_copy(nodes, nodes_h);
2216 }
2217 
2218 template <typename ArrayT>
2219 void STK_Interface::getElementNodes_FromCoordsNoResize(const std::vector<stk::mesh::Entity> & elements, ArrayT & nodes) const
2220 {
2221  // nothing to do! silently return
2222  if(elements.size()==0) {
2223  return;
2224  }
2225 
2226  //
2227  // gather from the intrinsic mesh coordinates (non-lagrangian)
2228  //
2229 
2230  // get *master* cell toplogy...(belongs to first element)
2231  unsigned masterNodeCount
2232  = stk::mesh::get_cell_topology(bulkData_->bucket(elements[0]).topology()).getCellTopologyData()->node_count;
2233 
2234  // loop over each requested element
2235  unsigned dim = getDimension();
2236  auto nodes_h = Kokkos::create_mirror_view(nodes);
2237  for(std::size_t cell=0;cell<elements.size();cell++) {
2238  stk::mesh::Entity element = elements[cell];
2239  TEUCHOS_ASSERT(element!=0);
2240 
2241  unsigned nodeCount
2242  = stk::mesh::get_cell_topology(bulkData_->bucket(element).topology()).getCellTopologyData()->node_count;
2243  TEUCHOS_TEST_FOR_EXCEPTION(nodeCount!=masterNodeCount,std::runtime_error,
2244  "In call to STK_Interface::getElementNodes all elements "
2245  "must have the same node count!");
2246 
2247  // loop over all element nodes
2248  const size_t num_nodes = bulkData_->num_nodes(element);
2249  stk::mesh::Entity const* elem_nodes = bulkData_->begin_nodes(element);
2250  TEUCHOS_TEST_FOR_EXCEPTION(num_nodes!=masterNodeCount,std::runtime_error,
2251  "In call to STK_Interface::getElementNodes cardinality of "
2252  "element node relations must be the node count!");
2253  for(std::size_t node=0; node<num_nodes; ++node) {
2254  const double * coord = getNodeCoordinates(elem_nodes[node]);
2255 
2256  // set each dimension of the coordinate
2257  for(unsigned d=0;d<dim;d++)
2258  nodes_h(cell,node,d) = coord[d];
2259  }
2260  }
2261  Kokkos::deep_copy(nodes, nodes_h);
2262 }
2263 
2264 template <typename ArrayT>
2265 void STK_Interface::getElementNodes_FromField(const std::vector<stk::mesh::Entity> & elements,const std::string & eBlock, ArrayT & nodes) const
2266 {
2268 
2269  // nothing to do! silently return
2270  if(elements.size()==0) {
2271  nodes = Kokkos::createDynRankView(nodes,"nodes",0,0,0);
2272  return;
2273  }
2274 
2275  // get *master* cell toplogy...(belongs to first element)
2276  unsigned masterNodeCount
2277  = stk::mesh::get_cell_topology(bulkData_->bucket(elements[0]).topology()).getCellTopologyData()->node_count;
2278 
2279  // allocate space
2280  nodes = Kokkos::createDynRankView(nodes,"nodes",elements.size(),masterNodeCount,getDimension());
2281  auto nodes_h = Kokkos::create_mirror_view(nodes);
2282  std::map<std::string,std::vector<std::string> >::const_iterator itr = meshCoordFields_.find(eBlock);
2283  if(itr==meshCoordFields_.end()) {
2284  // no coordinate field set for this element block
2285  TEUCHOS_ASSERT(false);
2286  }
2287 
2288  const std::vector<std::string> & coordField = itr->second;
2289  std::vector<SolutionFieldType*> fields(getDimension());
2290  for(std::size_t d=0;d<fields.size();d++) {
2291  fields[d] = this->getSolutionField(coordField[d],eBlock);
2292  }
2293 
2294  // loop over elements
2295  for(std::size_t cell=0;cell<elements.size();cell++) {
2296  stk::mesh::Entity element = elements[cell];
2297 
2298  // loop over nodes set solution values
2299  const size_t num_nodes = bulkData_->num_nodes(element);
2300  stk::mesh::Entity const* elem_nodes = bulkData_->begin_nodes(element);
2301  for(std::size_t i=0; i<num_nodes; ++i) {
2302  stk::mesh::Entity node = elem_nodes[i];
2303 
2304  const double * coord = getNodeCoordinates(node);
2305 
2306  for(unsigned d=0;d<getDimension();d++) {
2307  double * solnData = stk::mesh::field_data(*fields[d],node);
2308 
2309  // recall mesh field coordinates are stored as displacements
2310  // from the mesh coordinates, make sure to add them together
2311  nodes_h(cell,i,d) = solnData[0]+coord[d];
2312  }
2313  }
2314  }
2315  Kokkos::deep_copy(nodes, nodes_h);
2316 }
2317 
2318 template <typename ArrayT>
2319 void STK_Interface::getElementNodes_FromFieldNoResize(const std::vector<stk::mesh::Entity> & elements,
2320  const std::string & eBlock, ArrayT & nodes) const
2321 {
2323 
2324  // nothing to do! silently return
2325  if(elements.size()==0) {
2326  return;
2327  }
2328 
2329  std::map<std::string,std::vector<std::string> >::const_iterator itr = meshCoordFields_.find(eBlock);
2330  if(itr==meshCoordFields_.end()) {
2331  // no coordinate field set for this element block
2332  TEUCHOS_ASSERT(false);
2333  }
2334 
2335  const std::vector<std::string> & coordField = itr->second;
2336  std::vector<SolutionFieldType*> fields(getDimension());
2337  for(std::size_t d=0;d<fields.size();d++) {
2338  fields[d] = this->getSolutionField(coordField[d],eBlock);
2339  }
2340 
2341  // loop over elements
2342  for(std::size_t cell=0;cell<elements.size();cell++) {
2343  stk::mesh::Entity element = elements[cell];
2344 
2345  // loop over nodes set solution values
2346  const size_t num_nodes = bulkData_->num_nodes(element);
2347  stk::mesh::Entity const* elem_nodes = bulkData_->begin_nodes(element);
2348  for(std::size_t i=0; i<num_nodes; ++i) {
2349  stk::mesh::Entity node = elem_nodes[i];
2350 
2351  const double * coord = getNodeCoordinates(node);
2352 
2353  for(unsigned d=0;d<getDimension();d++) {
2354  double * solnData = stk::mesh::field_data(*fields[d],node);
2355 
2356  // recall mesh field coordinates are stored as displacements
2357  // from the mesh coordinates, make sure to add them together
2358  nodes(cell,i,d) = solnData[0]+coord[d];
2359  }
2360  }
2361  }
2362 }
2363 
2364 }
2365 
2366 #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
const bool & useBoundingBoxSearch() const
Teuchos::RCP< std::vector< stk::mesh::Entity > > orderedEdgeVector_
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 getFaceBlockNames(std::vector< std::string > &names) const
void print(std::ostream &os) const
void getElementNodesNoResize(const std::vector< std::size_t > &localIds, ArrayT &nodes) const
unsigned entityOwnerRank(stk::mesh::Entity entity) const
void getMyEdges(std::vector< stk::mesh::Entity > &edges) const
void setEdgeFieldData(const std::string &fieldName, const std::string &blockId, const std::vector< std::size_t > &localEdgeIds, const ArrayT &edgeValues, double scaleValue=1.0)
stk::mesh::Field< double > VectorFieldType
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 getElementNodes_FromCoordsNoResize(const std::vector< stk::mesh::Entity > &elements, ArrayT &nodes) const
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 part
std::vector< stk::mesh::Part * > edgesPartVec_
bool nonnull(const std::shared_ptr< T > &p)
bool isFaceLocal(stk::mesh::Entity face) const
std::map< std::string, double > blockWeights_
Teuchos::RCP< const std::vector< stk::mesh::Entity > > getEdgesOrderedByLID() const
#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 getElementNodes_FromFieldNoResize(const std::vector< stk::mesh::Entity > &elements, const std::string &eBlock, ArrayT &nodes) const
void getElementVertices(const std::vector< std::size_t > &localIds, ArrayT &vertices) const
void addEntityToSideset(stk::mesh::Entity entity, stk::mesh::Part *sideset)
std::map< std::string, stk::mesh::Part * > faceBlocks_
stk::mesh::EntityId edgeGlobalId(std::size_t lid) const
void getElementVertices_FromCoords(const std::vector< stk::mesh::Entity > &elements, ArrayT &vertices) const
void getElementNodes(const std::vector< std::size_t > &localIds, ArrayT &nodes) const
void addSolutionField(const std::string &fieldName, const std::string &blockId)
std::map< std::string, Teuchos::RCP< shards::CellTopology > > elementBlockCT_
void addEntityToFaceBlock(stk::mesh::Entity entity, stk::mesh::Part *faceblock)
stk::mesh::EntityId getMaxEntityId(unsigned entityRank) const
get max entity ID of type entityRank
std::map< std::pair< std::string, std::string >, SolutionFieldType * > fieldNameToFaceField_
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
void getEdgeBlockNames(std::vector< std::string > &names) const
stk::mesh::Part * getFaceBlock(const std::string &name) const
get the block part
void addInformationRecords(const std::vector< std::string > &info_records)
stk::mesh::EntityId faceGlobalId(std::size_t lid) const
void getAllFaces(const std::string &faceBlockName, std::vector< stk::mesh::Entity > &faces) const
stk::mesh::Field< double > * getSolutionField(const std::string &fieldName, const std::string &blockId) const
void addEdgeBlock(const std::string &elemBlockName, const std::string &edgeBlockName, const stk::topology &topology)
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
std::set< std::string > informationRecords_
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
stk::mesh::Field< double > * getEdgeField(const std::string &fieldName, const std::string &blockId) const
std::size_t elementLocalId(stk::mesh::Entity elmt) const
stk::mesh::EntityId faceGlobalId(stk::mesh::Entity face) 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 setupExodusFile(const std::string &filename, const bool append=false, const bool append_after_restart_time=false, const double restart_time=0.0)
Set up an output Exodus file for writing results.
void setUseFieldCoordinates(bool useFieldCoordinates)
stk::mesh::EntityRank getSideRank() const
void getAllEdges(const std::string &edgeBlockName, std::vector< stk::mesh::Entity > &edges) const
void initialize(stk::ParallelMachine parallelMach, bool setupIO=true, const bool buildRefinementSupport=false)
static const std::string coordsString
void setFaceFieldData(const std::string &fieldName, const std::string &blockId, const std::vector< std::size_t > &localFaceIds, const ArrayT &faceValues, double scaleValue=1.0)
std::size_t getNumElementBlocks() const
get the block count
void getElementNodes_FromField(const std::vector< stk::mesh::Entity > &elements, const std::string &eBlock, ArrayT &nodes) const
Teuchos::RCP< stk::mesh::BulkData > getBulkData() const
std::map< std::string, stk::mesh::Part * > edgeBlocks_
void getElementVertices_FromField(const std::vector< stk::mesh::Entity > &elements, const std::string &eBlock, ArrayT &vertices) const
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 edgeBlockString
static const std::string edgesString
std::string containingBlockId(stk::mesh::Entity elmt) const
stk::mesh::EntityRank getFaceRank() const
std::size_t edgeLocalId(stk::mesh::Entity elmt) 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 addEntitiesToEdgeBlock(std::vector< stk::mesh::Entity > entities, stk::mesh::Part *edgeblock)
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 writeToExodus(const std::string &filename, const bool append=false)
Write this mesh and associated fields to the given output file.
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 addEdgeField(const std::string &fieldName, const std::string &blockId)
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)
static const std::string faceBlockString
std::unordered_map< stk::mesh::EntityId, std::size_t > localFaceIDHash_
void addEntityToEdgeBlock(stk::mesh::Entity entity, stk::mesh::Part *edgeblock)
std::map< std::pair< std::string, std::string >, SolutionFieldType * > fieldNameToEdgeField_
stk::mesh::Part * getEdgeBlock(const std::string &name) const
get the block part
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
std::unordered_map< stk::mesh::EntityId, std::size_t > localEdgeIDHash_
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 addEntitiesToFaceBlock(std::vector< stk::mesh::Entity > entities, stk::mesh::Part *faceblock)
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)
stk::mesh::Field< double > * getFaceField(const std::string &fieldName, const std::string &blockId) const
std::size_t faceLocalId(stk::mesh::Entity elmt) const
void addFaceBlock(const std::string &elemBlockName, const std::string &faceBlockName, const stk::topology &topology)
Teuchos::RCP< std::vector< stk::mesh::Entity > > orderedFaceVector_
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
void addFaceField(const std::string &fieldName, const std::string &blockId)
#define TEUCHOS_ASSERT(assertion_test)
void getMyElements(std::vector< stk::mesh::Entity > &elements) const
stk::mesh::EntityRank getNodeRank() const
void getMyFaces(std::vector< stk::mesh::Entity > &faces) 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
Teuchos::RCP< const std::vector< stk::mesh::Entity > > getFacesOrderedByLID() const
bool operator()(stk::mesh::Entity a, stk::mesh::Entity b)
stk::mesh::EntityRank getElementRank() const
void getElementNodes_FromCoords(const std::vector< stk::mesh::Entity > &elements, ArrayT &nodes) const
stk::mesh::EntityId edgeGlobalId(stk::mesh::Entity edge) 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 setBoundingBoxSearchFlag(const bool &searchFlag)
bool isEdgeLocal(stk::mesh::Entity edge) 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_