Panzer  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Panzer_STK_SquareTriMeshFactory.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> SquareTriMeshFactory::buildMesh(stk::ParallelMachine parallelMach) const
32 {
33  PANZER_FUNC_TIME_MONITOR("panzer::SquareTriMeshFactory::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::SquareTriMeshFactory::buildUncomittedMesh()");
50 
51  RCP<STK_Interface> mesh = rcp(new STK_Interface(2));
52 
53  machRank_ = stk::parallel_machine_rank(parallelMach);
54  machSize_ = stk::parallel_machine_size(parallelMach);
55 
56  if (xProcs_ == -1 && yProcs_ == -1) {
57  // copied from galeri
58  xProcs_ = yProcs_ = Teuchos::as<int>(pow(Teuchos::as<double>(machSize_), 0.5));
59 
60  if (xProcs_ * yProcs_ != Teuchos::as<int>(machSize_)) {
61  // Simple method to find a set of processor assignments
62  xProcs_ = yProcs_ = 1;
63 
64  // This means that this works correctly up to about maxFactor^2
65  // processors.
66  const int maxFactor = 100;
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_ = xProcs_*jj;
88  else yProcs_ = yProcs_*jj;
89  factors[jj]--;
90  }
91  }
92  }
93 
94  } else if(xProcs_==-1) {
95  // default x only decomposition
97  yProcs_ = 1;
98  }
99  TEUCHOS_TEST_FOR_EXCEPTION(int(machSize_)!=xProcs_*yProcs_,std::logic_error,
100  "Cannot build SquareTriMeshFactory, the product of \"X Procs\" and \"Y Procs\""
101  " must equal the number of processors.");
103 
104  // build meta information: blocks and side set setups
105  buildMetaData(parallelMach,*mesh);
106 
107  mesh->addPeriodicBCs(periodicBCVec_);
108  mesh->setBoundingBoxSearchFlag(useBBoxSearch_);
109 
110  return mesh;
111 }
112 
113 void SquareTriMeshFactory::completeMeshConstruction(STK_Interface & mesh,stk::ParallelMachine parallelMach) const
114 {
115  PANZER_FUNC_TIME_MONITOR("panzer::SquareTriMeshFactory::completeMeshConstruction()");
116 
117  if(not mesh.isInitialized())
118  mesh.initialize(parallelMach);
119 
120  // add node and element information
121  buildElements(parallelMach,mesh);
122 
123  // finish up the edges
124  mesh.buildSubcells();
125  mesh.buildLocalElementIDs();
126  if(createEdgeBlocks_) {
127  mesh.buildLocalEdgeIDs();
128  }
129 
130  // now that edges are built, sidets can be added
131  addSideSets(mesh);
132 
133  // add nodesets
134  addNodeSets(mesh);
135 
136  if(createEdgeBlocks_) {
137  addEdgeBlocks(mesh);
138  }
139 
140  // calls Stk_MeshFactory::rebalance
141  this->rebalance(mesh);
142 }
143 
146 {
148 
149  setMyParamList(paramList);
150 
151  x0_ = paramList->get<double>("X0");
152  y0_ = paramList->get<double>("Y0");
153 
154  xf_ = paramList->get<double>("Xf");
155  yf_ = paramList->get<double>("Yf");
156 
157  xBlocks_ = paramList->get<int>("X Blocks");
158  yBlocks_ = paramList->get<int>("Y Blocks");
159 
160  nXElems_ = paramList->get<int>("X Elements");
161  nYElems_ = paramList->get<int>("Y Elements");
162 
163  xProcs_ = paramList->get<int>("X Procs");
164  yProcs_ = paramList->get<int>("Y Procs");
165 
166  createEdgeBlocks_ = paramList->get<bool>("Create Edge Blocks");
167 
168  // read in periodic boundary conditions
169  parsePeriodicBCList(Teuchos::rcpFromRef(paramList->sublist("Periodic BCs")),periodicBCVec_,useBBoxSearch_);
170 }
171 
174 {
175  static RCP<Teuchos::ParameterList> defaultParams;
176 
177  // fill with default values
178  if(defaultParams == Teuchos::null) {
179  defaultParams = rcp(new Teuchos::ParameterList);
180 
181  defaultParams->set<double>("X0",0.0);
182  defaultParams->set<double>("Y0",0.0);
183 
184  defaultParams->set<double>("Xf",1.0);
185  defaultParams->set<double>("Yf",1.0);
186 
187  defaultParams->set<int>("X Blocks",1);
188  defaultParams->set<int>("Y Blocks",1);
189 
190  defaultParams->set<int>("X Procs",-1);
191  defaultParams->set<int>("Y Procs",1);
192 
193  defaultParams->set<int>("X Elements",5);
194  defaultParams->set<int>("Y Elements",5);
195 
196  // default to false for backward compatibility
197  defaultParams->set<bool>("Create Edge Blocks",false,"Create edge blocks in the mesh");
198 
199  Teuchos::ParameterList & bcs = defaultParams->sublist("Periodic BCs");
200  bcs.set<int>("Count",0); // no default periodic boundary conditions
201  }
202 
203  return defaultParams;
204 }
205 
207 {
208  // get valid parameters
210 
211  // set that parameter list
212  setParameterList(validParams);
213 
214  /* This is a tri mesh factory so all elements in all element blocks
215  * will be tri3. This means that all the edges will be line2.
216  * The edge block name is hard coded to reflect this.
217  */
219 }
220 
221 void SquareTriMeshFactory::buildMetaData(stk::ParallelMachine /* parallelMach */, STK_Interface & mesh) const
222 {
223  typedef shards::Triangle<> TriTopo;
224  const CellTopologyData * ctd = shards::getCellTopologyData<TriTopo>();
225  const CellTopologyData * side_ctd = shards::CellTopology(ctd).getBaseCellTopologyData(1,0);
226  const CellTopologyData * edge_ctd = shards::CellTopology(ctd).getBaseCellTopologyData(1,0);
227 
228  // build meta data
229  //mesh.setDimension(2);
230  for(int bx=0;bx<xBlocks_;bx++) {
231  for(int by=0;by<yBlocks_;by++) {
232 
233  // add this element block
234  {
235  std::stringstream ebPostfix;
236  ebPostfix << "-" << bx << "_" << by;
237 
238  // add element blocks
239  mesh.addElementBlock("eblock"+ebPostfix.str(),ctd);
240  if(createEdgeBlocks_) {
241  mesh.addEdgeBlock("eblock"+ebPostfix.str(),
243  edge_ctd);
244  }
245  }
246 
247  }
248  }
249 
250  // add sidesets
251  mesh.addSideset("left",side_ctd);
252  mesh.addSideset("right",side_ctd);
253  mesh.addSideset("top",side_ctd);
254  mesh.addSideset("bottom",side_ctd);
255 
256  // add nodesets
257  mesh.addNodeset("origin");
258 }
259 
260 void SquareTriMeshFactory::buildElements(stk::ParallelMachine parallelMach,STK_Interface & mesh) const
261 {
262  mesh.beginModification();
263  // build each block
264  for(int xBlock=0;xBlock<xBlocks_;xBlock++) {
265  for(int yBlock=0;yBlock<yBlocks_;yBlock++) {
266  buildBlock(parallelMach,xBlock,yBlock,mesh);
267  }
268  }
269  mesh.endModification();
270 }
271 
272 void SquareTriMeshFactory::buildBlock(stk::ParallelMachine /* parallelMach */, int xBlock, int yBlock, STK_Interface& mesh) const
273 {
274  // grab this processors rank and machine size
275  std::pair<int,int> sizeAndStartX = determineXElemSizeAndStart(xBlock,xProcs_,machRank_);
276  std::pair<int,int> sizeAndStartY = determineYElemSizeAndStart(yBlock,yProcs_,machRank_);
277 
278  int myXElems_start = sizeAndStartX.first;
279  int myXElems_end = myXElems_start+sizeAndStartX.second;
280  int myYElems_start = sizeAndStartY.first;
281  int myYElems_end = myYElems_start+sizeAndStartY.second;
282  int totalXElems = nXElems_*xBlocks_;
283  int totalYElems = nYElems_*yBlocks_;
284 
285  double deltaX = (xf_-x0_)/double(totalXElems);
286  double deltaY = (yf_-y0_)/double(totalYElems);
287 
288  std::vector<double> coord(2,0.0);
289 
290  // build the nodes
291  for(int nx=myXElems_start;nx<myXElems_end+1;++nx) {
292  coord[0] = this->getMeshCoord(nx, deltaX, x0_);
293  for(int ny=myYElems_start;ny<myYElems_end+1;++ny) {
294  coord[1] = this->getMeshCoord(ny, deltaY, y0_);
295 
296  mesh.addNode(ny*(totalXElems+1)+nx+1,coord);
297  }
298  }
299 
300  std::stringstream blockName;
301  blockName << "eblock-" << xBlock << "_" << yBlock;
302  stk::mesh::Part * block = mesh.getElementBlockPart(blockName.str());
303 
304  // build the elements
305  for(int nx=myXElems_start;nx<myXElems_end;++nx) {
306  for(int ny=myYElems_start;ny<myYElems_end;++ny) {
307  stk::mesh::EntityId gid_a = 2*(totalXElems*ny+nx+1)-1;
308  stk::mesh::EntityId gid_b = gid_a+1;
309  std::vector<stk::mesh::EntityId> nodes(3);
310  stk::mesh::EntityId sw,se,ne,nw;
311  sw = nx+1+ny*(totalXElems+1);
312  se = sw+1;
313  ne = se+(totalXElems+1);
314  nw = ne-1;
315 
316  nodes[0] = sw;
317  nodes[1] = se;
318  nodes[2] = ne;
319  mesh.addElement(rcp(new ElementDescriptor(gid_a,nodes)),block);
320 
321  nodes[0] = sw;
322  nodes[1] = ne;
323  nodes[2] = nw;
324  mesh.addElement(rcp(new ElementDescriptor(gid_b,nodes)),block);
325  }
326  }
327 }
328 
329 std::pair<int,int> SquareTriMeshFactory::determineXElemSizeAndStart(int xBlock,unsigned int size,unsigned int /* rank */) const
330 {
331  std::size_t xProcLoc = procTuple_[0];
332  unsigned int minElements = nXElems_/size;
333  unsigned int extra = nXElems_ - minElements*size;
334 
335  TEUCHOS_ASSERT(minElements>0);
336 
337  // first "extra" elements get an extra column of elements
338  // this determines the starting X index and number of elements
339  int nume=0, start=0;
340  if(xProcLoc<extra) {
341  nume = minElements+1;
342  start = xProcLoc*(minElements+1);
343  }
344  else {
345  nume = minElements;
346  start = extra*(minElements+1)+(xProcLoc-extra)*minElements;
347  }
348 
349  return std::make_pair(start+nXElems_*xBlock,nume);
350 }
351 
352 std::pair<int,int> SquareTriMeshFactory::determineYElemSizeAndStart(int yBlock,unsigned int size,unsigned int /* rank */) const
353 {
354  std::size_t yProcLoc = procTuple_[1];
355  unsigned int minElements = nYElems_/size;
356  unsigned int extra = nYElems_ - minElements*size;
357 
358  TEUCHOS_ASSERT(minElements>0);
359 
360  // first "extra" elements get an extra column of elements
361  // this determines the starting X index and number of elements
362  int nume=0, start=0;
363  if(yProcLoc<extra) {
364  nume = minElements+1;
365  start = yProcLoc*(minElements+1);
366  }
367  else {
368  nume = minElements;
369  start = extra*(minElements+1)+(yProcLoc-extra)*minElements;
370  }
371 
372  return std::make_pair(start+nYElems_*yBlock,nume);
373 }
374 
376 {
377  mesh.beginModification();
378 
379  std::size_t totalXElems = nXElems_*xBlocks_;
380  std::size_t totalYElems = nYElems_*yBlocks_;
381 
382  // get all part vectors
383  stk::mesh::Part * left = mesh.getSideset("left");
384  stk::mesh::Part * right = mesh.getSideset("right");
385  stk::mesh::Part * top = mesh.getSideset("top");
386  stk::mesh::Part * bottom = mesh.getSideset("bottom");
387 
388  std::vector<stk::mesh::Entity> localElmts;
389  mesh.getMyElements(localElmts);
390 
391  // loop over elements adding edges to sidesets
392  std::vector<stk::mesh::Entity>::const_iterator itr;
393  for(itr=localElmts.begin();itr!=localElmts.end();++itr) {
394  stk::mesh::Entity element = (*itr);
395  stk::mesh::EntityId gid = mesh.elementGlobalId(element);
396 
397  bool lower = (gid%2 != 0);
398  std::size_t block = lower ? (gid+1)/2 : gid/2;
399  std::size_t nx,ny;
400  ny = (block-1) / totalXElems;
401  nx = block-ny*totalXElems-1;
402 
403  // vertical boundaries
405 
406  if(nx+1==totalXElems && lower) {
407  stk::mesh::Entity edge = mesh.findConnectivityById(element, stk::topology::EDGE_RANK, 1);
408 
409  // on the right
410  if(mesh.entityOwnerRank(edge)==machRank_)
411  mesh.addEntityToSideset(edge,right);
412  }
413 
414  if(nx==0 && !lower) {
415  stk::mesh::Entity edge = mesh.findConnectivityById(element, stk::topology::EDGE_RANK, 2);
416 
417  // on the left
418  if(mesh.entityOwnerRank(edge)==machRank_)
419  mesh.addEntityToSideset(edge,left);
420  }
421 
422  // horizontal boundaries
424 
425  if(ny==0 && lower) {
426  stk::mesh::Entity edge = mesh.findConnectivityById(element, stk::topology::EDGE_RANK, 0);
427 
428  // on the bottom
429  if(mesh.entityOwnerRank(edge)==machRank_)
430  mesh.addEntityToSideset(edge,bottom);
431  }
432 
433  if(ny+1==totalYElems && !lower) {
434  stk::mesh::Entity edge = mesh.findConnectivityById(element, stk::topology::EDGE_RANK, 1);
435 
436  // on the top
437  if(mesh.entityOwnerRank(edge)==machRank_)
438  mesh.addEntityToSideset(edge,top);
439  }
440  }
441 
442  mesh.endModification();
443 }
444 
446 {
447  mesh.beginModification();
448 
449  // get all part vectors
450  stk::mesh::Part * origin = mesh.getNodeset("origin");
451 
453  if(machRank_==0)
454  {
455  stk::mesh::Entity node = bulkData->get_entity(mesh.getNodeRank(),1);
456 
457  // add zero node to origin node set
458  mesh.addEntityToNodeset(node,origin);
459  }
460 
461  mesh.endModification();
462 }
463 
465 {
466  mesh.beginModification();
467 
470 
471  stk::mesh::Part * edge_block = mesh.getEdgeBlock(edgeBlockName_);
472 
473  stk::mesh::Selector owned_block = metaData->locally_owned_part();
474 
475  std::vector<stk::mesh::Entity> edges;
476  bulkData->get_entities(mesh.getEdgeRank(), owned_block, edges);
477  mesh.addEntitiesToEdgeBlock(edges, edge_block);
478 
479  mesh.endModification();
480 }
481 
484 {
485  std::size_t i=0,j=0;
486 
487  j = procRank/xProcs_;
488  procRank = procRank % xProcs_;
489  i = procRank;
490 
491  return Teuchos::tuple(i,j);
492 }
493 
494 } // end panzer_stk
virtual Teuchos::RCP< STK_Interface > buildUncommitedMesh(stk::ParallelMachine parallelMach) const
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
Teuchos::Tuple< std::size_t, 2 > procTuple_
void addEntityToNodeset(stk::mesh::Entity entity, stk::mesh::Part *nodeset)
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const
From ParameterListAcceptor.
static void parsePeriodicBCList(const Teuchos::RCP< Teuchos::ParameterList > &pl, std::vector< Teuchos::RCP< const PeriodicBC_MatcherBase > > &periodicBC, bool &useBBoxSearch)
T & get(const std::string &name, T def_value)
stk::mesh::Part * getElementBlockPart(const std::string &name) const
get the block part
void setParameterList(const Teuchos::RCP< Teuchos::ParameterList > &paramList)
From ParameterListAcceptor.
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
void addEntityToSideset(stk::mesh::Entity entity, stk::mesh::Part *sideset)
void addNodeSets(STK_Interface &mesh) const
ParameterList & set(std::string const &name, T &&value, std::string const &docString="", RCP< const ParameterEntryValidator > const &validator=null)
stk::mesh::Part * getNodeset(const std::string &name) const
void addEdgeBlock(const std::string &elemBlockName, const std::string &edgeBlockName, const stk::topology &topology)
Teuchos::RCP< STK_Interface > buildMesh(stk::ParallelMachine parallelMach) const
Build the mesh object.
void buildBlock(stk::ParallelMachine machRank, int xBlock, int yBlock, STK_Interface &mesh) const
stk::mesh::Part * getSideset(const std::string &name) const
void buildMetaData(stk::ParallelMachine parallelMach, STK_Interface &mesh) const
void addElement(const Teuchos::RCP< ElementDescriptor > &ed, stk::mesh::Part *block)
void addEdgeBlocks(STK_Interface &mesh) const
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
void initialize(stk::ParallelMachine parallelMach, bool setupIO=true, const bool buildRefinementSupport=false)
void rebalance(STK_Interface &mesh) const
void addSideSets(STK_Interface &mesh) const
Teuchos::RCP< stk::mesh::BulkData > getBulkData() const
void buildElements(stk::ParallelMachine parallelMach, STK_Interface &mesh) 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
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
std::pair< int, int > determineXElemSizeAndStart(int xBlock, unsigned int size, unsigned int rank) const
bool isInitialized() const
Has initialize been called on this mesh object?
Teuchos::Tuple< std::size_t, 2 > procRankToProcTuple(std::size_t procRank) const
what is the 2D tuple describe this processor distribution
stk::mesh::Entity findConnectivityById(stk::mesh::Entity src, stk::mesh::EntityRank tgt_rank, unsigned rel_id) const
virtual void completeMeshConstruction(STK_Interface &mesh, stk::ParallelMachine parallelMach) const
void addSideset(const std::string &name, const CellTopologyData *ctData)
double getMeshCoord(const int nx, const double deltaX, const double x0) const
stk::mesh::Part * getEdgeBlock(const std::string &name) const
get the block part
std::vector< Teuchos::RCP< const PeriodicBC_MatcherBase > > periodicBCVec_
stk::mesh::EntityRank getEdgeRank() const
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
void addElementBlock(const std::string &name, const CellTopologyData *ctData)