Panzer  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
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/CreateAdjacentEntities.hpp>
56 
57 // #include <stk_rebalance/Rebalance.hpp>
58 // #include <stk_rebalance/Partition.hpp>
59 // #include <stk_rebalance/ZoltanPartition.hpp>
60 // #include <stk_rebalance_utils/RebalanceUtils.hpp>
61 
62 #include <stk_util/parallel/ParallelReduce.hpp>
63 #include <stk_util/parallel/CommSparse.hpp>
64 
65 #ifdef PANZER_HAVE_IOSS
66 #include <Ionit_Initializer.h>
67 #include <stk_io/IossBridge.hpp>
68 #endif
69 
70 #ifdef PANZER_HAVE_PERCEPT
71 #include <percept/PerceptMesh.hpp>
72 #include <adapt/UniformRefinerPattern.hpp>
73 #include <adapt/UniformRefiner.hpp>
74 #endif
75 
77 
78 #include <set>
79 
80 using Teuchos::RCP;
81 using Teuchos::rcp;
82 
83 namespace panzer_stk {
84 
86 ElementDescriptor::ElementDescriptor(stk::mesh::EntityId gid,const std::vector<stk::mesh::EntityId> & nodes)
87  : gid_(gid), nodes_(nodes) {}
89 
93 buildElementDescriptor(stk::mesh::EntityId elmtId,std::vector<stk::mesh::EntityId> & nodes)
94 {
95  return Teuchos::rcp(new ElementDescriptor(elmtId,nodes));
96 }
97 
98 const std::string STK_Interface::coordsString = "coordinates";
99 const std::string STK_Interface::nodesString = "nodes";
100 const std::string STK_Interface::edgesString = "edges";
101 const std::string STK_Interface::facesString = "faces";
102 
104  : dimension_(0), initialized_(false), currentLocalId_(0), initialStateTime_(0.0), currentStateTime_(0.0), useFieldCoordinates_(false)
105 {
106  metaData_ = rcp(new stk::mesh::MetaData());
107 }
108 
110  : dimension_(0), initialized_(false), currentLocalId_(0), initialStateTime_(0.0), currentStateTime_(0.0), useFieldCoordinates_(false)
111 {
112  metaData_ = metaData;
113 }
114 
116  : dimension_(dim), initialized_(false), currentLocalId_(0), useFieldCoordinates_(false)
117 {
118  std::vector<std::string> entity_rank_names = stk::mesh::entity_rank_names();
119  entity_rank_names.push_back("FAMILY_TREE");
120 
121  metaData_ = rcp(new stk::mesh::MetaData(dimension_,entity_rank_names));
122 
124 }
125 
126 void STK_Interface::addSideset(const std::string & name,const CellTopologyData * ctData)
127 {
130 
131  stk::mesh::Part * sideset = metaData_->get_part(name);
132  if(sideset==nullptr)
133  sideset = &metaData_->declare_part_with_topology(name,
134  stk::mesh::get_topology(shards::CellTopology(ctData), dimension_));
135  sidesets_.insert(std::make_pair(name,sideset));
136 }
137 
138 void STK_Interface::addNodeset(const std::string & name)
139 {
142 
143  stk::mesh::Part * nodeset = metaData_->get_part(name);
144  if(nodeset==nullptr) {
145  const CellTopologyData * ctData = shards::getCellTopologyData<shards::Node>();
146  nodeset = &metaData_->declare_part_with_topology(name,
147  stk::mesh::get_topology(shards::CellTopology(ctData), dimension_));
148  }
149  nodesets_.insert(std::make_pair(name,nodeset));
150 }
151 
152 void STK_Interface::addSolutionField(const std::string & fieldName,const std::string & blockId)
153 {
155  "Unknown element block \"" << blockId << "\"");
156  std::pair<std::string,std::string> key = std::make_pair(fieldName,blockId);
157 
158  // add & declare field if not already added...currently assuming linears
159  if(fieldNameToSolution_.find(key)==fieldNameToSolution_.end()) {
160  SolutionFieldType * field = metaData_->get_field<SolutionFieldType>(stk::topology::NODE_RANK, fieldName);
161  if(field==0)
162  field = &metaData_->declare_field<SolutionFieldType>(stk::topology::NODE_RANK, fieldName);
163  if ( initialized_ ) {
164  metaData_->enable_late_fields();
165  stk::mesh::FieldTraits<SolutionFieldType>::data_type* init_sol = nullptr;
166  stk::mesh::put_field_on_mesh(*field, metaData_->universal_part(),init_sol );
167  }
169  }
170 }
171 
172 void STK_Interface::addCellField(const std::string & fieldName,const std::string & blockId)
173 {
175  "Unknown element block \"" << blockId << "\"");
176  std::pair<std::string,std::string> key = std::make_pair(fieldName,blockId);
177 
178  // add & declare field if not already added...currently assuming linears
179  if(fieldNameToCellField_.find(key)==fieldNameToCellField_.end()) {
180  SolutionFieldType * field = metaData_->get_field<SolutionFieldType>(stk::topology::ELEMENT_RANK, fieldName);
181  if(field==0)
182  field = &metaData_->declare_field<SolutionFieldType>(stk::topology::ELEMENT_RANK, fieldName);
183 
184  if ( initialized_ ) {
185  metaData_->enable_late_fields();
186  stk::mesh::FieldTraits<SolutionFieldType>::data_type* init_sol = nullptr;
187  stk::mesh::put_field_on_mesh(*field, metaData_->universal_part(),init_sol );
188  }
190  }
191 }
192 
193 void STK_Interface::addMeshCoordFields(const std::string & blockId,
194  const std::vector<std::string> & coordNames,
195  const std::string & dispPrefix)
196 {
198  TEUCHOS_ASSERT(dimension_==coordNames.size());
201  "Unknown element block \"" << blockId << "\"");
202 
203  // we only allow one alternative coordinate field
204  TEUCHOS_TEST_FOR_EXCEPTION(meshCoordFields_.find(blockId)!=meshCoordFields_.end(),std::invalid_argument,
205  "STK_Interface::addMeshCoordFields: Can't set more than one set of coordinate "
206  "fields for element block \""+blockId+"\".");
207 
208  // Note that there is a distinction between the key which is used for lookups
209  // and the field that lives on the mesh, which is used for printing the displacement.
210 
211  // just copy the coordinate names
212  meshCoordFields_[blockId] = coordNames;
213 
214  // must fill in the displacement fields
215  std::vector<std::string> & dispFields = meshDispFields_[blockId];
216  dispFields.resize(dimension_);
217 
218  for(unsigned i=0;i<dimension_;i++) {
219  std::pair<std::string,std::string> key = std::make_pair(coordNames[i],blockId);
220  std::string dispName = dispPrefix+coordNames[i];
221 
222  dispFields[i] = dispName; // record this field as a
223  // displacement field
224 
225  // add & declare field if not already added...currently assuming linears
226  if(fieldNameToSolution_.find(key)==fieldNameToSolution_.end()) {
227 
228  SolutionFieldType * field = metaData_->get_field<SolutionFieldType>(stk::topology::NODE_RANK, dispName);
229  if(field==0) {
230  field = &metaData_->declare_field<SolutionFieldType>(stk::topology::NODE_RANK, dispName);
231  }
233  }
234  }
235 }
236 
237 void STK_Interface::initialize(stk::ParallelMachine parallelMach,bool setupIO,
238  const bool buildRefinementSupport)
239 {
241  TEUCHOS_ASSERT(dimension_!=0); // no zero dimensional meshes!
242 
243  if(mpiComm_==Teuchos::null)
244  mpiComm_ = getSafeCommunicator(parallelMach);
245 
246  procRank_ = stk::parallel_machine_rank(*mpiComm_->getRawMpiComm());
247 
248  // associating the field with a part: universal part!
249  stk::mesh::FieldTraits<VectorFieldType>::data_type* init_vf = nullptr; // gcc 4.8 hack
250  stk::mesh::FieldTraits<ProcIdFieldType>::data_type* init_pid = nullptr; // gcc 4.8 hack
251  stk::mesh::FieldTraits<SolutionFieldType>::data_type* init_sol = nullptr; // gcc 4.8 hack
252  stk::mesh::put_field_on_mesh( *coordinatesField_ , metaData_->universal_part(), getDimension(),init_vf);
253  stk::mesh::put_field_on_mesh( *edgesField_ , metaData_->universal_part(), getDimension(),init_vf);
254  if (dimension_ > 2)
255  stk::mesh::put_field_on_mesh( *facesField_ , metaData_->universal_part(), getDimension(),init_vf);
256  stk::mesh::put_field_on_mesh( *processorIdField_ , metaData_->universal_part(),init_pid);
257  stk::mesh::put_field_on_mesh( *loadBalField_ , metaData_->universal_part(),init_sol);
258 
261 
262 #ifdef PANZER_HAVE_IOSS
263  if(setupIO) {
264  // setup Exodus file IO
266 
267  // add element blocks
268  {
269  std::map<std::string, stk::mesh::Part*>::iterator itr;
270  for(itr=elementBlocks_.begin();itr!=elementBlocks_.end();++itr)
271  if(!stk::io::is_part_io_part(*itr->second))
272  stk::io::put_io_part_attribute(*itr->second); // this can only be called once per part
273  }
274 
275  // add side sets
276  {
277  std::map<std::string, stk::mesh::Part*>::iterator itr;
278  for(itr=sidesets_.begin();itr!=sidesets_.end();++itr)
279  if(!stk::io::is_part_io_part(*itr->second))
280  stk::io::put_io_part_attribute(*itr->second); // this can only be called once per part
281  }
282 
283  // add node sets
284  {
285  std::map<std::string, stk::mesh::Part*>::iterator itr;
286  for(itr=nodesets_.begin();itr!=nodesets_.end();++itr)
287  if(!stk::io::is_part_io_part(*itr->second))
288  stk::io::put_io_part_attribute(*itr->second); // this can only be called once per part
289  }
290 
291  // add nodes
292  if(!stk::io::is_part_io_part(*nodesPart_))
293  stk::io::put_io_part_attribute(*nodesPart_);
294 
295  stk::io::set_field_role(*coordinatesField_, Ioss::Field::MESH);
296  stk::io::set_field_role(*edgesField_, Ioss::Field::MESH);
297  if (dimension_ > 2)
298  stk::io::set_field_role(*facesField_, Ioss::Field::MESH);
299  stk::io::set_field_role(*processorIdField_, Ioss::Field::TRANSIENT);
300  // stk::io::set_field_role(*loadBalField_, Ioss::Field::TRANSIENT);
301  }
302 #endif
303 
304  if (buildRefinementSupport) {
305 #ifdef PANZER_HAVE_PERCEPT
306  refinedMesh_ = Teuchos::rcp(new percept::PerceptMesh(this->getMetaData().get(),
307  this->getBulkData().get(),
308  true));
309 
310  breakPattern_ = Teuchos::rcp(new percept::URP_Heterogeneous_3D(*refinedMesh_));
311 #else
312  TEUCHOS_TEST_FOR_EXCEPTION(true,std::logic_error,
313  "ERROR: Uniform refinement requested. This requires the Percept package to be enabled in Trilinos!");
314 #endif
315  }
316 
317  if(bulkData_==Teuchos::null)
318  instantiateBulkData(*mpiComm_->getRawMpiComm());
319 
320  metaData_->commit();
321 
322  initialized_ = true;
323 }
324 
325 void STK_Interface::initializeFieldsInSTK(const std::map<std::pair<std::string,std::string>,SolutionFieldType*> & nameToField,
326  bool setupIO)
327 {
328  std::set<SolutionFieldType*> uniqueFields;
329  std::map<std::pair<std::string,std::string>,SolutionFieldType*>::const_iterator fieldIter;
330  for(fieldIter=nameToField.begin();fieldIter!=nameToField.end();++fieldIter)
331  uniqueFields.insert(fieldIter->second); // this makes setting up IO easier!
332 
333  {
334  std::set<SolutionFieldType*>::const_iterator uniqueFieldIter;
335  stk::mesh::FieldTraits<SolutionFieldType>::data_type* init_sol = nullptr; // gcc 4.8 hack
336  for(uniqueFieldIter=uniqueFields.begin();uniqueFieldIter!=uniqueFields.end();++uniqueFieldIter)
337  stk::mesh::put_field_on_mesh(*(*uniqueFieldIter),metaData_->universal_part(),init_sol);
338  }
339 
340 #ifdef PANZER_HAVE_IOSS
341  if(setupIO) {
342  // add solution fields
343  std::set<SolutionFieldType*>::const_iterator uniqueFieldIter;
344  for(uniqueFieldIter=uniqueFields.begin();uniqueFieldIter!=uniqueFields.end();++uniqueFieldIter)
345  stk::io::set_field_role(*(*uniqueFieldIter), Ioss::Field::TRANSIENT);
346  }
347 #endif
348 }
349 
350 void STK_Interface::instantiateBulkData(stk::ParallelMachine parallelMach)
351 {
352  TEUCHOS_ASSERT(bulkData_==Teuchos::null);
353  if(mpiComm_==Teuchos::null)
354  mpiComm_ = getSafeCommunicator(parallelMach);
355 
356  bulkData_ = rcp(new stk::mesh::BulkData(*metaData_, *mpiComm_->getRawMpiComm()));
357 }
358 
360 {
361  TEUCHOS_TEST_FOR_EXCEPTION(bulkData_==Teuchos::null,std::logic_error,
362  "STK_Interface: Must call \"initialized\" or \"instantiateBulkData\" before \"beginModification\"");
363 
364  bulkData_->modification_begin();
365 }
366 
368 {
369  TEUCHOS_TEST_FOR_EXCEPTION(bulkData_==Teuchos::null,std::logic_error,
370  "STK_Interface: Must call \"initialized\" or \"instantiateBulkData\" before \"endModification\"");
371 
372  // TODO: Resolving sharing this way comes at a cost in performance. The STK
373  // team has decided that users need to declare their own sharing. We should
374  // find where shared entities are being created in Panzer and declare it.
375  // Once this is done, the extra code below can be deleted.
376 
377  stk::CommSparse comm(bulkData_->parallel());
378 
379  for (int phase=0;phase<2;++phase) {
380  for (int i=0;i<bulkData_->parallel_size();++i) {
381  if ( i != bulkData_->parallel_rank() ) {
382  const stk::mesh::BucketVector& buckets = bulkData_->buckets(stk::topology::NODE_RANK);
383  for (size_t j=0;j<buckets.size();++j) {
384  const stk::mesh::Bucket& bucket = *buckets[j];
385  if ( bucket.owned() ) {
386  for (size_t k=0;k<bucket.size();++k) {
387  stk::mesh::EntityKey key = bulkData_->entity_key(bucket[k]);
388  comm.send_buffer(i).pack<stk::mesh::EntityKey>(key);
389  }
390  }
391  }
392  }
393  }
394 
395  if (phase == 0 ) {
396  comm.allocate_buffers();
397  }
398  else {
399  comm.communicate();
400  }
401  }
402 
403  for (int i=0;i<bulkData_->parallel_size();++i) {
404  if ( i != bulkData_->parallel_rank() ) {
405  while(comm.recv_buffer(i).remaining()) {
406  stk::mesh::EntityKey key;
407  comm.recv_buffer(i).unpack<stk::mesh::EntityKey>(key);
408  stk::mesh::Entity node = bulkData_->get_entity(key);
409  if ( bulkData_->is_valid(node) ) {
410  bulkData_->add_node_sharing(node, i);
411  }
412  }
413  }
414  }
415 
416 
417  bulkData_->modification_end();
418 
421 }
422 
423 void STK_Interface::addNode(stk::mesh::EntityId gid, const std::vector<double> & coord)
424 {
425  TEUCHOS_TEST_FOR_EXCEPTION(not isModifiable(),std::logic_error,
426  "STK_Interface::addNode: STK_Interface must be modifiable to add a node");
427  TEUCHOS_TEST_FOR_EXCEPTION(coord.size()!=getDimension(),std::logic_error,
428  "STK_Interface::addNode: number of coordinates in vector must mation dimension");
429  TEUCHOS_TEST_FOR_EXCEPTION(gid==0,std::logic_error,
430  "STK_Interface::addNode: STK has STUPID restriction of no zero GIDs, pick something else");
431  stk::mesh::EntityRank nodeRank = getNodeRank();
432 
433  stk::mesh::Entity node = bulkData_->declare_entity(nodeRank,gid,nodesPartVec_);
434 
435  // set coordinate vector
436  double * fieldCoords = stk::mesh::field_data(*coordinatesField_,node);
437  for(std::size_t i=0;i<coord.size();++i)
438  fieldCoords[i] = coord[i];
439 }
440 
441 void STK_Interface::addEntityToSideset(stk::mesh::Entity entity,stk::mesh::Part * sideset)
442 {
443  std::vector<stk::mesh::Part*> sidesetV;
444  sidesetV.push_back(sideset);
445 
446  bulkData_->change_entity_parts(entity,sidesetV);
447 }
448 
449 void STK_Interface::addEntityToNodeset(stk::mesh::Entity entity,stk::mesh::Part * nodeset)
450 {
451  std::vector<stk::mesh::Part*> nodesetV;
452  nodesetV.push_back(nodeset);
453 
454  bulkData_->change_entity_parts(entity,nodesetV);
455 }
456 
457 void STK_Interface::addElement(const Teuchos::RCP<ElementDescriptor> & ed,stk::mesh::Part * block)
458 {
459  std::vector<stk::mesh::Part*> blockVec;
460  blockVec.push_back(block);
461 
462  stk::mesh::EntityRank elementRank = getElementRank();
463  stk::mesh::EntityRank nodeRank = getNodeRank();
464  stk::mesh::Entity element = bulkData_->declare_entity(elementRank,ed->getGID(),blockVec);
465 
466  // build relations that give the mesh structure
467  const std::vector<stk::mesh::EntityId> & nodes = ed->getNodes();
468  for(std::size_t i=0;i<nodes.size();++i) {
469  // add element->node relation
470  stk::mesh::Entity node = bulkData_->get_entity(nodeRank,nodes[i]);
471  TEUCHOS_ASSERT(bulkData_->is_valid(node));
472  bulkData_->declare_relation(element,node,i);
473  }
474 
475  ProcIdData * procId = stk::mesh::field_data(*processorIdField_,element);
476  procId[0] = Teuchos::as<ProcIdData>(procRank_);
477 }
478 
479 
481 {
482  // loop over elements
483  stk::mesh::EntityRank edgeRank = getEdgeRank();
484  std::vector<stk::mesh::Entity> localElmts;
485  getMyElements(localElmts);
486  std::vector<stk::mesh::Entity>::const_iterator itr;
487  for(itr=localElmts.begin();itr!=localElmts.end();++itr) {
488  stk::mesh::Entity element = (*itr);
489  stk::mesh::EntityId gid = bulkData_->identifier(element);
490  std::vector<stk::mesh::EntityId> subcellIds;
491  getSubcellIndices(edgeRank,gid,subcellIds);
492 
493  for(std::size_t i=0;i<subcellIds.size();++i) {
494  stk::mesh::Entity edge = bulkData_->get_entity(edgeRank,subcellIds[i]);
495  stk::mesh::Entity const* relations = bulkData_->begin_nodes(edge);
496 
497  double * node_coord_1 = stk::mesh::field_data(*coordinatesField_,relations[0]);
498  double * node_coord_2 = stk::mesh::field_data(*coordinatesField_,relations[1]);
499 
500  // set coordinate vector
501  double * edgeCoords = stk::mesh::field_data(*edgesField_,edge);
502  for(std::size_t j=0;j<getDimension();++j)
503  edgeCoords[j] = (node_coord_1[j]+node_coord_2[j])/2.0;
504  }
505  }
506 }
507 
509 {
510  // loop over elements
511  stk::mesh::EntityRank faceRank = getFaceRank();
512  std::vector<stk::mesh::Entity> localElmts;
513  getMyElements(localElmts);
514  std::vector<stk::mesh::Entity>::const_iterator itr;
515  for(itr=localElmts.begin();itr!=localElmts.end();++itr) {
516  stk::mesh::Entity element = (*itr);
517  stk::mesh::EntityId gid = elementGlobalId(element);
518  std::vector<stk::mesh::EntityId> subcellIds;
519  getSubcellIndices(faceRank,gid,subcellIds);
520 
521  for(std::size_t i=0;i<subcellIds.size();++i) {
522  stk::mesh::Entity face = bulkData_->get_entity(faceRank,subcellIds[i]);
523  stk::mesh::Entity const* relations = bulkData_->begin_nodes(face);
524  const size_t num_relations = bulkData_->num_nodes(face);
525 
526  // set coordinate vector
527  double * faceCoords = stk::mesh::field_data(*facesField_,face);
528  for(std::size_t j=0;j<getDimension();++j){
529  faceCoords[j] = 0.0;
530  for(std::size_t k=0;k<num_relations;++k)
531  faceCoords[j] += stk::mesh::field_data(*coordinatesField_,relations[k])[j];
532  faceCoords[j] /= double(num_relations);
533  }
534  }
535  }
536 }
537 
538 stk::mesh::Entity STK_Interface::findConnectivityById(stk::mesh::Entity src, stk::mesh::EntityRank tgt_rank, unsigned rel_id) const
539 {
540  const size_t num_rels = bulkData_->num_connectivity(src, tgt_rank);
541  stk::mesh::Entity const* relations = bulkData_->begin(src, tgt_rank);
542  stk::mesh::ConnectivityOrdinal const* ordinals = bulkData_->begin_ordinals(src, tgt_rank);
543  for (size_t i = 0; i < num_rels; ++i) {
544  if (ordinals[i] == static_cast<stk::mesh::ConnectivityOrdinal>(rel_id)) {
545  return relations[i];
546  }
547  }
548 
549  return stk::mesh::Entity();
550 }
551 
553 //
554 // writeToExodus()
555 //
557 void
560  const std::string& filename)
561 {
562  setupExodusFile(filename);
563  writeToExodus(0.0);
564 } // end of writeToExodus()
565 
567 //
568 // setupExodusFile()
569 //
571 void
574  const std::string& filename)
575 {
576  using std::runtime_error;
577  using stk::io::StkMeshIoBroker;
578  using stk::mesh::FieldVector;
579  using stk::ParallelMachine;
580  using Teuchos::rcp;
581  PANZER_FUNC_TIME_MONITOR("STK_Interface::setupExodusFile(filename)");
582 #ifdef PANZER_HAVE_IOSS
584  ParallelMachine comm = *mpiComm_->getRawMpiComm();
585  meshData_ = rcp(new StkMeshIoBroker(comm));
586  meshData_->set_bulk_data(bulkData_);
587  meshIndex_ = meshData_->create_output_mesh(filename, stk::io::WRITE_RESULTS);
588  const FieldVector& fields = metaData_->get_fields();
589  for (size_t i(0); i < fields.size(); ++i) {
590  // Do NOT add MESH type stk fields to exodus io, but do add everything
591  // else. This allows us to avoid having to catch a throw for
592  // re-registering coordinates, sidesets, etc... Note that some
593  // fields like LOAD_BAL don't always have a role assigned, so for
594  // roles that point to null, we need to register them as well.
595  auto role = stk::io::get_field_role(*fields[i]);
596  if (role != nullptr) {
597  if (*role != Ioss::Field::MESH)
598  meshData_->add_field(meshIndex_, *fields[i]);
599  } else {
600  meshData_->add_field(meshIndex_, *fields[i]);
601  }
602  }
603 #else
604  TEUCHOS_ASSERT(false)
605 #endif
606 } // end of setupExodusFile()
607 
609 //
610 // writeToExodus()
611 //
613 void
616  double timestep)
617 {
618  using std::cout;
619  using std::endl;
620  using Teuchos::FancyOStream;
621  using Teuchos::rcpFromRef;
622  PANZER_FUNC_TIME_MONITOR("STK_Interface::writeToExodus(timestep)");
623 #ifdef PANZER_HAVE_IOSS
624  if (not meshData_.is_null())
625  {
626  currentStateTime_ = timestep;
627  globalToExodus(GlobalVariable::ADD);
628  meshData_->begin_output_step(meshIndex_, currentStateTime_);
629  meshData_->write_defined_output_fields(meshIndex_);
630  globalToExodus(GlobalVariable::WRITE);
631  meshData_->end_output_step(meshIndex_);
632  }
633  else // if (meshData_.is_null())
634  {
635  FancyOStream out(rcpFromRef(cout));
636  out.setOutputToRootOnly(0);
637  out << "WARNING: Exodus I/O has been disabled or not setup properly; "
638  << "not writing to Exodus." << endl;
639  } // end if (meshData_.is_null()) or not
640 #else
641  TEUCHOS_ASSERT(false);
642 #endif
643 } // end of writeToExodus()
644 
646 //
647 // globalToExodus()
648 //
650 void
651 STK_Interface::
652 globalToExodus(
653  const GlobalVariable& flag)
654 {
655  using std::invalid_argument;
656  using std::string;
657  using Teuchos::Array;
658 
659  // Loop over all the global variables to be added to the Exodus output file.
660  // For each global variable, we determine the data type, and then add or
661  // write it accordingly, depending on the value of flag.
662  for (auto i = globalData_.begin(); i != globalData_.end(); ++i)
663  {
664  const string& name = globalData_.name(i);
665 
666  // Integers.
667  if (globalData_.isType<int>(name))
668  {
669  const auto& value = globalData_.get<int>(name);
670  if (flag == GlobalVariable::ADD)
671  {
672  try
673  {
674  meshData_->add_global(meshIndex_, name, value,
675  stk::util::ParameterType::INTEGER);
676  }
677  catch (...)
678  {
679  return;
680  }
681  }
682  else // if (flag == GlobalVariable::WRITE)
683  meshData_->write_global(meshIndex_, name, value,
684  stk::util::ParameterType::INTEGER);
685  }
686 
687  // Doubles.
688  else if (globalData_.isType<double>(name))
689  {
690  const auto& value = globalData_.get<double>(name);
691  if (flag == GlobalVariable::ADD)
692  {
693  try
694  {
695  meshData_->add_global(meshIndex_, name, value,
696  stk::util::ParameterType::DOUBLE);
697  }
698  catch (...)
699  {
700  return;
701  }
702  }
703  else // if (flag == GlobalVariable::WRITE)
704  meshData_->write_global(meshIndex_, name, value,
705  stk::util::ParameterType::DOUBLE);
706  }
707 
708  // Vectors of integers.
709  else if (globalData_.isType<Array<int>>(name))
710  {
711  const auto& value = globalData_.get<Array<int>>(name).toVector();
712  if (flag == GlobalVariable::ADD)
713  {
714  try
715  {
716  meshData_->add_global(meshIndex_, name, value,
717  stk::util::ParameterType::INTEGERVECTOR);
718  }
719  catch (...)
720  {
721  return;
722  }
723  }
724  else // if (flag == GlobalVariable::WRITE)
725  meshData_->write_global(meshIndex_, name, value,
726  stk::util::ParameterType::INTEGERVECTOR);
727  }
728 
729  // Vectors of doubles.
730  else if (globalData_.isType<Array<double>>(name))
731  {
732  const auto& value = globalData_.get<Array<double>>(name).toVector();
733  if (flag == GlobalVariable::ADD)
734  {
735  try
736  {
737  meshData_->add_global(meshIndex_, name, value,
738  stk::util::ParameterType::DOUBLEVECTOR);
739  }
740  catch (...)
741  {
742  return;
743  }
744  }
745  else // if (flag == GlobalVariable::WRITE)
746  meshData_->write_global(meshIndex_, name, value,
747  stk::util::ParameterType::DOUBLEVECTOR);
748  }
749 
750  // If the data type is something else, throw an exception.
751  else
752  TEUCHOS_TEST_FOR_EXCEPTION(true, invalid_argument,
753  "STK_Interface::globalToExodus(): The global variable to be added " \
754  "to the Exodus output file is of an invalid type. Valid types are " \
755  "int and double, along with std::vectors of those types.")
756  } // end loop over globalData_
757 } // end of globalToExodus()
758 
760 //
761 // addGlobalToExodus()
762 //
764 void
767  const std::string& key,
768  const int& value)
769 {
770  globalData_.set(key, value);
771 } // end of addGlobalToExodus()
772 
774 //
775 // addGlobalToExodus()
776 //
778 void
781  const std::string& key,
782  const double& value)
783 {
784  globalData_.set(key, value);
785 } // end of addGlobalToExodus()
786 
788 //
789 // addGlobalToExodus()
790 //
792 void
795  const std::string& key,
796  const std::vector<int>& value)
797 {
798  using Teuchos::Array;
799  globalData_.set(key, Array<int>(value));
800 } // end of addGlobalToExodus()
801 
803 //
804 // addGlobalToExodus()
805 //
807 void
810  const std::string& key,
811  const std::vector<double>& value)
812 {
813  using Teuchos::Array;
814  globalData_.set(key, Array<double>(value));
815 } // end of addGlobalToExodus()
816 
818 {
819  #ifdef PANZER_HAVE_IOSS
820  return true;
821  #else
822  return false;
823  #endif
824 }
825 
826 void STK_Interface::getElementsSharingNode(stk::mesh::EntityId nodeId,std::vector<stk::mesh::Entity> & elements) const
827 {
828  stk::mesh::EntityRank nodeRank = getNodeRank();
829 
830  // get all relations for node
831  stk::mesh::Entity node = bulkData_->get_entity(nodeRank,nodeId);
832  const size_t numElements = bulkData_->num_elements(node);
833  stk::mesh::Entity const* relations = bulkData_->begin_elements(node);
834 
835  // extract elements sharing nodes
836  elements.insert(elements.end(), relations, relations + numElements);
837 }
838 
839 void STK_Interface::getOwnedElementsSharingNode(stk::mesh::Entity node,std::vector<stk::mesh::Entity> & elements,
840  std::vector<int> & relIds) const
841 {
842  // get all relations for node
843  const size_t numElements = bulkData_->num_elements(node);
844  stk::mesh::Entity const* relations = bulkData_->begin_elements(node);
845  stk::mesh::ConnectivityOrdinal const* rel_ids = bulkData_->begin_element_ordinals(node);
846 
847  // extract elements sharing nodes
848  for (size_t i = 0; i < numElements; ++i) {
849  stk::mesh::Entity element = relations[i];
850 
851  // if owned by this processor
852  if(bulkData_->parallel_owner_rank(element) == static_cast<int>(procRank_)) {
853  elements.push_back(element);
854  relIds.push_back(rel_ids[i]);
855  }
856  }
857 }
858 
859 void STK_Interface::getOwnedElementsSharingNode(stk::mesh::EntityId nodeId,std::vector<stk::mesh::Entity> & elements,
860  std::vector<int> & relIds, unsigned int matchType) const
861 {
862  stk::mesh::EntityRank rank;
863  if(matchType == 0)
864  rank = getNodeRank();
865  else if(matchType == 1)
866  rank = getEdgeRank();
867  else if(matchType == 2)
868  rank = getFaceRank();
869  else
870  TEUCHOS_ASSERT(false);
871 
872  stk::mesh::Entity node = bulkData_->get_entity(rank,nodeId);
873 
874  getOwnedElementsSharingNode(node,elements,relIds);
875 }
876 
877 void STK_Interface::getElementsSharingNodes(const std::vector<stk::mesh::EntityId> nodeIds,std::vector<stk::mesh::Entity> & elements) const
878 {
879  std::vector<stk::mesh::Entity> current;
880 
881  getElementsSharingNode(nodeIds[0],current); // fill it with elements touching first node
882  std::sort(current.begin(),current.end()); // sort for intersection on the pointer
883 
884  // find intersection with remaining nodes
885  for(std::size_t n=1;n<nodeIds.size();++n) {
886  // get elements associated with next node
887  std::vector<stk::mesh::Entity> nextNode;
888  getElementsSharingNode(nodeIds[n],nextNode); // fill it with elements touching first node
889  std::sort(nextNode.begin(),nextNode.end()); // sort for intersection on the pointer ID
890 
891  // intersect next node elements with current element list
892  std::vector<stk::mesh::Entity> intersection(std::min(nextNode.size(),current.size()));
893  std::vector<stk::mesh::Entity>::const_iterator endItr
894  = std::set_intersection(current.begin(),current.end(),
895  nextNode.begin(),nextNode.end(),
896  intersection.begin());
897  std::size_t newLength = endItr-intersection.begin();
898  intersection.resize(newLength);
899 
900  // store intersection
901  current.clear();
902  current = intersection;
903  }
904 
905  // return the elements computed
906  elements = current;
907 }
908 
909 void STK_Interface::getNodeIdsForElement(stk::mesh::Entity element,std::vector<stk::mesh::EntityId> & nodeIds) const
910 {
911  stk::mesh::Entity const* nodeRel = getBulkData()->begin_nodes(element);
912  const size_t numNodes = getBulkData()->num_nodes(element);
913 
914  nodeIds.reserve(numNodes);
915  for(size_t i = 0; i < numNodes; ++i) {
916  nodeIds.push_back(elementGlobalId(nodeRel[i]));
917  }
918 }
919 
921 {
922  entityCounts_.clear();
923  stk::mesh::comm_mesh_counts(*bulkData_,entityCounts_);
924 }
925 
927 {
928  // developed to mirror "comm_mesh_counts" in stk_mesh/base/Comm.cpp
929 
930  const auto entityRankCount = metaData_->entity_rank_count();
931  const size_t commCount = 10; // entityRankCount
932 
933  TEUCHOS_ASSERT(entityRankCount<10);
934 
935  // stk::ParallelMachine mach = bulkData_->parallel();
936  stk::ParallelMachine mach = *mpiComm_->getRawMpiComm();
937 
938  std::vector<stk::mesh::EntityId> local(commCount,0);
939 
940  // determine maximum ID for this processor for each entity type
941  stk::mesh::Selector ownedPart = metaData_->locally_owned_part();
942  for(stk::mesh::EntityRank i=stk::topology::NODE_RANK;
943  i < static_cast<stk::mesh::EntityRank>(entityRankCount); ++i) {
944  std::vector<stk::mesh::Entity> entities;
945 
946  stk::mesh::get_selected_entities(ownedPart,bulkData_->buckets(i),entities);
947 
948  // determine maximum ID for this processor
949  std::vector<stk::mesh::Entity>::const_iterator itr;
950  for(itr=entities.begin();itr!=entities.end();++itr) {
951  stk::mesh::EntityId id = bulkData_->identifier(*itr);
952  if(id>local[i])
953  local[i] = id;
954  }
955  }
956 
957  // get largest IDs across processors
958  stk::all_reduce(mach,stk::ReduceMax<10>(&local[0]));
959  maxEntityId_.assign(local.begin(),local.begin()+entityRankCount+1);
960 }
961 
962 std::size_t STK_Interface::getEntityCounts(unsigned entityRank) const
963 {
964  TEUCHOS_TEST_FOR_EXCEPTION(entityRank>=entityCounts_.size(),std::logic_error,
965  "STK_Interface::getEntityCounts: Entity counts do not include rank: " << entityRank);
966 
967  return entityCounts_[entityRank];
968 }
969 
970 stk::mesh::EntityId STK_Interface::getMaxEntityId(unsigned entityRank) const
971 {
972  TEUCHOS_TEST_FOR_EXCEPTION(entityRank>=maxEntityId_.size(),std::logic_error,
973  "STK_Interface::getMaxEntityId: Max entity ids do not include rank: " << entityRank);
974 
975  return maxEntityId_[entityRank];
976 }
977 
979 {
980  stk::mesh::PartVector emptyPartVector;
981  stk::mesh::create_adjacent_entities(*bulkData_,emptyPartVector);
982 
985 
986  addEdges();
987  if (dimension_ > 2)
988  addFaces();
989 }
990 
991 const double * STK_Interface::getNodeCoordinates(stk::mesh::EntityId nodeId) const
992 {
993  stk::mesh::Entity node = bulkData_->get_entity(getNodeRank(),nodeId);
994  return stk::mesh::field_data(*coordinatesField_,node);
995 }
996 
997 const double * STK_Interface::getNodeCoordinates(stk::mesh::Entity node) const
998 {
999  return stk::mesh::field_data(*coordinatesField_,node);
1000 }
1001 
1002 void STK_Interface::getSubcellIndices(unsigned entityRank,stk::mesh::EntityId elementId,
1003  std::vector<stk::mesh::EntityId> & subcellIds) const
1004 {
1005  stk::mesh::EntityRank elementRank = getElementRank();
1006  stk::mesh::Entity cell = bulkData_->get_entity(elementRank,elementId);
1007 
1008  TEUCHOS_TEST_FOR_EXCEPTION(!bulkData_->is_valid(cell),std::logic_error,
1009  "STK_Interface::getSubcellIndices: could not find element requested (GID = " << elementId << ")");
1010 
1011  const size_t numSubcells = bulkData_->num_connectivity(cell, static_cast<stk::mesh::EntityRank>(entityRank));
1012  stk::mesh::Entity const* subcells = bulkData_->begin(cell, static_cast<stk::mesh::EntityRank>(entityRank));
1013  subcellIds.clear();
1014  subcellIds.resize(numSubcells,0);
1015 
1016  // loop over relations and fill subcell vector
1017  for(size_t i = 0; i < numSubcells; ++i) {
1018  stk::mesh::Entity subcell = subcells[i];
1019  subcellIds[i] = bulkData_->identifier(subcell);
1020  }
1021 }
1022 
1023 void STK_Interface::getMyElements(std::vector<stk::mesh::Entity> & elements) const
1024 {
1025  // setup local ownership
1026  stk::mesh::Selector ownedPart = metaData_->locally_owned_part();
1027 
1028  // grab elements
1029  stk::mesh::EntityRank elementRank = getElementRank();
1030  stk::mesh::get_selected_entities(ownedPart,bulkData_->buckets(elementRank),elements);
1031 }
1032 
1033 void STK_Interface::getMyElements(const std::string & blockID,std::vector<stk::mesh::Entity> & elements) const
1034 {
1035  stk::mesh::Part * elementBlock = getElementBlockPart(blockID);
1036 
1037  TEUCHOS_TEST_FOR_EXCEPTION(elementBlock==0,std::logic_error,"Could not find element block \"" << blockID << "\"");
1038 
1039  // setup local ownership
1040  // stk::mesh::Selector block = *elementBlock;
1041  stk::mesh::Selector ownedBlock = metaData_->locally_owned_part() & (*elementBlock);
1042 
1043  // grab elements
1044  stk::mesh::EntityRank elementRank = getElementRank();
1045  stk::mesh::get_selected_entities(ownedBlock,bulkData_->buckets(elementRank),elements);
1046 }
1047 
1048 void STK_Interface::getNeighborElements(std::vector<stk::mesh::Entity> & elements) const
1049 {
1050  // setup local ownership
1051  stk::mesh::Selector neighborBlock = (!metaData_->locally_owned_part());
1052 
1053  // grab elements
1054  stk::mesh::EntityRank elementRank = getElementRank();
1055  stk::mesh::get_selected_entities(neighborBlock,bulkData_->buckets(elementRank),elements);
1056 }
1057 
1058 void STK_Interface::getNeighborElements(const std::string & blockID,std::vector<stk::mesh::Entity> & elements) const
1059 {
1060  stk::mesh::Part * elementBlock = getElementBlockPart(blockID);
1061 
1062  TEUCHOS_TEST_FOR_EXCEPTION(elementBlock==0,std::logic_error,"Could not find element block \"" << blockID << "\"");
1063 
1064  // setup local ownership
1065  stk::mesh::Selector neighborBlock = (!metaData_->locally_owned_part()) & (*elementBlock);
1066 
1067  // grab elements
1068  stk::mesh::EntityRank elementRank = getElementRank();
1069  stk::mesh::get_selected_entities(neighborBlock,bulkData_->buckets(elementRank),elements);
1070 }
1071 
1072 void STK_Interface::getMySides(const std::string & sideName,std::vector<stk::mesh::Entity> & sides) const
1073 {
1074  stk::mesh::Part * sidePart = getSideset(sideName);
1075  TEUCHOS_TEST_FOR_EXCEPTION(sidePart==0,std::logic_error,
1076  "Unknown side set \"" << sideName << "\"");
1077 
1078  stk::mesh::Selector side = *sidePart;
1079  stk::mesh::Selector ownedBlock = metaData_->locally_owned_part() & side;
1080 
1081  // grab elements
1082  stk::mesh::get_selected_entities(ownedBlock,bulkData_->buckets(getSideRank()),sides);
1083 }
1084 
1085 void STK_Interface::getMySides(const std::string & sideName,const std::string & blockName,std::vector<stk::mesh::Entity> & sides) const
1086 {
1087  stk::mesh::Part * sidePart = getSideset(sideName);
1088  stk::mesh::Part * elmtPart = getElementBlockPart(blockName);
1090  "Unknown side set \"" << sideName << "\"");
1092  "Unknown element block \"" << blockName << "\"");
1093 
1094  stk::mesh::Selector side = *sidePart;
1095  stk::mesh::Selector block = *elmtPart;
1096  stk::mesh::Selector ownedBlock = metaData_->locally_owned_part() & block & side;
1097 
1098  // grab elements
1099  stk::mesh::get_selected_entities(ownedBlock,bulkData_->buckets(getSideRank()),sides);
1100 }
1101 
1102 void STK_Interface::getAllSides(const std::string & sideName,std::vector<stk::mesh::Entity> & sides) const
1103 {
1104  stk::mesh::Part * sidePart = getSideset(sideName);
1105  TEUCHOS_TEST_FOR_EXCEPTION(sidePart==0,std::logic_error,
1106  "Unknown side set \"" << sideName << "\"");
1107 
1108  stk::mesh::Selector side = *sidePart;
1109 
1110  // grab elements
1111  stk::mesh::get_selected_entities(side,bulkData_->buckets(getSideRank()),sides);
1112 }
1113 
1114 void STK_Interface::getAllSides(const std::string & sideName,const std::string & blockName,std::vector<stk::mesh::Entity> & sides) const
1115 {
1116  stk::mesh::Part * sidePart = getSideset(sideName);
1117  stk::mesh::Part * elmtPart = getElementBlockPart(blockName);
1119  "Unknown side set \"" << sideName << "\"");
1121  "Unknown element block \"" << blockName << "\"");
1122 
1123  stk::mesh::Selector side = *sidePart;
1124  stk::mesh::Selector block = *elmtPart;
1125  stk::mesh::Selector sideBlock = block & side;
1126 
1127  // grab elements
1128  stk::mesh::get_selected_entities(sideBlock,bulkData_->buckets(getSideRank()),sides);
1129 }
1130 
1131 void STK_Interface::getMyNodes(const std::string & nodesetName,const std::string & blockName,std::vector<stk::mesh::Entity> & nodes) const
1132 {
1133  stk::mesh::Part * nodePart = getNodeset(nodesetName);
1134  stk::mesh::Part * elmtPart = getElementBlockPart(blockName);
1136  "Unknown node set \"" << nodesetName << "\"");
1138  "Unknown element block \"" << blockName << "\"");
1139 
1140  stk::mesh::Selector nodeset = *nodePart;
1141  stk::mesh::Selector block = *elmtPart;
1142  stk::mesh::Selector ownedBlock = metaData_->locally_owned_part() & block & nodeset;
1143 
1144  // grab elements
1145  stk::mesh::get_selected_entities(ownedBlock,bulkData_->buckets(getNodeRank()),nodes);
1146 }
1147 
1148 void STK_Interface::getElementBlockNames(std::vector<std::string> & names) const
1149 {
1150  // TEUCHOS_ASSERT(initialized_); // all blocks must have been added
1151 
1152  names.clear();
1153 
1154  // fill vector with automagically ordered string values
1155  std::map<std::string, stk::mesh::Part*>::const_iterator blkItr; // Element blocks
1156  for(blkItr=elementBlocks_.begin();blkItr!=elementBlocks_.end();++blkItr)
1157  names.push_back(blkItr->first);
1158 }
1159 
1160 void STK_Interface::getSidesetNames(std::vector<std::string> & names) const
1161 {
1162  // TEUCHOS_ASSERT(initialized_); // all blocks must have been added
1163 
1164  names.clear();
1165 
1166  // fill vector with automagically ordered string values
1167  std::map<std::string, stk::mesh::Part*>::const_iterator sideItr; // Element blocks
1168  for(sideItr=sidesets_.begin();sideItr!=sidesets_.end();++sideItr)
1169  names.push_back(sideItr->first);
1170 }
1171 
1172 void STK_Interface::getNodesetNames(std::vector<std::string> & names) const
1173 {
1174  names.clear();
1175 
1176  // fill vector with automagically ordered string values
1177  std::map<std::string, stk::mesh::Part*>::const_iterator nodeItr; // Element blocks
1178  for(nodeItr=nodesets_.begin();nodeItr!=nodesets_.end();++nodeItr)
1179  names.push_back(nodeItr->first);
1180 }
1181 
1182 std::size_t STK_Interface::elementLocalId(stk::mesh::Entity elmt) const
1183 {
1184  return elementLocalId(bulkData_->identifier(elmt));
1185  // const std::size_t * fieldCoords = stk::mesh::field_data(*localIdField_,*elmt);
1186  // return fieldCoords[0];
1187 }
1188 
1189 std::size_t STK_Interface::elementLocalId(stk::mesh::EntityId gid) const
1190 {
1191  // stk::mesh::EntityRank elementRank = getElementRank();
1192  // stk::mesh::Entity elmt = bulkData_->get_entity(elementRank,gid);
1193  // TEUCHOS_ASSERT(elmt->owner_rank()==procRank_);
1194  // return elementLocalId(elmt);
1195  std::unordered_map<stk::mesh::EntityId,std::size_t>::const_iterator itr = localIDHash_.find(gid);
1196  TEUCHOS_ASSERT(itr!=localIDHash_.end());
1197  return itr->second;
1198 }
1199 
1200 
1201 std::string STK_Interface::containingBlockId(stk::mesh::Entity elmt) const
1202 {
1203  for(const auto & eb_pair : elementBlocks_)
1204  if(bulkData_->bucket(elmt).member(*(eb_pair.second)))
1205  return eb_pair.first;
1206  return "";
1207 }
1208 
1209 stk::mesh::Field<double> * STK_Interface::getSolutionField(const std::string & fieldName,
1210  const std::string & blockId) const
1211 {
1212  // look up field in map
1213  std::map<std::pair<std::string,std::string>, SolutionFieldType*>::const_iterator
1214  iter = fieldNameToSolution_.find(std::make_pair(fieldName,blockId));
1215 
1216  // check to make sure field was actually found
1217  TEUCHOS_TEST_FOR_EXCEPTION(iter==fieldNameToSolution_.end(),std::runtime_error,
1218  "Solution field name \"" << fieldName << "\" in block ID \"" << blockId << "\" was not found");
1219 
1220  return iter->second;
1221 }
1222 
1223 stk::mesh::Field<double> * STK_Interface::getCellField(const std::string & fieldName,
1224  const std::string & blockId) const
1225 {
1226  // look up field in map
1227  std::map<std::pair<std::string,std::string>, SolutionFieldType*>::const_iterator
1228  iter = fieldNameToCellField_.find(std::make_pair(fieldName,blockId));
1229 
1230  // check to make sure field was actually found
1231  TEUCHOS_TEST_FOR_EXCEPTION(iter==fieldNameToCellField_.end(),std::runtime_error,
1232  "Cell field named \"" << fieldName << "\" in block ID \"" << blockId << "\" was not found");
1233 
1234  return iter->second;
1235 }
1236 
1238 {
1239  using Teuchos::RCP;
1240  using Teuchos::rcp;
1241 
1242  if(orderedElementVector_==Teuchos::null) {
1243  // safe because essentially this is a call to modify a mutable object
1244  const_cast<STK_Interface*>(this)->buildLocalElementIDs();
1245  }
1246 
1248 }
1249 
1250 void STK_Interface::addElementBlock(const std::string & name,const CellTopologyData * ctData)
1251 {
1253 
1254  stk::mesh::Part * block = metaData_->get_part(name);
1255  if(block==0) {
1256  block = &metaData_->declare_part_with_topology(name, stk::mesh::get_topology(shards::CellTopology(ctData), dimension_));
1257  }
1258 
1259  // construct cell topology object for this block
1261  = Teuchos::rcp(new shards::CellTopology(ctData));
1262 
1263  // add element block part and cell topology
1264  elementBlocks_.insert(std::make_pair(name,block));
1265  elementBlockCT_.insert(std::make_pair(name,ct));
1266 }
1267 
1269 {
1270  dimension_ = metaData_->spatial_dimension();
1271 
1272  // declare coordinates and node parts
1273  coordinatesField_ = &metaData_->declare_field<VectorFieldType>(stk::topology::NODE_RANK, coordsString);
1274  edgesField_ = &metaData_->declare_field<VectorFieldType>(stk::topology::EDGE_RANK, edgesString);
1275  if (dimension_ > 2)
1276  facesField_ = &metaData_->declare_field<VectorFieldType>(stk::topology::FACE_RANK, facesString);
1277  processorIdField_ = &metaData_->declare_field<ProcIdFieldType>(stk::topology::ELEMENT_RANK, "PROC_ID");
1278  loadBalField_ = &metaData_->declare_field<SolutionFieldType>(stk::topology::ELEMENT_RANK, "LOAD_BAL");
1279 
1280  // stk::mesh::put_field( *coordinatesField_ , getNodeRank() , metaData_->universal_part() );
1281 
1282  nodesPart_ = &metaData_->declare_part(nodesString,getNodeRank());
1283  nodesPartVec_.push_back(nodesPart_);
1284  edgesPart_ = &metaData_->declare_part(edgesString,getEdgeRank());
1285  edgesPartVec_.push_back(edgesPart_);
1286  if (dimension_ > 2) {
1287  facesPart_ = &metaData_->declare_part(facesString,getFaceRank());
1288  facesPartVec_.push_back(facesPart_);
1289  }
1290 }
1291 
1293 {
1294  currentLocalId_ = 0;
1295 
1296  orderedElementVector_ = Teuchos::null; // forces rebuild of ordered lists
1297 
1298  // might be better (faster) to do this by buckets
1299  std::vector<stk::mesh::Entity> elements;
1300  getMyElements(elements);
1301 
1302  for(std::size_t index=0;index<elements.size();++index) {
1303  stk::mesh::Entity element = elements[index];
1304 
1305  // set processor rank
1306  ProcIdData * procId = stk::mesh::field_data(*processorIdField_,element);
1307  procId[0] = Teuchos::as<ProcIdData>(procRank_);
1308 
1309  localIDHash_[bulkData_->identifier(element)] = currentLocalId_;
1310 
1311  currentLocalId_++;
1312  }
1313 
1314  // copy elements into the ordered element vector
1315  orderedElementVector_ = Teuchos::rcp(new std::vector<stk::mesh::Entity>(elements));
1316 
1317  elements.clear();
1318  getNeighborElements(elements);
1319 
1320  for(std::size_t index=0;index<elements.size();++index) {
1321  stk::mesh::Entity element = elements[index];
1322 
1323  // set processor rank
1324  ProcIdData * procId = stk::mesh::field_data(*processorIdField_,element);
1325  procId[0] = Teuchos::as<ProcIdData>(procRank_);
1326 
1327  localIDHash_[bulkData_->identifier(element)] = currentLocalId_;
1328 
1329  currentLocalId_++;
1330  }
1331 
1332  orderedElementVector_->insert(orderedElementVector_->end(),elements.begin(),elements.end());
1333 }
1334 
1336 {
1337  std::vector<std::string> names;
1338  getElementBlockNames(names);
1339 
1340  for(std::size_t b=0;b<names.size();b++) {
1341  // find user specified block weight, otherwise use 1.0
1342  std::map<std::string,double>::const_iterator bw_itr = blockWeights_.find(names[b]);
1343  double blockWeight = (bw_itr!=blockWeights_.end()) ? bw_itr->second : 1.0;
1344 
1345  std::vector<stk::mesh::Entity> elements;
1346  getMyElements(names[b],elements);
1347 
1348  for(std::size_t index=0;index<elements.size();++index) {
1349  // set local element ID
1350  double * loadBal = stk::mesh::field_data(*loadBalField_,elements[index]);
1351  loadBal[0] = blockWeight;
1352  }
1353  }
1354 }
1355 
1356 bool
1357 STK_Interface::isMeshCoordField(const std::string & eBlock,
1358  const std::string & fieldName,
1359  int & axis) const
1360 {
1361  std::map<std::string,std::vector<std::string> >::const_iterator blkItr = meshCoordFields_.find(eBlock);
1362  if(blkItr==meshCoordFields_.end()) {
1363  return false;
1364  }
1365 
1366  axis = 0;
1367  for(axis=0;axis<Teuchos::as<int>(blkItr->second.size());axis++) {
1368  if(blkItr->second[axis]==fieldName)
1369  break; // found axis, break
1370  }
1371 
1372  if(axis>=Teuchos::as<int>(blkItr->second.size()))
1373  return false;
1374 
1375  return true;
1376 }
1377 
1378 std::pair<Teuchos::RCP<std::vector<std::pair<std::size_t,std::size_t> > >, Teuchos::RCP<std::vector<unsigned int> > >
1380 {
1382  Teuchos::RCP<std::vector<unsigned int > > type_vec = rcp(new std::vector<unsigned int>);
1383  const std::vector<Teuchos::RCP<const PeriodicBC_MatcherBase> > & matchers = getPeriodicBCVector();
1384 
1385  // build up the vectors by looping over the matched pair
1386  for(std::size_t m=0;m<matchers.size();m++){
1387  vec = matchers[m]->getMatchedPair(*this,vec);
1388  unsigned int type;
1389  if(matchers[m]->getType() == "coord")
1390  type = 0;
1391  else if(matchers[m]->getType() == "edge")
1392  type = 1;
1393  else if(matchers[m]->getType() == "face")
1394  type = 2;
1395  else
1396  TEUCHOS_ASSERT(false);
1397  type_vec->insert(type_vec->begin(),vec->size()-type_vec->size(),type);
1398  }
1399 
1400  return std::make_pair(vec,type_vec);
1401 }
1402 
1403 bool STK_Interface::validBlockId(const std::string & blockId) const
1404 {
1405  std::map<std::string, stk::mesh::Part*>::const_iterator blkItr = elementBlocks_.find(blockId);
1406 
1407  return blkItr!=elementBlocks_.end();
1408 }
1409 
1410 void STK_Interface::print(std::ostream & os) const
1411 {
1412  std::vector<std::string> blockNames, sidesetNames, nodesetNames;
1413 
1414  getElementBlockNames(blockNames);
1415  getSidesetNames(sidesetNames);
1416  getNodesetNames(nodesetNames);
1417 
1418  os << "STK Mesh data:\n";
1419  os << " Spatial dim = " << getDimension() << "\n";
1420  if(getDimension()==2)
1421  os << " Entity counts (Nodes, Edges, Cells) = ( "
1422  << getEntityCounts(getNodeRank()) << ", "
1423  << getEntityCounts(getEdgeRank()) << ", "
1424  << getEntityCounts(getElementRank()) << " )\n";
1425  else if(getDimension()==3)
1426  os << " Entity counts (Nodes, Edges, Faces, Cells) = ( "
1427  << getEntityCounts(getNodeRank()) << ", "
1428  << getEntityCounts(getEdgeRank()) << ", "
1429  << getEntityCounts(getSideRank()) << ", "
1430  << getEntityCounts(getElementRank()) << " )\n";
1431  else
1432  os << " Entity counts (Nodes, Cells) = ( "
1433  << getEntityCounts(getNodeRank()) << ", "
1434  << getEntityCounts(getElementRank()) << " )\n";
1435 
1436  os << " Element blocks = ";
1437  for(std::size_t i=0;i<blockNames.size();i++)
1438  os << "\"" << blockNames[i] << "\" ";
1439  os << "\n";
1440  os << " Sidesets = ";
1441  for(std::size_t i=0;i<sidesetNames.size();i++)
1442  os << "\"" << sidesetNames[i] << "\" ";
1443  os << "\n";
1444  os << " Nodesets = ";
1445  for(std::size_t i=0;i<nodesetNames.size();i++)
1446  os << "\"" << nodesetNames[i] << "\" ";
1447  os << std::endl;
1448 
1449  // print out periodic boundary conditions
1450  const std::vector<Teuchos::RCP<const PeriodicBC_MatcherBase> > & bcVector
1451  = getPeriodicBCVector();
1452  if(bcVector.size()!=0) {
1453  os << " Periodic BCs:\n";
1454  for(std::size_t i=0;i<bcVector.size();i++)
1455  os << " " << bcVector[i]->getString() << "\n";
1456  os << std::endl;
1457  }
1458 }
1459 
1460 void STK_Interface::printMetaData(std::ostream & os) const
1461 {
1462  std::vector<std::string> blockNames, sidesetNames, nodesetNames;
1463 
1464  getElementBlockNames(blockNames);
1465  getSidesetNames(sidesetNames);
1466  getNodesetNames(nodesetNames);
1467 
1468  os << "STK Meta data:\n";
1469  os << " Element blocks = ";
1470  for(std::size_t i=0;i<blockNames.size();i++)
1471  os << "\"" << blockNames[i] << "\" ";
1472  os << "\n";
1473  os << " Sidesets = ";
1474  for(std::size_t i=0;i<sidesetNames.size();i++)
1475  os << "\"" << sidesetNames[i] << "\" ";
1476  os << "\n";
1477  os << " Nodesets = ";
1478  for(std::size_t i=0;i<nodesetNames.size();i++)
1479  os << "\"" << nodesetNames[i] << "\" ";
1480  os << std::endl;
1481 
1482  // print out periodic boundary conditions
1483  const std::vector<Teuchos::RCP<const PeriodicBC_MatcherBase> > & bcVector
1484  = getPeriodicBCVector();
1485  if(bcVector.size()!=0) {
1486  os << " Periodic BCs:\n";
1487  for(std::size_t i=0;i<bcVector.size();i++)
1488  os << " " << bcVector[i]->getString() << "\n";
1489  os << std::endl;
1490  }
1491 
1492  // print all available fields in meta data
1493  os << " Fields = ";
1494  const stk::mesh::FieldVector & fv = metaData_->get_fields();
1495  for(std::size_t i=0;i<fv.size();i++)
1496  os << "\"" << fv[i]->name() << "\" ";
1497  os << std::endl;
1498 }
1499 
1501 {
1502  std::map<std::string, Teuchos::RCP<shards::CellTopology> >::const_iterator itr;
1503  itr = elementBlockCT_.find(eBlock);
1504 
1505  if(itr==elementBlockCT_.end()) {
1506  std::stringstream ss;
1507  printMetaData(ss);
1508  TEUCHOS_TEST_FOR_EXCEPTION(itr==elementBlockCT_.end(),std::logic_error,
1509  "STK_Interface::getCellTopology: No such element block \"" +eBlock +"\" available.\n\n"
1510  << "STK Meta Data follows: \n" << ss.str());
1511  }
1512 
1513  return itr->second;
1514 }
1515 
1517 {
1518  MPI_Comm newComm;
1519  const int err = MPI_Comm_dup (parallelMach, &newComm);
1520  TEUCHOS_TEST_FOR_EXCEPTION(err != MPI_SUCCESS, std::runtime_error,
1521  "panzer::STK_Interface: MPI_Comm_dup failed with error \""
1522  << Teuchos::mpiErrorCodeToString (err) << "\".");
1523 
1524  return Teuchos::rcp(new Teuchos::MpiComm<int>(Teuchos::opaqueWrapper (newComm,MPI_Comm_free)));
1525 }
1526 
1528 {
1529  TEUCHOS_TEST_FOR_EXCEPTION(true, std::runtime_error, "Rebalance not currently supported");
1530 #if 0
1531  // make sure weights have been set (a local operation)
1533 
1534  stk::mesh::Selector selector(getMetaData()->universal_part());
1535  stk::mesh::Selector owned_selector(getMetaData()->locally_owned_part());
1536 
1537  Teuchos::FancyOStream out(Teuchos::rcpFromRef(std::cout));
1538  out.setOutputToRootOnly(0);
1539  out << "Load balance before: " << stk::rebalance::check_balance(*getBulkData(), loadBalField_, getElementRank(), &selector) << std::endl;
1540 
1541  // perform reblance
1542  Teuchos::ParameterList graph;
1543  if(params.begin()!=params.end())
1544  graph.sublist(stk::rebalance::Zoltan::default_parameters_name()) = params;
1545  stk::rebalance::Zoltan zoltan_partition(*mpiComm_->getRawMpiComm(), getDimension(), graph);
1546  stk::rebalance::rebalance(*getBulkData(), owned_selector, &getCoordinatesField(), loadBalField_, zoltan_partition);
1547 
1548  out << "Load balance after: " << stk::rebalance::check_balance(*getBulkData(), loadBalField_, getElementRank(), &selector) << std::endl;
1549 
1550  currentLocalId_ = 0;
1551  orderedElementVector_ = Teuchos::null; // forces rebuild of ordered lists
1552 #endif
1553 }
1554 
1557 {
1558  TEUCHOS_ASSERT(this->isInitialized());
1559  return mpiComm_;
1560 }
1561 
1562 void STK_Interface::refineMesh(const int numberOfLevels, const bool deleteParentElements) {
1563 #ifdef PANZER_HAVE_PERCEPT
1564  TEUCHOS_TEST_FOR_EXCEPTION(numberOfLevels < 1,std::runtime_error,
1565  "ERROR: Number of levels for uniform refinement must be greater than 0");
1566  TEUCHOS_ASSERT(nonnull(refinedMesh_));
1567  TEUCHOS_ASSERT(nonnull(breakPattern_));
1568 
1569  refinedMesh_->setCoordinatesField();
1570 
1571  percept::UniformRefiner breaker(*refinedMesh_,*breakPattern_);
1572 
1573  for (int i=0; i < numberOfLevels; ++i)
1574  breaker.doBreak();
1575 
1576  if (deleteParentElements)
1577  breaker.deleteParentElements();
1578 
1579 #else
1580  TEUCHOS_TEST_FOR_EXCEPTION(true,std::logic_error,
1581  "ERROR: Uniform refinement requested. This requires the Percept package to be enabled in Trilinos!");
1582 #endif
1583 }
1584 
1585 
1586 } // end namespace panzer_stk
Teuchos::RCP< ElementDescriptor > buildElementDescriptor(stk::mesh::EntityId elmtId, std::vector< stk::mesh::EntityId > &nodes)
Teuchos::RCP< Teuchos::MpiComm< int > > getSafeCommunicator(stk::ParallelMachine parallelMach) const
void getSidesetNames(std::vector< std::string > &name) const
RCP< const T > getConst() const
std::pair< Teuchos::RCP< std::vector< std::pair< std::size_t, std::size_t > > >, Teuchos::RCP< std::vector< unsigned int > > > getPeriodicNodePairing() const
void addNodeset(const std::string &name)
Teuchos::RCP< stk::mesh::BulkData > bulkData_
void print(std::ostream &os) const
Teuchos::RCP< const std::vector< stk::mesh::Entity > > getElementsOrderedByLID() const
std::vector< stk::mesh::Part * > facesPartVec_
std::map< std::pair< std::string, std::string >, SolutionFieldType * > fieldNameToSolution_
void addEntityToNodeset(stk::mesh::Entity entity, stk::mesh::Part *nodeset)
void addGlobalToExodus(const std::string &key, const int &value)
Add an int global variable to the information to be written to the Exodus output file.
void getElementBlockNames(std::vector< std::string > &names) const
void getSubcellIndices(unsigned entityRank, stk::mesh::EntityId elementId, std::vector< stk::mesh::EntityId > &subcellIds) 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
stk::mesh::Field< double > SolutionFieldType
stk::mesh::Part * getElementBlockPart(const std::string &name) const
get the block count
std::vector< stk::mesh::Part * > edgesPartVec_
std::map< std::string, double > blockWeights_
basic_FancyOStream< char > FancyOStream
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
std::size_t getEntityCounts(unsigned entityRank) const
get the global counts for the entity of specified rank
void printMetaData(std::ostream &os) const
void addEntityToSideset(stk::mesh::Entity entity, stk::mesh::Part *sideset)
void addSolutionField(const std::string &fieldName, const std::string &blockId)
std::map< std::string, Teuchos::RCP< shards::CellTopology > > elementBlockCT_
stk::mesh::EntityId getMaxEntityId(unsigned entityRank) const
get max entity ID of type entityRank
std::map< std::string, stk::mesh::Part * > sidesets_
std::vector< std::size_t > entityCounts_
bool validBlockId(const std::string &blockId) const
stk::mesh::Part * getNodeset(const std::string &name) const
stk::mesh::Field< double > * getSolutionField(const std::string &fieldName, const std::string &blockId) const
const double * getNodeCoordinates(stk::mesh::EntityId nodeId) const
stk::mesh::Part * getSideset(const std::string &name) const
static const std::string nodesString
void instantiateBulkData(stk::ParallelMachine parallelMach)
void getMySides(const std::string &sideName, std::vector< stk::mesh::Entity > &sides) const
std::size_t elementLocalId(stk::mesh::Entity elmt) 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_
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
stk::mesh::EntityRank getSideRank() const
void initialize(stk::ParallelMachine parallelMach, bool setupIO=true, const bool buildRefinementSupport=false)
static const std::string coordsString
stk::mesh::Field< double, stk::mesh::Cartesian > VectorFieldType
Teuchos::RCP< stk::mesh::BulkData > getBulkData() const
void writeToExodus(const std::string &filename)
Write this mesh and associated fields to the given output file.
stk::mesh::EntityId elementGlobalId(std::size_t lid) const
std::map< std::string, stk::mesh::Part * > elementBlocks_
std::vector< stk::mesh::Part * > nodesPartVec_
static const std::string edgesString
basic_FancyOStream & setOutputToRootOnly(const int rootRank)
std::string containingBlockId(stk::mesh::Entity elmt) const
stk::mesh::EntityRank getFaceRank() const
void getOwnedElementsSharingNode(stk::mesh::Entity node, std::vector< stk::mesh::Entity > &elements, std::vector< int > &relIds) const
void addNode(stk::mesh::EntityId gid, const std::vector< double > &coord)
void buildSubcells()
force the mesh to build subcells: edges and faces
Teuchos::RCP< Teuchos::MpiComm< int > > mpiComm_
Teuchos::RCP< stk::mesh::MetaData > metaData_
Teuchos::RCP< stk::mesh::MetaData > getMetaData() const
Teuchos::RCP< std::vector< stk::mesh::Entity > > orderedElementVector_
bool isInitialized() const
Has initialize been called on this mesh object?
TransListIter iter
void getMyNodes(const std::string &sideName, const std::string &blockName, std::vector< stk::mesh::Entity > &nodes) const
bool nonnull(const boost::shared_ptr< T > &p)
void addMeshCoordFields(const std::string &blockId, const std::vector< std::string > &coordField, const std::string &dispPrefix)
stk::mesh::Entity findConnectivityById(stk::mesh::Entity src, stk::mesh::EntityRank tgt_rank, unsigned rel_id) const
bool isMeshCoordField(const std::string &eBlock, const std::string &fieldName, int &axis) const
void getAllSides(const std::string &sideName, std::vector< stk::mesh::Entity > &sides) const
void addSideset(const std::string &name, const CellTopologyData *ctData)
PHX::MDField< ScalarT, panzer::Cell, panzer::BASIS > field
A field to which we&#39;ll contribute, or in which we&#39;ll store, the result of computing this integral...
void rebalance(const Teuchos::ParameterList &params)
std::map< std::string, std::vector< std::string > > meshCoordFields_
stk::mesh::EntityRank getEdgeRank() const
static const std::string facesString
Teuchos::RCP< const Teuchos::Comm< int > > getComm() const
get the comm associated with this mesh
void initializeFieldsInSTK(const std::map< std::pair< std::string, std::string >, SolutionFieldType * > &nameToField, bool setupIO)
const std::vector< Teuchos::RCP< const PeriodicBC_MatcherBase > > & getPeriodicBCVector() const
void refineMesh(const int numberOfLevels, const bool deleteParentElements)
void setupExodusFile(const std::string &filename)
Set up an output Exodus file for writing results.
std::unordered_map< stk::mesh::EntityId, std::size_t > localIDHash_
ParameterList & sublist(const std::string &name, bool mustAlreadyExist=false, const std::string &docString="")
std::vector< stk::mesh::EntityId > maxEntityId_
stk::mesh::Field< ProcIdData > ProcIdFieldType
void getElementsSharingNode(stk::mesh::EntityId nodeId, std::vector< stk::mesh::Entity > &elements) const
get a set of elements sharing a single node
void getMyElements(std::vector< stk::mesh::Entity > &elements) const
#define TEUCHOS_ASSERT(assertion_test)
stk::mesh::EntityRank getNodeRank() const
void addCellField(const std::string &fieldName, const std::string &blockId)
void getNodeIdsForElement(stk::mesh::Entity element, std::vector< stk::mesh::EntityId > &nodeIds) const
get a list of node ids for nodes connected to an element
stk::mesh::EntityRank getElementRank() const
void addElementBlock(const std::string &name, const CellTopologyData *ctData)
stk::mesh::Field< double > * getCellField(const std::string &fieldName, const std::string &blockId) const
const VectorFieldType & getCoordinatesField() const
void getNeighborElements(std::vector< stk::mesh::Entity > &elements) const
void getNodesetNames(std::vector< std::string > &name) const
Teuchos::RCP< const shards::CellTopology > getCellTopology(const std::string &eBlock) const
std::map< std::pair< std::string, std::string >, SolutionFieldType * > fieldNameToCellField_
std::map< std::string, stk::mesh::Part * > nodesets_
bool is_null() const