Panzer  Version of the Day
Panzer_STK_Interface.cpp
Go to the documentation of this file.
1 // @HEADER
2 // ***********************************************************************
3 //
4 // Panzer: A partial differential equation assembly
5 // engine for strongly coupled complex multiphysics systems
6 // Copyright (2011) Sandia Corporation
7 //
8 // Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
9 // the U.S. Government retains certain rights in this software.
10 //
11 // Redistribution and use in source and binary forms, with or without
12 // modification, are permitted provided that the following conditions are
13 // met:
14 //
15 // 1. Redistributions of source code must retain the above copyright
16 // notice, this list of conditions and the following disclaimer.
17 //
18 // 2. Redistributions in binary form must reproduce the above copyright
19 // notice, this list of conditions and the following disclaimer in the
20 // documentation and/or other materials provided with the distribution.
21 //
22 // 3. Neither the name of the Corporation nor the names of the
23 // contributors may be used to endorse or promote products derived from
24 // this software without specific prior written permission.
25 //
26 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
27 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
28 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
29 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
30 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
31 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
32 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
33 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
34 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
35 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
36 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
37 //
38 // Questions? Contact Roger P. Pawlowski (rppawlo@sandia.gov) and
39 // Eric C. Cyr (eccyr@sandia.gov)
40 // ***********************************************************************
41 // @HEADER
42 
43 #include <PanzerAdaptersSTK_config.hpp>
44 #include <Panzer_STK_Interface.hpp>
45 
46 #include <Teuchos_as.hpp>
47 
48 #include <limits>
49 
50 #include <stk_mesh/base/FieldBase.hpp>
51 #include <stk_mesh/base/Comm.hpp>
52 #include <stk_mesh/base/Selector.hpp>
53 #include <stk_mesh/base/GetEntities.hpp>
54 #include <stk_mesh/base/GetBuckets.hpp>
55 #include <stk_mesh/base/MeshBuilder.hpp>
56 #include <stk_mesh/base/CreateAdjacentEntities.hpp>
57 
58 // #include <stk_rebalance/Rebalance.hpp>
59 // #include <stk_rebalance/Partition.hpp>
60 // #include <stk_rebalance/ZoltanPartition.hpp>
61 // #include <stk_rebalance_utils/RebalanceUtils.hpp>
62 
63 #include <stk_util/parallel/ParallelReduce.hpp>
64 #include <stk_util/parallel/CommSparse.hpp>
65 
66 #ifdef PANZER_HAVE_IOSS
67 #include <Ionit_Initializer.h>
68 #include <stk_io/IossBridge.hpp>
69 #include <stk_io/WriteMesh.hpp>
70 #endif
71 
72 #ifdef PANZER_HAVE_PERCEPT
73 #include <percept/PerceptMesh.hpp>
74 #include <adapt/UniformRefinerPattern.hpp>
75 #include <adapt/UniformRefiner.hpp>
76 #endif
77 
79 
80 #include <set>
81 
82 using Teuchos::RCP;
83 using Teuchos::rcp;
84 
85 namespace panzer_stk {
86 
88 ElementDescriptor::ElementDescriptor(stk::mesh::EntityId gid,const std::vector<stk::mesh::EntityId> & nodes)
89  : gid_(gid), nodes_(nodes) {}
91 
94 Teuchos::RCP<ElementDescriptor>
95 buildElementDescriptor(stk::mesh::EntityId elmtId,std::vector<stk::mesh::EntityId> & nodes)
96 {
97  return Teuchos::rcp(new ElementDescriptor(elmtId,nodes));
98 }
99 
100 const std::string STK_Interface::coordsString = "coordinates";
101 const std::string STK_Interface::nodesString = "nodes";
102 const std::string STK_Interface::edgesString = "edges";
103 const std::string STK_Interface::facesString = "faces";
104 const std::string STK_Interface::edgeBlockString = "edge_block";
105 const std::string STK_Interface::faceBlockString = "face_block";
106 
108  : dimension_(0), initialized_(false), currentLocalId_(0), initialStateTime_(0.0), currentStateTime_(0.0), useFieldCoordinates_(false)
109 {
110  metaData_ = rcp(new stk::mesh::MetaData());
111 }
112 
113 STK_Interface::STK_Interface(Teuchos::RCP<stk::mesh::MetaData> metaData)
114  : dimension_(0), initialized_(false), currentLocalId_(0), initialStateTime_(0.0), currentStateTime_(0.0), useFieldCoordinates_(false)
115 {
116  metaData_ = metaData;
117 }
118 
120  : dimension_(dim), initialized_(false), currentLocalId_(0), useFieldCoordinates_(false)
121 {
122  std::vector<std::string> entity_rank_names = stk::mesh::entity_rank_names();
123  entity_rank_names.push_back("FAMILY_TREE");
124 
125  metaData_ = rcp(new stk::mesh::MetaData(dimension_,entity_rank_names));
126 
128 }
129 
130 void STK_Interface::addSideset(const std::string & name,const CellTopologyData * ctData)
131 {
132  TEUCHOS_ASSERT(not initialized_);
133  TEUCHOS_ASSERT(dimension_!=0);
134 
135  stk::mesh::Part * sideset = metaData_->get_part(name);
136  if(sideset==nullptr)
137  sideset = &metaData_->declare_part_with_topology(name,
138  stk::mesh::get_topology(shards::CellTopology(ctData), dimension_));
139  sidesets_.insert(std::make_pair(name,sideset));
140 }
141 
142 void STK_Interface::addNodeset(const std::string & name)
143 {
144  TEUCHOS_ASSERT(not initialized_);
145  TEUCHOS_ASSERT(dimension_!=0);
146 
147  stk::mesh::Part * nodeset = metaData_->get_part(name);
148  if(nodeset==nullptr) {
149  const CellTopologyData * ctData = shards::getCellTopologyData<shards::Node>();
150  nodeset = &metaData_->declare_part_with_topology(name,
151  stk::mesh::get_topology(shards::CellTopology(ctData), dimension_));
152  }
153  nodesets_.insert(std::make_pair(name,nodeset));
154 }
155 
156 void STK_Interface::addSolutionField(const std::string & fieldName,const std::string & blockId)
157 {
158  TEUCHOS_TEST_FOR_EXCEPTION(!validBlockId(blockId),ElementBlockException,
159  "Unknown element block \"" << blockId << "\"");
160  std::pair<std::string,std::string> key = std::make_pair(fieldName,blockId);
161 
162  // add & declare field if not already added...currently assuming linears
163  if(fieldNameToSolution_.find(key)==fieldNameToSolution_.end()) {
164  SolutionFieldType * field = metaData_->get_field<SolutionFieldType>(stk::topology::NODE_RANK, fieldName);
165  if(field==0)
166  field = &metaData_->declare_field<SolutionFieldType>(stk::topology::NODE_RANK, fieldName);
167  if ( initialized_ ) {
168  metaData_->enable_late_fields();
169  stk::mesh::FieldTraits<SolutionFieldType>::data_type* init_sol = nullptr;
170  stk::mesh::put_field_on_mesh(*field, metaData_->universal_part(),init_sol );
171  }
173  }
174 }
175 
176 void STK_Interface::addCellField(const std::string & fieldName,const std::string & blockId)
177 {
178  TEUCHOS_TEST_FOR_EXCEPTION(!validBlockId(blockId),ElementBlockException,
179  "Unknown element block \"" << blockId << "\"");
180  std::pair<std::string,std::string> key = std::make_pair(fieldName,blockId);
181 
182  // add & declare field if not already added...currently assuming linears
183  if(fieldNameToCellField_.find(key)==fieldNameToCellField_.end()) {
184  SolutionFieldType * field = metaData_->get_field<SolutionFieldType>(stk::topology::ELEMENT_RANK, fieldName);
185  if(field==0)
186  field = &metaData_->declare_field<SolutionFieldType>(stk::topology::ELEMENT_RANK, fieldName);
187 
188  if ( initialized_ ) {
189  metaData_->enable_late_fields();
190  stk::mesh::FieldTraits<SolutionFieldType>::data_type* init_sol = nullptr;
191  stk::mesh::put_field_on_mesh(*field, metaData_->universal_part(),init_sol );
192  }
194  }
195 }
196 
197 void STK_Interface::addEdgeField(const std::string & fieldName,const std::string & blockId)
198 {
199  TEUCHOS_TEST_FOR_EXCEPTION(!validBlockId(blockId),ElementBlockException,
200  "Unknown element block \"" << blockId << "\"");
201  std::pair<std::string,std::string> key = std::make_pair(fieldName,blockId);
202 
203  // add & declare field if not already added...currently assuming linears
204  if(fieldNameToEdgeField_.find(key)==fieldNameToEdgeField_.end()) {
205  SolutionFieldType * field = metaData_->get_field<SolutionFieldType>(stk::topology::EDGE_RANK, fieldName);
206  if(field==0) {
207  field = &metaData_->declare_field<SolutionFieldType>(stk::topology::EDGE_RANK, fieldName);
208  }
209 
210  if ( initialized_ ) {
211  metaData_->enable_late_fields();
212  stk::mesh::FieldTraits<SolutionFieldType>::data_type* init_sol = nullptr;
213  stk::mesh::put_field_on_mesh(*field, metaData_->universal_part(),init_sol );
214  }
216  }
217 }
218 
219 void STK_Interface::addFaceField(const std::string & fieldName,const std::string & blockId)
220 {
221  TEUCHOS_TEST_FOR_EXCEPTION(!validBlockId(blockId),ElementBlockException,
222  "Unknown element block \"" << blockId << "\"");
223  std::pair<std::string,std::string> key = std::make_pair(fieldName,blockId);
224 
225  // add & declare field if not already added...currently assuming linears
226  if(fieldNameToFaceField_.find(key)==fieldNameToFaceField_.end()) {
227  SolutionFieldType * field = metaData_->get_field<SolutionFieldType>(stk::topology::FACE_RANK, fieldName);
228  if(field==0) {
229  field = &metaData_->declare_field<SolutionFieldType>(stk::topology::FACE_RANK, fieldName);
230  }
231 
232  if ( initialized_ ) {
233  metaData_->enable_late_fields();
234  stk::mesh::FieldTraits<SolutionFieldType>::data_type* init_sol = nullptr;
235  stk::mesh::put_field_on_mesh(*field, metaData_->universal_part(),init_sol );
236  }
238  }
239 }
240 
241 void STK_Interface::addMeshCoordFields(const std::string & blockId,
242  const std::vector<std::string> & coordNames,
243  const std::string & dispPrefix)
244 {
245  TEUCHOS_ASSERT(dimension_!=0);
246  TEUCHOS_ASSERT(dimension_==coordNames.size());
247  TEUCHOS_ASSERT(not initialized_);
248  TEUCHOS_TEST_FOR_EXCEPTION(!validBlockId(blockId),ElementBlockException,
249  "Unknown element block \"" << blockId << "\"");
250 
251  // we only allow one alternative coordinate field
252  TEUCHOS_TEST_FOR_EXCEPTION(meshCoordFields_.find(blockId)!=meshCoordFields_.end(),std::invalid_argument,
253  "STK_Interface::addMeshCoordFields: Can't set more than one set of coordinate "
254  "fields for element block \""+blockId+"\".");
255 
256  // Note that there is a distinction between the key which is used for lookups
257  // and the field that lives on the mesh, which is used for printing the displacement.
258 
259  // just copy the coordinate names
260  meshCoordFields_[blockId] = coordNames;
261 
262  // must fill in the displacement fields
263  std::vector<std::string> & dispFields = meshDispFields_[blockId];
264  dispFields.resize(dimension_);
265 
266  for(unsigned i=0;i<dimension_;i++) {
267  std::pair<std::string,std::string> key = std::make_pair(coordNames[i],blockId);
268  std::string dispName = dispPrefix+coordNames[i];
269 
270  dispFields[i] = dispName; // record this field as a
271  // displacement field
272 
273  // add & declare field if not already added...currently assuming linears
274  if(fieldNameToSolution_.find(key)==fieldNameToSolution_.end()) {
275 
276  SolutionFieldType * field = metaData_->get_field<SolutionFieldType>(stk::topology::NODE_RANK, dispName);
277  if(field==0) {
278  field = &metaData_->declare_field<SolutionFieldType>(stk::topology::NODE_RANK, dispName);
279  }
281  }
282  }
283 }
284 
285 void STK_Interface::addInformationRecords(const std::vector<std::string> & info_records)
286 {
287  informationRecords_.insert(info_records.begin(), info_records.end());
288 }
289 
290 void STK_Interface::initialize(stk::ParallelMachine parallelMach,bool setupIO,
291  const bool buildRefinementSupport)
292 {
293  TEUCHOS_ASSERT(not initialized_);
294  TEUCHOS_ASSERT(dimension_!=0); // no zero dimensional meshes!
295 
296  if(mpiComm_==Teuchos::null)
297  mpiComm_ = getSafeCommunicator(parallelMach);
298 
299  procRank_ = stk::parallel_machine_rank(*mpiComm_->getRawMpiComm());
300 
301  // associating the field with a part: universal part!
302  stk::mesh::FieldTraits<VectorFieldType>::data_type* init_vf = nullptr; // gcc 4.8 hack
303  stk::mesh::FieldTraits<ProcIdFieldType>::data_type* init_pid = nullptr; // gcc 4.8 hack
304  stk::mesh::FieldTraits<SolutionFieldType>::data_type* init_sol = nullptr; // gcc 4.8 hack
305  stk::mesh::put_field_on_mesh( *coordinatesField_ , metaData_->universal_part(), getDimension(),init_vf);
306  stk::mesh::put_field_on_mesh( *edgesField_ , metaData_->universal_part(), getDimension(),init_vf);
307  if (dimension_ > 2)
308  stk::mesh::put_field_on_mesh( *facesField_ , metaData_->universal_part(), getDimension(),init_vf);
309  stk::mesh::put_field_on_mesh( *processorIdField_ , metaData_->universal_part(),init_pid);
310  stk::mesh::put_field_on_mesh( *loadBalField_ , metaData_->universal_part(),init_sol);
311 
316 
317 #ifdef PANZER_HAVE_IOSS
318  if(setupIO) {
319  // setup Exodus file IO
321 
322  // add element blocks
323  {
324  std::map<std::string, stk::mesh::Part*>::iterator itr;
325  for(itr=elementBlocks_.begin();itr!=elementBlocks_.end();++itr)
326  if(!stk::io::is_part_io_part(*itr->second))
327  stk::io::put_io_part_attribute(*itr->second); // this can only be called once per part
328  }
329 
330  // add edge blocks
331  {
332  std::map<std::string, stk::mesh::Part*>::iterator itr;
333  for(itr=edgeBlocks_.begin();itr!=edgeBlocks_.end();++itr)
334  if(!stk::io::is_part_edge_block_io_part(*itr->second))
335  stk::io::put_edge_block_io_part_attribute(*itr->second); // this can only be called once per part
336  }
337 
338  // add face blocks
339  {
340  std::map<std::string, stk::mesh::Part*>::iterator itr;
341  for(itr=faceBlocks_.begin();itr!=faceBlocks_.end();++itr)
342  if(!stk::io::is_part_face_block_io_part(*itr->second))
343  stk::io::put_face_block_io_part_attribute(*itr->second); // this can only be called once per part
344  }
345 
346  // add side sets
347  {
348  std::map<std::string, stk::mesh::Part*>::iterator itr;
349  for(itr=sidesets_.begin();itr!=sidesets_.end();++itr)
350  if(!stk::io::is_part_io_part(*itr->second))
351  stk::io::put_io_part_attribute(*itr->second); // this can only be called once per part
352  }
353 
354  // add node sets
355  {
356  std::map<std::string, stk::mesh::Part*>::iterator itr;
357  for(itr=nodesets_.begin();itr!=nodesets_.end();++itr)
358  if(!stk::io::is_part_io_part(*itr->second))
359  stk::io::put_io_part_attribute(*itr->second); // this can only be called once per part
360  }
361 
362  // add nodes
363  if(!stk::io::is_part_io_part(*nodesPart_))
364  stk::io::put_io_part_attribute(*nodesPart_);
365 
366  stk::io::set_field_role(*coordinatesField_, Ioss::Field::MESH);
367  stk::io::set_field_role(*edgesField_, Ioss::Field::MESH);
368  if (dimension_ > 2)
369  stk::io::set_field_role(*facesField_, Ioss::Field::MESH);
370  stk::io::set_field_role(*processorIdField_, Ioss::Field::TRANSIENT);
371  // stk::io::set_field_role(*loadBalField_, Ioss::Field::TRANSIENT);
372  }
373 #endif
374 
375  if (buildRefinementSupport) {
376 #ifdef PANZER_HAVE_PERCEPT
377  refinedMesh_ = Teuchos::rcp(new percept::PerceptMesh(this->getMetaData().get(),
378  this->getBulkData().get(),
379  true));
380 
381  breakPattern_ = Teuchos::rcp(new percept::URP_Heterogeneous_3D(*refinedMesh_));
382 #else
383  TEUCHOS_TEST_FOR_EXCEPTION(true,std::logic_error,
384  "ERROR: Uniform refinement requested. This requires the Percept package to be enabled in Trilinos!");
385 #endif
386  }
387 
388  if(bulkData_==Teuchos::null)
389  instantiateBulkData(*mpiComm_->getRawMpiComm());
390 
391  metaData_->commit();
392 
393  initialized_ = true;
394 }
395 
396 void STK_Interface::initializeFieldsInSTK(const std::map<std::pair<std::string,std::string>,SolutionFieldType*> & nameToField,
397  bool setupIO)
398 {
399  std::set<SolutionFieldType*> uniqueFields;
400  std::map<std::pair<std::string,std::string>,SolutionFieldType*>::const_iterator fieldIter;
401  for(fieldIter=nameToField.begin();fieldIter!=nameToField.end();++fieldIter)
402  uniqueFields.insert(fieldIter->second); // this makes setting up IO easier!
403 
404  {
405  std::set<SolutionFieldType*>::const_iterator uniqueFieldIter;
406  stk::mesh::FieldTraits<SolutionFieldType>::data_type* init_sol = nullptr; // gcc 4.8 hack
407  for(uniqueFieldIter=uniqueFields.begin();uniqueFieldIter!=uniqueFields.end();++uniqueFieldIter)
408  stk::mesh::put_field_on_mesh(*(*uniqueFieldIter),metaData_->universal_part(),init_sol);
409  }
410 
411 #ifdef PANZER_HAVE_IOSS
412  if(setupIO) {
413  // add solution fields
414  std::set<SolutionFieldType*>::const_iterator uniqueFieldIter;
415  for(uniqueFieldIter=uniqueFields.begin();uniqueFieldIter!=uniqueFields.end();++uniqueFieldIter)
416  stk::io::set_field_role(*(*uniqueFieldIter), Ioss::Field::TRANSIENT);
417  }
418 #endif
419 }
420 
421 void STK_Interface::instantiateBulkData(stk::ParallelMachine parallelMach)
422 {
423  TEUCHOS_ASSERT(bulkData_==Teuchos::null);
424  if(mpiComm_==Teuchos::null)
425  mpiComm_ = getSafeCommunicator(parallelMach);
426 
427  std::unique_ptr<stk::mesh::BulkData> bulkUPtr = stk::mesh::MeshBuilder(*mpiComm_->getRawMpiComm()).create(Teuchos::get_shared_ptr(metaData_));
428  bulkData_ = rcp(bulkUPtr.release());
429 }
430 
432 {
433  TEUCHOS_TEST_FOR_EXCEPTION(bulkData_==Teuchos::null,std::logic_error,
434  "STK_Interface: Must call \"initialized\" or \"instantiateBulkData\" before \"beginModification\"");
435 
436  bulkData_->modification_begin();
437 }
438 
440 {
441  TEUCHOS_TEST_FOR_EXCEPTION(bulkData_==Teuchos::null,std::logic_error,
442  "STK_Interface: Must call \"initialized\" or \"instantiateBulkData\" before \"endModification\"");
443 
444  // TODO: Resolving sharing this way comes at a cost in performance. The STK
445  // team has decided that users need to declare their own sharing. We should
446  // find where shared entities are being created in Panzer and declare it.
447  // Once this is done, the extra code below can be deleted.
448 
449  stk::CommSparse comm(bulkData_->parallel());
450 
451  for (int phase=0;phase<2;++phase) {
452  for (int i=0;i<bulkData_->parallel_size();++i) {
453  if ( i != bulkData_->parallel_rank() ) {
454  const stk::mesh::BucketVector& buckets = bulkData_->buckets(stk::topology::NODE_RANK);
455  for (size_t j=0;j<buckets.size();++j) {
456  const stk::mesh::Bucket& bucket = *buckets[j];
457  if ( bucket.owned() ) {
458  for (size_t k=0;k<bucket.size();++k) {
459  stk::mesh::EntityKey key = bulkData_->entity_key(bucket[k]);
460  comm.send_buffer(i).pack<stk::mesh::EntityKey>(key);
461  }
462  }
463  }
464  }
465  }
466 
467  if (phase == 0 ) {
468  comm.allocate_buffers();
469  }
470  else {
471  comm.communicate();
472  }
473  }
474 
475  for (int i=0;i<bulkData_->parallel_size();++i) {
476  if ( i != bulkData_->parallel_rank() ) {
477  while(comm.recv_buffer(i).remaining()) {
478  stk::mesh::EntityKey key;
479  comm.recv_buffer(i).unpack<stk::mesh::EntityKey>(key);
480  stk::mesh::Entity node = bulkData_->get_entity(key);
481  if ( bulkData_->is_valid(node) ) {
482  bulkData_->add_node_sharing(node, i);
483  }
484  }
485  }
486  }
487 
488 
489  bulkData_->modification_end();
490 
493 }
494 
495 void STK_Interface::addNode(stk::mesh::EntityId gid, const std::vector<double> & coord)
496 {
497  TEUCHOS_TEST_FOR_EXCEPTION(not isModifiable(),std::logic_error,
498  "STK_Interface::addNode: STK_Interface must be modifiable to add a node");
499  TEUCHOS_TEST_FOR_EXCEPTION(coord.size()!=getDimension(),std::logic_error,
500  "STK_Interface::addNode: number of coordinates in vector must mation dimension");
501  TEUCHOS_TEST_FOR_EXCEPTION(gid==0,std::logic_error,
502  "STK_Interface::addNode: STK has STUPID restriction of no zero GIDs, pick something else");
503  stk::mesh::EntityRank nodeRank = getNodeRank();
504 
505  stk::mesh::Entity node = bulkData_->declare_entity(nodeRank,gid,nodesPartVec_);
506 
507  // set coordinate vector
508  double * fieldCoords = stk::mesh::field_data(*coordinatesField_,node);
509  for(std::size_t i=0;i<coord.size();++i)
510  fieldCoords[i] = coord[i];
511 }
512 
513 void STK_Interface::addEntityToSideset(stk::mesh::Entity entity,stk::mesh::Part * sideset)
514 {
515  std::vector<stk::mesh::Part*> sidesetV;
516  sidesetV.push_back(sideset);
517 
518  bulkData_->change_entity_parts(entity,sidesetV);
519 }
520 
521 void STK_Interface::addEntityToNodeset(stk::mesh::Entity entity,stk::mesh::Part * nodeset)
522 {
523  std::vector<stk::mesh::Part*> nodesetV;
524  nodesetV.push_back(nodeset);
525 
526  bulkData_->change_entity_parts(entity,nodesetV);
527 }
528 
529 void STK_Interface::addEntityToEdgeBlock(stk::mesh::Entity entity,stk::mesh::Part * edgeblock)
530 {
531  std::vector<stk::mesh::Part*> edgeblockV;
532  edgeblockV.push_back(edgeblock);
533 
534  bulkData_->change_entity_parts(entity,edgeblockV);
535 }
536 void STK_Interface::addEntitiesToEdgeBlock(std::vector<stk::mesh::Entity> entities,stk::mesh::Part * edgeblock)
537 {
538  if (entities.size() > 0) {
539  std::vector<stk::mesh::Part*> edgeblockV;
540  edgeblockV.push_back(edgeblock);
541 
542  bulkData_->change_entity_parts(entities,edgeblockV);
543  }
544 }
545 
546 void STK_Interface::addEntityToFaceBlock(stk::mesh::Entity entity,stk::mesh::Part * faceblock)
547 {
548  std::vector<stk::mesh::Part*> faceblockV;
549  faceblockV.push_back(faceblock);
550 
551  bulkData_->change_entity_parts(entity,faceblockV);
552 }
553 void STK_Interface::addEntitiesToFaceBlock(std::vector<stk::mesh::Entity> entities,stk::mesh::Part * faceblock)
554 {
555  if (entities.size() > 0) {
556  std::vector<stk::mesh::Part*> faceblockV;
557  faceblockV.push_back(faceblock);
558 
559  bulkData_->change_entity_parts(entities,faceblockV);
560  }
561 }
562 
563 void STK_Interface::addElement(const Teuchos::RCP<ElementDescriptor> & ed,stk::mesh::Part * block)
564 {
565  std::vector<stk::mesh::Part*> blockVec;
566  blockVec.push_back(block);
567 
568  stk::mesh::EntityRank elementRank = getElementRank();
569  stk::mesh::EntityRank nodeRank = getNodeRank();
570  stk::mesh::Entity element = bulkData_->declare_entity(elementRank,ed->getGID(),blockVec);
571 
572  // build relations that give the mesh structure
573  const std::vector<stk::mesh::EntityId> & nodes = ed->getNodes();
574  for(std::size_t i=0;i<nodes.size();++i) {
575  // add element->node relation
576  stk::mesh::Entity node = bulkData_->get_entity(nodeRank,nodes[i]);
577  TEUCHOS_ASSERT(bulkData_->is_valid(node));
578  bulkData_->declare_relation(element,node,i);
579  }
580 
581  ProcIdData * procId = stk::mesh::field_data(*processorIdField_,element);
582  procId[0] = Teuchos::as<ProcIdData>(procRank_);
583 }
584 
585 
587 {
588  // loop over elements
589  stk::mesh::EntityRank edgeRank = getEdgeRank();
590  std::vector<stk::mesh::Entity> localElmts;
591  getMyElements(localElmts);
592  std::vector<stk::mesh::Entity>::const_iterator itr;
593  for(itr=localElmts.begin();itr!=localElmts.end();++itr) {
594  stk::mesh::Entity element = (*itr);
595  stk::mesh::EntityId gid = bulkData_->identifier(element);
596  std::vector<stk::mesh::EntityId> subcellIds;
597  getSubcellIndices(edgeRank,gid,subcellIds);
598 
599  for(std::size_t i=0;i<subcellIds.size();++i) {
600  stk::mesh::Entity edge = bulkData_->get_entity(edgeRank,subcellIds[i]);
601  stk::mesh::Entity const* relations = bulkData_->begin_nodes(edge);
602 
603  double * node_coord_1 = stk::mesh::field_data(*coordinatesField_,relations[0]);
604  double * node_coord_2 = stk::mesh::field_data(*coordinatesField_,relations[1]);
605 
606  // set coordinate vector
607  double * edgeCoords = stk::mesh::field_data(*edgesField_,edge);
608  for(std::size_t j=0;j<getDimension();++j)
609  edgeCoords[j] = (node_coord_1[j]+node_coord_2[j])/2.0;
610  }
611  }
612 }
613 
615 {
616  // loop over elements
617  stk::mesh::EntityRank faceRank = getFaceRank();
618  std::vector<stk::mesh::Entity> localElmts;
619  getMyElements(localElmts);
620  std::vector<stk::mesh::Entity>::const_iterator itr;
621  for(itr=localElmts.begin();itr!=localElmts.end();++itr) {
622  stk::mesh::Entity element = (*itr);
623  stk::mesh::EntityId gid = elementGlobalId(element);
624  std::vector<stk::mesh::EntityId> subcellIds;
625  getSubcellIndices(faceRank,gid,subcellIds);
626 
627  for(std::size_t i=0;i<subcellIds.size();++i) {
628  stk::mesh::Entity face = bulkData_->get_entity(faceRank,subcellIds[i]);
629  stk::mesh::Entity const* relations = bulkData_->begin_nodes(face);
630  const size_t num_relations = bulkData_->num_nodes(face);
631 
632  // set coordinate vector
633  double * faceCoords = stk::mesh::field_data(*facesField_,face);
634  for(std::size_t j=0;j<getDimension();++j){
635  faceCoords[j] = 0.0;
636  for(std::size_t k=0;k<num_relations;++k)
637  faceCoords[j] += stk::mesh::field_data(*coordinatesField_,relations[k])[j];
638  faceCoords[j] /= double(num_relations);
639  }
640  }
641  }
642 }
643 
644 stk::mesh::Entity STK_Interface::findConnectivityById(stk::mesh::Entity src, stk::mesh::EntityRank tgt_rank, unsigned rel_id) const
645 {
646  const size_t num_rels = bulkData_->num_connectivity(src, tgt_rank);
647  stk::mesh::Entity const* relations = bulkData_->begin(src, tgt_rank);
648  stk::mesh::ConnectivityOrdinal const* ordinals = bulkData_->begin_ordinals(src, tgt_rank);
649  for (size_t i = 0; i < num_rels; ++i) {
650  if (ordinals[i] == static_cast<stk::mesh::ConnectivityOrdinal>(rel_id)) {
651  return relations[i];
652  }
653  }
654 
655  return stk::mesh::Entity();
656 }
657 
659 //
660 // writeToExodus()
661 //
663 void
665 writeToExodus(const std::string& filename,
666  const bool append)
667 {
668  setupExodusFile(filename,append);
669  writeToExodus(0.0);
670 } // end of writeToExodus()
671 
673 //
674 // setupExodusFile()
675 //
677 void
679 setupExodusFile(const std::string& filename,
680  const bool append,
681  const bool append_after_restart_time,
682  const double restart_time)
683 {
684  using std::runtime_error;
685  using stk::io::StkMeshIoBroker;
686  using stk::mesh::FieldVector;
687  using stk::ParallelMachine;
688  using Teuchos::rcp;
689  PANZER_FUNC_TIME_MONITOR("STK_Interface::setupExodusFile(filename)");
690 #ifdef PANZER_HAVE_IOSS
691  TEUCHOS_ASSERT(not mpiComm_.is_null())
692 
693  ParallelMachine comm = *mpiComm_->getRawMpiComm();
694  meshData_ = rcp(new StkMeshIoBroker(comm));
695  meshData_->set_bulk_data(Teuchos::get_shared_ptr(bulkData_));
696  Ioss::PropertyManager props;
697  props.add(Ioss::Property("LOWER_CASE_VARIABLE_NAMES", "FALSE"));
698  if (append) {
699  if (append_after_restart_time) {
700  meshIndex_ = meshData_->create_output_mesh(filename, stk::io::APPEND_RESULTS,
701  props, restart_time);
702  }
703  else // Append results to the end of the file
704  meshIndex_ = meshData_->create_output_mesh(filename, stk::io::APPEND_RESULTS, props);
705  }
706  else
707  meshIndex_ = meshData_->create_output_mesh(filename, stk::io::WRITE_RESULTS, props);
708  const FieldVector& fields = metaData_->get_fields();
709  for (size_t i(0); i < fields.size(); ++i) {
710  // Do NOT add MESH type stk fields to exodus io, but do add everything
711  // else. This allows us to avoid having to catch a throw for
712  // re-registering coordinates, sidesets, etc... Note that some
713  // fields like LOAD_BAL don't always have a role assigned, so for
714  // roles that point to null, we need to register them as well.
715  auto role = stk::io::get_field_role(*fields[i]);
716  if (role != nullptr) {
717  if (*role != Ioss::Field::MESH)
718  meshData_->add_field(meshIndex_, *fields[i]);
719  } else {
720  meshData_->add_field(meshIndex_, *fields[i]);
721  }
722  }
723 
724  // convert the set to a vector
725  std::vector<std::string> deduped_info_records(informationRecords_.begin(),informationRecords_.end());
726  meshData_->add_info_records(meshIndex_, deduped_info_records);
727 #else
728  TEUCHOS_ASSERT(false)
729 #endif
730 } // end of setupExodusFile()
731 
733 //
734 // writeToExodus()
735 //
737 void
740  double timestep)
741 {
742  using std::cout;
743  using std::endl;
744  using Teuchos::FancyOStream;
745  using Teuchos::rcpFromRef;
746  PANZER_FUNC_TIME_MONITOR("STK_Interface::writeToExodus(timestep)");
747 #ifdef PANZER_HAVE_IOSS
748  if (not meshData_.is_null())
749  {
750  currentStateTime_ = timestep;
751  globalToExodus(GlobalVariable::ADD);
752  meshData_->begin_output_step(meshIndex_, currentStateTime_);
753  meshData_->write_defined_output_fields(meshIndex_);
754  globalToExodus(GlobalVariable::WRITE);
755  meshData_->end_output_step(meshIndex_);
756  }
757  else // if (meshData_.is_null())
758  {
759  FancyOStream out(rcpFromRef(cout));
760  out.setOutputToRootOnly(0);
761  out << "WARNING: Exodus I/O has been disabled or not setup properly; "
762  << "not writing to Exodus." << endl;
763  } // end if (meshData_.is_null()) or not
764 #else
765  TEUCHOS_ASSERT(false);
766 #endif
767 } // end of writeToExodus()
768 
770 //
771 // globalToExodus()
772 //
774 void
775 STK_Interface::
776 globalToExodus(
777  const GlobalVariable& flag)
778 {
779  using std::invalid_argument;
780  using std::string;
781  using Teuchos::Array;
782 
783  // Loop over all the global variables to be added to the Exodus output file.
784  // For each global variable, we determine the data type, and then add or
785  // write it accordingly, depending on the value of flag.
786  for (auto i = globalData_.begin(); i != globalData_.end(); ++i)
787  {
788  const string& name = globalData_.name(i);
789 
790  // Integers.
791  if (globalData_.isType<int>(name))
792  {
793  const auto& value = globalData_.get<int>(name);
794  if (flag == GlobalVariable::ADD)
795  {
796  try
797  {
798  meshData_->add_global(meshIndex_, name, value,
799  stk::util::ParameterType::INTEGER);
800  }
801  catch (...)
802  {
803  return;
804  }
805  }
806  else // if (flag == GlobalVariable::WRITE)
807  meshData_->write_global(meshIndex_, name, value,
808  stk::util::ParameterType::INTEGER);
809  }
810 
811  // Doubles.
812  else if (globalData_.isType<double>(name))
813  {
814  const auto& value = globalData_.get<double>(name);
815  if (flag == GlobalVariable::ADD)
816  {
817  try
818  {
819  meshData_->add_global(meshIndex_, name, value,
820  stk::util::ParameterType::DOUBLE);
821  }
822  catch (...)
823  {
824  return;
825  }
826  }
827  else // if (flag == GlobalVariable::WRITE)
828  meshData_->write_global(meshIndex_, name, value,
829  stk::util::ParameterType::DOUBLE);
830  }
831 
832  // Vectors of integers.
833  else if (globalData_.isType<Array<int>>(name))
834  {
835  const auto& value = globalData_.get<Array<int>>(name).toVector();
836  if (flag == GlobalVariable::ADD)
837  {
838  try
839  {
840  meshData_->add_global(meshIndex_, name, value,
841  stk::util::ParameterType::INTEGERVECTOR);
842  }
843  catch (...)
844  {
845  return;
846  }
847  }
848  else // if (flag == GlobalVariable::WRITE)
849  meshData_->write_global(meshIndex_, name, value,
850  stk::util::ParameterType::INTEGERVECTOR);
851  }
852 
853  // Vectors of doubles.
854  else if (globalData_.isType<Array<double>>(name))
855  {
856  const auto& value = globalData_.get<Array<double>>(name).toVector();
857  if (flag == GlobalVariable::ADD)
858  {
859  try
860  {
861  meshData_->add_global(meshIndex_, name, value,
862  stk::util::ParameterType::DOUBLEVECTOR);
863  }
864  catch (...)
865  {
866  return;
867  }
868  }
869  else // if (flag == GlobalVariable::WRITE)
870  meshData_->write_global(meshIndex_, name, value,
871  stk::util::ParameterType::DOUBLEVECTOR);
872  }
873 
874  // If the data type is something else, throw an exception.
875  else
876  TEUCHOS_TEST_FOR_EXCEPTION(true, invalid_argument,
877  "STK_Interface::globalToExodus(): The global variable to be added " \
878  "to the Exodus output file is of an invalid type. Valid types are " \
879  "int and double, along with std::vectors of those types.")
880  } // end loop over globalData_
881 } // end of globalToExodus()
882 
884 //
885 // addGlobalToExodus()
886 //
888 void
891  const std::string& key,
892  const int& value)
893 {
894  globalData_.set(key, value);
895 } // end of addGlobalToExodus()
896 
898 //
899 // addGlobalToExodus()
900 //
902 void
905  const std::string& key,
906  const double& value)
907 {
908  globalData_.set(key, value);
909 } // end of addGlobalToExodus()
910 
912 //
913 // addGlobalToExodus()
914 //
916 void
919  const std::string& key,
920  const std::vector<int>& value)
921 {
922  using Teuchos::Array;
923  globalData_.set(key, Array<int>(value));
924 } // end of addGlobalToExodus()
925 
927 //
928 // addGlobalToExodus()
929 //
931 void
934  const std::string& key,
935  const std::vector<double>& value)
936 {
937  using Teuchos::Array;
938  globalData_.set(key, Array<double>(value));
939 } // end of addGlobalToExodus()
940 
942 {
943  #ifdef PANZER_HAVE_IOSS
944  return true;
945  #else
946  return false;
947  #endif
948 }
949 
950 void STK_Interface::getElementsSharingNode(stk::mesh::EntityId nodeId,std::vector<stk::mesh::Entity> & elements) const
951 {
952  stk::mesh::EntityRank nodeRank = getNodeRank();
953 
954  // get all relations for node
955  stk::mesh::Entity node = bulkData_->get_entity(nodeRank,nodeId);
956  const size_t numElements = bulkData_->num_elements(node);
957  stk::mesh::Entity const* relations = bulkData_->begin_elements(node);
958 
959  // extract elements sharing nodes
960  elements.insert(elements.end(), relations, relations + numElements);
961 }
962 
963 void STK_Interface::getOwnedElementsSharingNode(stk::mesh::Entity node,std::vector<stk::mesh::Entity> & elements,
964  std::vector<int> & relIds) const
965 {
966  // get all relations for node
967  const size_t numElements = bulkData_->num_elements(node);
968  stk::mesh::Entity const* relations = bulkData_->begin_elements(node);
969  stk::mesh::ConnectivityOrdinal const* rel_ids = bulkData_->begin_element_ordinals(node);
970 
971  // extract elements sharing nodes
972  for (size_t i = 0; i < numElements; ++i) {
973  stk::mesh::Entity element = relations[i];
974 
975  // if owned by this processor
976  if(bulkData_->parallel_owner_rank(element) == static_cast<int>(procRank_)) {
977  elements.push_back(element);
978  relIds.push_back(rel_ids[i]);
979  }
980  }
981 }
982 
983 void STK_Interface::getOwnedElementsSharingNode(stk::mesh::EntityId nodeId,std::vector<stk::mesh::Entity> & elements,
984  std::vector<int> & relIds, unsigned int matchType) const
985 {
986  stk::mesh::EntityRank rank;
987  if(matchType == 0)
988  rank = getNodeRank();
989  else if(matchType == 1)
990  rank = getEdgeRank();
991  else if(matchType == 2)
992  rank = getFaceRank();
993  else
994  TEUCHOS_ASSERT(false);
995 
996  stk::mesh::Entity node = bulkData_->get_entity(rank,nodeId);
997 
998  getOwnedElementsSharingNode(node,elements,relIds);
999 }
1000 
1001 void STK_Interface::getElementsSharingNodes(const std::vector<stk::mesh::EntityId> nodeIds,std::vector<stk::mesh::Entity> & elements) const
1002 {
1003  std::vector<stk::mesh::Entity> current;
1004 
1005  getElementsSharingNode(nodeIds[0],current); // fill it with elements touching first node
1006  std::sort(current.begin(),current.end()); // sort for intersection on the pointer
1007 
1008  // find intersection with remaining nodes
1009  for(std::size_t n=1;n<nodeIds.size();++n) {
1010  // get elements associated with next node
1011  std::vector<stk::mesh::Entity> nextNode;
1012  getElementsSharingNode(nodeIds[n],nextNode); // fill it with elements touching first node
1013  std::sort(nextNode.begin(),nextNode.end()); // sort for intersection on the pointer ID
1014 
1015  // intersect next node elements with current element list
1016  std::vector<stk::mesh::Entity> intersection(std::min(nextNode.size(),current.size()));
1017  std::vector<stk::mesh::Entity>::const_iterator endItr
1018  = std::set_intersection(current.begin(),current.end(),
1019  nextNode.begin(),nextNode.end(),
1020  intersection.begin());
1021  std::size_t newLength = endItr-intersection.begin();
1022  intersection.resize(newLength);
1023 
1024  // store intersection
1025  current.clear();
1026  current = intersection;
1027  }
1028 
1029  // return the elements computed
1030  elements = current;
1031 }
1032 
1033 void STK_Interface::getNodeIdsForElement(stk::mesh::Entity element,std::vector<stk::mesh::EntityId> & nodeIds) const
1034 {
1035  stk::mesh::Entity const* nodeRel = getBulkData()->begin_nodes(element);
1036  const size_t numNodes = getBulkData()->num_nodes(element);
1037 
1038  nodeIds.reserve(numNodes);
1039  for(size_t i = 0; i < numNodes; ++i) {
1040  nodeIds.push_back(elementGlobalId(nodeRel[i]));
1041  }
1042 }
1043 
1045 {
1046  entityCounts_.clear();
1047  stk::mesh::comm_mesh_counts(*bulkData_,entityCounts_);
1048 }
1049 
1051 {
1052  // developed to mirror "comm_mesh_counts" in stk_mesh/base/Comm.cpp
1053 
1054  const auto entityRankCount = metaData_->entity_rank_count();
1055  const size_t commCount = 10; // entityRankCount
1056 
1057  TEUCHOS_ASSERT(entityRankCount<10);
1058 
1059  // stk::ParallelMachine mach = bulkData_->parallel();
1060  stk::ParallelMachine mach = *mpiComm_->getRawMpiComm();
1061 
1062  std::vector<stk::mesh::EntityId> local(commCount,0);
1063 
1064  // determine maximum ID for this processor for each entity type
1065  stk::mesh::Selector ownedPart = metaData_->locally_owned_part();
1066  for(stk::mesh::EntityRank i=stk::topology::NODE_RANK;
1067  i < static_cast<stk::mesh::EntityRank>(entityRankCount); ++i) {
1068  std::vector<stk::mesh::Entity> entities;
1069 
1070  stk::mesh::get_selected_entities(ownedPart,bulkData_->buckets(i),entities);
1071 
1072  // determine maximum ID for this processor
1073  std::vector<stk::mesh::Entity>::const_iterator itr;
1074  for(itr=entities.begin();itr!=entities.end();++itr) {
1075  stk::mesh::EntityId id = bulkData_->identifier(*itr);
1076  if(id>local[i])
1077  local[i] = id;
1078  }
1079  }
1080 
1081  // get largest IDs across processors
1082  stk::all_reduce(mach,stk::ReduceMax<10>(&local[0]));
1083  maxEntityId_.assign(local.begin(),local.begin()+entityRankCount+1);
1084 }
1085 
1086 std::size_t STK_Interface::getEntityCounts(unsigned entityRank) const
1087 {
1088  TEUCHOS_TEST_FOR_EXCEPTION(entityRank>=entityCounts_.size(),std::logic_error,
1089  "STK_Interface::getEntityCounts: Entity counts do not include rank: " << entityRank);
1090 
1091  return entityCounts_[entityRank];
1092 }
1093 
1094 stk::mesh::EntityId STK_Interface::getMaxEntityId(unsigned entityRank) const
1095 {
1096  TEUCHOS_TEST_FOR_EXCEPTION(entityRank>=maxEntityId_.size(),std::logic_error,
1097  "STK_Interface::getMaxEntityId: Max entity ids do not include rank: " << entityRank);
1098 
1099  return maxEntityId_[entityRank];
1100 }
1101 
1103 {
1104  stk::mesh::PartVector emptyPartVector;
1105  stk::mesh::create_adjacent_entities(*bulkData_,emptyPartVector);
1106 
1109 
1110  addEdges();
1111  if (dimension_ > 2)
1112  addFaces();
1113 }
1114 
1115 const double * STK_Interface::getNodeCoordinates(stk::mesh::EntityId nodeId) const
1116 {
1117  stk::mesh::Entity node = bulkData_->get_entity(getNodeRank(),nodeId);
1118  return stk::mesh::field_data(*coordinatesField_,node);
1119 }
1120 
1121 const double * STK_Interface::getNodeCoordinates(stk::mesh::Entity node) const
1122 {
1123  return stk::mesh::field_data(*coordinatesField_,node);
1124 }
1125 
1126 void STK_Interface::getSubcellIndices(unsigned entityRank,stk::mesh::EntityId elementId,
1127  std::vector<stk::mesh::EntityId> & subcellIds) const
1128 {
1129  stk::mesh::EntityRank elementRank = getElementRank();
1130  stk::mesh::Entity cell = bulkData_->get_entity(elementRank,elementId);
1131 
1132  TEUCHOS_TEST_FOR_EXCEPTION(!bulkData_->is_valid(cell),std::logic_error,
1133  "STK_Interface::getSubcellIndices: could not find element requested (GID = " << elementId << ")");
1134 
1135  const size_t numSubcells = bulkData_->num_connectivity(cell, static_cast<stk::mesh::EntityRank>(entityRank));
1136  stk::mesh::Entity const* subcells = bulkData_->begin(cell, static_cast<stk::mesh::EntityRank>(entityRank));
1137  subcellIds.clear();
1138  subcellIds.resize(numSubcells,0);
1139 
1140  // loop over relations and fill subcell vector
1141  for(size_t i = 0; i < numSubcells; ++i) {
1142  stk::mesh::Entity subcell = subcells[i];
1143  subcellIds[i] = bulkData_->identifier(subcell);
1144  }
1145 }
1146 
1147 void STK_Interface::getMyElements(std::vector<stk::mesh::Entity> & elements) const
1148 {
1149  // setup local ownership
1150  stk::mesh::Selector ownedPart = metaData_->locally_owned_part();
1151 
1152  // grab elements
1153  stk::mesh::EntityRank elementRank = getElementRank();
1154  stk::mesh::get_selected_entities(ownedPart,bulkData_->buckets(elementRank),elements);
1155 }
1156 
1157 void STK_Interface::getMyElements(const std::string & blockID,std::vector<stk::mesh::Entity> & elements) const
1158 {
1159  stk::mesh::Part * elementBlock = getElementBlockPart(blockID);
1160 
1161  TEUCHOS_TEST_FOR_EXCEPTION(elementBlock==0,std::logic_error,"Could not find element block \"" << blockID << "\"");
1162 
1163  // setup local ownership
1164  // stk::mesh::Selector block = *elementBlock;
1165  stk::mesh::Selector ownedBlock = metaData_->locally_owned_part() & (*elementBlock);
1166 
1167  // grab elements
1168  stk::mesh::EntityRank elementRank = getElementRank();
1169  stk::mesh::get_selected_entities(ownedBlock,bulkData_->buckets(elementRank),elements);
1170 }
1171 
1172 void STK_Interface::getNeighborElements(std::vector<stk::mesh::Entity> & elements) const
1173 {
1174  // setup local ownership
1175  stk::mesh::Selector neighborBlock = (!metaData_->locally_owned_part());
1176 
1177  // grab elements
1178  stk::mesh::EntityRank elementRank = getElementRank();
1179  stk::mesh::get_selected_entities(neighborBlock,bulkData_->buckets(elementRank),elements);
1180 }
1181 
1182 void STK_Interface::getNeighborElements(const std::string & blockID,std::vector<stk::mesh::Entity> & elements) const
1183 {
1184  stk::mesh::Part * elementBlock = getElementBlockPart(blockID);
1185 
1186  TEUCHOS_TEST_FOR_EXCEPTION(elementBlock==0,std::logic_error,"Could not find element block \"" << blockID << "\"");
1187 
1188  // setup local ownership
1189  stk::mesh::Selector neighborBlock = (!metaData_->locally_owned_part()) & (*elementBlock);
1190 
1191  // grab elements
1192  stk::mesh::EntityRank elementRank = getElementRank();
1193  stk::mesh::get_selected_entities(neighborBlock,bulkData_->buckets(elementRank),elements);
1194 }
1195 
1196 void STK_Interface::getMyEdges(std::vector<stk::mesh::Entity> & edges) const
1197 {
1198  // setup local ownership
1199  stk::mesh::Selector ownedPart = metaData_->locally_owned_part();
1200 
1201  // grab elements
1202  stk::mesh::get_selected_entities(ownedPart,bulkData_->buckets(getEdgeRank()),edges);
1203 }
1204 
1205 void STK_Interface::getMyEdges(const std::string & edgeBlockName,std::vector<stk::mesh::Entity> & edges) const
1206 {
1207  stk::mesh::Part * edgeBlockPart = getEdgeBlock(edgeBlockName);
1208  TEUCHOS_TEST_FOR_EXCEPTION(edgeBlockPart==0,std::logic_error,
1209  "Unknown edge block \"" << edgeBlockName << "\"");
1210 
1211  stk::mesh::Selector edge_block = *edgeBlockPart;
1212  stk::mesh::Selector owned_block = metaData_->locally_owned_part() & edge_block;
1213 
1214  // grab elements
1215  stk::mesh::get_selected_entities(owned_block,bulkData_->buckets(getEdgeRank()),edges);
1216 }
1217 
1218 void STK_Interface::getMyEdges(const std::string & edgeBlockName,const std::string & blockName,std::vector<stk::mesh::Entity> & edges) const
1219 {
1220  stk::mesh::Part * edgeBlockPart = getEdgeBlock(edgeBlockName);
1221  stk::mesh::Part * elmtPart = getElementBlockPart(blockName);
1222  TEUCHOS_TEST_FOR_EXCEPTION(edgeBlockPart==0,EdgeBlockException,
1223  "Unknown edge block \"" << edgeBlockName << "\"");
1224  TEUCHOS_TEST_FOR_EXCEPTION(elmtPart==0,ElementBlockException,
1225  "Unknown element block \"" << blockName << "\"");
1226 
1227  stk::mesh::Selector edge_block = *edgeBlockPart;
1228  stk::mesh::Selector element_block = *elmtPart;
1229  stk::mesh::Selector owned_block = metaData_->locally_owned_part() & element_block & edge_block;
1230 
1231  // grab elements
1232  stk::mesh::get_selected_entities(owned_block,bulkData_->buckets(getEdgeRank()),edges);
1233 }
1234 
1235 void STK_Interface::getAllEdges(const std::string & edgeBlockName,std::vector<stk::mesh::Entity> & edges) const
1236 {
1237  stk::mesh::Part * edgeBlockPart = getEdgeBlock(edgeBlockName);
1238  TEUCHOS_TEST_FOR_EXCEPTION(edgeBlockPart==0,std::logic_error,
1239  "Unknown edge block \"" << edgeBlockName << "\"");
1240 
1241  stk::mesh::Selector edge_block = *edgeBlockPart;
1242 
1243  // grab elements
1244  stk::mesh::get_selected_entities(edge_block,bulkData_->buckets(getEdgeRank()),edges);
1245 }
1246 
1247 void STK_Interface::getAllEdges(const std::string & edgeBlockName,const std::string & blockName,std::vector<stk::mesh::Entity> & edges) const
1248 {
1249  stk::mesh::Part * edgeBlockPart = getEdgeBlock(edgeBlockName);
1250  stk::mesh::Part * elmtPart = getElementBlockPart(blockName);
1251  TEUCHOS_TEST_FOR_EXCEPTION(edgeBlockPart==0,EdgeBlockException,
1252  "Unknown edge block \"" << edgeBlockName << "\"");
1253  TEUCHOS_TEST_FOR_EXCEPTION(elmtPart==0,ElementBlockException,
1254  "Unknown element block \"" << blockName << "\"");
1255 
1256  stk::mesh::Selector edge_block = *edgeBlockPart;
1257  stk::mesh::Selector element_block = *elmtPart;
1258  stk::mesh::Selector element_edge_block = element_block & edge_block;
1259 
1260  // grab elements
1261  stk::mesh::get_selected_entities(element_edge_block,bulkData_->buckets(getEdgeRank()),edges);
1262 }
1263 
1264 void STK_Interface::getMyFaces(std::vector<stk::mesh::Entity> & faces) const
1265 {
1266  // setup local ownership
1267  stk::mesh::Selector ownedPart = metaData_->locally_owned_part();
1268 
1269  // grab elements
1270  stk::mesh::EntityRank faceRank = getFaceRank();
1271  stk::mesh::get_selected_entities(ownedPart,bulkData_->buckets(faceRank),faces);
1272 }
1273 
1274 void STK_Interface::getMyFaces(const std::string & faceBlockName,std::vector<stk::mesh::Entity> & faces) const
1275 {
1276  stk::mesh::Part * faceBlockPart = getFaceBlock(faceBlockName);
1277  TEUCHOS_TEST_FOR_EXCEPTION(faceBlockPart==0,std::logic_error,
1278  "Unknown face block \"" << faceBlockName << "\"");
1279 
1280  stk::mesh::Selector face_block = *faceBlockPart;
1281  stk::mesh::Selector owned_block = metaData_->locally_owned_part() & face_block;
1282 
1283  // grab elements
1284  stk::mesh::get_selected_entities(owned_block,bulkData_->buckets(getFaceRank()),faces);
1285 }
1286 
1287 void STK_Interface::getMyFaces(const std::string & faceBlockName,const std::string & blockName,std::vector<stk::mesh::Entity> & faces) const
1288 {
1289  stk::mesh::Part * faceBlockPart = getFaceBlock(faceBlockName);
1290  stk::mesh::Part * elmtPart = getElementBlockPart(blockName);
1291  TEUCHOS_TEST_FOR_EXCEPTION(faceBlockPart==0,FaceBlockException,
1292  "Unknown face block \"" << faceBlockName << "\"");
1293  TEUCHOS_TEST_FOR_EXCEPTION(elmtPart==0,ElementBlockException,
1294  "Unknown element block \"" << blockName << "\"");
1295 
1296  stk::mesh::Selector face_block = *faceBlockPart;
1297  stk::mesh::Selector element_block = *elmtPart;
1298  stk::mesh::Selector owned_block = metaData_->locally_owned_part() & element_block & face_block;
1299 
1300  // grab elements
1301  stk::mesh::get_selected_entities(owned_block,bulkData_->buckets(getFaceRank()),faces);
1302 }
1303 
1304 void STK_Interface::getAllFaces(const std::string & faceBlockName,std::vector<stk::mesh::Entity> & faces) const
1305 {
1306  stk::mesh::Part * faceBlockPart = getFaceBlock(faceBlockName);
1307  TEUCHOS_TEST_FOR_EXCEPTION(faceBlockPart==0,std::logic_error,
1308  "Unknown face block \"" << faceBlockName << "\"");
1309 
1310  stk::mesh::Selector face_block = *faceBlockPart;
1311 
1312  // grab elements
1313  stk::mesh::get_selected_entities(face_block,bulkData_->buckets(getFaceRank()),faces);
1314 }
1315 
1316 void STK_Interface::getAllFaces(const std::string & faceBlockName,const std::string & blockName,std::vector<stk::mesh::Entity> & faces) const
1317 {
1318  stk::mesh::Part * faceBlockPart = getFaceBlock(faceBlockName);
1319  stk::mesh::Part * elmtPart = getElementBlockPart(blockName);
1320  TEUCHOS_TEST_FOR_EXCEPTION(faceBlockPart==0,FaceBlockException,
1321  "Unknown face block \"" << faceBlockName << "\"");
1322  TEUCHOS_TEST_FOR_EXCEPTION(elmtPart==0,ElementBlockException,
1323  "Unknown element block \"" << blockName << "\"");
1324 
1325  stk::mesh::Selector face_block = *faceBlockPart;
1326  stk::mesh::Selector element_block = *elmtPart;
1327  stk::mesh::Selector element_face_block = element_block & face_block;
1328 
1329  // grab elements
1330  stk::mesh::get_selected_entities(element_face_block,bulkData_->buckets(getFaceRank()),faces);
1331 }
1332 
1333 void STK_Interface::getMySides(const std::string & sideName,std::vector<stk::mesh::Entity> & sides) const
1334 {
1335  stk::mesh::Part * sidePart = getSideset(sideName);
1336  TEUCHOS_TEST_FOR_EXCEPTION(sidePart==0,std::logic_error,
1337  "Unknown side set \"" << sideName << "\"");
1338 
1339  stk::mesh::Selector side = *sidePart;
1340  stk::mesh::Selector ownedBlock = metaData_->locally_owned_part() & side;
1341 
1342  // grab elements
1343  stk::mesh::get_selected_entities(ownedBlock,bulkData_->buckets(getSideRank()),sides);
1344 }
1345 
1346 void STK_Interface::getMySides(const std::string & sideName,const std::string & blockName,std::vector<stk::mesh::Entity> & sides) const
1347 {
1348  stk::mesh::Part * sidePart = getSideset(sideName);
1349  stk::mesh::Part * elmtPart = getElementBlockPart(blockName);
1350  TEUCHOS_TEST_FOR_EXCEPTION(sidePart==0,SidesetException,
1351  "Unknown side set \"" << sideName << "\"");
1352  TEUCHOS_TEST_FOR_EXCEPTION(elmtPart==0,ElementBlockException,
1353  "Unknown element block \"" << blockName << "\"");
1354 
1355  stk::mesh::Selector side = *sidePart;
1356  stk::mesh::Selector block = *elmtPart;
1357  stk::mesh::Selector ownedBlock = metaData_->locally_owned_part() & block & side;
1358 
1359  // grab elements
1360  stk::mesh::get_selected_entities(ownedBlock,bulkData_->buckets(getSideRank()),sides);
1361 }
1362 
1363 void STK_Interface::getAllSides(const std::string & sideName,std::vector<stk::mesh::Entity> & sides) const
1364 {
1365  stk::mesh::Part * sidePart = getSideset(sideName);
1366  TEUCHOS_TEST_FOR_EXCEPTION(sidePart==0,std::logic_error,
1367  "Unknown side set \"" << sideName << "\"");
1368 
1369  stk::mesh::Selector side = *sidePart;
1370 
1371  // grab elements
1372  stk::mesh::get_selected_entities(side,bulkData_->buckets(getSideRank()),sides);
1373 }
1374 
1375 void STK_Interface::getAllSides(const std::string & sideName,const std::string & blockName,std::vector<stk::mesh::Entity> & sides) const
1376 {
1377  stk::mesh::Part * sidePart = getSideset(sideName);
1378  stk::mesh::Part * elmtPart = getElementBlockPart(blockName);
1379  TEUCHOS_TEST_FOR_EXCEPTION(sidePart==0,SidesetException,
1380  "Unknown side set \"" << sideName << "\"");
1381  TEUCHOS_TEST_FOR_EXCEPTION(elmtPart==0,ElementBlockException,
1382  "Unknown element block \"" << blockName << "\"");
1383 
1384  stk::mesh::Selector side = *sidePart;
1385  stk::mesh::Selector block = *elmtPart;
1386  stk::mesh::Selector sideBlock = block & side;
1387 
1388  // grab elements
1389  stk::mesh::get_selected_entities(sideBlock,bulkData_->buckets(getSideRank()),sides);
1390 }
1391 
1392 void STK_Interface::getMyNodes(const std::string & nodesetName,const std::string & blockName,std::vector<stk::mesh::Entity> & nodes) const
1393 {
1394  stk::mesh::Part * nodePart = getNodeset(nodesetName);
1395  stk::mesh::Part * elmtPart = getElementBlockPart(blockName);
1396  TEUCHOS_TEST_FOR_EXCEPTION(nodePart==0,SidesetException,
1397  "Unknown node set \"" << nodesetName << "\"");
1398  TEUCHOS_TEST_FOR_EXCEPTION(elmtPart==0,ElementBlockException,
1399  "Unknown element block \"" << blockName << "\"");
1400 
1401  stk::mesh::Selector nodeset = *nodePart;
1402  stk::mesh::Selector block = *elmtPart;
1403  stk::mesh::Selector ownedBlock = metaData_->locally_owned_part() & block & nodeset;
1404 
1405  // grab elements
1406  stk::mesh::get_selected_entities(ownedBlock,bulkData_->buckets(getNodeRank()),nodes);
1407 }
1408 
1409 void STK_Interface::getElementBlockNames(std::vector<std::string> & names) const
1410 {
1411  // TEUCHOS_ASSERT(initialized_); // all blocks must have been added
1412 
1413  names.clear();
1414 
1415  // fill vector with automagically ordered string values
1416  std::map<std::string, stk::mesh::Part*>::const_iterator blkItr; // Element blocks
1417  for(blkItr=elementBlocks_.begin();blkItr!=elementBlocks_.end();++blkItr)
1418  names.push_back(blkItr->first);
1419 }
1420 
1421 void STK_Interface::getSidesetNames(std::vector<std::string> & names) const
1422 {
1423  // TEUCHOS_ASSERT(initialized_); // all blocks must have been added
1424 
1425  names.clear();
1426 
1427  // fill vector with automagically ordered string values
1428  std::map<std::string, stk::mesh::Part*>::const_iterator sideItr; // Side sets
1429  for(sideItr=sidesets_.begin();sideItr!=sidesets_.end();++sideItr)
1430  names.push_back(sideItr->first);
1431 }
1432 
1433 void STK_Interface::getNodesetNames(std::vector<std::string> & names) const
1434 {
1435  names.clear();
1436 
1437  // fill vector with automagically ordered string values
1438  std::map<std::string, stk::mesh::Part*>::const_iterator nodeItr; // Node sets
1439  for(nodeItr=nodesets_.begin();nodeItr!=nodesets_.end();++nodeItr)
1440  names.push_back(nodeItr->first);
1441 }
1442 
1443 void STK_Interface::getEdgeBlockNames(std::vector<std::string> & names) const
1444 {
1445  names.clear();
1446 
1447  // fill vector with automagically ordered string values
1448  std::map<std::string, stk::mesh::Part*>::const_iterator edgeBlockItr; // Edge blocks
1449  for(edgeBlockItr=edgeBlocks_.begin();edgeBlockItr!=edgeBlocks_.end();++edgeBlockItr)
1450  names.push_back(edgeBlockItr->first);
1451 }
1452 
1453 void STK_Interface::getFaceBlockNames(std::vector<std::string> & names) const
1454 {
1455  names.clear();
1456 
1457  // fill vector with automagically ordered string values
1458  std::map<std::string, stk::mesh::Part*>::const_iterator faceBlockItr; // Face blocks
1459  for(faceBlockItr=faceBlocks_.begin();faceBlockItr!=faceBlocks_.end();++faceBlockItr)
1460  names.push_back(faceBlockItr->first);
1461 }
1462 
1463 std::size_t STK_Interface::elementLocalId(stk::mesh::Entity elmt) const
1464 {
1465  return elementLocalId(bulkData_->identifier(elmt));
1466  // const std::size_t * fieldCoords = stk::mesh::field_data(*localIdField_,*elmt);
1467  // return fieldCoords[0];
1468 }
1469 
1470 std::size_t STK_Interface::elementLocalId(stk::mesh::EntityId gid) const
1471 {
1472  // stk::mesh::EntityRank elementRank = getElementRank();
1473  // stk::mesh::Entity elmt = bulkData_->get_entity(elementRank,gid);
1474  // TEUCHOS_ASSERT(elmt->owner_rank()==procRank_);
1475  // return elementLocalId(elmt);
1476  std::unordered_map<stk::mesh::EntityId,std::size_t>::const_iterator itr = localIDHash_.find(gid);
1477  TEUCHOS_ASSERT(itr!=localIDHash_.end());
1478  return itr->second;
1479 }
1480 
1481 bool STK_Interface::isEdgeLocal(stk::mesh::Entity edge) const
1482 {
1483  return isEdgeLocal(bulkData_->identifier(edge));
1484 }
1485 
1486 bool STK_Interface::isEdgeLocal(stk::mesh::EntityId gid) const
1487 {
1488  std::unordered_map<stk::mesh::EntityId,std::size_t>::const_iterator itr = localEdgeIDHash_.find(gid);
1489  if (itr==localEdgeIDHash_.end()) {
1490  return false;
1491  }
1492  return true;
1493 }
1494 
1495 std::size_t STK_Interface::edgeLocalId(stk::mesh::Entity edge) const
1496 {
1497  return edgeLocalId(bulkData_->identifier(edge));
1498 }
1499 
1500 std::size_t STK_Interface::edgeLocalId(stk::mesh::EntityId gid) const
1501 {
1502  std::unordered_map<stk::mesh::EntityId,std::size_t>::const_iterator itr = localEdgeIDHash_.find(gid);
1503  TEUCHOS_ASSERT(itr!=localEdgeIDHash_.end());
1504  return itr->second;
1505 }
1506 
1507 bool STK_Interface::isFaceLocal(stk::mesh::Entity face) const
1508 {
1509  return isFaceLocal(bulkData_->identifier(face));
1510 }
1511 
1512 bool STK_Interface::isFaceLocal(stk::mesh::EntityId gid) const
1513 {
1514  std::unordered_map<stk::mesh::EntityId,std::size_t>::const_iterator itr = localFaceIDHash_.find(gid);
1515  if (itr==localFaceIDHash_.end()) {
1516  return false;
1517  }
1518  return true;
1519 }
1520 
1521 std::size_t STK_Interface::faceLocalId(stk::mesh::Entity face) const
1522 {
1523  return faceLocalId(bulkData_->identifier(face));
1524 }
1525 
1526 std::size_t STK_Interface::faceLocalId(stk::mesh::EntityId gid) const
1527 {
1528  std::unordered_map<stk::mesh::EntityId,std::size_t>::const_iterator itr = localFaceIDHash_.find(gid);
1529  TEUCHOS_ASSERT(itr!=localFaceIDHash_.end());
1530  return itr->second;
1531 }
1532 
1533 
1534 std::string STK_Interface::containingBlockId(stk::mesh::Entity elmt) const
1535 {
1536  for(const auto & eb_pair : elementBlocks_)
1537  if(bulkData_->bucket(elmt).member(*(eb_pair.second)))
1538  return eb_pair.first;
1539  return "";
1540 }
1541 
1542 stk::mesh::Field<double> * STK_Interface::getSolutionField(const std::string & fieldName,
1543  const std::string & blockId) const
1544 {
1545  // look up field in map
1546  std::map<std::pair<std::string,std::string>, SolutionFieldType*>::const_iterator
1547  iter = fieldNameToSolution_.find(std::make_pair(fieldName,blockId));
1548 
1549  // check to make sure field was actually found
1550  TEUCHOS_TEST_FOR_EXCEPTION(iter==fieldNameToSolution_.end(),std::runtime_error,
1551  "Solution field name \"" << fieldName << "\" in block ID \"" << blockId << "\" was not found");
1552 
1553  return iter->second;
1554 }
1555 
1556 stk::mesh::Field<double> * STK_Interface::getCellField(const std::string & fieldName,
1557  const std::string & blockId) const
1558 {
1559  // look up field in map
1560  std::map<std::pair<std::string,std::string>, SolutionFieldType*>::const_iterator
1561  iter = fieldNameToCellField_.find(std::make_pair(fieldName,blockId));
1562 
1563  // check to make sure field was actually found
1564  TEUCHOS_TEST_FOR_EXCEPTION(iter==fieldNameToCellField_.end(),std::runtime_error,
1565  "Cell field named \"" << fieldName << "\" in block ID \"" << blockId << "\" was not found");
1566 
1567  return iter->second;
1568 }
1569 
1570 stk::mesh::Field<double> * STK_Interface::getEdgeField(const std::string & fieldName,
1571  const std::string & blockId) const
1572 {
1573  // look up field in map
1574  std::map<std::pair<std::string,std::string>, SolutionFieldType*>::const_iterator
1575  iter = fieldNameToEdgeField_.find(std::make_pair(fieldName,blockId));
1576 
1577  // check to make sure field was actually found
1578  TEUCHOS_TEST_FOR_EXCEPTION(iter==fieldNameToEdgeField_.end(),std::runtime_error,
1579  "Edge field named \"" << fieldName << "\" in block ID \"" << blockId << "\" was not found");
1580 
1581  return iter->second;
1582 }
1583 
1584 stk::mesh::Field<double> * STK_Interface::getFaceField(const std::string & fieldName,
1585  const std::string & blockId) const
1586 {
1587  // look up field in map
1588  std::map<std::pair<std::string,std::string>, SolutionFieldType*>::const_iterator
1589  iter = fieldNameToFaceField_.find(std::make_pair(fieldName,blockId));
1590 
1591  // check to make sure field was actually found
1592  TEUCHOS_TEST_FOR_EXCEPTION(iter==fieldNameToFaceField_.end(),std::runtime_error,
1593  "Face field named \"" << fieldName << "\" in block ID \"" << blockId << "\" was not found");
1594 
1595  return iter->second;
1596 }
1597 
1598 Teuchos::RCP<const std::vector<stk::mesh::Entity> > STK_Interface::getElementsOrderedByLID() const
1599 {
1600  using Teuchos::RCP;
1601  using Teuchos::rcp;
1602 
1603  if(orderedElementVector_==Teuchos::null) {
1604  // safe because essentially this is a call to modify a mutable object
1605  const_cast<STK_Interface*>(this)->buildLocalElementIDs();
1606  }
1607 
1608  return orderedElementVector_.getConst();
1609 }
1610 
1611 void STK_Interface::addElementBlock(const std::string & name,const CellTopologyData * ctData)
1612 {
1613  TEUCHOS_ASSERT(not initialized_);
1614 
1615  stk::mesh::Part * block = metaData_->get_part(name);
1616  if(block==0) {
1617  block = &metaData_->declare_part_with_topology(name, stk::mesh::get_topology(shards::CellTopology(ctData), dimension_));
1618  }
1619 
1620  // construct cell topology object for this block
1621  Teuchos::RCP<shards::CellTopology> ct
1622  = Teuchos::rcp(new shards::CellTopology(ctData));
1623 
1624  // add element block part and cell topology
1625  elementBlocks_.insert(std::make_pair(name,block));
1626  elementBlockCT_.insert(std::make_pair(name,ct));
1627 }
1628 
1629 Teuchos::RCP<const std::vector<stk::mesh::Entity> > STK_Interface::getEdgesOrderedByLID() const
1630 {
1631  using Teuchos::RCP;
1632  using Teuchos::rcp;
1633 
1634  if(orderedEdgeVector_==Teuchos::null) {
1635  // safe because essentially this is a call to modify a mutable object
1636  const_cast<STK_Interface*>(this)->buildLocalEdgeIDs();
1637  }
1638 
1639  return orderedEdgeVector_.getConst();
1640 }
1641 
1642 void STK_Interface::addEdgeBlock(const std::string & elemBlockName,
1643  const std::string & edgeBlockName,
1644  const stk::topology & topology)
1645 {
1646  TEUCHOS_ASSERT(not initialized_);
1647 
1648  stk::mesh::Part * edge_block = metaData_->get_part(edgeBlockName);
1649  if(edge_block==0) {
1650  edge_block = &metaData_->declare_part_with_topology(edgeBlockName, topology);
1651  }
1652 
1653  /* There is only one edge block for each edge topology, so declare it
1654  * as a subset of the element block even if it wasn't just created.
1655  */
1656  stk::mesh::Part * elem_block = metaData_->get_part(elemBlockName);
1657  metaData_->declare_part_subset(*elem_block, *edge_block);
1658 
1659  // add edge block part
1660  edgeBlocks_.insert(std::make_pair(edgeBlockName,edge_block));
1661 }
1662 
1663 void STK_Interface::addEdgeBlock(const std::string & elemBlockName,
1664  const std::string & edgeBlockName,
1665  const CellTopologyData * ctData)
1666 {
1667  TEUCHOS_ASSERT(not initialized_);
1668 
1669  addEdgeBlock(elemBlockName,
1670  edgeBlockName,
1671  stk::mesh::get_topology(shards::CellTopology(ctData), dimension_));
1672 }
1673 
1674 Teuchos::RCP<const std::vector<stk::mesh::Entity> > STK_Interface::getFacesOrderedByLID() const
1675 {
1676  using Teuchos::RCP;
1677  using Teuchos::rcp;
1678 
1679  if(orderedFaceVector_==Teuchos::null) {
1680  // safe because essentially this is a call to modify a mutable object
1681  const_cast<STK_Interface*>(this)->buildLocalFaceIDs();
1682  }
1683 
1684  return orderedFaceVector_.getConst();
1685 }
1686 
1687 void STK_Interface::addFaceBlock(const std::string & elemBlockName,
1688  const std::string & faceBlockName,
1689  const stk::topology & topology)
1690 {
1691  TEUCHOS_ASSERT(not initialized_);
1692 
1693  stk::mesh::Part * face_block = metaData_->get_part(faceBlockName);
1694  if(face_block==0) {
1695  face_block = &metaData_->declare_part_with_topology(faceBlockName, topology);
1696  }
1697 
1698  /* There is only one face block for each edge topology, so declare it
1699  * as a subset of the element block even if it wasn't just created.
1700  */
1701  stk::mesh::Part * elem_block = metaData_->get_part(elemBlockName);
1702  metaData_->declare_part_subset(*elem_block, *face_block);
1703 
1704  // add face block part
1705  faceBlocks_.insert(std::make_pair(faceBlockName,face_block));
1706 }
1707 
1708 void STK_Interface::addFaceBlock(const std::string & elemBlockName,
1709  const std::string & faceBlockName,
1710  const CellTopologyData * ctData)
1711 {
1712  TEUCHOS_ASSERT(not initialized_);
1713 
1714  addFaceBlock(elemBlockName,
1715  faceBlockName,
1716  stk::mesh::get_topology(shards::CellTopology(ctData), dimension_));
1717 }
1718 
1720 {
1721  dimension_ = metaData_->spatial_dimension();
1722 
1723  // declare coordinates and node parts
1724  coordinatesField_ = &metaData_->declare_field<VectorFieldType>(stk::topology::NODE_RANK, coordsString);
1725  edgesField_ = &metaData_->declare_field<VectorFieldType>(stk::topology::EDGE_RANK, edgesString);
1726  if (dimension_ > 2)
1727  facesField_ = &metaData_->declare_field<VectorFieldType>(stk::topology::FACE_RANK, facesString);
1728  processorIdField_ = &metaData_->declare_field<ProcIdFieldType>(stk::topology::ELEMENT_RANK, "PROC_ID");
1729  loadBalField_ = &metaData_->declare_field<SolutionFieldType>(stk::topology::ELEMENT_RANK, "LOAD_BAL");
1730 
1731  // stk::mesh::put_field( *coordinatesField_ , getNodeRank() , metaData_->universal_part() );
1732 
1733  nodesPart_ = &metaData_->declare_part(nodesString,getNodeRank());
1734  nodesPartVec_.push_back(nodesPart_);
1735  edgesPart_ = &metaData_->declare_part(edgesString,getEdgeRank());
1736  edgesPartVec_.push_back(edgesPart_);
1737  if (dimension_ > 2) {
1738  facesPart_ = &metaData_->declare_part(facesString,getFaceRank());
1739  facesPartVec_.push_back(facesPart_);
1740  }
1741 }
1742 
1744 {
1745  currentLocalId_ = 0;
1746 
1747  orderedElementVector_ = Teuchos::null; // forces rebuild of ordered lists
1748 
1749  // might be better (faster) to do this by buckets
1750  std::vector<stk::mesh::Entity> elements;
1751  getMyElements(elements);
1752 
1753  for(std::size_t index=0;index<elements.size();++index) {
1754  stk::mesh::Entity element = elements[index];
1755 
1756  // set processor rank
1757  ProcIdData * procId = stk::mesh::field_data(*processorIdField_,element);
1758  procId[0] = Teuchos::as<ProcIdData>(procRank_);
1759 
1760  localIDHash_[bulkData_->identifier(element)] = currentLocalId_;
1761 
1762  currentLocalId_++;
1763  }
1764 
1765  // copy elements into the ordered element vector
1766  orderedElementVector_ = Teuchos::rcp(new std::vector<stk::mesh::Entity>(elements));
1767 
1768  elements.clear();
1769  getNeighborElements(elements);
1770 
1771  for(std::size_t index=0;index<elements.size();++index) {
1772  stk::mesh::Entity element = elements[index];
1773 
1774  // set processor rank
1775  ProcIdData * procId = stk::mesh::field_data(*processorIdField_,element);
1776  procId[0] = Teuchos::as<ProcIdData>(procRank_);
1777 
1778  localIDHash_[bulkData_->identifier(element)] = currentLocalId_;
1779 
1780  currentLocalId_++;
1781  }
1782 
1783  orderedElementVector_->insert(orderedElementVector_->end(),elements.begin(),elements.end());
1784 }
1785 
1787 {
1788  std::vector<std::string> names;
1789  getElementBlockNames(names);
1790 
1791  for(std::size_t b=0;b<names.size();b++) {
1792  // find user specified block weight, otherwise use 1.0
1793  std::map<std::string,double>::const_iterator bw_itr = blockWeights_.find(names[b]);
1794  double blockWeight = (bw_itr!=blockWeights_.end()) ? bw_itr->second : 1.0;
1795 
1796  std::vector<stk::mesh::Entity> elements;
1797  getMyElements(names[b],elements);
1798 
1799  for(std::size_t index=0;index<elements.size();++index) {
1800  // set local element ID
1801  double * loadBal = stk::mesh::field_data(*loadBalField_,elements[index]);
1802  loadBal[0] = blockWeight;
1803  }
1804  }
1805 }
1806 
1808 {
1809  currentLocalId_ = 0;
1810 
1811  orderedEdgeVector_ = Teuchos::null; // forces rebuild of ordered lists
1812 
1813  // might be better (faster) to do this by buckets
1814  std::vector<stk::mesh::Entity> edges;
1815  getMyEdges(edges);
1816 
1817  for(std::size_t index=0;index<edges.size();++index) {
1818  stk::mesh::Entity edge = edges[index];
1819  localEdgeIDHash_[bulkData_->identifier(edge)] = currentLocalId_;
1820  currentLocalId_++;
1821  }
1822 
1823  // copy edges into the ordered edge vector
1824  orderedEdgeVector_ = Teuchos::rcp(new std::vector<stk::mesh::Entity>(edges));
1825 }
1826 
1828 {
1829  currentLocalId_ = 0;
1830 
1831  orderedFaceVector_ = Teuchos::null; // forces rebuild of ordered lists
1832 
1833  // might be better (faster) to do this by buckets
1834  std::vector<stk::mesh::Entity> faces;
1835  getMyFaces(faces);
1836 
1837  for(std::size_t index=0;index<faces.size();++index) {
1838  stk::mesh::Entity face = faces[index];
1839  localFaceIDHash_[bulkData_->identifier(face)] = currentLocalId_;
1840  currentLocalId_++;
1841  }
1842 
1843  // copy faces into the ordered face vector
1844  orderedFaceVector_ = Teuchos::rcp(new std::vector<stk::mesh::Entity>(faces));
1845 }
1846 
1847 bool
1848 STK_Interface::isMeshCoordField(const std::string & eBlock,
1849  const std::string & fieldName,
1850  int & axis) const
1851 {
1852  std::map<std::string,std::vector<std::string> >::const_iterator blkItr = meshCoordFields_.find(eBlock);
1853  if(blkItr==meshCoordFields_.end()) {
1854  return false;
1855  }
1856 
1857  axis = 0;
1858  for(axis=0;axis<Teuchos::as<int>(blkItr->second.size());axis++) {
1859  if(blkItr->second[axis]==fieldName)
1860  break; // found axis, break
1861  }
1862 
1863  if(axis>=Teuchos::as<int>(blkItr->second.size()))
1864  return false;
1865 
1866  return true;
1867 }
1868 
1869 std::pair<Teuchos::RCP<std::vector<std::pair<std::size_t,std::size_t> > >, Teuchos::RCP<std::vector<unsigned int> > >
1871 {
1872  Teuchos::RCP<std::vector<std::pair<std::size_t,std::size_t> > > vec;
1873  Teuchos::RCP<std::vector<unsigned int > > type_vec = rcp(new std::vector<unsigned int>);
1874  const std::vector<Teuchos::RCP<const PeriodicBC_MatcherBase> > & matchers = getPeriodicBCVector();
1875  const bool & useBBoxSearch = useBoundingBoxSearch();
1876  std::vector<std::vector<std::string> > matchedSides(3); // (coord,edge,face)
1877 
1878  // build up the vectors by looping over the matched pair
1879  for(std::size_t m=0;m<matchers.size();m++){
1880  unsigned int type;
1881  if(matchers[m]->getType() == "coord")
1882  type = 0;
1883  else if(matchers[m]->getType() == "edge")
1884  type = 1;
1885  else if(matchers[m]->getType() == "face")
1886  type = 2;
1887  else
1888  TEUCHOS_ASSERT(false);
1889 #ifdef PANZER_HAVE_STKSEARCH
1890 
1891  if (useBBoxSearch) {
1892  vec = matchers[m]->getMatchedPair(*this,matchedSides[type],vec);
1893  } else {
1894  vec = matchers[m]->getMatchedPair(*this,vec);
1895  }
1896 #else
1897  TEUCHOS_TEST_FOR_EXCEPTION(useBBoxSearch,std::logic_error,
1898  "panzer::STK_Interface::getPeriodicNodePairing(): Requested bounding box search, but "
1899  "did not compile with STK_SEARCH enabled.");
1900  vec = matchers[m]->getMatchedPair(*this,vec);
1901 #endif
1902  type_vec->insert(type_vec->begin(),vec->size()-type_vec->size(),type);
1903  matchedSides[type].push_back(matchers[m]->getLeftSidesetName());
1904  }
1905 
1906  return std::make_pair(vec,type_vec);
1907 }
1908 
1909 bool STK_Interface::validBlockId(const std::string & blockId) const
1910 {
1911  std::map<std::string, stk::mesh::Part*>::const_iterator blkItr = elementBlocks_.find(blockId);
1912 
1913  return blkItr!=elementBlocks_.end();
1914 }
1915 
1916 void STK_Interface::print(std::ostream & os) const
1917 {
1918  std::vector<std::string> blockNames, sidesetNames, nodesetNames;
1919 
1920  getElementBlockNames(blockNames);
1921  getSidesetNames(sidesetNames);
1922  getNodesetNames(nodesetNames);
1923 
1924  os << "STK Mesh data:\n";
1925  os << " Spatial dim = " << getDimension() << "\n";
1926  if(getDimension()==2)
1927  os << " Entity counts (Nodes, Edges, Cells) = ( "
1928  << getEntityCounts(getNodeRank()) << ", "
1929  << getEntityCounts(getEdgeRank()) << ", "
1930  << getEntityCounts(getElementRank()) << " )\n";
1931  else if(getDimension()==3)
1932  os << " Entity counts (Nodes, Edges, Faces, Cells) = ( "
1933  << getEntityCounts(getNodeRank()) << ", "
1934  << getEntityCounts(getEdgeRank()) << ", "
1935  << getEntityCounts(getSideRank()) << ", "
1936  << getEntityCounts(getElementRank()) << " )\n";
1937  else
1938  os << " Entity counts (Nodes, Cells) = ( "
1939  << getEntityCounts(getNodeRank()) << ", "
1940  << getEntityCounts(getElementRank()) << " )\n";
1941 
1942  os << " Element blocks = ";
1943  for(std::size_t i=0;i<blockNames.size();i++)
1944  os << "\"" << blockNames[i] << "\" ";
1945  os << "\n";
1946  os << " Sidesets = ";
1947  for(std::size_t i=0;i<sidesetNames.size();i++)
1948  os << "\"" << sidesetNames[i] << "\" ";
1949  os << "\n";
1950  os << " Nodesets = ";
1951  for(std::size_t i=0;i<nodesetNames.size();i++)
1952  os << "\"" << nodesetNames[i] << "\" ";
1953  os << std::endl;
1954 
1955  // print out periodic boundary conditions
1956  const std::vector<Teuchos::RCP<const PeriodicBC_MatcherBase> > & bcVector
1957  = getPeriodicBCVector();
1958  if(bcVector.size()!=0) {
1959  os << " Periodic BCs:\n";
1960  for(std::size_t i=0;i<bcVector.size();i++)
1961  os << " " << bcVector[i]->getString() << "\n";
1962  os << std::endl;
1963  }
1964 }
1965 
1966 void STK_Interface::printMetaData(std::ostream & os) const
1967 {
1968  std::vector<std::string> blockNames, sidesetNames, nodesetNames;
1969 
1970  getElementBlockNames(blockNames);
1971  getSidesetNames(sidesetNames);
1972  getNodesetNames(nodesetNames);
1973 
1974  os << "STK Meta data:\n";
1975  os << " Element blocks = ";
1976  for(std::size_t i=0;i<blockNames.size();i++)
1977  os << "\"" << blockNames[i] << "\" ";
1978  os << "\n";
1979  os << " Sidesets = ";
1980  for(std::size_t i=0;i<sidesetNames.size();i++)
1981  os << "\"" << sidesetNames[i] << "\" ";
1982  os << "\n";
1983  os << " Nodesets = ";
1984  for(std::size_t i=0;i<nodesetNames.size();i++)
1985  os << "\"" << nodesetNames[i] << "\" ";
1986  os << std::endl;
1987 
1988  // print out periodic boundary conditions
1989  const std::vector<Teuchos::RCP<const PeriodicBC_MatcherBase> > & bcVector
1990  = getPeriodicBCVector();
1991  if(bcVector.size()!=0) {
1992  os << " Periodic BCs:\n";
1993  for(std::size_t i=0;i<bcVector.size();i++)
1994  os << " " << bcVector[i]->getString() << "\n";
1995  os << std::endl;
1996  }
1997 
1998  // print all available fields in meta data
1999  os << " Fields = ";
2000  const stk::mesh::FieldVector & fv = metaData_->get_fields();
2001  for(std::size_t i=0;i<fv.size();i++)
2002  os << "\"" << fv[i]->name() << "\" ";
2003  os << std::endl;
2004 }
2005 
2006 Teuchos::RCP<const shards::CellTopology> STK_Interface::getCellTopology(const std::string & eBlock) const
2007 {
2008  std::map<std::string, Teuchos::RCP<shards::CellTopology> >::const_iterator itr;
2009  itr = elementBlockCT_.find(eBlock);
2010 
2011  if(itr==elementBlockCT_.end()) {
2012  std::stringstream ss;
2013  printMetaData(ss);
2014  TEUCHOS_TEST_FOR_EXCEPTION(itr==elementBlockCT_.end(),std::logic_error,
2015  "STK_Interface::getCellTopology: No such element block \"" +eBlock +"\" available.\n\n"
2016  << "STK Meta Data follows: \n" << ss.str());
2017  }
2018 
2019  return itr->second;
2020 }
2021 
2022 Teuchos::RCP<Teuchos::MpiComm<int> > STK_Interface::getSafeCommunicator(stk::ParallelMachine parallelMach) const
2023 {
2024  MPI_Comm newComm;
2025  const int err = MPI_Comm_dup (parallelMach, &newComm);
2026  TEUCHOS_TEST_FOR_EXCEPTION(err != MPI_SUCCESS, std::runtime_error,
2027  "panzer::STK_Interface: MPI_Comm_dup failed with error \""
2028  << Teuchos::mpiErrorCodeToString (err) << "\".");
2029 
2030  return Teuchos::rcp(new Teuchos::MpiComm<int>(Teuchos::opaqueWrapper (newComm,MPI_Comm_free)));
2031 }
2032 
2033 void STK_Interface::rebalance(const Teuchos::ParameterList & /* params */)
2034 {
2035  TEUCHOS_TEST_FOR_EXCEPTION(true, std::runtime_error, "Rebalance not currently supported");
2036 #if 0
2037  // make sure weights have been set (a local operation)
2039 
2040  stk::mesh::Selector selector(getMetaData()->universal_part());
2041  stk::mesh::Selector owned_selector(getMetaData()->locally_owned_part());
2042 
2043  Teuchos::FancyOStream out(Teuchos::rcpFromRef(std::cout));
2044  out.setOutputToRootOnly(0);
2045  out << "Load balance before: " << stk::rebalance::check_balance(*getBulkData(), loadBalField_, getElementRank(), &selector) << std::endl;
2046 
2047  // perform reblance
2048  Teuchos::ParameterList graph;
2049  if(params.begin()!=params.end())
2050  graph.sublist(stk::rebalance::Zoltan::default_parameters_name()) = params;
2051  stk::rebalance::Zoltan zoltan_partition(*mpiComm_->getRawMpiComm(), getDimension(), graph);
2052  stk::rebalance::rebalance(*getBulkData(), owned_selector, &getCoordinatesField(), loadBalField_, zoltan_partition);
2053 
2054  out << "Load balance after: " << stk::rebalance::check_balance(*getBulkData(), loadBalField_, getElementRank(), &selector) << std::endl;
2055 
2056  currentLocalId_ = 0;
2057  orderedElementVector_ = Teuchos::null; // forces rebuild of ordered lists
2058 #endif
2059 }
2060 
2061 Teuchos::RCP<const Teuchos::Comm<int> >
2063 {
2064  TEUCHOS_ASSERT(this->isInitialized());
2065  return mpiComm_;
2066 }
2067 
2068 void STK_Interface::refineMesh(const int numberOfLevels, const bool deleteParentElements) {
2069 #ifdef PANZER_HAVE_PERCEPT
2070  TEUCHOS_TEST_FOR_EXCEPTION(numberOfLevels < 1,std::runtime_error,
2071  "ERROR: Number of levels for uniform refinement must be greater than 0");
2072  TEUCHOS_ASSERT(nonnull(refinedMesh_));
2073  TEUCHOS_ASSERT(nonnull(breakPattern_));
2074 
2075  refinedMesh_->setCoordinatesField();
2076 
2077  percept::UniformRefiner breaker(*refinedMesh_,*breakPattern_);
2078  breaker.setRemoveOldElements(deleteParentElements);
2079 
2080  for (int i=0; i < numberOfLevels; ++i)
2081  breaker.doBreak();
2082 
2083 #else
2084  TEUCHOS_TEST_FOR_EXCEPTION(true,std::logic_error,
2085  "ERROR: Uniform refinement requested. This requires the Percept package to be enabled in Trilinos!");
2086 #endif
2087 }
2088 
2089 
2090 } // end namespace panzer_stk
std::size_t edgeLocalId(stk::mesh::Entity elmt) const
Teuchos::RCP< ElementDescriptor > buildElementDescriptor(stk::mesh::EntityId elmtId, std::vector< stk::mesh::EntityId > &nodes)
Teuchos::RCP< stk::mesh::MetaData > getMetaData() const
stk::mesh::Part * getElementBlockPart(const std::string &name) const
get the block part
Teuchos::RCP< std::vector< stk::mesh::Entity > > orderedEdgeVector_
bool isMeshCoordField(const std::string &eBlock, const std::string &fieldName, int &axis) const
void addNodeset(const std::string &name)
Teuchos::RCP< stk::mesh::BulkData > bulkData_
void getMySides(const std::string &sideName, std::vector< stk::mesh::Entity > &sides) const
std::vector< stk::mesh::Part * > facesPartVec_
std::map< std::pair< std::string, std::string >, SolutionFieldType * > fieldNameToSolution_
std::string containingBlockId(stk::mesh::Entity elmt) const
void addEntityToNodeset(stk::mesh::Entity entity, stk::mesh::Part *nodeset)
bool isFaceLocal(stk::mesh::Entity face) const
stk::mesh::EntityRank getFaceRank() const
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.
std::size_t elementLocalId(stk::mesh::Entity elmt) const
stk::mesh::Field< double > SolutionFieldType
std::vector< stk::mesh::Part * > edgesPartVec_
std::map< std::string, double > blockWeights_
void getMyElements(std::vector< stk::mesh::Entity > &elements) const
void addEntityToSideset(stk::mesh::Entity entity, stk::mesh::Part *sideset)
std::map< std::string, stk::mesh::Part * > faceBlocks_
stk::mesh::EntityRank getNodeRank() 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::Field< double > * getEdgeField(const std::string &fieldName, const std::string &blockId) const
Teuchos::RCP< stk::mesh::BulkData > getBulkData() const
void getFaceBlockNames(std::vector< std::string > &names) const
std::map< std::pair< std::string, std::string >, SolutionFieldType * > fieldNameToFaceField_
std::map< std::string, stk::mesh::Part * > sidesets_
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 getSubcellIndices(unsigned entityRank, stk::mesh::EntityId elementId, std::vector< stk::mesh::EntityId > &subcellIds) const
void addInformationRecords(const std::vector< std::string > &info_records)
void getMyFaces(std::vector< stk::mesh::Entity > &faces) const
stk::mesh::EntityRank getSideRank() const
stk::mesh::Entity findConnectivityById(stk::mesh::Entity src, stk::mesh::EntityRank tgt_rank, unsigned rel_id) const
void addEdgeBlock(const std::string &elemBlockName, const std::string &edgeBlockName, const stk::topology &topology)
static const std::string nodesString
const double * getNodeCoordinates(stk::mesh::EntityId nodeId) const
void instantiateBulkData(stk::ParallelMachine parallelMach)
stk::mesh::Part * getEdgeBlock(const std::string &name) const
get the block part
void getElementBlockNames(std::vector< std::string > &names) const
std::set< std::string > informationRecords_
stk::mesh::EntityRank getEdgeRank() const
stk::mesh::Part * getSideset(const std::string &name) 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.
stk::mesh::EntityId elementGlobalId(std::size_t lid) 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 initialize(stk::ParallelMachine parallelMach, bool setupIO=true, const bool buildRefinementSupport=false)
static const std::string coordsString
stk::mesh::Field< double, stk::mesh::Cartesian > VectorFieldType
stk::mesh::Field< double > * getCellField(const std::string &fieldName, const std::string &blockId) const
std::size_t getEntityCounts(unsigned entityRank) const
get the global counts for the entity of specified rank
std::map< std::string, stk::mesh::Part * > edgeBlocks_
stk::mesh::Part * getFaceBlock(const std::string &name) const
get the block part
Teuchos::RCP< const std::vector< stk::mesh::Entity > > getFacesOrderedByLID() const
std::map< std::string, stk::mesh::Part * > elementBlocks_
bool isInitialized() const
Has initialize been called on this mesh object?
std::vector< stk::mesh::Part * > nodesPartVec_
static const std::string edgeBlockString
static const std::string edgesString
void printMetaData(std::ostream &os) const
void initializeFieldsInSTK(const std::map< std::pair< std::string, std::string >, SolutionFieldType *> &nameToField, bool setupIO)
Teuchos::RCP< const std::vector< stk::mesh::Entity > > getElementsOrderedByLID() const
void addNode(stk::mesh::EntityId gid, const std::vector< double > &coord)
void getAllFaces(const std::string &faceBlockName, std::vector< stk::mesh::Entity > &faces) const
void addEntitiesToEdgeBlock(std::vector< stk::mesh::Entity > entities, stk::mesh::Part *edgeblock)
const std::vector< Teuchos::RCP< const PeriodicBC_MatcherBase > > & getPeriodicBCVector() const
void buildSubcells()
force the mesh to build subcells: edges and faces
Teuchos::RCP< Teuchos::MpiComm< int > > mpiComm_
void getAllEdges(const std::string &edgeBlockName, std::vector< stk::mesh::Entity > &edges) const
Teuchos::RCP< stk::mesh::MetaData > metaData_
Teuchos::RCP< std::vector< stk::mesh::Entity > > orderedElementVector_
void getEdgeBlockNames(std::vector< std::string > &names) const
Teuchos::RCP< const Teuchos::Comm< int > > getComm() const
get the comm associated with this mesh
void writeToExodus(const std::string &filename, const bool append=false)
Write this mesh and associated fields to the given output file.
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)
bool validBlockId(const std::string &blockId) 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_
std::unordered_map< stk::mesh::EntityId, std::size_t > localEdgeIDHash_
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)
Teuchos::RCP< const shards::CellTopology > getCellTopology(const std::string &eBlock) const
void getAllSides(const std::string &sideName, std::vector< stk::mesh::Entity > &sides) const
Teuchos::RCP< const std::vector< stk::mesh::Entity > > getEdgesOrderedByLID() const
void getSidesetNames(std::vector< std::string > &name) const
const bool & useBoundingBoxSearch() const
std::map< std::string, std::vector< std::string > > meshCoordFields_
static const std::string facesString
Teuchos::RCP< Teuchos::MpiComm< int > > getSafeCommunicator(stk::ParallelMachine parallelMach) const
void addEntitiesToFaceBlock(std::vector< stk::mesh::Entity > entities, stk::mesh::Part *faceblock)
std::pair< Teuchos::RCP< std::vector< std::pair< std::size_t, std::size_t > > >, Teuchos::RCP< std::vector< unsigned int > > > getPeriodicNodePairing() const
bool isEdgeLocal(stk::mesh::Entity edge) const
void refineMesh(const int numberOfLevels, const bool deleteParentElements)
void getNodesetNames(std::vector< std::string > &name) 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_
std::vector< stk::mesh::EntityId > maxEntityId_
const VectorFieldType & getCoordinatesField() const
void print(std::ostream &os) const
stk::mesh::EntityRank getElementRank() const
void getMyEdges(std::vector< stk::mesh::Entity > &edges) const
stk::mesh::Field< ProcIdData > ProcIdFieldType
stk::mesh::EntityId getMaxEntityId(unsigned entityRank) const
get max entity ID of type entityRank
void addFaceField(const std::string &fieldName, const std::string &blockId)
void addCellField(const std::string &fieldName, const std::string &blockId)
void getMyNodes(const std::string &sideName, const std::string &blockName, std::vector< stk::mesh::Entity > &nodes) const
void getElementsSharingNodes(const std::vector< stk::mesh::EntityId > nodeId, std::vector< stk::mesh::Entity > &elements) const
get a set of elements sharing multiple nodes
void addElementBlock(const std::string &name, const CellTopologyData *ctData)
void getOwnedElementsSharingNode(stk::mesh::Entity node, std::vector< stk::mesh::Entity > &elements, std::vector< int > &relIds) const
std::size_t faceLocalId(stk::mesh::Entity elmt) const
stk::mesh::Part * getNodeset(const std::string &name) const
stk::mesh::Field< double > * getFaceField(const std::string &fieldName, const std::string &blockId) const
stk::mesh::Field< double > * getSolutionField(const std::string &fieldName, const std::string &blockId) const
std::map< std::pair< std::string, std::string >, SolutionFieldType * > fieldNameToCellField_
std::map< std::string, stk::mesh::Part * > nodesets_
void getNeighborElements(std::vector< stk::mesh::Entity > &elements) const