Panzer  Version of the Day
Panzer_STK_SquareQuadMeshFactory.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 #include "Teuchos_StandardParameterEntryValidators.hpp" // for plist validation
47 
48 // #define ENABLE_UNIFORM
49 
50 using Teuchos::RCP;
51 using Teuchos::rcp;
52 
53 namespace panzer_stk {
54 
56 {
58 }
59 
62 {
63 }
64 
66 Teuchos::RCP<STK_Interface> SquareQuadMeshFactory::buildMesh(stk::ParallelMachine parallelMach) const
67 {
68  PANZER_FUNC_TIME_MONITOR("panzer::SquareQuadMeshFactory::buildMesh()");
69 
70  // build all meta data
71  RCP<STK_Interface> mesh = buildUncommitedMesh(parallelMach);
72 
73  // commit meta data
74  mesh->initialize(parallelMach);
75 
76  // build bulk data
77  completeMeshConstruction(*mesh,parallelMach);
78 
79  return mesh;
80 }
81 
82 Teuchos::RCP<STK_Interface> SquareQuadMeshFactory::buildUncommitedMesh(stk::ParallelMachine parallelMach) const
83 {
84  PANZER_FUNC_TIME_MONITOR("panzer::SquareQuadMeshFactory::buildUncomittedMesh()");
85 
86  RCP<STK_Interface> mesh = rcp(new STK_Interface(2));
87 
88  machRank_ = stk::parallel_machine_rank(parallelMach);
89  machSize_ = stk::parallel_machine_size(parallelMach);
90 
91  if (xProcs_ == -1 && yProcs_ == -1) {
92  // copied from galeri
93  xProcs_ = yProcs_ = Teuchos::as<int>(pow(Teuchos::as<double>(machSize_), 0.5));
94 
95  if (xProcs_ * yProcs_ != Teuchos::as<int>(machSize_)) {
96  // Simple method to find a set of processor assignments
97  xProcs_ = yProcs_ = 1;
98 
99  // This means that this works correctly up to about maxFactor^2
100  // processors.
101  const int maxFactor = 100;
102 
103  int ProcTemp = machSize_;
104  int factors[maxFactor];
105  for (int jj = 0; jj < maxFactor; jj++) factors[jj] = 0;
106  for (int jj = 2; jj < maxFactor; jj++) {
107  bool flag = true;
108  while (flag) {
109  int temp = ProcTemp/jj;
110  if (temp*jj == ProcTemp) {
111  factors[jj]++;
112  ProcTemp = temp;
113 
114  } else {
115  flag = false;
116  }
117  }
118  }
119  xProcs_ = ProcTemp;
120  for (int jj = maxFactor-1; jj > 0; jj--) {
121  while (factors[jj] != 0) {
122  if (xProcs_ <= yProcs_) xProcs_ = xProcs_*jj;
123  else yProcs_ = yProcs_*jj;
124  factors[jj]--;
125  }
126  }
127  }
128 
129  } else if(xProcs_==-1) {
130  // default x only decomposition
131  xProcs_ = machSize_;
132  yProcs_ = 1;
133  }
134  TEUCHOS_TEST_FOR_EXCEPTION(int(machSize_) != xProcs_ * yProcs_, std::logic_error,
135  "Cannot build SquareQuadMeshFactory. The product of 'X Procs * Y Procs = " << xProcs_ << "*" << yProcs_ << " = " << xProcs_*yProcs_
136  << "' must equal the number of processors = " << machSize_
137  << "\n\n\t==> Run the simulation with an appropriate number of processors, i.e. #procs = " << xProcs_*yProcs_ << ".\n");
139 
140  // build meta information: blocks and side set setups
141  buildMetaData(parallelMach,*mesh);
142 
143  mesh->addPeriodicBCs(periodicBCVec_);
144  mesh->setBoundingBoxSearchFlag(useBBoxSearch_);
145 
146  return mesh;
147 }
148 
149 void SquareQuadMeshFactory::completeMeshConstruction(STK_Interface & mesh,stk::ParallelMachine parallelMach) const
150 {
151  PANZER_FUNC_TIME_MONITOR("panzer::SquareQuadMeshFactory::completeMeshConstruction()");
152 
153  if(not mesh.isInitialized())
154  mesh.initialize(parallelMach);
155 
156  // add node and element information
157  buildElements(parallelMach,mesh);
158 
159  // finish up the edges
160 #ifndef ENABLE_UNIFORM
161  mesh.buildSubcells();
162 #endif
163  mesh.buildLocalElementIDs();
164  if(createEdgeBlocks_) {
165  mesh.buildLocalEdgeIDs();
166  }
167 
168  // now that edges are built, sidsets can be added
169 #ifndef ENABLE_UNIFORM
170  addSideSets(mesh);
171 #endif
172 
173  // add nodesets
174  addNodeSets(mesh);
175 
176  if(createEdgeBlocks_) {
177  addEdgeBlocks(mesh);
178  }
179 
180  // calls Stk_MeshFactory::rebalance
181  this->rebalance(mesh);
182 }
183 
185 void SquareQuadMeshFactory::setParameterList(const Teuchos::RCP<Teuchos::ParameterList> & paramList)
186 {
187  paramList->validateParametersAndSetDefaults(*getValidParameters(),0);
188 
189  setMyParamList(paramList);
190 
191  x0_ = paramList->get<double>("X0");
192  y0_ = paramList->get<double>("Y0");
193 
194  xf_ = paramList->get<double>("Xf");
195  yf_ = paramList->get<double>("Yf");
196 
197  xBlocks_ = paramList->get<int>("X Blocks");
198  yBlocks_ = paramList->get<int>("Y Blocks");
199 
200  nXElems_ = paramList->get<int>("X Elements");
201  nYElems_ = paramList->get<int>("Y Elements");
202 
203  xProcs_ = paramList->get<int>("X Procs");
204  yProcs_ = paramList->get<int>("Y Procs");
205 
206  offsetGIDs_ = (paramList->get<std::string>("Offset mesh GIDs above 32-bit int limit") == "ON") ? true : false;
207 
208  createEdgeBlocks_ = paramList->get<bool>("Create Edge Blocks");
209 
210  // read in periodic boundary conditions
211  parsePeriodicBCList(Teuchos::rcpFromRef(paramList->sublist("Periodic BCs")),periodicBCVec_,useBBoxSearch_);
212 }
213 
215 Teuchos::RCP<const Teuchos::ParameterList> SquareQuadMeshFactory::getValidParameters() const
216 {
217  static RCP<Teuchos::ParameterList> defaultParams;
218 
219  // fill with default values
220  if(defaultParams == Teuchos::null) {
221  defaultParams = rcp(new Teuchos::ParameterList);
222 
223  defaultParams->set<double>("X0",0.0);
224  defaultParams->set<double>("Y0",0.0);
225 
226  defaultParams->set<double>("Xf",1.0);
227  defaultParams->set<double>("Yf",1.0);
228 
229  defaultParams->set<int>("X Blocks",1);
230  defaultParams->set<int>("Y Blocks",1);
231 
232  defaultParams->set<int>("X Procs",-1);
233  defaultParams->set<int>("Y Procs",1);
234 
235  defaultParams->set<int>("X Elements",5);
236  defaultParams->set<int>("Y Elements",5);
237 
238  // default to false for backward compatibility
239  defaultParams->set<bool>("Create Edge Blocks",false,"Create edge blocks in the mesh");
240 
241  Teuchos::setStringToIntegralParameter<int>(
242  "Offset mesh GIDs above 32-bit int limit",
243  "OFF",
244  "If 64-bit GIDs are supported, the mesh element and node global indices will start at a value greater than 32-bit limit.",
245  Teuchos::tuple<std::string>("OFF", "ON"),
246  defaultParams.get());
247 
248  Teuchos::ParameterList & bcs = defaultParams->sublist("Periodic BCs");
249  bcs.set<int>("Count",0); // no default periodic boundary conditions
250  }
251 
252  return defaultParams;
253 }
254 
256 {
257  // get valid parameters
258  RCP<Teuchos::ParameterList> validParams = rcp(new Teuchos::ParameterList(*getValidParameters()));
259 
260  // set that parameter list
261  setParameterList(validParams);
262 
263  /* This is a quad mesh factory so all elements in all element blocks
264  * will be quad4. This means that all the edges will be line2.
265  * The edge block name is hard coded to reflect this.
266  */
268 }
269 
270 void SquareQuadMeshFactory::buildMetaData(stk::ParallelMachine /* parallelMach */, STK_Interface & mesh) const
271 {
272  typedef shards::Quadrilateral<4> QuadTopo;
273  const CellTopologyData * ctd = shards::getCellTopologyData<QuadTopo>();
274  const CellTopologyData * side_ctd = shards::CellTopology(ctd).getBaseCellTopologyData(1,0);
275  const CellTopologyData * edge_ctd = shards::CellTopology(ctd).getBaseCellTopologyData(1,0);
276 
277  // build meta data
278  //mesh.setDimension(2);
279  for(int bx=0;bx<xBlocks_;bx++) {
280  for(int by=0;by<yBlocks_;by++) {
281 
282  // add this element block
283  {
284  std::stringstream ebPostfix;
285  ebPostfix << "-" << bx << "_" << by;
286 
287  // add element blocks
288  mesh.addElementBlock("eblock"+ebPostfix.str(),ctd);
289  if(createEdgeBlocks_) {
290  mesh.addEdgeBlock("eblock"+ebPostfix.str(),
292  edge_ctd);
293  }
294  }
295 
296  }
297  }
298 
299  // add sidesets
300 #ifndef ENABLE_UNIFORM
301  mesh.addSideset("left",side_ctd);
302  mesh.addSideset("right",side_ctd);
303  mesh.addSideset("top",side_ctd);
304  mesh.addSideset("bottom",side_ctd);
305 
306  for(int bx=1;bx<xBlocks_;bx++) {
307  std::stringstream ss;
308  ss << "vertical_" << bx-1;
309  mesh.addSideset(ss.str(),side_ctd);
310  }
311  for(int by=1;by<yBlocks_;by++) {
312  std::stringstream ss;
313  ss << "horizontal_" << by-1;
314  mesh.addSideset(ss.str(),side_ctd);
315  }
316 #endif
317 
318  // add nodesets
319  mesh.addNodeset("lower_left");
320  mesh.addNodeset("origin");
321 }
322 
323 void SquareQuadMeshFactory::buildElements(stk::ParallelMachine parallelMach,STK_Interface & mesh) const
324 {
325  mesh.beginModification();
326  // build each block
327  for(int xBlock=0;xBlock<xBlocks_;xBlock++) {
328  for(int yBlock=0;yBlock<yBlocks_;yBlock++) {
329  buildBlock(parallelMach,xBlock,yBlock,mesh);
330  }
331  }
332  mesh.endModification();
333 }
334 
335 void SquareQuadMeshFactory::buildBlock(stk::ParallelMachine /* parallelMach */,int xBlock,int yBlock,STK_Interface & mesh) const
336 {
337  // grab this processors rank and machine size
338  std::pair<int,int> sizeAndStartX = determineXElemSizeAndStart(xBlock,xProcs_,machRank_);
339  std::pair<int,int> sizeAndStartY = determineYElemSizeAndStart(yBlock,yProcs_,machRank_);
340 
341  int myXElems_start = sizeAndStartX.first;
342  int myXElems_end = myXElems_start+sizeAndStartX.second;
343  int myYElems_start = sizeAndStartY.first;
344  int myYElems_end = myYElems_start+sizeAndStartY.second;
345  int totalXElems = nXElems_*xBlocks_;
346  int totalYElems = nYElems_*yBlocks_;
347 
348  double deltaX = (xf_-x0_)/double(totalXElems);
349  double deltaY = (yf_-y0_)/double(totalYElems);
350 
351  std::vector<double> coord(2,0.0);
352 
353  offset_ = 0;
354  if (offsetGIDs_) {
355  if (std::numeric_limits<panzer::GlobalOrdinal>::max() > std::numeric_limits<unsigned int>::max())
356  offset_ = panzer::GlobalOrdinal(std::numeric_limits<unsigned int>::max()) + 1;
357  }
358 
359  // build the nodes
360  for(int nx=myXElems_start;nx<myXElems_end+1;++nx) {
361  coord[0] = this->getMeshCoord(nx, deltaX, x0_);
362  for(int ny=myYElems_start;ny<myYElems_end+1;++ny) {
363  coord[1] = this->getMeshCoord(ny, deltaY, y0_);
364 
365  mesh.addNode(ny*(totalXElems+1)+nx+1+offset_,coord);
366  }
367  }
368 
369  std::stringstream blockName;
370  blockName << "eblock-" << xBlock << "_" << yBlock;
371  stk::mesh::Part * block = mesh.getElementBlockPart(blockName.str());
372 
373  // build the elements
374  for(int nx=myXElems_start;nx<myXElems_end;++nx) {
375  for(int ny=myYElems_start;ny<myYElems_end;++ny) {
376  stk::mesh::EntityId gid = totalXElems*ny+nx+1 + offset_;
377  std::vector<stk::mesh::EntityId> nodes(4);
378  nodes[0] = nx+1+ny*(totalXElems+1);
379  nodes[1] = nodes[0]+1;
380  nodes[2] = nodes[1]+(totalXElems+1);
381  nodes[3] = nodes[2]-1;
382 
383  for (int i=0; i < 4; ++i)
384  nodes[i] += offset_;
385 
386  RCP<ElementDescriptor> ed = rcp(new ElementDescriptor(gid,nodes));
387  mesh.addElement(ed,block);
388  }
389  }
390 }
391 
392 std::pair<int,int> SquareQuadMeshFactory::determineXElemSizeAndStart(int xBlock,unsigned int size,unsigned int /* rank */) const
393 {
394  std::size_t xProcLoc = procTuple_[0];
395  unsigned int minElements = nXElems_/size;
396  unsigned int extra = nXElems_ - minElements*size;
397 
398  TEUCHOS_ASSERT(minElements>0);
399 
400  // first "extra" elements get an extra column of elements
401  // this determines the starting X index and number of elements
402  int nume=0, start=0;
403  if(xProcLoc<extra) {
404  nume = minElements+1;
405  start = xProcLoc*(minElements+1);
406  }
407  else {
408  nume = minElements;
409  start = extra*(minElements+1)+(xProcLoc-extra)*minElements;
410  }
411 
412  return std::make_pair(start+nXElems_*xBlock,nume);
413 }
414 
415 std::pair<int,int> SquareQuadMeshFactory::determineYElemSizeAndStart(int yBlock,unsigned int size,unsigned int /* rank */) const
416 {
417  std::size_t yProcLoc = procTuple_[1];
418  unsigned int minElements = nYElems_/size;
419  unsigned int extra = nYElems_ - minElements*size;
420 
421  TEUCHOS_ASSERT(minElements>0);
422 
423  // first "extra" elements get an extra column of elements
424  // this determines the starting X index and number of elements
425  int nume=0, start=0;
426  if(yProcLoc<extra) {
427  nume = minElements+1;
428  start = yProcLoc*(minElements+1);
429  }
430  else {
431  nume = minElements;
432  start = extra*(minElements+1)+(yProcLoc-extra)*minElements;
433  }
434 
435  return std::make_pair(start+nYElems_*yBlock,nume);
436 }
437 
439 {
440  mesh.beginModification();
441 
442  std::size_t totalXElems = nXElems_*xBlocks_;
443  std::size_t totalYElems = nYElems_*yBlocks_;
444 
445  // get all part vectors
446  stk::mesh::Part * left = mesh.getSideset("left");
447  stk::mesh::Part * right = mesh.getSideset("right");
448  stk::mesh::Part * top = mesh.getSideset("top");
449  stk::mesh::Part * bottom = mesh.getSideset("bottom");
450 
451  std::vector<stk::mesh::Part*> vertical;
452  std::vector<stk::mesh::Part*> horizontal;
453  for(int bx=1;bx<xBlocks_;bx++) {
454  std::stringstream ss;
455  ss << "vertical_" << bx-1;
456  vertical.push_back(mesh.getSideset(ss.str()));
457  }
458  for(int by=1;by<yBlocks_;by++) {
459  std::stringstream ss;
460  ss << "horizontal_" << by-1;
461  horizontal.push_back(mesh.getSideset(ss.str()));
462  }
463 
464  std::vector<stk::mesh::Entity> localElmts;
465  mesh.getMyElements(localElmts);
466 
467  // loop over elements adding edges to sidesets
468  std::vector<stk::mesh::Entity>::const_iterator itr;
469  for(itr=localElmts.begin();itr!=localElmts.end();++itr) {
470  stk::mesh::Entity element = (*itr);
471  stk::mesh::EntityId gid = mesh.elementGlobalId(element);
472 
473  // reverse the offset for local gid numbering scheme
474  gid -= offset_;
475 
476  std::size_t nx,ny;
477  ny = (gid-1) / totalXElems;
478  nx = gid-ny*totalXElems-1;
479 
480  // vertical boundaries
482 
483  if(nx+1==totalXElems) {
484  stk::mesh::Entity edge = mesh.findConnectivityById(element, stk::topology::EDGE_RANK, 1);
485 
486  // on the right
487  if(mesh.entityOwnerRank(edge)==machRank_)
488  mesh.addEntityToSideset(edge,right);
489  }
490 
491  if(nx==0) {
492  stk::mesh::Entity edge = mesh.findConnectivityById(element, stk::topology::EDGE_RANK, 3);
493 
494  // on the left
495  if(mesh.entityOwnerRank(edge)==machRank_)
496  mesh.addEntityToSideset(edge,left);
497  }
498 
499  if(nx+1!=totalXElems && ((nx+1) % nXElems_==0)) {
500  stk::mesh::Entity edge = mesh.findConnectivityById(element, stk::topology::EDGE_RANK, 1);
501 
502  // on the right
503  if(mesh.entityOwnerRank(edge)==machRank_) {
504  int index = (nx+1)/nXElems_-1;
505  mesh.addEntityToSideset(edge,vertical[index]);
506  }
507  }
508 
509  if(nx!=0 && (nx % nXElems_==0)) {
510  stk::mesh::Entity edge = mesh.findConnectivityById(element, stk::topology::EDGE_RANK, 3);
511 
512  // on the left
513  if(mesh.entityOwnerRank(edge)==machRank_) {
514  int index = nx/nXElems_-1;
515  mesh.addEntityToSideset(edge,vertical[index]);
516  }
517  }
518 
519  // horizontal boundaries
521 
522  if(ny==0) {
523  stk::mesh::Entity edge = mesh.findConnectivityById(element, stk::topology::EDGE_RANK, 0);
524 
525  // on the bottom
526  if(mesh.entityOwnerRank(edge)==machRank_)
527  mesh.addEntityToSideset(edge,bottom);
528  }
529 
530  if(ny+1==totalYElems) {
531  stk::mesh::Entity edge = mesh.findConnectivityById(element, stk::topology::EDGE_RANK, 2);
532 
533  // on the top
534  if(mesh.entityOwnerRank(edge)==machRank_)
535  mesh.addEntityToSideset(edge,top);
536  }
537 
538  if(ny!=0 && (ny % nYElems_==0)) {
539  stk::mesh::Entity edge = mesh.findConnectivityById(element, stk::topology::EDGE_RANK, 0);
540 
541  // on the bottom
542  if(mesh.entityOwnerRank(edge)==machRank_) {
543  int index = ny/nYElems_-1;
544  mesh.addEntityToSideset(edge,horizontal[index]);
545  }
546  }
547 
548  if(ny+1!=totalYElems && ((ny+1) % nYElems_==0)) {
549  stk::mesh::Entity edge = mesh.findConnectivityById(element, stk::topology::EDGE_RANK, 2);
550 
551  // on the top
552  if(mesh.entityOwnerRank(edge)==machRank_) {
553  int index = (ny+1)/nYElems_-1;
554  mesh.addEntityToSideset(edge,horizontal[index]);
555  }
556  }
557  }
558 
559  mesh.endModification();
560 }
561 
563 {
564  mesh.beginModification();
565 
566  // get all part vectors
567  stk::mesh::Part * lower_left = mesh.getNodeset("lower_left");
568  stk::mesh::Part * origin = mesh.getNodeset("origin");
569 
570  // std::vector<stk::mesh::Entity> localElmts;
571  // mesh.getMyElements(localElmts);
572 
573  Teuchos::RCP<stk::mesh::BulkData> bulkData = mesh.getBulkData();
574  if(machRank_==0)
575  {
576  // add zero node to lower_left node set
577  stk::mesh::Entity node = bulkData->get_entity(mesh.getNodeRank(),1 + offset_);
578  mesh.addEntityToNodeset(node,lower_left);
579 
580  // add zero node to origin node set
581  mesh.addEntityToNodeset(node,origin);
582  }
583 
584  mesh.endModification();
585 }
586 
588 {
589  mesh.beginModification();
590 
591  Teuchos::RCP<stk::mesh::BulkData> bulkData = mesh.getBulkData();
592  Teuchos::RCP<stk::mesh::MetaData> metaData = mesh.getMetaData();
593 
594  stk::mesh::Part * edge_block = mesh.getEdgeBlock(edgeBlockName_);
595 
596  stk::mesh::Selector owned_block = metaData->locally_owned_part();
597 
598  std::vector<stk::mesh::Entity> edges;
599  bulkData->get_entities(mesh.getEdgeRank(), owned_block, edges);
600  mesh.addEntitiesToEdgeBlock(edges, edge_block);
601 
602  mesh.endModification();
603 }
604 
606 Teuchos::Tuple<std::size_t,2> SquareQuadMeshFactory::procRankToProcTuple(std::size_t procRank) const
607 {
608  std::size_t i=0,j=0;
609 
610  j = procRank/xProcs_;
611  procRank = procRank % xProcs_;
612  i = procRank;
613 
614  return Teuchos::tuple(i,j);
615 }
616 
617 } // end panzer_stk
Teuchos::RCP< stk::mesh::MetaData > getMetaData() const
stk::mesh::Part * getElementBlockPart(const std::string &name) const
get the block part
void buildMetaData(stk::ParallelMachine parallelMach, STK_Interface &mesh) const
void addNodeset(const std::string &name)
void buildBlock(stk::ParallelMachine machRank, int xBlock, int yBlock, STK_Interface &mesh) 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)
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::mesh::BulkData > getBulkData() const
virtual Teuchos::RCP< STK_Interface > buildUncommitedMesh(stk::ParallelMachine parallelMach) const
stk::mesh::Entity findConnectivityById(stk::mesh::Entity src, stk::mesh::EntityRank tgt_rank, unsigned rel_id) const
void addEdgeBlock(const std::string &elemBlockName, const std::string &edgeBlockName, const stk::topology &topology)
bool offsetGIDs_
If true, offset mesh GIDs to exercise 32-bit limits.
stk::mesh::Part * getEdgeBlock(const std::string &name) const
get the block part
stk::mesh::EntityRank getEdgeRank() 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 initialize(stk::ParallelMachine parallelMach, bool setupIO=true, const bool buildRefinementSupport=false)
SquareQuadMeshFactory(bool enableRebalance=false)
Constructor.
virtual void completeMeshConstruction(STK_Interface &mesh, stk::ParallelMachine parallelMach) const
double getMeshCoord(const int nx, const double deltaX, const double x0) const
bool isInitialized() const
Has initialize been called on this mesh object?
static const std::string edgeBlockString
std::pair< int, int > determineYElemSizeAndStart(int yBlock, unsigned int size, unsigned int rank) 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
void buildElements(stk::ParallelMachine parallelMach, STK_Interface &mesh) const
void addSideset(const std::string &name, const CellTopologyData *ctData)
std::pair< int, int > determineXElemSizeAndStart(int xBlock, unsigned int size, unsigned int rank) const
std::vector< Teuchos::RCP< const PeriodicBC_MatcherBase > > periodicBCVec_
Teuchos::Tuple< std::size_t, 2 > procRankToProcTuple(std::size_t procRank) const
what is the 2D tuple describe this processor distribution
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const
From ParameterListAcceptor.
void setParameterList(const Teuchos::RCP< Teuchos::ParameterList > &paramList)
From ParameterListAcceptor.
Teuchos::RCP< STK_Interface > buildMesh(stk::ParallelMachine parallelMach) const
Build the mesh object.
unsigned entityOwnerRank(stk::mesh::Entity entity) const
void rebalance(STK_Interface &mesh) const
void addElementBlock(const std::string &name, const CellTopologyData *ctData)
stk::mesh::Part * getNodeset(const std::string &name) const