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 //
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> SquareTriMeshFactory::buildMesh(stk::ParallelMachine parallelMach) const
64 {
65  PANZER_FUNC_TIME_MONITOR("panzer::SquareTriMeshFactory::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 
80 {
81  PANZER_FUNC_TIME_MONITOR("panzer::SquareTriMeshFactory::buildUncomittedMesh()");
82 
83  RCP<STK_Interface> mesh = rcp(new STK_Interface(2));
84 
85  machRank_ = stk::parallel_machine_rank(parallelMach);
86  machSize_ = stk::parallel_machine_size(parallelMach);
87 
88  if (xProcs_ == -1 && yProcs_ == -1) {
89  // copied from galeri
90  xProcs_ = yProcs_ = Teuchos::as<int>(pow(Teuchos::as<double>(machSize_), 0.5));
91 
92  if (xProcs_ * yProcs_ != Teuchos::as<int>(machSize_)) {
93  // Simple method to find a set of processor assignments
94  xProcs_ = yProcs_ = 1;
95 
96  // This means that this works correctly up to about maxFactor^2
97  // processors.
98  const int maxFactor = 100;
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_ = xProcs_*jj;
120  else yProcs_ = yProcs_*jj;
121  factors[jj]--;
122  }
123  }
124  }
125 
126  } else if(xProcs_==-1) {
127  // default x only decomposition
128  xProcs_ = machSize_;
129  yProcs_ = 1;
130  }
131  TEUCHOS_TEST_FOR_EXCEPTION(int(machSize_)!=xProcs_*yProcs_,std::logic_error,
132  "Cannot build SquareTriMeshFactory, the product of \"X Procs\" and \"Y Procs\""
133  " must equal the number of processors.");
135 
136  // build meta information: blocks and side set setups
137  buildMetaData(parallelMach,*mesh);
138 
139  mesh->addPeriodicBCs(periodicBCVec_);
140 
141  return mesh;
142 }
143 
144 void SquareTriMeshFactory::completeMeshConstruction(STK_Interface & mesh,stk::ParallelMachine parallelMach) const
145 {
146  PANZER_FUNC_TIME_MONITOR("panzer::SquareTriMeshFactory::completeMeshConstruction()");
147 
148  if(not mesh.isInitialized())
149  mesh.initialize(parallelMach);
150 
151  // add node and element information
152  buildElements(parallelMach,mesh);
153 
154  // finish up the edges
155  mesh.buildSubcells();
156  mesh.buildLocalElementIDs();
157 
158  // now that edges are built, sidets can be added
159  addSideSets(mesh);
160 
161  // add nodesets
162  addNodeSets(mesh);
163 
164  // calls Stk_MeshFactory::rebalance
165  this->rebalance(mesh);
166 }
167 
170 {
172 
173  setMyParamList(paramList);
174 
175  x0_ = paramList->get<double>("X0");
176  y0_ = paramList->get<double>("Y0");
177 
178  xf_ = paramList->get<double>("Xf");
179  yf_ = paramList->get<double>("Yf");
180 
181  xBlocks_ = paramList->get<int>("X Blocks");
182  yBlocks_ = paramList->get<int>("Y Blocks");
183 
184  nXElems_ = paramList->get<int>("X Elements");
185  nYElems_ = paramList->get<int>("Y Elements");
186 
187  xProcs_ = paramList->get<int>("X Procs");
188  yProcs_ = paramList->get<int>("Y Procs");
189 
190  // read in periodic boundary conditions
191  parsePeriodicBCList(Teuchos::rcpFromRef(paramList->sublist("Periodic BCs")),periodicBCVec_);
192 }
193 
196 {
197  static RCP<Teuchos::ParameterList> defaultParams;
198 
199  // fill with default values
200  if(defaultParams == Teuchos::null) {
201  defaultParams = rcp(new Teuchos::ParameterList);
202 
203  defaultParams->set<double>("X0",0.0);
204  defaultParams->set<double>("Y0",0.0);
205 
206  defaultParams->set<double>("Xf",1.0);
207  defaultParams->set<double>("Yf",1.0);
208 
209  defaultParams->set<int>("X Blocks",1);
210  defaultParams->set<int>("Y Blocks",1);
211 
212  defaultParams->set<int>("X Procs",-1);
213  defaultParams->set<int>("Y Procs",1);
214 
215  defaultParams->set<int>("X Elements",5);
216  defaultParams->set<int>("Y Elements",5);
217 
218  Teuchos::ParameterList & bcs = defaultParams->sublist("Periodic BCs");
219  bcs.set<int>("Count",0); // no default periodic boundary conditions
220  }
221 
222  return defaultParams;
223 }
224 
226 {
227  // get valid parameters
229 
230  // set that parameter list
231  setParameterList(validParams);
232 }
233 
234 void SquareTriMeshFactory::buildMetaData(stk::ParallelMachine /* parallelMach */, STK_Interface & mesh) const
235 {
236  typedef shards::Triangle<> TriTopo;
237  const CellTopologyData * ctd = shards::getCellTopologyData<TriTopo>();
238  const CellTopologyData * side_ctd = shards::CellTopology(ctd).getBaseCellTopologyData(1,0);
239 
240  // build meta data
241  //mesh.setDimension(2);
242  for(int bx=0;bx<xBlocks_;bx++) {
243  for(int by=0;by<yBlocks_;by++) {
244 
245  // add this element block
246  {
247  std::stringstream ebPostfix;
248  ebPostfix << "-" << bx << "_" << by;
249 
250  // add element blocks
251  mesh.addElementBlock("eblock"+ebPostfix.str(),ctd);
252  }
253 
254  }
255  }
256 
257  // add sidesets
258  mesh.addSideset("left",side_ctd);
259  mesh.addSideset("right",side_ctd);
260  mesh.addSideset("top",side_ctd);
261  mesh.addSideset("bottom",side_ctd);
262 
263  // add nodesets
264  mesh.addNodeset("origin");
265 }
266 
267 void SquareTriMeshFactory::buildElements(stk::ParallelMachine parallelMach,STK_Interface & mesh) const
268 {
269  mesh.beginModification();
270  // build each block
271  for(int xBlock=0;xBlock<xBlocks_;xBlock++) {
272  for(int yBlock=0;yBlock<yBlocks_;yBlock++) {
273  buildBlock(parallelMach,xBlock,yBlock,mesh);
274  }
275  }
276  mesh.endModification();
277 }
278 
279 void SquareTriMeshFactory::buildBlock(stk::ParallelMachine /* parallelMach */, int xBlock, int yBlock, STK_Interface& mesh) const
280 {
281  // grab this processors rank and machine size
282  std::pair<int,int> sizeAndStartX = determineXElemSizeAndStart(xBlock,xProcs_,machRank_);
283  std::pair<int,int> sizeAndStartY = determineYElemSizeAndStart(yBlock,yProcs_,machRank_);
284 
285  int myXElems_start = sizeAndStartX.first;
286  int myXElems_end = myXElems_start+sizeAndStartX.second;
287  int myYElems_start = sizeAndStartY.first;
288  int myYElems_end = myYElems_start+sizeAndStartY.second;
289  int totalXElems = nXElems_*xBlocks_;
290  int totalYElems = nYElems_*yBlocks_;
291 
292  double deltaX = (xf_-x0_)/double(totalXElems);
293  double deltaY = (yf_-y0_)/double(totalYElems);
294 
295  std::vector<double> coord(2,0.0);
296 
297  // build the nodes
298  for(int nx=myXElems_start;nx<myXElems_end+1;++nx) {
299  coord[0] = double(nx)*deltaX+x0_;
300  for(int ny=myYElems_start;ny<myYElems_end+1;++ny) {
301  coord[1] = double(ny)*deltaY+y0_;
302 
303  mesh.addNode(ny*(totalXElems+1)+nx+1,coord);
304  }
305  }
306 
307  std::stringstream blockName;
308  blockName << "eblock-" << xBlock << "_" << yBlock;
309  stk::mesh::Part * block = mesh.getElementBlockPart(blockName.str());
310 
311  // build the elements
312  for(int nx=myXElems_start;nx<myXElems_end;++nx) {
313  for(int ny=myYElems_start;ny<myYElems_end;++ny) {
314  stk::mesh::EntityId gid_a = 2*(totalXElems*ny+nx+1)-1;
315  stk::mesh::EntityId gid_b = gid_a+1;
316  std::vector<stk::mesh::EntityId> nodes(3);
317  stk::mesh::EntityId sw,se,ne,nw;
318  sw = nx+1+ny*(totalXElems+1);
319  se = sw+1;
320  ne = se+(totalXElems+1);
321  nw = ne-1;
322 
323  nodes[0] = sw;
324  nodes[1] = se;
325  nodes[2] = ne;
326  mesh.addElement(rcp(new ElementDescriptor(gid_a,nodes)),block);
327 
328  nodes[0] = sw;
329  nodes[1] = ne;
330  nodes[2] = nw;
331  mesh.addElement(rcp(new ElementDescriptor(gid_b,nodes)),block);
332  }
333  }
334 }
335 
336 std::pair<int,int> SquareTriMeshFactory::determineXElemSizeAndStart(int xBlock,unsigned int size,unsigned int /* rank */) const
337 {
338  std::size_t xProcLoc = procTuple_[0];
339  unsigned int minElements = nXElems_/size;
340  unsigned int extra = nXElems_ - minElements*size;
341 
342  TEUCHOS_ASSERT(minElements>0);
343 
344  // first "extra" elements get an extra column of elements
345  // this determines the starting X index and number of elements
346  int nume=0, start=0;
347  if(xProcLoc<extra) {
348  nume = minElements+1;
349  start = xProcLoc*(minElements+1);
350  }
351  else {
352  nume = minElements;
353  start = extra*(minElements+1)+(xProcLoc-extra)*minElements;
354  }
355 
356  return std::make_pair(start+nXElems_*xBlock,nume);
357 }
358 
359 std::pair<int,int> SquareTriMeshFactory::determineYElemSizeAndStart(int yBlock,unsigned int size,unsigned int /* rank */) const
360 {
361  std::size_t yProcLoc = procTuple_[1];
362  unsigned int minElements = nYElems_/size;
363  unsigned int extra = nYElems_ - minElements*size;
364 
365  TEUCHOS_ASSERT(minElements>0);
366 
367  // first "extra" elements get an extra column of elements
368  // this determines the starting X index and number of elements
369  int nume=0, start=0;
370  if(yProcLoc<extra) {
371  nume = minElements+1;
372  start = yProcLoc*(minElements+1);
373  }
374  else {
375  nume = minElements;
376  start = extra*(minElements+1)+(yProcLoc-extra)*minElements;
377  }
378 
379  return std::make_pair(start+nYElems_*yBlock,nume);
380 }
381 
383 {
384  mesh.beginModification();
385 
386  std::size_t totalXElems = nXElems_*xBlocks_;
387  std::size_t totalYElems = nYElems_*yBlocks_;
388 
389  // get all part vectors
390  stk::mesh::Part * left = mesh.getSideset("left");
391  stk::mesh::Part * right = mesh.getSideset("right");
392  stk::mesh::Part * top = mesh.getSideset("top");
393  stk::mesh::Part * bottom = mesh.getSideset("bottom");
394 
395  std::vector<stk::mesh::Entity> localElmts;
396  mesh.getMyElements(localElmts);
397 
398  // loop over elements adding edges to sidesets
399  std::vector<stk::mesh::Entity>::const_iterator itr;
400  for(itr=localElmts.begin();itr!=localElmts.end();++itr) {
401  stk::mesh::Entity element = (*itr);
402  stk::mesh::EntityId gid = mesh.elementGlobalId(element);
403 
404  bool lower = (gid%2 != 0);
405  std::size_t block = lower ? (gid+1)/2 : gid/2;
406  std::size_t nx,ny;
407  ny = (block-1) / totalXElems;
408  nx = block-ny*totalXElems-1;
409 
410  // vertical boundaries
412 
413  if(nx+1==totalXElems && lower) {
414  stk::mesh::Entity edge = mesh.findConnectivityById(element, stk::topology::EDGE_RANK, 1);
415 
416  // on the right
417  if(mesh.entityOwnerRank(edge)==machRank_)
418  mesh.addEntityToSideset(edge,right);
419  }
420 
421  if(nx==0 && !lower) {
422  stk::mesh::Entity edge = mesh.findConnectivityById(element, stk::topology::EDGE_RANK, 2);
423 
424  // on the left
425  if(mesh.entityOwnerRank(edge)==machRank_)
426  mesh.addEntityToSideset(edge,left);
427  }
428 
429  // horizontal boundaries
431 
432  if(ny==0 && lower) {
433  stk::mesh::Entity edge = mesh.findConnectivityById(element, stk::topology::EDGE_RANK, 0);
434 
435  // on the bottom
436  if(mesh.entityOwnerRank(edge)==machRank_)
437  mesh.addEntityToSideset(edge,bottom);
438  }
439 
440  if(ny+1==totalYElems && !lower) {
441  stk::mesh::Entity edge = mesh.findConnectivityById(element, stk::topology::EDGE_RANK, 1);
442 
443  // on the top
444  if(mesh.entityOwnerRank(edge)==machRank_)
445  mesh.addEntityToSideset(edge,top);
446  }
447  }
448 
449  mesh.endModification();
450 }
451 
453 {
454  mesh.beginModification();
455 
456  // get all part vectors
457  stk::mesh::Part * origin = mesh.getNodeset("origin");
458 
460  if(machRank_==0)
461  {
462  stk::mesh::Entity node = bulkData->get_entity(mesh.getNodeRank(),1);
463 
464  // add zero node to origin node set
465  mesh.addEntityToNodeset(node,origin);
466  }
467 
468  mesh.endModification();
469 }
470 
473 {
474  std::size_t i=0,j=0;
475 
476  j = procRank/xProcs_;
477  procRank = procRank % xProcs_;
478  i = procRank;
479 
480  return Teuchos::tuple(i,j);
481 }
482 
483 } // 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
void addEntityToNodeset(stk::mesh::Entity entity, stk::mesh::Part *nodeset)
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const
From ParameterListAcceptor.
T & get(const std::string &name, T def_value)
stk::mesh::Part * getElementBlockPart(const std::string &name) const
get the block count
ParameterList & set(std::string const &name, T const &value, std::string const &docString="", RCP< const ParameterEntryValidator > const &validator=null)
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
stk::mesh::Part * getNodeset(const std::string &name) const
void initialize(stk::ParallelMachine parallelMach, bool setupIO=true)
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)
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
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)
void addNode(stk::mesh::EntityId gid, const std::vector< double > &coord)
void buildSubcells()
force the mesh to build subcells: edges and faces
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)
std::vector< Teuchos::RCP< const PeriodicBC_MatcherBase > > periodicBCVec_
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)
static void parsePeriodicBCList(const Teuchos::RCP< Teuchos::ParameterList > &pl, std::vector< Teuchos::RCP< const PeriodicBC_MatcherBase > > &periodicBC)