Panzer  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Panzer_STK_CubeTetMeshFactory.cpp
Go to the documentation of this file.
1 // @HEADER
2 // *****************************************************************************
3 // Panzer: A partial differential equation assembly
4 // engine for strongly coupled complex multiphysics systems
5 //
6 // Copyright 2011 NTESS and the Panzer contributors.
7 // SPDX-License-Identifier: BSD-3-Clause
8 // *****************************************************************************
9 // @HEADER
10 
12 #include <Teuchos_TimeMonitor.hpp>
13 #include <PanzerAdaptersSTK_config.hpp>
14 
15 using Teuchos::RCP;
16 using Teuchos::rcp;
17 
18 namespace panzer_stk {
19 
21 {
23 }
24 
27 {
28 }
29 
31 Teuchos::RCP<STK_Interface> CubeTetMeshFactory::buildMesh(stk::ParallelMachine parallelMach) const
32 {
33  PANZER_FUNC_TIME_MONITOR("panzer::CubeTetMeshFactory::buildMesh()");
34 
35  // build all meta data
36  RCP<STK_Interface> mesh = buildUncommitedMesh(parallelMach);
37 
38  // commit meta data
39  mesh->initialize(parallelMach);
40 
41  // build bulk data
42  completeMeshConstruction(*mesh,parallelMach);
43 
44  return mesh;
45 }
46 
48 {
49  PANZER_FUNC_TIME_MONITOR("panzer::CubeTetMeshFactory::buildUncomittedMesh()");
50 
51  RCP<STK_Interface> mesh = rcp(new STK_Interface(3));
52 
53  machRank_ = stk::parallel_machine_rank(parallelMach);
54  machSize_ = stk::parallel_machine_size(parallelMach);
55 
56  if (xProcs_ == -1 && yProcs_ == -1 && zProcs_ == -1) {
57  // copied from galeri
58  xProcs_ = yProcs_ = zProcs_ = Teuchos::as<int>(pow(Teuchos::as<double>(machSize_), 0.333334));
59 
60  if (xProcs_ * yProcs_ * zProcs_ != Teuchos::as<int>(machSize_)) {
61  // Simple method to find a set of processor assignments
62  xProcs_ = yProcs_ = zProcs_ = 1;
63 
64  // This means that this works correctly up to about maxFactor^3
65  // processors.
66  const int maxFactor = 50;
67 
68  int ProcTemp = machSize_;
69  int factors[maxFactor];
70  for (int jj = 0; jj < maxFactor; jj++) factors[jj] = 0;
71  for (int jj = 2; jj < maxFactor; jj++) {
72  bool flag = true;
73  while (flag) {
74  int temp = ProcTemp/jj;
75  if (temp*jj == ProcTemp) {
76  factors[jj]++;
77  ProcTemp = temp;
78 
79  } else {
80  flag = false;
81  }
82  }
83  }
84  xProcs_ = ProcTemp;
85  for (int jj = maxFactor-1; jj > 0; jj--) {
86  while (factors[jj] != 0) {
87  if ((xProcs_ <= yProcs_) && (xProcs_ <= zProcs_)) xProcs_ = xProcs_*jj;
88  else if ((yProcs_ <= xProcs_) && (yProcs_ <= zProcs_)) yProcs_ = yProcs_*jj;
89  else zProcs_ = zProcs_*jj;
90  factors[jj]--;
91  }
92  }
93  }
94 
95  } else if(xProcs_==-1) {
96  // default x only decomposition
98  yProcs_ = 1;
99  zProcs_ = 1;
100  }
102  "Cannot build CubeTetMeshFactory, the product of \"X Procs\", \"Y Procs\", and \"Z Procs\""
103  " must equal the number of processors.");
105 
106  // build meta information: blocks and side set setups
107  buildMetaData(parallelMach,*mesh);
108 
109  mesh->addPeriodicBCs(periodicBCVec_);
110  mesh->setBoundingBoxSearchFlag(useBBoxSearch_);
111 
112  return mesh;
113 }
114 
115 void CubeTetMeshFactory::completeMeshConstruction(STK_Interface & mesh,stk::ParallelMachine parallelMach) const
116 {
117  PANZER_FUNC_TIME_MONITOR("panzer::CubeTetMeshFactory::completeMeshConstruction()");
118 
119  if(not mesh.isInitialized())
120  mesh.initialize(parallelMach);
121 
122  // add node and element information
123  buildElements(parallelMach,mesh);
124 
125  // finish up the edges and faces
126  mesh.buildSubcells();
127  mesh.buildLocalElementIDs();
128  if(createEdgeBlocks_) {
129  mesh.buildLocalEdgeIDs();
130  }
131  if(createFaceBlocks_) {
132  mesh.buildLocalFaceIDs();
133  }
134 
135  // now that edges are built, sidets can be added
136  addSideSets(mesh);
137  addNodeSets(mesh);
138 
139  mesh.beginModification();
140  if(createEdgeBlocks_) {
141  addEdgeBlocks(mesh);
142  }
143  if(createFaceBlocks_) {
144  addFaceBlocks(mesh);
145  }
146  mesh.endModification();
147 
148  // calls Stk_MeshFactory::rebalance
149  this->rebalance(mesh);
150 }
151 
154 {
156 
157  setMyParamList(paramList);
158 
159  x0_ = paramList->get<double>("X0");
160  y0_ = paramList->get<double>("Y0");
161  z0_ = paramList->get<double>("Z0");
162 
163  xf_ = paramList->get<double>("Xf");
164  yf_ = paramList->get<double>("Yf");
165  zf_ = paramList->get<double>("Zf");
166 
167  xBlocks_ = paramList->get<int>("X Blocks");
168  yBlocks_ = paramList->get<int>("Y Blocks");
169  zBlocks_ = paramList->get<int>("Z Blocks");
170 
171  xProcs_ = paramList->get<int>("X Procs");
172  yProcs_ = paramList->get<int>("Y Procs");
173  zProcs_ = paramList->get<int>("Z Procs");
174 
175  nXElems_ = paramList->get<int>("X Elements");
176  nYElems_ = paramList->get<int>("Y Elements");
177  nZElems_ = paramList->get<int>("Z Elements");
178 
179  createEdgeBlocks_ = paramList->get<bool>("Create Edge Blocks");
180  createFaceBlocks_ = paramList->get<bool>("Create Face Blocks");
181 
182  // read in periodic boundary conditions
183  parsePeriodicBCList(Teuchos::rcpFromRef(paramList->sublist("Periodic BCs")),periodicBCVec_,useBBoxSearch_);
184 }
185 
188 {
189  static RCP<Teuchos::ParameterList> defaultParams;
190 
191  // fill with default values
192  if(defaultParams == Teuchos::null) {
193  defaultParams = rcp(new Teuchos::ParameterList);
194 
195  defaultParams->set<double>("X0",0.0);
196  defaultParams->set<double>("Y0",0.0);
197  defaultParams->set<double>("Z0",0.0);
198 
199  defaultParams->set<double>("Xf",1.0);
200  defaultParams->set<double>("Yf",1.0);
201  defaultParams->set<double>("Zf",1.0);
202 
203  defaultParams->set<int>("X Blocks",1);
204  defaultParams->set<int>("Y Blocks",1);
205  defaultParams->set<int>("Z Blocks",1);
206 
207  defaultParams->set<int>("X Procs",-1);
208  defaultParams->set<int>("Y Procs",1);
209  defaultParams->set<int>("Z Procs",1);
210 
211  defaultParams->set<int>("X Elements",5);
212  defaultParams->set<int>("Y Elements",5);
213  defaultParams->set<int>("Z Elements",5);
214 
215  // default to false for backward compatibility
216  defaultParams->set<bool>("Create Edge Blocks",false,"Create edge blocks in the mesh");
217  defaultParams->set<bool>("Create Face Blocks",false,"Create face blocks in the mesh");
218 
219  Teuchos::ParameterList & bcs = defaultParams->sublist("Periodic BCs");
220  bcs.set<int>("Count",0); // no default periodic boundary conditions
221  }
222 
223  return defaultParams;
224 }
225 
227 {
228  // get valid parameters
230 
231  // set that parameter list
232  setParameterList(validParams);
233 
234  /* This is a tet mesh factory so all elements in all element blocks
235  * will be tet4. This means that all the edges will be line2 and
236  * all the faces will be tri3. The edge and face block names are
237  * hard coded to reflect this.
238  */
241 }
242 
243 void CubeTetMeshFactory::buildMetaData(stk::ParallelMachine /* parallelMach */, STK_Interface & mesh) const
244 {
245  typedef shards::Tetrahedron<4> TetTopo;
246  const CellTopologyData * ctd = shards::getCellTopologyData<TetTopo>();
247  const CellTopologyData * side_ctd = shards::CellTopology(ctd).getBaseCellTopologyData(2,0);
248  const CellTopologyData * edge_ctd = shards::CellTopology(ctd).getBaseCellTopologyData(1,0);
249  const CellTopologyData * face_ctd = shards::CellTopology(ctd).getBaseCellTopologyData(2,0);
250 
251  // build meta data
252  //mesh.setDimension(2);
253  for(int bx=0;bx<xBlocks_;bx++) {
254  for(int by=0;by<yBlocks_;by++) {
255  for(int bz=0;bz<zBlocks_;bz++) {
256 
257  std::stringstream ebPostfix;
258  ebPostfix << "-" << bx << "_" << by << "_" << bz;
259 
260  // add element blocks
261  mesh.addElementBlock("eblock"+ebPostfix.str(),ctd);
262  if(createEdgeBlocks_) {
263  mesh.addEdgeBlock("eblock"+ebPostfix.str(),
265  edge_ctd);
266  }
267  if(createFaceBlocks_) {
268  mesh.addFaceBlock("eblock"+ebPostfix.str(),
270  face_ctd);
271  }
272  }
273  }
274  }
275 
276  // add sidesets
277  mesh.addSideset("left",side_ctd);
278  mesh.addSideset("right",side_ctd);
279  mesh.addSideset("top",side_ctd);
280  mesh.addSideset("bottom",side_ctd);
281  mesh.addSideset("front",side_ctd);
282  mesh.addSideset("back",side_ctd);
283 
284  mesh.addNodeset("origin");
285 }
286 
287 void CubeTetMeshFactory::buildElements(stk::ParallelMachine parallelMach,STK_Interface & mesh) const
288 {
289  mesh.beginModification();
290  // build each block
291  for(int xBlock=0;xBlock<xBlocks_;xBlock++) {
292  for(int yBlock=0;yBlock<yBlocks_;yBlock++) {
293  for(int zBlock=0;zBlock<zBlocks_;zBlock++) {
294  buildBlock(parallelMach,xBlock,yBlock,zBlock,mesh);
295  }
296  }
297  }
298  mesh.endModification();
299 }
300 
301 void CubeTetMeshFactory::buildBlock(stk::ParallelMachine /* parallelMach */,int xBlock,int yBlock,int zBlock,STK_Interface & mesh) const
302 {
303  // grab this processors rank and machine size
304  std::pair<int,int> sizeAndStartX = determineXElemSizeAndStart(xBlock,xProcs_,machRank_);
305  std::pair<int,int> sizeAndStartY = determineYElemSizeAndStart(yBlock,yProcs_,machRank_);
306  std::pair<int,int> sizeAndStartZ = determineZElemSizeAndStart(zBlock,zProcs_,machRank_);
307 
308  int myXElems_start = sizeAndStartX.first;
309  int myXElems_end = myXElems_start+sizeAndStartX.second;
310  int myYElems_start = sizeAndStartY.first;
311  int myYElems_end = myYElems_start+sizeAndStartY.second;
312  int myZElems_start = sizeAndStartZ.first;
313  int myZElems_end = myZElems_start+sizeAndStartZ.second;
314 
315  int totalXElems = nXElems_*xBlocks_;
316  int totalYElems = nYElems_*yBlocks_;
317  int totalZElems = nZElems_*zBlocks_;
318 
319  double deltaX = (xf_-x0_)/double(totalXElems);
320  double deltaY = (yf_-y0_)/double(totalYElems);
321  double deltaZ = (zf_-z0_)/double(totalZElems);
322 
323  std::vector<double> coord(3,0.0);
324 
325  // build the nodes
326  for(int nx=myXElems_start;nx<myXElems_end+1;++nx) {
327  coord[0] = this->getMeshCoord(nx, deltaX, x0_);
328  for(int ny=myYElems_start;ny<myYElems_end+1;++ny) {
329  coord[1] = this->getMeshCoord(ny, deltaY, y0_);
330  for(int nz=myZElems_start;nz<myZElems_end+1;++nz) {
331  coord[2] = this->getMeshCoord(nz, deltaZ, z0_);
332 
333  mesh.addNode(nz*(totalYElems+1)*(totalXElems+1)+ny*(totalXElems+1)+nx+1,coord);
334  }
335  }
336  }
337 
338  std::stringstream blockName;
339  blockName << "eblock-" << xBlock << "_" << yBlock << "_" << zBlock;
340  stk::mesh::Part * block = mesh.getElementBlockPart(blockName.str());
341 
342  // build the elements
343  for(int nx=myXElems_start;nx<myXElems_end;++nx) {
344  for(int ny=myYElems_start;ny<myYElems_end;++ny) {
345  for(int nz=myZElems_start;nz<myZElems_end;++nz) {
346 
347  std::vector<stk::mesh::EntityId> nodes(8);
348  nodes[0] = nx+1+ny*(totalXElems+1) +nz*(totalYElems+1)*(totalXElems+1);
349  nodes[1] = nodes[0]+1;
350  nodes[2] = nodes[1]+(totalXElems+1);
351  nodes[3] = nodes[2]-1;
352  nodes[4] = nodes[0]+(totalYElems+1)*(totalXElems+1);
353  nodes[5] = nodes[1]+(totalYElems+1)*(totalXElems+1);
354  nodes[6] = nodes[2]+(totalYElems+1)*(totalXElems+1);
355  nodes[7] = nodes[3]+(totalYElems+1)*(totalXElems+1);
356 
357  buildTetsOnHex(Teuchos::tuple(totalXElems,totalYElems,totalZElems),
358  Teuchos::tuple(nx,ny,nz),
359  block,nodes,mesh);
360  }
361  }
362  }
363 }
364 
366  const Teuchos::Tuple<int,3> & element,
367  stk::mesh::Part * block,
368  const std::vector<stk::mesh::EntityId> & h_nodes,
369  STK_Interface & mesh) const
370 {
371  Teuchos::FancyOStream out(Teuchos::rcpFromRef(std::cout));
372  out.setShowProcRank(true);
373  out.setOutputToRootOnly(-1);
374 
375  int totalXElems = meshDesc[0]; int totalYElems = meshDesc[1]; int totalZElems = meshDesc[2];
376  int nx = element[0]; int ny = element[1]; int nz = element[2];
377 
378  stk::mesh::EntityId hex_id = totalXElems*totalYElems*nz+totalXElems*ny+nx+1;
379  stk::mesh::EntityId gid_0 = 12*(hex_id-1)+1;
380  std::vector<stk::mesh::EntityId> nodes(4);
381 
382  // add centroid node
383  stk::mesh::EntityId centroid = 0;
384  {
385  stk::mesh::EntityId largestNode = (totalXElems+1)*(totalYElems+1)*(totalZElems+1);
386  centroid = hex_id+largestNode;
387 
388  // compute average of coordinates
389  std::vector<double> coord(3,0.0);
390  for(std::size_t i=0;i<h_nodes.size();i++) {
391  const double * node_coord = mesh.getNodeCoordinates(h_nodes[i]);
392  coord[0] += node_coord[0];
393  coord[1] += node_coord[1];
394  coord[2] += node_coord[2];
395  }
396  coord[0] /= 8.0;
397  coord[1] /= 8.0;
398  coord[2] /= 8.0;
399 
400  mesh.addNode(centroid,coord);
401  }
402 
403  //
404  int idSet[][3] = { { 0, 1, 2}, // back
405  { 0, 2, 3},
406  { 0, 5, 1}, // bottom
407  { 0, 4, 5},
408  { 0, 7, 4}, // left
409  { 0, 3, 7},
410  { 6, 1, 5}, // right
411  { 6, 2, 1},
412  { 6, 3, 2}, // top
413  { 6, 7, 3},
414  { 6, 4, 7}, // front
415  { 6, 5, 4} };
416 
417  for(int i=0;i<12;i++) {
418  nodes[0] = h_nodes[idSet[i][0]];
419  nodes[1] = h_nodes[idSet[i][1]];
420  nodes[2] = h_nodes[idSet[i][2]];
421  nodes[3] = centroid;
422 
423  // add element to mesh
424  mesh.addElement(rcp(new ElementDescriptor(gid_0+i,nodes)),block);
425  }
426 }
427 
428 std::pair<int,int> CubeTetMeshFactory::determineXElemSizeAndStart(int xBlock,unsigned int size,unsigned int /* rank */) const
429 {
430  std::size_t xProcLoc = procTuple_[0];
431  unsigned int minElements = nXElems_/size;
432  unsigned int extra = nXElems_ - minElements*size;
433 
434  TEUCHOS_ASSERT(minElements>0);
435 
436  // first "extra" elements get an extra column of elements
437  // this determines the starting X index and number of elements
438  int nume=0, start=0;
439  if(xProcLoc<extra) {
440  nume = minElements+1;
441  start = xProcLoc*(minElements+1);
442  }
443  else {
444  nume = minElements;
445  start = extra*(minElements+1)+(xProcLoc-extra)*minElements;
446  }
447 
448  return std::make_pair(start+nXElems_*xBlock,nume);
449 }
450 
451 std::pair<int,int> CubeTetMeshFactory::determineYElemSizeAndStart(int yBlock,unsigned int size,unsigned int /* rank */) const
452 {
453  // int start = yBlock*nYElems_;
454  // return std::make_pair(start,nYElems_);
455 
456  std::size_t yProcLoc = procTuple_[1];
457  unsigned int minElements = nYElems_/size;
458  unsigned int extra = nYElems_ - minElements*size;
459 
460  TEUCHOS_ASSERT(minElements>0);
461 
462  // first "extra" elements get an extra column of elements
463  // this determines the starting X index and number of elements
464  int nume=0, start=0;
465  if(yProcLoc<extra) {
466  nume = minElements+1;
467  start = yProcLoc*(minElements+1);
468  }
469  else {
470  nume = minElements;
471  start = extra*(minElements+1)+(yProcLoc-extra)*minElements;
472  }
473 
474  return std::make_pair(start+nYElems_*yBlock,nume);
475 }
476 
477 std::pair<int,int> CubeTetMeshFactory::determineZElemSizeAndStart(int zBlock,unsigned int size,unsigned int /* rank */) const
478 {
479  // int start = zBlock*nZElems_;
480  // return std::make_pair(start,nZElems_);
481  std::size_t zProcLoc = procTuple_[2];
482  unsigned int minElements = nZElems_/size;
483  unsigned int extra = nZElems_ - minElements*size;
484 
485  TEUCHOS_ASSERT(minElements>0);
486 
487  // first "extra" elements get an extra column of elements
488  // this determines the starting X index and number of elements
489  int nume=0, start=0;
490  if(zProcLoc<extra) {
491  nume = minElements+1;
492  start = zProcLoc*(minElements+1);
493  }
494  else {
495  nume = minElements;
496  start = extra*(minElements+1)+(zProcLoc-extra)*minElements;
497  }
498 
499  return std::make_pair(start+nZElems_*zBlock,nume);
500 }
501 
503 {
504  mesh.beginModification();
505  const stk::mesh::EntityRank side_rank = mesh.getSideRank();
506 
507  std::size_t totalXElems = nXElems_*xBlocks_;
508  std::size_t totalYElems = nYElems_*yBlocks_;
509  std::size_t totalZElems = nZElems_*zBlocks_;
510 
511  // get all part vectors
512  stk::mesh::Part * left = mesh.getSideset("left");
513  stk::mesh::Part * right = mesh.getSideset("right");
514  stk::mesh::Part * top = mesh.getSideset("top");
515  stk::mesh::Part * bottom = mesh.getSideset("bottom");
516  stk::mesh::Part * front = mesh.getSideset("front");
517  stk::mesh::Part * back = mesh.getSideset("back");
518 
519  std::vector<stk::mesh::Entity> localElmts;
520  mesh.getMyElements(localElmts);
521 
522  // gid = totalXElems*totalYElems*nz+totalXElems*ny+nx+1
523 
524  // loop over elements adding sides to sidesets
525  std::vector<stk::mesh::Entity>::const_iterator itr;
526  for(itr=localElmts.begin();itr!=localElmts.end();++itr) {
527  stk::mesh::Entity element = (*itr);
528  stk::mesh::EntityId gid = mesh.elementGlobalId(element);
529 
530  // get hex global id
531  stk::mesh::EntityId h_gid = (gid-1)/12+1;
532  stk::mesh::EntityId t_offset = gid - (12*(h_gid-1)+1);
533 
534  std::size_t nx,ny,nz;
535  nz = (h_gid-1) / (totalXElems*totalYElems);
536  h_gid = (h_gid-1)-nz*(totalXElems*totalYElems);
537  ny = h_gid / totalXElems;
538  nx = h_gid-ny*totalXElems;
539 
540  if(nz==0 && (t_offset==0 || t_offset==1)) {
541  stk::mesh::Entity side = mesh.findConnectivityById(element, side_rank, 3);
542 
543  // on the back
544  if(mesh.entityOwnerRank(side)==machRank_)
545  mesh.addEntityToSideset(side,back);
546  }
547  if(nz+1==totalZElems && (t_offset==10 || t_offset==11)) {
548  stk::mesh::Entity side = mesh.findConnectivityById(element, side_rank, 3);
549 
550  // on the front
551  if(mesh.entityOwnerRank(side)==machRank_)
552  mesh.addEntityToSideset(side,front);
553  }
554 
555  if(ny==0 && (t_offset==2 || t_offset==3)) {
556  stk::mesh::Entity side = mesh.findConnectivityById(element, side_rank, 3);
557 
558  // on the bottom
559  if(mesh.entityOwnerRank(side)==machRank_)
560  mesh.addEntityToSideset(side,bottom);
561  }
562  if(ny+1==totalYElems && (t_offset==8 || t_offset==9)) {
563  stk::mesh::Entity side = mesh.findConnectivityById(element, side_rank, 3);
564 
565  // on the top
566  if(mesh.entityOwnerRank(side)==machRank_)
567  mesh.addEntityToSideset(side,top);
568  }
569 
570  if(nx==0 && (t_offset==4 || t_offset==5)) {
571  stk::mesh::Entity side = mesh.findConnectivityById(element, side_rank, 3);
572 
573  // on the left
574  if(mesh.entityOwnerRank(side)==machRank_)
575  mesh.addEntityToSideset(side,left);
576  }
577  if(nx+1==totalXElems && (t_offset==6 || t_offset==7)) {
578  stk::mesh::Entity side = mesh.findConnectivityById(element, side_rank, 3);
579 
580  // on the right
581  if(mesh.entityOwnerRank(side)==machRank_)
582  mesh.addEntityToSideset(side,right);
583  }
584  }
585 
586  mesh.endModification();
587 }
588 
590 {
591  mesh.beginModification();
592 
593  // get all part vectors
594  stk::mesh::Part * origin = mesh.getNodeset("origin");
595 
597  if(machRank_==0)
598  {
599  // add zero node to origin node set
600  stk::mesh::Entity node = bulkData->get_entity(mesh.getNodeRank(),1);
601  mesh.addEntityToNodeset(node,origin);
602  }
603 
604  mesh.endModification();
605 }
606 
607 // Pre-Condition: call beginModification() before entry
608 // Post-Condition: call endModification() after exit
610 {
613 
614  stk::mesh::Part * edge_block = mesh.getEdgeBlock(edgeBlockName_);
615 
616  stk::mesh::Selector owned_block = metaData->locally_owned_part();
617 
618  std::vector<stk::mesh::Entity> edges;
619  bulkData->get_entities(mesh.getEdgeRank(), owned_block, edges);
620  mesh.addEntitiesToEdgeBlock(edges, edge_block);
621 }
622 
623 // Pre-Condition: call beginModification() before entry
624 // Post-Condition: call endModification() after exit
626 {
629 
630  stk::mesh::Part * face_block = mesh.getFaceBlock(faceBlockName_);
631 
632  stk::mesh::Selector owned_block = metaData->locally_owned_part();
633 
634  std::vector<stk::mesh::Entity> faces;
635  bulkData->get_entities(mesh.getFaceRank(), owned_block, faces);
636  mesh.addEntitiesToFaceBlock(faces, face_block);
637 }
638 
641 {
642  std::size_t i=0,j=0,k=0;
643 
644  k = procRank/(xProcs_*yProcs_); procRank = procRank % (xProcs_*yProcs_);
645  j = procRank/xProcs_; procRank = procRank % xProcs_;
646  i = procRank;
647 
648  return Teuchos::tuple(i,j,k);
649 }
650 
651 } // end panzer_stk
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const
From ParameterListAcceptor.
std::pair< int, int > determineYElemSizeAndStart(int yBlock, unsigned int size, unsigned int rank) const
void addNodeset(const std::string &name)
unsigned entityOwnerRank(stk::mesh::Entity entity) const
void addEntityToNodeset(stk::mesh::Entity entity, stk::mesh::Part *nodeset)
static void parsePeriodicBCList(const Teuchos::RCP< Teuchos::ParameterList > &pl, std::vector< Teuchos::RCP< const PeriodicBC_MatcherBase > > &periodicBC, bool &useBBoxSearch)
basic_FancyOStream & setShowProcRank(const bool showProcRank)
T & get(const std::string &name, T def_value)
Teuchos::RCP< STK_Interface > buildMesh(stk::ParallelMachine parallelMach) const
Build the mesh object.
stk::mesh::Part * getElementBlockPart(const std::string &name) const
get the block part
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
void addEntityToSideset(stk::mesh::Entity entity, stk::mesh::Part *sideset)
void addFaceBlocks(STK_Interface &mesh) const
ParameterList & set(std::string const &name, T &&value, std::string const &docString="", RCP< const ParameterEntryValidator > const &validator=null)
std::pair< int, int > determineZElemSizeAndStart(int zBlock, unsigned int size, unsigned int rank) const
stk::mesh::Part * getNodeset(const std::string &name) const
stk::mesh::Part * getFaceBlock(const std::string &name) const
get the block part
void addEdgeBlock(const std::string &elemBlockName, const std::string &edgeBlockName, const stk::topology &topology)
const double * getNodeCoordinates(stk::mesh::EntityId nodeId) const
stk::mesh::Part * getSideset(const std::string &name) const
virtual Teuchos::RCP< STK_Interface > buildUncommitedMesh(stk::ParallelMachine parallelMach) const
void buildElements(stk::ParallelMachine parallelMach, STK_Interface &mesh) const
void addSideSets(STK_Interface &mesh) const
void addElement(const Teuchos::RCP< ElementDescriptor > &ed, stk::mesh::Part *block)
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)
void rebalance(STK_Interface &mesh) const
void buildBlock(stk::ParallelMachine machRank, int xBlock, int yBlock, int zBlock, STK_Interface &mesh) const
Teuchos::RCP< stk::mesh::BulkData > getBulkData() const
stk::mesh::EntityId elementGlobalId(std::size_t lid) const
void setMyParamList(const RCP< ParameterList > &paramList)
void validateParametersAndSetDefaults(ParameterList const &validParamList, int const depth=1000)
static const std::string edgeBlockString
basic_FancyOStream & setOutputToRootOnly(const int rootRank)
void setParameterList(const Teuchos::RCP< Teuchos::ParameterList > &paramList)
From ParameterListAcceptor.
std::pair< int, int > determineXElemSizeAndStart(int xBlock, unsigned int size, unsigned int rank) const
stk::mesh::EntityRank getFaceRank() const
void addEdgeBlocks(STK_Interface &mesh) const
void addNode(stk::mesh::EntityId gid, const std::vector< double > &coord)
void addEntitiesToEdgeBlock(std::vector< stk::mesh::Entity > entities, stk::mesh::Part *edgeblock)
void buildSubcells()
force the mesh to build subcells: edges and faces
Teuchos::RCP< stk::mesh::MetaData > getMetaData() const
bool isInitialized() const
Has initialize been called on this mesh object?
Teuchos::Tuple< std::size_t, 3 > procRankToProcTuple(std::size_t procRank) const
what is the 3D tuple describe this processor distribution
stk::mesh::Entity findConnectivityById(stk::mesh::Entity src, stk::mesh::EntityRank tgt_rank, unsigned rel_id) const
void addSideset(const std::string &name, const CellTopologyData *ctData)
static const std::string faceBlockString
double getMeshCoord(const int nx, const double deltaX, const double x0) const
void buildMetaData(stk::ParallelMachine parallelMach, STK_Interface &mesh) const
stk::mesh::Part * getEdgeBlock(const std::string &name) const
get the block part
void buildTetsOnHex(const Teuchos::Tuple< int, 3 > &meshDesc, const Teuchos::Tuple< int, 3 > &element, stk::mesh::Part *block, const std::vector< stk::mesh::EntityId > &h_nodes, STK_Interface &mesh) const
std::vector< Teuchos::RCP< const PeriodicBC_MatcherBase > > periodicBCVec_
stk::mesh::EntityRank getEdgeRank() const
void addEntitiesToFaceBlock(std::vector< stk::mesh::Entity > entities, stk::mesh::Part *faceblock)
void addFaceBlock(const std::string &elemBlockName, const std::string &faceBlockName, const stk::topology &topology)
ParameterList & sublist(const std::string &name, bool mustAlreadyExist=false, const std::string &docString="")
void getMyElements(std::vector< stk::mesh::Entity > &elements) const
#define TEUCHOS_ASSERT(assertion_test)
stk::mesh::EntityRank getNodeRank() const
virtual void completeMeshConstruction(STK_Interface &mesh, stk::ParallelMachine parallelMach) const
Teuchos::Tuple< std::size_t, 3 > procTuple_
void addElementBlock(const std::string &name, const CellTopologyData *ctData)
void addNodeSets(STK_Interface &mesh) const