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