Panzer Version of the Day
Loading...
Searching...
No Matches
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>
15#include <Teuchos_DefaultMpiComm.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
38namespace percept {
39 class PerceptMesh;
40 class URP_Heterogeneous_3D;
41}
42#endif
43
44namespace panzer_stk {
45
46class PeriodicBC_MatcherBase;
47
52public:
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_; }
58protected:
59 stk::mesh::EntityId gid_;
60 std::vector<stk::mesh::EntityId> nodes_;
61
63};
64
67Teuchos::RCP<ElementDescriptor>
68buildElementDescriptor(stk::mesh::EntityId elmtId,std::vector<stk::mesh::EntityId> & nodes);
69
71public:
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
91
94 STK_Interface(unsigned dim);
95
96 STK_Interface(Teuchos::RCP<stk::mesh::MetaData> metaData);
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
208 void endModification(const bool find_and_set_shared_nodes_in_stk=true);
209
215 void addNode(stk::mesh::EntityId gid, const std::vector<double> & coord);
216
217 void addElement(const Teuchos::RCP<ElementDescriptor> & ed,stk::mesh::Part * block);
218
219 void addEdges();
220
221 void addFaces();
222
225 void addEntityToSideset(stk::mesh::Entity entity,stk::mesh::Part * sideset);
226
229 void addEntityToNodeset(stk::mesh::Entity entity,stk::mesh::Part * nodeset);
230
233 void addEntityToEdgeBlock(stk::mesh::Entity entity,stk::mesh::Part * edgeblock);
236 void addEntitiesToEdgeBlock(std::vector<stk::mesh::Entity> entities,stk::mesh::Part * edgeblock);
237
240 void addEntityToFaceBlock(stk::mesh::Entity entity,stk::mesh::Part * faceblock);
243 void addEntitiesToFaceBlock(std::vector<stk::mesh::Entity> entities,stk::mesh::Part * faceblock);
244
245 // Methods to interrogate the mesh topology and structure
247
251 { return *coordinatesField_; }
252
256 { return *edgesField_; }
257
259 { return *facesField_; }
260
263 const double * getNodeCoordinates(stk::mesh::EntityId nodeId) const;
264
267 const double * getNodeCoordinates(stk::mesh::Entity node) const;
268
271 void getSubcellIndices(unsigned entityRank,stk::mesh::EntityId elementId,
272 std::vector<stk::mesh::EntityId> & subcellIds) const;
273
276 void getMyElements(std::vector<stk::mesh::Entity> & elements) const;
277
280 void getMyElements(const std::string & blockID,std::vector<stk::mesh::Entity> & elements) const;
281
285 void getNeighborElements(std::vector<stk::mesh::Entity> & elements) const;
286
289 void getNeighborElements(const std::string & blockID,std::vector<stk::mesh::Entity> & elements) const;
290
293 void getMyEdges(std::vector<stk::mesh::Entity> & edges) const;
294
302 void getMyEdges(const std::string & edgeBlockName,std::vector<stk::mesh::Entity> & edges) const;
303
312 void getMyEdges(const std::string & edgeBlockName,const std::string & blockName,std::vector<stk::mesh::Entity> & edges) const;
313
321 void getAllEdges(const std::string & edgeBlockName,std::vector<stk::mesh::Entity> & edges) const;
322
331 void getAllEdges(const std::string & edgeBlockName,const std::string & blockName,std::vector<stk::mesh::Entity> & edges) const;
332
335 void getMyFaces(std::vector<stk::mesh::Entity> & faces) const;
336
344 void getMyFaces(const std::string & faceBlockName,std::vector<stk::mesh::Entity> & faces) const;
345
354 void getMyFaces(const std::string & faceBlockName,const std::string & blockName,std::vector<stk::mesh::Entity> & faces) const;
355
363 void getAllFaces(const std::string & faceBlockName,std::vector<stk::mesh::Entity> & faces) const;
364
373 void getAllFaces(const std::string & faceBlockName,const std::string & blockName,std::vector<stk::mesh::Entity> & faces) const;
374
382 void getMySides(const std::string & sideName,std::vector<stk::mesh::Entity> & sides) const;
383
392 void getMySides(const std::string & sideName,const std::string & blockName,std::vector<stk::mesh::Entity> & sides) const;
393
401 void getAllSides(const std::string & sideName,std::vector<stk::mesh::Entity> & sides) const;
402
411 void getAllSides(const std::string & sideName,const std::string & blockName,std::vector<stk::mesh::Entity> & sides) const;
412
422 void getMyNodes(const std::string & sideName,const std::string & blockName,std::vector<stk::mesh::Entity> & nodes) const;
423
432 stk::mesh::Entity findConnectivityById(stk::mesh::Entity src, stk::mesh::EntityRank tgt_rank, unsigned rel_id) const;
433
434 // Utility functions
436
445 void
446 writeToExodus(const std::string& filename,
447 const bool append = false);
448
466 void
467 setupExodusFile(const std::string& filename,
468 const bool append = false,
469 const bool append_after_restart_time = false,
470 const double restart_time = 0.0);
471
472 void
473 setupExodusFile(const std::string& filename,
474 const std::vector<Ioss::Property>& ioss_properties,
475 const bool append = false,
476 const bool append_after_restart_time = false,
477 const double restart_time = 0.0);
478
491 void
493 double timestep);
494
514 void
516 const std::string& key,
517 const int& value);
518
538 void
540 const std::string& key,
541 const double& value);
542
562 void
564 const std::string& key,
565 const std::vector<int>& value);
566
586 void
588 const std::string& key,
589 const std::vector<double>& value);
590
591 // Accessor functions
593
595 Teuchos::RCP<const Teuchos::Comm<int> > getComm() const;
596
597 Teuchos::RCP<stk::mesh::BulkData> getBulkData() const { return bulkData_; }
598 Teuchos::RCP<stk::mesh::MetaData> getMetaData() const { return metaData_; }
599
600#ifdef PANZER_HAVE_PERCEPT
602 Teuchos::RCP<percept::PerceptMesh> getRefinedMesh() const
603 { TEUCHOS_ASSERT(Teuchos::nonnull(refinedMesh_)); return refinedMesh_; }
604#endif
605
606 bool isWritable() const;
607
608 bool isModifiable() const
609 { if(bulkData_==Teuchos::null) return false;
610 return bulkData_->in_modifiable_state(); }
611
613 unsigned getDimension() const
614 { return dimension_; }
615
617 std::size_t getNumElementBlocks() const
618 { return elementBlocks_.size(); }
619
627 void getElementBlockNames(std::vector<std::string> & names) const;
628
636 void getSidesetNames(std::vector<std::string> & name) const;
637
645 void getNodesetNames(std::vector<std::string> & name) const;
646
654 void getEdgeBlockNames(std::vector<std::string> & names) const;
655
663 void getFaceBlockNames(std::vector<std::string> & names) const;
664
666 stk::mesh::Part * getOwnedPart() const
667 { return &getMetaData()->locally_owned_part(); } // I don't like the pointer access here, but it will do for now!
668
670 stk::mesh::Part * getElementBlockPart(const std::string & name) const
671 {
672 std::map<std::string, stk::mesh::Part*>::const_iterator itr = elementBlocks_.find(name); // Element blocks
673 if(itr==elementBlocks_.end()) return 0;
674 return itr->second;
675 }
676
678 stk::mesh::Part * getEdgeBlock(const std::string & name) const
679 {
680 std::map<std::string, stk::mesh::Part*>::const_iterator itr = edgeBlocks_.find(name); // edge blocks
681 if(itr==edgeBlocks_.end()) return 0;
682 return itr->second;
683 }
684
686 stk::mesh::Part * getFaceBlock(const std::string & name) const
687 {
688 std::map<std::string, stk::mesh::Part*>::const_iterator itr = faceBlocks_.find(name); // face blocks
689 if(itr==faceBlocks_.end()) return 0;
690 return itr->second;
691 }
692
694 std::size_t getNumSidesets() const
695 { return sidesets_.size(); }
696
697 stk::mesh::Part * getSideset(const std::string & name) const
698 {
699 auto itr = sidesets_.find(name);
700 return (itr != sidesets_.end()) ? itr->second : nullptr;
701 }
702
704 std::size_t getNumNodesets() const
705 { return nodesets_.size(); }
706
707 stk::mesh::Part * getNodeset(const std::string & name) const
708 {
709 auto itr = nodesets_.find(name);
710 return (itr != nodesets_.end()) ? itr->second : nullptr;
711 }
712
714 std::size_t getEntityCounts(unsigned entityRank) const;
715
717 stk::mesh::EntityId getMaxEntityId(unsigned entityRank) const;
718
719 // Utilities
721
723 void getElementsSharingNode(stk::mesh::EntityId nodeId,std::vector<stk::mesh::Entity> & elements) const;
724
726 void getNodeIdsForElement(stk::mesh::Entity element,std::vector<stk::mesh::EntityId> & nodeIds) const;
727
730 void getOwnedElementsSharingNode(stk::mesh::Entity node,std::vector<stk::mesh::Entity> & elements,
731 std::vector<int> & relIds) const;
732
735 void getOwnedElementsSharingNode(stk::mesh::EntityId nodeId,std::vector<stk::mesh::Entity> & elements,
736 std::vector<int> & relIds, unsigned int matchType) const;
737
738
740 void getElementsSharingNodes(const std::vector<stk::mesh::EntityId> nodeId,std::vector<stk::mesh::Entity> & elements) const;
741
743 void buildSubcells();
744
747 std::size_t elementLocalId(stk::mesh::Entity elmt) const;
748
751 std::size_t elementLocalId(stk::mesh::EntityId gid) const;
752
755 inline stk::mesh::EntityId elementGlobalId(std::size_t lid) const
756 { return bulkData_->identifier((*orderedElementVector_)[lid]); }
757
760 inline stk::mesh::EntityId elementGlobalId(stk::mesh::Entity elmt) const
761 { return bulkData_->identifier(elmt); }
762
765 bool isEdgeLocal(stk::mesh::Entity edge) const;
766
769 bool isEdgeLocal(stk::mesh::EntityId gid) const;
770
773 std::size_t edgeLocalId(stk::mesh::Entity elmt) const;
774
777 std::size_t edgeLocalId(stk::mesh::EntityId gid) const;
778
781 inline stk::mesh::EntityId edgeGlobalId(std::size_t lid) const
782 { return bulkData_->identifier((*orderedEdgeVector_)[lid]); }
783
786 inline stk::mesh::EntityId edgeGlobalId(stk::mesh::Entity edge) const
787 { return bulkData_->identifier(edge); }
788
791 bool isFaceLocal(stk::mesh::Entity face) const;
792
795 bool isFaceLocal(stk::mesh::EntityId gid) const;
796
799 std::size_t faceLocalId(stk::mesh::Entity elmt) const;
800
803 std::size_t faceLocalId(stk::mesh::EntityId gid) const;
804
807 inline stk::mesh::EntityId faceGlobalId(std::size_t lid) const
808 { return bulkData_->identifier((*orderedFaceVector_)[lid]); }
809
812 inline stk::mesh::EntityId faceGlobalId(stk::mesh::Entity face) const
813 { return bulkData_->identifier(face); }
814
817 inline unsigned entityOwnerRank(stk::mesh::Entity entity) const
818 { return bulkData_->parallel_owner_rank(entity); }
819
822 inline bool isValid(stk::mesh::Entity entity) const
823 { return bulkData_->is_valid(entity); }
824
827 std::string containingBlockId(stk::mesh::Entity elmt) const;
828
833 stk::mesh::Field<double> * getSolutionField(const std::string & fieldName,
834 const std::string & blockId) const;
835
840 stk::mesh::Field<double> * getCellField(const std::string & fieldName,
841 const std::string & blockId) const;
842
847 stk::mesh::Field<double> * getEdgeField(const std::string & fieldName,
848 const std::string & blockId) const;
849
854 stk::mesh::Field<double> * getFaceField(const std::string & fieldName,
855 const std::string & blockId) const;
856
858
860 bool isInitialized() const { return initialized_; }
861
865 Teuchos::RCP<const std::vector<stk::mesh::Entity> > getElementsOrderedByLID() const;
866
881 template <typename ArrayT>
882 void setSolutionFieldData(const std::string & fieldName,const std::string & blockId,
883 const std::vector<std::size_t> & localElementIds,const ArrayT & solutionValues,double scaleValue=1.0);
884
899 template <typename ArrayT>
900 void getSolutionFieldData(const std::string & fieldName,const std::string & blockId,
901 const std::vector<std::size_t> & localElementIds,ArrayT & solutionValues) const;
902
917 template <typename ArrayT>
918 void setCellFieldData(const std::string & fieldName,const std::string & blockId,
919 const std::vector<std::size_t> & localElementIds,const ArrayT & solutionValues,double scaleValue=1.0);
920
924 Teuchos::RCP<const std::vector<stk::mesh::Entity> > getEdgesOrderedByLID() const;
925
929 Teuchos::RCP<const std::vector<stk::mesh::Entity> > getFacesOrderedByLID() const;
930
945 template <typename ArrayT>
946 void setEdgeFieldData(const std::string & fieldName,const std::string & blockId,
947 const std::vector<std::size_t> & localEdgeIds,const ArrayT & edgeValues,double scaleValue=1.0);
948
963 template <typename ArrayT>
964 void setFaceFieldData(const std::string & fieldName,const std::string & blockId,
965 const std::vector<std::size_t> & localFaceIds,const ArrayT & faceValues,double scaleValue=1.0);
966
968
977 template <typename ArrayT>
978 void getElementVertices(const std::vector<std::size_t> & localIds, ArrayT & vertices) const;
979
988 template <typename ArrayT>
989 void getElementVertices(const std::vector<stk::mesh::Entity> & elements, ArrayT & vertices) const;
990
1000 template <typename ArrayT>
1001 void getElementVertices(const std::vector<std::size_t> & localIds,const std::string & eBlock, ArrayT & vertices) const;
1002
1012 template <typename ArrayT>
1013 void getElementVertices(const std::vector<stk::mesh::Entity> & elements,const std::string & eBlock, ArrayT & vertices) const;
1014
1023 template <typename ArrayT>
1024 void getElementVerticesNoResize(const std::vector<std::size_t> & localIds, ArrayT & vertices) const;
1025
1034 template <typename ArrayT>
1035 void getElementVerticesNoResize(const std::vector<stk::mesh::Entity> & elements, ArrayT & vertices) const;
1036
1046 template <typename ArrayT>
1047 void getElementVerticesNoResize(const std::vector<std::size_t> & localIds,const std::string & eBlock, ArrayT & vertices) const;
1048
1058 template <typename ArrayT>
1059 void getElementVerticesNoResize(const std::vector<stk::mesh::Entity> & elements,const std::string & eBlock, ArrayT & vertices) const;
1060
1062
1072 template <typename ArrayT>
1073 void getElementNodes(const std::vector<std::size_t> & localIds, ArrayT & nodes) const;
1074
1083 template <typename ArrayT>
1084 void getElementNodes(const std::vector<stk::mesh::Entity> & elements, ArrayT & nodes) const;
1085
1095 template <typename ArrayT>
1096 void getElementNodes(const std::vector<std::size_t> & localIds,const std::string & eBlock, ArrayT & nodes) const;
1097
1107 template <typename ArrayT>
1108 void getElementNodes(const std::vector<stk::mesh::Entity> & elements,const std::string & eBlock, ArrayT & nodes) const;
1109
1118 template <typename ArrayT>
1119 void getElementNodesNoResize(const std::vector<std::size_t> & localIds, ArrayT & nodes) const;
1120
1129 template <typename ArrayT>
1130 void getElementNodesNoResize(const std::vector<stk::mesh::Entity> & elements, ArrayT & nodes) const;
1131
1141 template <typename ArrayT>
1142 void getElementNodesNoResize(const std::vector<std::size_t> & localIds,const std::string & eBlock, ArrayT & nodes) const;
1143
1153 template <typename ArrayT>
1154 void getElementNodesNoResize(const std::vector<stk::mesh::Entity> & elements,const std::string & eBlock, ArrayT & nodes) const;
1155
1156 // const stk::mesh::FEMInterface & getFEMInterface() const
1157 // { return *femPtr_; }
1158
1159 stk::mesh::EntityRank getElementRank() const { return stk::topology::ELEMENT_RANK; }
1160 stk::mesh::EntityRank getSideRank() const { return metaData_->side_rank(); }
1161 stk::mesh::EntityRank getFaceRank() const { return stk::topology::FACE_RANK; }
1162 stk::mesh::EntityRank getEdgeRank() const { return stk::topology::EDGE_RANK; }
1163 stk::mesh::EntityRank getNodeRank() const { return stk::topology::NODE_RANK; }
1164
1168
1171 void buildLocalElementIDs();
1172
1175 void buildLocalEdgeIDs();
1176
1179 void buildLocalFaceIDs();
1180
1183 const std::vector<Teuchos::RCP<const PeriodicBC_MatcherBase> > &
1185 { return periodicBCs_; }
1186
1189 std::vector<Teuchos::RCP<const PeriodicBC_MatcherBase> > &
1192
1195 const bool & useBoundingBoxSearch() const
1196 { return useBBoxSearch_; }
1197
1199 void setBoundingBoxSearchFlag(const bool & searchFlag)
1200 { useBBoxSearch_ = searchFlag; return; }
1201
1207 void addPeriodicBC(const Teuchos::RCP<const PeriodicBC_MatcherBase> & bc)
1208 { periodicBCs_.push_back(bc); }
1209
1215 void addPeriodicBCs(const std::vector<Teuchos::RCP<const PeriodicBC_MatcherBase> > & bc_vec)
1216 { periodicBCs_.insert(periodicBCs_.end(),bc_vec.begin(),bc_vec.end()); }
1217
1220 std::pair<Teuchos::RCP<std::vector<std::pair<std::size_t,std::size_t> > >, Teuchos::RCP<std::vector<unsigned int> > >
1221 getPeriodicNodePairing() const;
1222
1225 bool validBlockId(const std::string & blockId) const;
1226
1229 void print(std::ostream & os) const;
1230
1233 void printMetaData(std::ostream & os) const;
1234
1237 Teuchos::RCP<const shards::CellTopology> getCellTopology(const std::string & eBlock) const;
1238
1243 double getCurrentStateTime() const { return currentStateTime_; }
1244
1250 double getInitialStateTime() const { return initialStateTime_; }
1251
1256 void setInitialStateTime(double value) { initialStateTime_ = value; }
1257
1260 void rebalance(const Teuchos::ParameterList & params);
1261
1265 void setBlockWeight(const std::string & blockId,double weight)
1266 { blockWeights_[blockId] = weight; }
1267
1274 void setUseFieldCoordinates(bool useFieldCoordinates)
1275 { useFieldCoordinates_ = useFieldCoordinates; }
1276
1279 { return useFieldCoordinates_; }
1280
1282 void setUseLowerCaseForIO(bool useLowerCase)
1283 { useLowerCase_ = useLowerCase; }
1284
1287 { return useLowerCase_; }
1288
1290
1301 template <typename ArrayT>
1302 void getElementVertices_FromField(const std::vector<stk::mesh::Entity> & elements,const std::string & eBlock, ArrayT & vertices) const;
1303
1304 template <typename ArrayT>
1305 void getElementVertices_FromFieldNoResize(const std::vector<stk::mesh::Entity> & elements,
1306 const std::string & eBlock, ArrayT & vertices) const;
1307
1317 template <typename ArrayT>
1318 void getElementVertices_FromCoords(const std::vector<stk::mesh::Entity> & elements, ArrayT & vertices) const;
1319
1320 template <typename ArrayT>
1321 void getElementVertices_FromCoordsNoResize(const std::vector<stk::mesh::Entity> & elements, ArrayT & vertices) const;
1322
1324
1335 template <typename ArrayT>
1336 void getElementNodes_FromField(const std::vector<stk::mesh::Entity> & elements,const std::string & eBlock, ArrayT & nodes) const;
1337
1338 template <typename ArrayT>
1339 void getElementNodes_FromFieldNoResize(const std::vector<stk::mesh::Entity> & elements,
1340 const std::string & eBlock, ArrayT & nodes) const;
1341
1351 template <typename ArrayT>
1352 void getElementNodes_FromCoords(const std::vector<stk::mesh::Entity> & elements, ArrayT & nodes) const;
1353
1354 template <typename ArrayT>
1355 void getElementNodes_FromCoordsNoResize(const std::vector<stk::mesh::Entity> & elements, ArrayT & nodes) const;
1356
1362 void refineMesh(const int numberOfLevels, const bool deleteParentElements);
1363
1364public: // static operations
1365 static const std::string coordsString;
1366 static const std::string nodesString;
1367 static const std::string edgesString;
1368 static const std::string edgeBlockString;
1369 static const std::string faceBlockString;
1370 static const std::string facesString;
1371
1372protected:
1373
1376 void buildEntityCounts();
1377
1380 void buildMaxEntityIds();
1381
1385 void initializeFieldsInSTK(const std::map<std::pair<std::string,std::string>,SolutionFieldType*> & nameToField,
1386 bool setupIO);
1387
1392 Teuchos::RCP<Teuchos::MpiComm<int> > getSafeCommunicator(stk::ParallelMachine parallelMach) const;
1393
1400
1404 bool isMeshCoordField(const std::string & eBlock,const std::string & fieldName,int & axis) const;
1405
1421 template <typename ArrayT>
1422 void setDispFieldData(const std::string & fieldName,const std::string & blockId,int axis,
1423 const std::vector<std::size_t> & localElementIds,const ArrayT & solutionValues);
1424
1425 std::vector<Teuchos::RCP<const PeriodicBC_MatcherBase> > periodicBCs_;
1426 bool useBBoxSearch_ = false; // TODO swap this to change default periodic BC search (see also PeriodicBC_Parser.cpp)
1427
1428 Teuchos::RCP<stk::mesh::MetaData> metaData_;
1429 Teuchos::RCP<stk::mesh::BulkData> bulkData_;
1430#ifdef PANZER_HAVE_PERCEPT
1431 Teuchos::RCP<percept::PerceptMesh> refinedMesh_;
1432 Teuchos::RCP<percept::URP_Heterogeneous_3D> breakPattern_;
1433#endif
1434
1435 std::map<std::string, stk::mesh::Part*> elementBlocks_; // Element blocks
1436 std::map<std::string, stk::mesh::Part*> sidesets_; // Side sets
1437 std::map<std::string, stk::mesh::Part*> nodesets_; // Node sets
1438 std::map<std::string, stk::mesh::Part*> edgeBlocks_; // Edge blocks
1439 std::map<std::string, stk::mesh::Part*> faceBlocks_; // Face blocks
1440
1441 std::map<std::string, Teuchos::RCP<shards::CellTopology> > elementBlockCT_;
1442
1443 // for storing/accessing nodes
1444 stk::mesh::Part * nodesPart_;
1445 std::vector<stk::mesh::Part*> nodesPartVec_;
1446 stk::mesh::Part * edgesPart_;
1447 std::vector<stk::mesh::Part*> edgesPartVec_;
1448 stk::mesh::Part * facesPart_;
1449 std::vector<stk::mesh::Part*> facesPartVec_;
1450
1456
1457 // maps field names to solution field stk mesh handles
1458 std::map<std::pair<std::string,std::string>,SolutionFieldType*> fieldNameToSolution_;
1459 std::map<std::pair<std::string,std::string>,SolutionFieldType*> fieldNameToCellField_;
1460 std::map<std::pair<std::string,std::string>,SolutionFieldType*> fieldNameToEdgeField_;
1461 std::map<std::pair<std::string,std::string>,SolutionFieldType*> fieldNameToFaceField_;
1462
1463 // use a set to maintain a list of unique information records
1464 std::set<std::string> informationRecords_;
1465
1466 unsigned dimension_;
1467
1469
1470 // how many elements, faces, edges, and nodes are there globally
1471 std::vector<std::size_t> entityCounts_;
1472
1473 // what is maximum entity ID
1474 std::vector<stk::mesh::EntityId> maxEntityId_;
1475
1476 unsigned procRank_;
1477 std::size_t currentLocalId_;
1478
1479 Teuchos::RCP<Teuchos::MpiComm<int> > mpiComm_;
1480
1481 double initialStateTime_; // the time stamp at the time this object was constructed (default 0.0)
1482 double currentStateTime_; // the time stamp set by the user when writeToExodus is called (default 0.0)
1483
1484#ifdef PANZER_HAVE_IOSS
1485 // I/O support
1486 Teuchos::RCP<stk::io::StkMeshIoBroker> meshData_;
1487 int meshIndex_;
1488
1493 enum class GlobalVariable
1494 {
1495 ADD,
1496 WRITE
1497 }; // end of enum class GlobalVariable
1498
1515 void
1516 globalToExodus(
1517 const GlobalVariable& flag);
1518
1522 Teuchos::ParameterList globalData_;
1523#endif
1524
1525 // uses lazy evaluation
1526 mutable Teuchos::RCP<std::vector<stk::mesh::Entity> > orderedElementVector_;
1527
1528 // uses lazy evaluation
1529 mutable Teuchos::RCP<std::vector<stk::mesh::Entity> > orderedEdgeVector_;
1530
1531 // uses lazy evaluation
1532 mutable Teuchos::RCP<std::vector<stk::mesh::Entity> > orderedFaceVector_;
1533
1534 // for element block weights
1535 std::map<std::string,double> blockWeights_;
1536
1537 std::unordered_map<stk::mesh::EntityId,std::size_t> localIDHash_;
1538 std::unordered_map<stk::mesh::EntityId,std::size_t> localEdgeIDHash_;
1539 std::unordered_map<stk::mesh::EntityId,std::size_t> localFaceIDHash_;
1540
1541 // Store mesh displacement fields by element block. This map
1542 // goes like this meshCoordFields_[eBlock][axis_index] => coordinate FieldName
1543 // goes like this meshDispFields_[eBlock][axis_index] => displacement FieldName
1544 std::map<std::string,std::vector<std::string> > meshCoordFields_; // coordinate fields written by user
1545 std::map<std::string,std::vector<std::string> > meshDispFields_; // displacement fields, output to exodus
1546
1548
1550
1551 // Object describing how to sort a vector of elements using
1552 // local ID as the key, very short lived object
1554 public:
1555 LocalIdCompare(const STK_Interface * mesh) : mesh_(mesh) {}
1556
1557 // Compares two stk mesh entities based on local ID
1558 bool operator() (stk::mesh::Entity a,stk::mesh::Entity b)
1559 { return mesh_->elementLocalId(a) < mesh_->elementLocalId(b);}
1560
1561 private:
1563 };
1564};
1565
1566template <typename ArrayT>
1567void STK_Interface::setSolutionFieldData(const std::string & fieldName,const std::string & blockId,
1568 const std::vector<std::size_t> & localElementIds,const ArrayT & solutionValues,double scaleValue)
1569{
1570 const std::vector<stk::mesh::Entity> & elements = *(this->getElementsOrderedByLID());
1571 auto solutionValues_h = Kokkos::create_mirror_view(solutionValues);
1572 Kokkos::deep_copy(solutionValues_h, solutionValues);
1573
1574 int field_axis = -1;
1575 if(isMeshCoordField(blockId,fieldName,field_axis)) {
1576 setDispFieldData(fieldName,blockId,field_axis,localElementIds,solutionValues_h);
1577 return;
1578 }
1579
1580 SolutionFieldType * field = this->getSolutionField(fieldName,blockId);
1581
1582 for(std::size_t cell=0;cell<localElementIds.size();cell++) {
1583 std::size_t localId = localElementIds[cell];
1584 stk::mesh::Entity element = elements[localId];
1585
1586 // loop over nodes set solution values
1587 const size_t num_nodes = bulkData_->num_nodes(element);
1588 stk::mesh::Entity const* nodes = bulkData_->begin_nodes(element);
1589 for(std::size_t i=0; i<num_nodes; ++i) {
1590 stk::mesh::Entity node = nodes[i];
1591
1592 double * solnData = stk::mesh::field_data(*field,node);
1593 // TEUCHOS_ASSERT(solnData!=0); // only needed if blockId is not specified
1594 solnData[0] = scaleValue*solutionValues_h(cell,i);
1595 }
1596 }
1597}
1598
1599template <typename ArrayT>
1600void STK_Interface::setDispFieldData(const std::string & fieldName,const std::string & blockId,int axis,
1601 const std::vector<std::size_t> & localElementIds,const ArrayT & dispValues)
1602{
1603 TEUCHOS_ASSERT(axis>=0); // sanity check
1604
1605 const std::vector<stk::mesh::Entity> & elements = *(this->getElementsOrderedByLID());
1606
1607 SolutionFieldType * field = this->getSolutionField(fieldName,blockId);
1608 const VectorFieldType & coord_field = this->getCoordinatesField();
1609
1610 for(std::size_t cell=0;cell<localElementIds.size();cell++) {
1611 std::size_t localId = localElementIds[cell];
1612 stk::mesh::Entity element = elements[localId];
1613
1614 // loop over nodes set solution values
1615 const size_t num_nodes = bulkData_->num_nodes(element);
1616 stk::mesh::Entity const* nodes = bulkData_->begin_nodes(element);
1617 for(std::size_t i=0; i<num_nodes; ++i) {
1618 stk::mesh::Entity node = nodes[i];
1619
1620 double * solnData = stk::mesh::field_data(*field,node);
1621 double * coordData = stk::mesh::field_data(coord_field,node);
1622 // TEUCHOS_ASSERT(solnData!=0); // only needed if blockId is not specified
1623 solnData[0] = dispValues(cell,i)-coordData[axis];
1624 }
1625 }
1626}
1627
1628template <typename ArrayT>
1629void STK_Interface::getSolutionFieldData(const std::string & fieldName,const std::string & blockId,
1630 const std::vector<std::size_t> & localElementIds,ArrayT & solutionValues) const
1631{
1632 const std::vector<stk::mesh::Entity> & elements = *(this->getElementsOrderedByLID());
1633
1634 solutionValues = Kokkos::createDynRankView(solutionValues,
1635 "solutionValues",
1636 localElementIds.size(),
1637 bulkData_->num_nodes(elements[localElementIds[0]]));
1638
1639 SolutionFieldType * field = this->getSolutionField(fieldName,blockId);
1640
1641 for(std::size_t cell=0;cell<localElementIds.size();cell++) {
1642 std::size_t localId = localElementIds[cell];
1643 stk::mesh::Entity element = elements[localId];
1644
1645 // loop over nodes set solution values
1646 const size_t num_nodes = bulkData_->num_nodes(element);
1647 stk::mesh::Entity const* nodes = bulkData_->begin_nodes(element);
1648 for(std::size_t i=0; i<num_nodes; ++i) {
1649 stk::mesh::Entity node = nodes[i];
1650
1651 double * solnData = stk::mesh::field_data(*field,node);
1652 // TEUCHOS_ASSERT(solnData!=0); // only needed if blockId is not specified
1653 solutionValues(cell,i) = solnData[0];
1654 }
1655 }
1656}
1657
1658template <typename ArrayT>
1659void STK_Interface::setCellFieldData(const std::string & fieldName,const std::string & blockId,
1660 const std::vector<std::size_t> & localElementIds,const ArrayT & solutionValues,double scaleValue)
1661{
1662 const std::vector<stk::mesh::Entity> & elements = *(this->getElementsOrderedByLID());
1663
1664 SolutionFieldType * field = this->getCellField(fieldName,blockId);
1665
1666 auto solutionValues_h = Kokkos::create_mirror_view(solutionValues);
1667 Kokkos::deep_copy(solutionValues_h, solutionValues);
1668
1669 for(std::size_t cell=0;cell<localElementIds.size();cell++) {
1670 std::size_t localId = localElementIds[cell];
1671 stk::mesh::Entity element = elements[localId];
1672
1673 double * solnData = stk::mesh::field_data(*field,element);
1674 TEUCHOS_ASSERT(solnData!=0); // only needed if blockId is not specified
1675 solnData[0] = scaleValue*solutionValues_h.access(cell,0);
1676 }
1677}
1678
1679template <typename ArrayT>
1680void STK_Interface::setEdgeFieldData(const std::string & fieldName,const std::string & blockId,
1681 const std::vector<std::size_t> & localEdgeIds,const ArrayT & edgeValues,double scaleValue)
1682{
1683 const std::vector<stk::mesh::Entity> & edges = *(this->getEdgesOrderedByLID());
1684
1685 SolutionFieldType * field = this->getEdgeField(fieldName,blockId);
1686
1687 auto edgeValues_h = Kokkos::create_mirror_view(edgeValues);
1688 Kokkos::deep_copy(edgeValues_h, edgeValues);
1689
1690 for(std::size_t idx=0;idx<localEdgeIds.size();idx++) {
1691 std::size_t localId = localEdgeIds[idx];
1692 stk::mesh::Entity edge = edges[localId];
1693
1694 double * solnData = stk::mesh::field_data(*field,edge);
1695 TEUCHOS_ASSERT(solnData!=0); // only needed if blockId is not specified
1696 solnData[0] = scaleValue*edgeValues_h.access(idx,0);
1697 }
1698}
1699
1700template <typename ArrayT>
1701void STK_Interface::setFaceFieldData(const std::string & fieldName,const std::string & blockId,
1702 const std::vector<std::size_t> & localFaceIds,const ArrayT & faceValues,double scaleValue)
1703{
1704 const std::vector<stk::mesh::Entity> & faces = *(this->getFacesOrderedByLID());
1705
1706 SolutionFieldType * field = this->getFaceField(fieldName,blockId);
1707
1708 auto faceValues_h = Kokkos::create_mirror_view(faceValues);
1709 Kokkos::deep_copy(faceValues_h, faceValues);
1710
1711 for(std::size_t idx=0;idx<localFaceIds.size();idx++) {
1712 std::size_t localId = localFaceIds[idx];
1713 stk::mesh::Entity face = faces[localId];
1714
1715 double * solnData = stk::mesh::field_data(*field,face);
1716 TEUCHOS_ASSERT(solnData!=0); // only needed if blockId is not specified
1717 solnData[0] = scaleValue*faceValues_h.access(idx,0);
1718 }
1719}
1720
1722
1723template <typename ArrayT>
1724void STK_Interface::getElementVertices(const std::vector<std::size_t> & localElementIds, ArrayT & vertices) const
1725{
1727 //
1728 // gather from the intrinsic mesh coordinates (non-lagrangian)
1729 //
1730
1731 const std::vector<stk::mesh::Entity> & elements = *(this->getElementsOrderedByLID());
1732
1733 // convert to a vector of entity objects
1734 std::vector<stk::mesh::Entity> selected_elements;
1735 for(std::size_t cell=0;cell<localElementIds.size();cell++)
1736 selected_elements.push_back(elements[localElementIds[cell]]);
1737
1738 getElementVertices_FromCoords(selected_elements,vertices);
1739 }
1740 else {
1741 TEUCHOS_TEST_FOR_EXCEPTION(true,std::invalid_argument,
1742 "STK_Interface::getElementVertices: Cannot call this method when field coordinates are used "
1743 "without specifying an element block.");
1744 }
1745}
1746
1747template <typename ArrayT>
1748void STK_Interface::getElementVertices(const std::vector<stk::mesh::Entity> & elements, ArrayT & vertices) const
1749{
1751 getElementVertices_FromCoords(elements,vertices);
1752 }
1753 else {
1754 TEUCHOS_TEST_FOR_EXCEPTION(true,std::invalid_argument,
1755 "STK_Interface::getElementVertices: Cannot call this method when field coordinates are used "
1756 "without specifying an element block.");
1757 }
1758}
1759
1760template <typename ArrayT>
1761void STK_Interface::getElementVertices(const std::vector<stk::mesh::Entity> & elements,const std::string & eBlock, ArrayT & vertices) const
1762{
1764 getElementVertices_FromCoords(elements,vertices);
1765 }
1766 else {
1767 getElementVertices_FromField(elements,eBlock,vertices);
1768 }
1769}
1770
1771template <typename ArrayT>
1772void STK_Interface::getElementVertices(const std::vector<std::size_t> & localElementIds,const std::string & eBlock, ArrayT & vertices) const
1773{
1774 const std::vector<stk::mesh::Entity> & elements = *(this->getElementsOrderedByLID());
1775
1776 // convert to a vector of entity objects
1777 std::vector<stk::mesh::Entity> selected_elements;
1778 for(std::size_t cell=0;cell<localElementIds.size();cell++)
1779 selected_elements.push_back(elements[localElementIds[cell]]);
1780
1782 getElementVertices_FromCoords(selected_elements,vertices);
1783 }
1784 else {
1785 getElementVertices_FromField(selected_elements,eBlock,vertices);
1786 }
1787}
1788
1789template <typename ArrayT>
1790void STK_Interface::getElementVerticesNoResize(const std::vector<std::size_t> & localElementIds, ArrayT & vertices) const
1791{
1793 //
1794 // gather from the intrinsic mesh coordinates (non-lagrangian)
1795 //
1796
1797 const std::vector<stk::mesh::Entity> & elements = *(this->getElementsOrderedByLID());
1798
1799 // convert to a vector of entity objects
1800 std::vector<stk::mesh::Entity> selected_elements;
1801 for(std::size_t cell=0;cell<localElementIds.size();cell++)
1802 selected_elements.push_back(elements[localElementIds[cell]]);
1803
1804 getElementVertices_FromCoordsNoResize(selected_elements,vertices);
1805 }
1806 else {
1807 TEUCHOS_TEST_FOR_EXCEPTION(true,std::invalid_argument,
1808 "STK_Interface::getElementVerticesNoResize: Cannot call this method when field coordinates are used "
1809 "without specifying an element block.");
1810 }
1811}
1812
1813template <typename ArrayT>
1814void STK_Interface::getElementVerticesNoResize(const std::vector<stk::mesh::Entity> & elements, ArrayT & vertices) const
1815{
1817 getElementVertices_FromCoordsNoResize(elements,vertices);
1818 }
1819 else {
1820 TEUCHOS_TEST_FOR_EXCEPTION(true,std::invalid_argument,
1821 "STK_Interface::getElementVerticesNoResize: Cannot call this method when field coordinates are used "
1822 "without specifying an element block.");
1823 }
1824}
1825
1826template <typename ArrayT>
1827void STK_Interface::getElementVerticesNoResize(const std::vector<stk::mesh::Entity> & elements,const std::string & eBlock, ArrayT & vertices) const
1828{
1830 getElementVertices_FromCoordsNoResize(elements,vertices);
1831 }
1832 else {
1833 getElementVertices_FromFieldNoResize(elements,eBlock,vertices);
1834 }
1835}
1836
1837template <typename ArrayT>
1838void STK_Interface::getElementVerticesNoResize(const std::vector<std::size_t> & localElementIds,const std::string & eBlock, ArrayT & vertices) const
1839{
1840 const std::vector<stk::mesh::Entity> & elements = *(this->getElementsOrderedByLID());
1841
1842 // convert to a vector of entity objects
1843 std::vector<stk::mesh::Entity> selected_elements;
1844 for(std::size_t cell=0;cell<localElementIds.size();cell++)
1845 selected_elements.push_back(elements[localElementIds[cell]]);
1846
1848 getElementVertices_FromCoordsNoResize(selected_elements,vertices);
1849 }
1850 else {
1851 getElementVertices_FromFieldNoResize(selected_elements,eBlock,vertices);
1852 }
1853}
1854
1855template <typename ArrayT>
1856void STK_Interface::getElementVertices_FromCoords(const std::vector<stk::mesh::Entity> & elements, ArrayT & vertices) const
1857{
1858 // nothing to do! silently return
1859 if(elements.size() == 0) {
1860 vertices = Kokkos::createDynRankView(vertices, "vertices", 0, 0, 0);
1861 return;
1862 }
1863
1864 //
1865 // gather from the intrinsic mesh coordinates (non-lagrangian)
1866 //
1867
1868 // get *master* cell toplogy...(belongs to first element)
1869 const auto masterVertexCount
1870 = stk::mesh::get_cell_topology(bulkData_->bucket(elements[0]).topology()).getCellTopologyData()->vertex_count;
1871
1872 // allocate space
1873 vertices = Kokkos::createDynRankView(vertices, "vertices", elements.size(), masterVertexCount,getDimension());
1874 auto vertices_h = Kokkos::create_mirror_view(vertices);
1875 Kokkos::deep_copy(vertices_h, vertices);
1876
1877 // loop over each requested element
1878 const auto dim = getDimension();
1879 for(std::size_t cell = 0; cell < elements.size(); cell++) {
1880 const auto element = elements[cell];
1881 TEUCHOS_ASSERT(element != 0);
1882
1883 const auto vertexCount
1884 = stk::mesh::get_cell_topology(bulkData_->bucket(element).topology()).getCellTopologyData()->vertex_count;
1885 TEUCHOS_TEST_FOR_EXCEPTION(vertexCount != masterVertexCount, std::runtime_error,
1886 "In call to STK_Interface::getElementVertices all elements "
1887 "must have the same vertex count!");
1888
1889 // loop over all element nodes
1890 const size_t num_nodes = bulkData_->num_nodes(element);
1891 auto const* nodes = bulkData_->begin_nodes(element);
1892 TEUCHOS_TEST_FOR_EXCEPTION(num_nodes!=masterVertexCount,std::runtime_error,
1893 "In call to STK_Interface::getElementVertices cardinality of "
1894 "element node relations must be the vertex count!");
1895 for(std::size_t node = 0; node < num_nodes; ++node) {
1896 const double * coord = getNodeCoordinates(nodes[node]);
1897
1898 // set each dimension of the coordinate
1899 for(unsigned d=0;d<dim;d++)
1900 vertices_h(cell,node,d) = coord[d];
1901 }
1902 }
1903 Kokkos::deep_copy(vertices, vertices_h);
1904}
1905
1906template <typename ArrayT>
1907void STK_Interface::getElementVertices_FromCoordsNoResize(const std::vector<stk::mesh::Entity> & elements, ArrayT & vertices) const
1908{
1909 // nothing to do! silently return
1910 if(elements.size()==0) {
1911 return;
1912 }
1913
1914 //
1915 // gather from the intrinsic mesh coordinates (non-lagrangian)
1916 //
1917
1918 // get *master* cell toplogy...(belongs to first element)
1919 unsigned masterVertexCount
1920 = stk::mesh::get_cell_topology(bulkData_->bucket(elements[0]).topology()).getCellTopologyData()->vertex_count;
1921
1922 // loop over each requested element
1923 unsigned dim = getDimension();
1924 auto vertices_h = Kokkos::create_mirror_view(vertices);
1925 for(std::size_t cell=0;cell<elements.size();cell++) {
1926 stk::mesh::Entity element = elements[cell];
1927 TEUCHOS_ASSERT(element!=0);
1928
1929 unsigned vertexCount
1930 = stk::mesh::get_cell_topology(bulkData_->bucket(element).topology()).getCellTopologyData()->vertex_count;
1931 TEUCHOS_TEST_FOR_EXCEPTION(vertexCount!=masterVertexCount,std::runtime_error,
1932 "In call to STK_Interface::getElementVertices all elements "
1933 "must have the same vertex count!");
1934
1935 // loop over all element nodes
1936 const size_t num_nodes = bulkData_->num_nodes(element);
1937 stk::mesh::Entity const* nodes = bulkData_->begin_nodes(element);
1938 TEUCHOS_TEST_FOR_EXCEPTION(num_nodes!=masterVertexCount,std::runtime_error,
1939 "In call to STK_Interface::getElementVertices cardinality of "
1940 "element node relations must be the vertex count!");
1941 for(std::size_t node=0; node<num_nodes; ++node) {
1942 const double * coord = getNodeCoordinates(nodes[node]);
1943
1944 // set each dimension of the coordinate
1945 for(unsigned d=0;d<dim;d++)
1946 vertices_h(cell,node,d) = coord[d];
1947 }
1948 }
1949 Kokkos::deep_copy(vertices, vertices_h);
1950}
1951
1952template <typename ArrayT>
1953void STK_Interface::getElementVertices_FromField(const std::vector<stk::mesh::Entity> & elements,const std::string & eBlock, ArrayT & vertices) const
1954{
1955 TEUCHOS_ASSERT(useFieldCoordinates_);
1956
1957 // nothing to do! silently return
1958 if(elements.size()==0) {
1959 vertices = Kokkos::createDynRankView(vertices,"vertices",0,0,0);
1960 return;
1961 }
1962
1963 // get *master* cell toplogy...(belongs to first element)
1964 unsigned masterVertexCount
1965 = stk::mesh::get_cell_topology(bulkData_->bucket(elements[0]).topology()).getCellTopologyData()->vertex_count;
1966
1967 // allocate space
1968 vertices = Kokkos::createDynRankView(vertices,"vertices",elements.size(),masterVertexCount,getDimension());
1969 auto vertices_h = Kokkos::create_mirror_view(vertices);
1970 std::map<std::string,std::vector<std::string> >::const_iterator itr = meshCoordFields_.find(eBlock);
1971 if(itr==meshCoordFields_.end()) {
1972 // no coordinate field set for this element block
1973 TEUCHOS_ASSERT(false);
1974 }
1975
1976 const std::vector<std::string> & coordField = itr->second;
1977 std::vector<SolutionFieldType*> fields(getDimension());
1978 for(std::size_t d=0;d<fields.size();d++) {
1979 fields[d] = this->getSolutionField(coordField[d],eBlock);
1980 }
1981
1982 for(std::size_t cell=0;cell<elements.size();cell++) {
1983 stk::mesh::Entity element = elements[cell];
1984
1985 // loop over nodes set solution values
1986 const size_t num_nodes = bulkData_->num_nodes(element);
1987 stk::mesh::Entity const* nodes = bulkData_->begin_nodes(element);
1988 for(std::size_t i=0; i<num_nodes; ++i) {
1989 stk::mesh::Entity node = nodes[i];
1990
1991 const double * coord = getNodeCoordinates(node);
1992
1993 for(unsigned d=0;d<getDimension();d++) {
1994 double * solnData = stk::mesh::field_data(*fields[d],node);
1995
1996 // recall mesh field coordinates are stored as displacements
1997 // from the mesh coordinates, make sure to add them together
1998 vertices_h(cell,i,d) = solnData[0]+coord[d];
1999 }
2000 }
2001 }
2002 Kokkos::deep_copy(vertices, vertices_h);
2003}
2004
2005template <typename ArrayT>
2006void STK_Interface::getElementVertices_FromFieldNoResize(const std::vector<stk::mesh::Entity> & elements,
2007 const std::string & eBlock, ArrayT & vertices) const
2008{
2009 TEUCHOS_ASSERT(useFieldCoordinates_);
2010
2011 // nothing to do! silently return
2012 if(elements.size()==0) {
2013 return;
2014 }
2015
2016 std::map<std::string,std::vector<std::string> >::const_iterator itr = meshCoordFields_.find(eBlock);
2017 if(itr==meshCoordFields_.end()) {
2018 // no coordinate field set for this element block
2019 TEUCHOS_ASSERT(false);
2020 }
2021
2022 const std::vector<std::string> & coordField = itr->second;
2023 std::vector<SolutionFieldType*> fields(getDimension());
2024 for(std::size_t d=0;d<fields.size();d++) {
2025 fields[d] = this->getSolutionField(coordField[d],eBlock);
2026 }
2027
2028 for(std::size_t cell=0;cell<elements.size();cell++) {
2029 stk::mesh::Entity element = elements[cell];
2030
2031 // loop over nodes set solution values
2032 const size_t num_nodes = bulkData_->num_nodes(element);
2033 stk::mesh::Entity const* nodes = bulkData_->begin_nodes(element);
2034 for(std::size_t i=0; i<num_nodes; ++i) {
2035 stk::mesh::Entity node = nodes[i];
2036
2037 const double * coord = getNodeCoordinates(node);
2038
2039 for(unsigned d=0;d<getDimension();d++) {
2040 double * solnData = stk::mesh::field_data(*fields[d],node);
2041
2042 // recall mesh field coordinates are stored as displacements
2043 // from the mesh coordinates, make sure to add them together
2044 vertices(cell,i,d) = solnData[0]+coord[d];
2045 }
2046 }
2047 }
2048}
2049
2051
2052template <typename ArrayT>
2053void STK_Interface::getElementNodes(const std::vector<std::size_t> & localElementIds, ArrayT & nodes) const
2054{
2056 //
2057 // gather from the intrinsic mesh coordinates (non-lagrangian)
2058 //
2059
2060 const std::vector<stk::mesh::Entity> & elements = *(this->getElementsOrderedByLID());
2061
2062 // convert to a vector of entity objects
2063 std::vector<stk::mesh::Entity> selected_elements;
2064 for(std::size_t cell=0;cell<localElementIds.size();cell++)
2065 selected_elements.push_back(elements[localElementIds[cell]]);
2066
2067 getElementNodes_FromCoords(selected_elements,nodes);
2068 }
2069 else {
2070 TEUCHOS_TEST_FOR_EXCEPTION(true,std::invalid_argument,
2071 "STK_Interface::getElementNodes: Cannot call this method when field coordinates are used "
2072 "without specifying an element block.");
2073 }
2074}
2075
2076template <typename ArrayT>
2077void STK_Interface::getElementNodes(const std::vector<stk::mesh::Entity> & elements, ArrayT & nodes) const
2078{
2080 getElementNodes_FromCoords(elements,nodes);
2081 }
2082 else {
2083 TEUCHOS_TEST_FOR_EXCEPTION(true,std::invalid_argument,
2084 "STK_Interface::getElementNodes: Cannot call this method when field coordinates are used "
2085 "without specifying an element block.");
2086 }
2087}
2088
2089template <typename ArrayT>
2090void STK_Interface::getElementNodes(const std::vector<stk::mesh::Entity> & elements,const std::string & eBlock, ArrayT & nodes) const
2091{
2093 getElementNodes_FromCoords(elements,nodes);
2094 }
2095 else {
2096 getElementNodes_FromField(elements,eBlock,nodes);
2097 }
2098}
2099
2100template <typename ArrayT>
2101void STK_Interface::getElementNodes(const std::vector<std::size_t> & localElementIds,const std::string & eBlock, ArrayT & nodes) const
2102{
2103 const std::vector<stk::mesh::Entity> & elements = *(this->getElementsOrderedByLID());
2104
2105 // convert to a vector of entity objects
2106 std::vector<stk::mesh::Entity> selected_elements;
2107 for(std::size_t cell=0;cell<localElementIds.size();cell++)
2108 selected_elements.push_back(elements[localElementIds[cell]]);
2109
2111 getElementNodes_FromCoords(selected_elements,nodes);
2112 }
2113 else {
2114 getElementNodes_FromField(selected_elements,eBlock,nodes);
2115 }
2116}
2117
2118template <typename ArrayT>
2119void STK_Interface::getElementNodesNoResize(const std::vector<std::size_t> & localElementIds, ArrayT & nodes) const
2120{
2122 //
2123 // gather from the intrinsic mesh coordinates (non-lagrangian)
2124 //
2125
2126 const std::vector<stk::mesh::Entity> & elements = *(this->getElementsOrderedByLID());
2127
2128 // convert to a vector of entity objects
2129 std::vector<stk::mesh::Entity> selected_elements;
2130 for(std::size_t cell=0;cell<localElementIds.size();cell++)
2131 selected_elements.push_back(elements[localElementIds[cell]]);
2132
2133 getElementNodes_FromCoordsNoResize(selected_elements,nodes);
2134 }
2135 else {
2136 TEUCHOS_TEST_FOR_EXCEPTION(true,std::invalid_argument,
2137 "STK_Interface::getElementNodesNoResize: Cannot call this method when field coordinates are used "
2138 "without specifying an element block.");
2139 }
2140}
2141
2142template <typename ArrayT>
2143void STK_Interface::getElementNodesNoResize(const std::vector<stk::mesh::Entity> & elements, ArrayT & nodes) const
2144{
2146 getElementNodes_FromCoordsNoResize(elements,nodes);
2147 }
2148 else {
2149 TEUCHOS_TEST_FOR_EXCEPTION(true,std::invalid_argument,
2150 "STK_Interface::getElementNodesNoResize: Cannot call this method when field coordinates are used "
2151 "without specifying an element block.");
2152 }
2153}
2154
2155template <typename ArrayT>
2156void STK_Interface::getElementNodesNoResize(const std::vector<stk::mesh::Entity> & elements,const std::string & eBlock, ArrayT & nodes) const
2157{
2159 getElementNodes_FromCoordsNoResize(elements,nodes);
2160 }
2161 else {
2162 getElementNodes_FromFieldNoResize(elements,eBlock,nodes);
2163 }
2164}
2165
2166template <typename ArrayT>
2167void STK_Interface::getElementNodesNoResize(const std::vector<std::size_t> & localElementIds,const std::string & eBlock, ArrayT & nodes) const
2168{
2169 const std::vector<stk::mesh::Entity> & elements = *(this->getElementsOrderedByLID());
2170
2171 // convert to a vector of entity objects
2172 std::vector<stk::mesh::Entity> selected_elements;
2173 for(std::size_t cell=0;cell<localElementIds.size();cell++)
2174 selected_elements.push_back(elements[localElementIds[cell]]);
2175
2177 getElementNodes_FromCoordsNoResize(selected_elements,nodes);
2178 }
2179 else {
2180 getElementNodes_FromFieldNoResize(selected_elements,eBlock,nodes);
2181 }
2182}
2183
2184template <typename ArrayT>
2185void STK_Interface::getElementNodes_FromCoords(const std::vector<stk::mesh::Entity> & elements, ArrayT & nodes) const
2186{
2187 // nothing to do! silently return
2188 if(elements.size() == 0) {
2189 nodes = Kokkos::createDynRankView(nodes, "nodes", 0, 0, 0);
2190 return;
2191 }
2192
2193 //
2194 // gather from the intrinsic mesh coordinates (non-lagrangian)
2195 //
2196
2197 // get *master* cell toplogy...(belongs to first element)
2198 const auto masterNodeCount
2199 = stk::mesh::get_cell_topology(bulkData_->bucket(elements[0]).topology()).getCellTopologyData()->node_count;
2200
2201 // allocate space
2202 nodes = Kokkos::createDynRankView(nodes, "nodes", elements.size(), masterNodeCount,getDimension());
2203 auto nodes_h = Kokkos::create_mirror_view(nodes);
2204 Kokkos::deep_copy(nodes_h, nodes);
2205
2206 // loop over each requested element
2207 const auto dim = getDimension();
2208 for(std::size_t cell = 0; cell < elements.size(); cell++) {
2209 const auto element = elements[cell];
2210 TEUCHOS_ASSERT(element != 0);
2211
2212 const auto nodeCount
2213 = stk::mesh::get_cell_topology(bulkData_->bucket(element).topology()).getCellTopologyData()->node_count;
2214 TEUCHOS_TEST_FOR_EXCEPTION(nodeCount != masterNodeCount, std::runtime_error,
2215 "In call to STK_Interface::getElementNodes all elements "
2216 "must have the same node count!");
2217
2218 // loop over all element nodes
2219 const size_t num_nodes = bulkData_->num_nodes(element);
2220 auto const* elem_nodes = bulkData_->begin_nodes(element);
2221 TEUCHOS_TEST_FOR_EXCEPTION(num_nodes!=masterNodeCount,std::runtime_error,
2222 "In call to STK_Interface::getElementNodes cardinality of "
2223 "element node relations must be the node count!");
2224 for(std::size_t node = 0; node < num_nodes; ++node) {
2225 const double * coord = getNodeCoordinates(elem_nodes[node]);
2226
2227 // set each dimension of the coordinate
2228 for(unsigned d=0;d<dim;d++)
2229 nodes_h(cell,node,d) = coord[d];
2230 }
2231 }
2232 Kokkos::deep_copy(nodes, nodes_h);
2233}
2234
2235template <typename ArrayT>
2236void STK_Interface::getElementNodes_FromCoordsNoResize(const std::vector<stk::mesh::Entity> & elements, ArrayT & nodes) const
2237{
2238 // nothing to do! silently return
2239 if(elements.size()==0) {
2240 return;
2241 }
2242
2243 //
2244 // gather from the intrinsic mesh coordinates (non-lagrangian)
2245 //
2246
2247 // get *master* cell toplogy...(belongs to first element)
2248 unsigned masterNodeCount
2249 = stk::mesh::get_cell_topology(bulkData_->bucket(elements[0]).topology()).getCellTopologyData()->node_count;
2250
2251 // loop over each requested element
2252 unsigned dim = getDimension();
2253 auto nodes_h = Kokkos::create_mirror_view(nodes);
2254 for(std::size_t cell=0;cell<elements.size();cell++) {
2255 stk::mesh::Entity element = elements[cell];
2256 TEUCHOS_ASSERT(element!=0);
2257
2258 unsigned nodeCount
2259 = stk::mesh::get_cell_topology(bulkData_->bucket(element).topology()).getCellTopologyData()->node_count;
2260 TEUCHOS_TEST_FOR_EXCEPTION(nodeCount!=masterNodeCount,std::runtime_error,
2261 "In call to STK_Interface::getElementNodes all elements "
2262 "must have the same node count!");
2263
2264 // loop over all element nodes
2265 const size_t num_nodes = bulkData_->num_nodes(element);
2266 stk::mesh::Entity const* elem_nodes = bulkData_->begin_nodes(element);
2267 TEUCHOS_TEST_FOR_EXCEPTION(num_nodes!=masterNodeCount,std::runtime_error,
2268 "In call to STK_Interface::getElementNodes cardinality of "
2269 "element node relations must be the node count!");
2270 for(std::size_t node=0; node<num_nodes; ++node) {
2271 const double * coord = getNodeCoordinates(elem_nodes[node]);
2272
2273 // set each dimension of the coordinate
2274 for(unsigned d=0;d<dim;d++)
2275 nodes_h(cell,node,d) = coord[d];
2276 }
2277 }
2278 Kokkos::deep_copy(nodes, nodes_h);
2279}
2280
2281template <typename ArrayT>
2282void STK_Interface::getElementNodes_FromField(const std::vector<stk::mesh::Entity> & elements,const std::string & eBlock, ArrayT & nodes) const
2283{
2284 TEUCHOS_ASSERT(useFieldCoordinates_);
2285
2286 // nothing to do! silently return
2287 if(elements.size()==0) {
2288 nodes = Kokkos::createDynRankView(nodes,"nodes",0,0,0);
2289 return;
2290 }
2291
2292 // get *master* cell toplogy...(belongs to first element)
2293 unsigned masterNodeCount
2294 = stk::mesh::get_cell_topology(bulkData_->bucket(elements[0]).topology()).getCellTopologyData()->node_count;
2295
2296 // allocate space
2297 nodes = Kokkos::createDynRankView(nodes,"nodes",elements.size(),masterNodeCount,getDimension());
2298 auto nodes_h = Kokkos::create_mirror_view(nodes);
2299 std::map<std::string,std::vector<std::string> >::const_iterator itr = meshCoordFields_.find(eBlock);
2300 if(itr==meshCoordFields_.end()) {
2301 // no coordinate field set for this element block
2302 TEUCHOS_ASSERT(false);
2303 }
2304
2305 const std::vector<std::string> & coordField = itr->second;
2306 std::vector<SolutionFieldType*> fields(getDimension());
2307 for(std::size_t d=0;d<fields.size();d++) {
2308 fields[d] = this->getSolutionField(coordField[d],eBlock);
2309 }
2310
2311 // loop over elements
2312 for(std::size_t cell=0;cell<elements.size();cell++) {
2313 stk::mesh::Entity element = elements[cell];
2314
2315 // loop over nodes set solution values
2316 const size_t num_nodes = bulkData_->num_nodes(element);
2317 stk::mesh::Entity const* elem_nodes = bulkData_->begin_nodes(element);
2318 for(std::size_t i=0; i<num_nodes; ++i) {
2319 stk::mesh::Entity node = elem_nodes[i];
2320
2321 const double * coord = getNodeCoordinates(node);
2322
2323 for(unsigned d=0;d<getDimension();d++) {
2324 double * solnData = stk::mesh::field_data(*fields[d],node);
2325
2326 // recall mesh field coordinates are stored as displacements
2327 // from the mesh coordinates, make sure to add them together
2328 nodes_h(cell,i,d) = solnData[0]+coord[d];
2329 }
2330 }
2331 }
2332 Kokkos::deep_copy(nodes, nodes_h);
2333}
2334
2335template <typename ArrayT>
2336void STK_Interface::getElementNodes_FromFieldNoResize(const std::vector<stk::mesh::Entity> & elements,
2337 const std::string & eBlock, ArrayT & nodes) const
2338{
2339 TEUCHOS_ASSERT(useFieldCoordinates_);
2340
2341 // nothing to do! silently return
2342 if(elements.size()==0) {
2343 return;
2344 }
2345
2346 std::map<std::string,std::vector<std::string> >::const_iterator itr = meshCoordFields_.find(eBlock);
2347 if(itr==meshCoordFields_.end()) {
2348 // no coordinate field set for this element block
2349 TEUCHOS_ASSERT(false);
2350 }
2351
2352 const std::vector<std::string> & coordField = itr->second;
2353 std::vector<SolutionFieldType*> fields(getDimension());
2354 for(std::size_t d=0;d<fields.size();d++) {
2355 fields[d] = this->getSolutionField(coordField[d],eBlock);
2356 }
2357
2358 // loop over elements
2359 for(std::size_t cell=0;cell<elements.size();cell++) {
2360 stk::mesh::Entity element = elements[cell];
2361
2362 // loop over nodes set solution values
2363 const size_t num_nodes = bulkData_->num_nodes(element);
2364 stk::mesh::Entity const* elem_nodes = bulkData_->begin_nodes(element);
2365 for(std::size_t i=0; i<num_nodes; ++i) {
2366 stk::mesh::Entity node = elem_nodes[i];
2367
2368 const double * coord = getNodeCoordinates(node);
2369
2370 for(unsigned d=0;d<getDimension();d++) {
2371 double * solnData = stk::mesh::field_data(*fields[d],node);
2372
2373 // recall mesh field coordinates are stored as displacements
2374 // from the mesh coordinates, make sure to add them together
2375 nodes(cell,i,d) = solnData[0]+coord[d];
2376 }
2377 }
2378 }
2379}
2380
2381}
2382
2383#endif
PHX::MDField< ScalarT, panzer::Cell, panzer::BASIS > field
A field to which we'll contribute, or in which we'll store, the result of computing this integral.
const std::vector< stk::mesh::EntityId > & getNodes() const
std::vector< stk::mesh::EntityId > nodes_
stk::mesh::EntityId getGID() const
bool operator()(stk::mesh::Entity a, stk::mesh::Entity b)
void rebalance(const Teuchos::ParameterList &params)
const VectorFieldType & getEdgesField() const
stk::mesh::Entity findConnectivityById(stk::mesh::Entity src, stk::mesh::EntityRank tgt_rank, unsigned rel_id) const
Teuchos::RCP< stk::mesh::MetaData > metaData_
Teuchos::RCP< std::vector< stk::mesh::Entity > > orderedFaceVector_
std::pair< Teuchos::RCP< std::vector< std::pair< std::size_t, std::size_t > > >, Teuchos::RCP< std::vector< unsigned int > > > getPeriodicNodePairing() const
void setFaceFieldData(const std::string &fieldName, const std::string &blockId, const std::vector< std::size_t > &localFaceIds, const ArrayT &faceValues, double scaleValue=1.0)
bool isValid(stk::mesh::Entity entity) const
const std::vector< Teuchos::RCP< const PeriodicBC_MatcherBase > > & getPeriodicBCVector() const
void initialize(stk::ParallelMachine parallelMach, bool setupIO=true, const bool buildRefinementSupport=false)
stk::mesh::Field< double > SolutionFieldType
void getMyFaces(std::vector< stk::mesh::Entity > &faces) const
stk::mesh::Part * getElementBlockPart(const std::string &name) const
get the block part
static const std::string edgesString
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.
static const std::string edgeBlockString
void getElementVertices(const std::vector< std::size_t > &localIds, ArrayT &vertices) const
void getElementVertices_FromCoords(const std::vector< stk::mesh::Entity > &elements, ArrayT &vertices) const
std::size_t getNumSidesets() const
get the side set count
Teuchos::RCP< std::vector< stk::mesh::Entity > > orderedElementVector_
std::map< std::string, stk::mesh::Part * > sidesets_
void addPeriodicBC(const Teuchos::RCP< const PeriodicBC_MatcherBase > &bc)
stk::mesh::EntityId elementGlobalId(std::size_t lid) const
stk::mesh::Part * getEdgeBlock(const std::string &name) const
get the block part
std::unordered_map< stk::mesh::EntityId, std::size_t > localFaceIDHash_
void getElementsSharingNodes(const std::vector< stk::mesh::EntityId > nodeId, std::vector< stk::mesh::Entity > &elements) const
get a set of elements sharing multiple nodes
std::size_t getEntityCounts(unsigned entityRank) const
get the global counts for the entity of specified rank
void getElementNodes(const std::vector< std::size_t > &localIds, ArrayT &nodes) const
void getFaceBlockNames(std::vector< std::string > &names) const
bool isFaceLocal(stk::mesh::Entity face) const
void addEntityToEdgeBlock(stk::mesh::Entity entity, stk::mesh::Part *edgeblock)
Teuchos::RCP< const std::vector< stk::mesh::Entity > > getElementsOrderedByLID() const
bool isInitialized() const
Has initialize been called on this mesh object?
std::map< std::string, stk::mesh::Part * > faceBlocks_
void addEntitiesToFaceBlock(std::vector< stk::mesh::Entity > entities, stk::mesh::Part *faceblock)
stk::mesh::EntityRank getNodeRank() const
void addInformationRecords(const std::vector< std::string > &info_records)
std::string containingBlockId(stk::mesh::Entity elmt) const
void getEdgeBlockNames(std::vector< std::string > &names) const
void addEdgeBlock(const std::string &elemBlockName, const std::string &edgeBlockName, const stk::topology &topology)
void getElementVertices_FromField(const std::vector< stk::mesh::Entity > &elements, const std::string &eBlock, ArrayT &vertices) const
const VectorFieldType & getFacesField() const
std::vector< std::size_t > entityCounts_
void getElementsSharingNode(stk::mesh::EntityId nodeId, std::vector< stk::mesh::Entity > &elements) const
get a set of elements sharing a single node
void getElementNodes_FromCoords(const std::vector< stk::mesh::Entity > &elements, ArrayT &nodes) const
Teuchos::RCP< const std::vector< stk::mesh::Entity > > getEdgesOrderedByLID() const
void buildSubcells()
force the mesh to build subcells: edges and faces
void setCellFieldData(const std::string &fieldName, const std::string &blockId, const std::vector< std::size_t > &localElementIds, const ArrayT &solutionValues, double scaleValue=1.0)
stk::mesh::EntityId getMaxEntityId(unsigned entityRank) const
get max entity ID of type entityRank
std::vector< Teuchos::RCP< const PeriodicBC_MatcherBase > > & getPeriodicBCVector()
stk::mesh::EntityRank getElementRank() const
void getOwnedElementsSharingNode(stk::mesh::Entity node, std::vector< stk::mesh::Entity > &elements, std::vector< int > &relIds) const
const double * getNodeCoordinates(stk::mesh::EntityId nodeId) const
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
void setBlockWeight(const std::string &blockId, double weight)
void addEdgeField(const std::string &fieldName, const std::string &blockId)
stk::mesh::EntityId faceGlobalId(stk::mesh::Entity face) const
std::unordered_map< stk::mesh::EntityId, std::size_t > localEdgeIDHash_
void addEntityToFaceBlock(stk::mesh::Entity entity, stk::mesh::Part *faceblock)
void getNeighborElements(std::vector< stk::mesh::Entity > &elements) const
void getMyNodes(const std::string &sideName, const std::string &blockName, std::vector< stk::mesh::Entity > &nodes) const
std::vector< stk::mesh::Part * > edgesPartVec_
const bool & useBoundingBoxSearch() 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)
void addEntityToNodeset(stk::mesh::Entity entity, stk::mesh::Part *nodeset)
void getElementVerticesNoResize(const std::vector< std::size_t > &localIds, ArrayT &vertices) const
void addElement(const Teuchos::RCP< ElementDescriptor > &ed, stk::mesh::Part *block)
bool isEdgeLocal(stk::mesh::Entity edge) const
void addPeriodicBCs(const std::vector< Teuchos::RCP< const PeriodicBC_MatcherBase > > &bc_vec)
void refineMesh(const int numberOfLevels, const bool deleteParentElements)
std::map< std::string, double > blockWeights_
std::map< std::pair< std::string, std::string >, SolutionFieldType * > fieldNameToEdgeField_
void addNode(stk::mesh::EntityId gid, const std::vector< double > &coord)
Teuchos::RCP< const Teuchos::Comm< int > > getComm() const
get the comm associated with this mesh
void getElementNodesNoResize(const std::vector< std::size_t > &localIds, ArrayT &nodes) const
std::map< std::pair< std::string, std::string >, SolutionFieldType * > fieldNameToFaceField_
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::size_t elementLocalId(stk::mesh::Entity elmt) const
stk::mesh::Field< double > VectorFieldType
void addCellField(const std::string &fieldName, const std::string &blockId)
stk::mesh::Field< double > * getSolutionField(const std::string &fieldName, const std::string &blockId) const
std::map< std::pair< std::string, std::string >, SolutionFieldType * > fieldNameToSolution_
void getAllSides(const std::string &sideName, std::vector< stk::mesh::Entity > &sides) const
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
void addMeshCoordFields(const std::string &blockId, const std::vector< std::string > &coordField, const std::string &dispPrefix)
stk::mesh::Field< double > * getEdgeField(const std::string &fieldName, const std::string &blockId) const
std::vector< stk::mesh::EntityId > maxEntityId_
void getElementNodes_FromField(const std::vector< stk::mesh::Entity > &elements, const std::string &eBlock, ArrayT &nodes) const
void getElementBlockNames(std::vector< std::string > &names) const
Teuchos::RCP< stk::mesh::BulkData > bulkData_
void addNodeset(const std::string &name)
void getSidesetNames(std::vector< std::string > &name) const
void addSideset(const std::string &name, const CellTopologyData *ctData)
stk::mesh::EntityRank getFaceRank() const
void getMyEdges(std::vector< stk::mesh::Entity > &edges) const
std::size_t getNumNodesets() const
get the side set count
std::set< std::string > informationRecords_
void getAllEdges(const std::string &edgeBlockName, std::vector< stk::mesh::Entity > &edges) const
std::vector< stk::mesh::Part * > nodesPartVec_
void getElementNodes_FromFieldNoResize(const std::vector< stk::mesh::Entity > &elements, const std::string &eBlock, ArrayT &nodes) const
void initializeFieldsInSTK(const std::map< std::pair< std::string, std::string >, SolutionFieldType * > &nameToField, bool setupIO)
static const std::string nodesString
const VectorFieldType & getCoordinatesField() const
void setUseFieldCoordinates(bool useFieldCoordinates)
Teuchos::RCP< stk::mesh::MetaData > getMetaData() const
stk::mesh::Field< double > * getFaceField(const std::string &fieldName, const std::string &blockId) const
static const std::string facesString
void writeToExodus(const std::string &filename, const bool append=false)
Write this mesh and associated fields to the given output file.
unsigned getDimension() const
get the dimension
bool isMeshCoordField(const std::string &eBlock, const std::string &fieldName, int &axis) const
std::map< std::string, stk::mesh::Part * > nodesets_
void getElementVertices_FromCoordsNoResize(const std::vector< stk::mesh::Entity > &elements, ArrayT &vertices) const
std::map< std::string, std::vector< std::string > > meshCoordFields_
stk::mesh::Field< double > * getCellField(const std::string &fieldName, const std::string &blockId) const
std::size_t faceLocalId(stk::mesh::Entity elmt) const
void getAllFaces(const std::string &faceBlockName, std::vector< stk::mesh::Entity > &faces) const
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.
Teuchos::RCP< std::vector< stk::mesh::Entity > > orderedEdgeVector_
stk::mesh::Part * getNodeset(const std::string &name) const
Teuchos::RCP< const shards::CellTopology > getCellTopology(const std::string &eBlock) const
void endModification(const bool find_and_set_shared_nodes_in_stk=true)
stk::mesh::EntityId edgeGlobalId(std::size_t lid) const
unsigned entityOwnerRank(stk::mesh::Entity entity) const
void getNodesetNames(std::vector< std::string > &name) const
void addEntitiesToEdgeBlock(std::vector< stk::mesh::Entity > entities, stk::mesh::Part *edgeblock)
stk::mesh::Field< ProcIdData > ProcIdFieldType
std::unordered_map< stk::mesh::EntityId, std::size_t > localIDHash_
void addSolutionField(const std::string &fieldName, const std::string &blockId)
void getSubcellIndices(unsigned entityRank, stk::mesh::EntityId elementId, std::vector< stk::mesh::EntityId > &subcellIds) const
std::vector< Teuchos::RCP< const PeriodicBC_MatcherBase > > periodicBCs_
void addEntityToSideset(stk::mesh::Entity entity, stk::mesh::Part *sideset)
Teuchos::RCP< const std::vector< stk::mesh::Entity > > getFacesOrderedByLID() const
void setDispFieldData(const std::string &fieldName, const std::string &blockId, int axis, const std::vector< std::size_t > &localElementIds, const ArrayT &solutionValues)
stk::mesh::EntityId elementGlobalId(stk::mesh::Entity elmt) const
std::map< std::string, Teuchos::RCP< shards::CellTopology > > elementBlockCT_
stk::mesh::EntityRank getSideRank() const
Teuchos::RCP< Teuchos::MpiComm< int > > mpiComm_
bool validBlockId(const std::string &blockId) const
stk::mesh::EntityId edgeGlobalId(stk::mesh::Entity edge) const
void addFaceField(const std::string &fieldName, const std::string &blockId)
std::map< std::string, stk::mesh::Part * > edgeBlocks_
void print(std::ostream &os) const
void getMyElements(std::vector< stk::mesh::Entity > &elements) const
stk::mesh::Part * getOwnedPart() const
Get a pointer to the locally owned part.
void addElementBlock(const std::string &name, const CellTopologyData *ctData)
ProcIdFieldType * getProcessorIdField()
void addFaceBlock(const std::string &elemBlockName, const std::string &faceBlockName, const stk::topology &topology)
Teuchos::RCP< stk::mesh::BulkData > getBulkData() const
std::map< std::string, std::vector< std::string > > meshDispFields_
static const std::string faceBlockString
std::vector< stk::mesh::Part * > facesPartVec_
void printMetaData(std::ostream &os) const
void instantiateBulkData(stk::ParallelMachine parallelMach)
void getElementNodes_FromCoordsNoResize(const std::vector< stk::mesh::Entity > &elements, ArrayT &nodes) const
stk::mesh::EntityId faceGlobalId(std::size_t lid) const
std::map< std::pair< std::string, std::string >, SolutionFieldType * > fieldNameToCellField_
std::size_t edgeLocalId(stk::mesh::Entity elmt) const
std::map< std::string, stk::mesh::Part * > elementBlocks_
static const std::string coordsString
void getSolutionFieldData(const std::string &fieldName, const std::string &blockId, const std::vector< std::size_t > &localElementIds, ArrayT &solutionValues) const
void setUseLowerCaseForIO(bool useLowerCase)
Teuchos::RCP< Teuchos::MpiComm< int > > getSafeCommunicator(stk::ParallelMachine parallelMach) const
std::size_t getNumElementBlocks() const
get the block count
void setBoundingBoxSearchFlag(const bool &searchFlag)
stk::mesh::Part * getFaceBlock(const std::string &name) const
get the block part
stk::mesh::EntityRank getEdgeRank() const
stk::mesh::Part * getSideset(const std::string &name) const
Teuchos::RCP< ElementDescriptor > buildElementDescriptor(stk::mesh::EntityId elmtId, std::vector< stk::mesh::EntityId > &nodes)