Panzer  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Panzer_STK_MultiBlockMeshFactory.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 
13 using Teuchos::RCP;
14 using Teuchos::rcp;
15 
16 namespace panzer_stk {
17 
19 {
21 }
22 
25 {
26 }
27 
29 Teuchos::RCP<STK_Interface> MultiBlockMeshFactory::buildMesh(stk::ParallelMachine parallelMach) const
30 {
31  // build all meta data
32  RCP<STK_Interface> mesh = buildUncommitedMesh(parallelMach);
33 
34  // commit meta data
35  mesh->initialize(parallelMach);
36 
37  // build bulk data
38  completeMeshConstruction(*mesh,parallelMach);
39 
40  return mesh;
41 }
42 
44 {
45  RCP<STK_Interface> mesh = rcp(new STK_Interface(2));
46 
47  machRank_ = stk::parallel_machine_rank(parallelMach);
48  machSize_ = stk::parallel_machine_size(parallelMach);
49 
50  // build meta information: blocks and side set setups
51  buildMetaData(parallelMach,*mesh);
52 
53  return mesh;
54 }
55 
56 void MultiBlockMeshFactory::completeMeshConstruction(STK_Interface & mesh,stk::ParallelMachine parallelMach) const
57 {
58  if(not mesh.isInitialized())
59  mesh.initialize(parallelMach);
60 
61  // add node and element information
62  buildElements(parallelMach,mesh);
63 
64  // finish up the edges
65  mesh.buildSubcells();
66  mesh.buildLocalElementIDs();
67 
68  // now that edges are built, sidets can be added
69  addSideSets(mesh);
70 
71  // calls Stk_MeshFactory::rebalance
72  this->rebalance(mesh);
73 }
74 
77 {
78 }
79 
82 {
83  static RCP<Teuchos::ParameterList> defaultParams;
84 
85  // fill with default values
86  if(defaultParams == Teuchos::null) {
87  defaultParams = rcp(new Teuchos::ParameterList);
88 
89  defaultParams->set<int>("X Elements",2);
90  defaultParams->set<int>("Y Elements",2);
91  }
92 
93  return defaultParams;
94 }
95 
97 {
98  nXElems_ = 2;
99  nYElems_ = 2;
100 }
101 
102 void MultiBlockMeshFactory::buildMetaData(stk::ParallelMachine /* parallelMach */, STK_Interface& mesh) const
103 {
104  typedef shards::Quadrilateral<4> QuadTopo;
105  const CellTopologyData * ctd = shards::getCellTopologyData<QuadTopo>();
106  const CellTopologyData * side_ctd = shards::CellTopology(ctd).getBaseCellTopologyData(1,0);
107 
108  // build meta data
109  //mesh.setDimension(2);
110  for(int bx=0;bx<2;bx++) {
111  // add this element block
112  std::stringstream ebPostfix;
113  ebPostfix << "-" << "0_" << bx;
114 
115  // add element blocks
116  mesh.addElementBlock("eblock"+ebPostfix.str(),ctd);
117  }
118 
119  // add sidesets
120  mesh.addSideset("left",side_ctd);
121  mesh.addSideset("right",side_ctd);
122  mesh.addSideset("top",side_ctd);
123  mesh.addSideset("bottom",side_ctd);
124 }
125 
126 void MultiBlockMeshFactory::buildElements(stk::ParallelMachine parallelMach,STK_Interface & mesh) const
127 {
128  mesh.beginModification();
129 
130  if(machRank_==0) {
131  buildBlock(parallelMach,0,0,mesh);
132  }
133  else if(machRank_==1) {
134  buildBlock(parallelMach,0,1,mesh);
135  }
136  else TEUCHOS_ASSERT(false);
137 
138  mesh.endModification();
139 }
140 
141 void MultiBlockMeshFactory::buildBlock(stk::ParallelMachine /* parallelMach */, int xBlock, int yBlock, STK_Interface& mesh) const
142 {
143  int myXElems_start = (machRank_==0 ? 0 : 2);
144  int myXElems_end = (machRank_==0 ? 1 : 3);
145  int myYElems_start = 0;
146  int myYElems_end = 1;
147  int totalXElems = nXElems_*2;
148  int totalYElems = nYElems_*1;
149 
150  double x0_ = 0.0;
151  double xf_ = 1.0;
152  double y0_ = 0.0;
153  double yf_ = 1.0;
154  double deltaX = (xf_-x0_)/double(totalXElems);
155  double deltaY = (yf_-y0_)/double(totalYElems);
156 
157  std::vector<double> coord(2,0.0);
158 
159  // build the nodes
160  for(int nx=myXElems_start;nx<myXElems_end+1;++nx) {
161  coord[0] = double(nx)*deltaX+x0_;
162  for(int ny=myYElems_start;ny<myYElems_end+1;++ny) {
163  coord[1] = double(ny)*deltaY+y0_;
164 
165  mesh.addNode(ny*(totalXElems+1)+nx+1,coord);
166  }
167  }
168 
169  std::stringstream blockName;
170  blockName << "eblock-" << xBlock << "_" << yBlock;
171  stk::mesh::Part * block = mesh.getElementBlockPart(blockName.str());
172 
173  // build the elements
174  for(int nx=myXElems_start;nx<myXElems_end;++nx) {
175  for(int ny=myYElems_start;ny<myYElems_end;++ny) {
176  stk::mesh::EntityId gid = totalXElems*ny+nx+1;
177  std::vector<stk::mesh::EntityId> nodes(4);
178  nodes[0] = nx+1+ny*(totalXElems+1);
179  nodes[1] = nodes[0]+1;
180  nodes[2] = nodes[1]+(totalXElems+1);
181  nodes[3] = nodes[2]-1;
182 
183  RCP<ElementDescriptor> ed = rcp(new ElementDescriptor(gid,nodes));
184  mesh.addElement(ed,block);
185  }
186  }
187 }
188 
189 std::pair<int,int> MultiBlockMeshFactory::determineXElemSizeAndStart(int xBlock,unsigned int size,unsigned int rank) const
190 {
191  unsigned int minElements = nXElems_/size;
192  unsigned int extra = nXElems_ - minElements*size;
193 
194  TEUCHOS_ASSERT(minElements>0);
195 
196  // first "extra" elements get an extra column of elements
197  // this determines the starting X index and number of elements
198  int nume=0, start=0;
199  if(rank<extra) {
200  nume = minElements+1;
201  start = rank*(minElements+1);
202  }
203  else {
204  nume = minElements;
205  start = extra*(minElements+1)+(rank-extra)*minElements;
206  }
207 
208  return std::make_pair(start+nXElems_*xBlock,nume);
209 }
210 
211 std::pair<int,int> MultiBlockMeshFactory::determineYElemSizeAndStart(int yBlock, unsigned int /* size */, unsigned int /* rank */) const
212 {
213  int start = yBlock*nYElems_;
214 
215  return std::make_pair(start,nYElems_);
216 }
217 
219 {
220  mesh.beginModification();
221 
222  std::size_t totalXElems = nXElems_*2;
223  std::size_t totalYElems = nYElems_*1;
224 
225  // get all part vectors
226  stk::mesh::Part * left = mesh.getSideset("left");
227  stk::mesh::Part * right = mesh.getSideset("right");
228  stk::mesh::Part * top = mesh.getSideset("top");
229  stk::mesh::Part * bottom = mesh.getSideset("bottom");
230 
231  std::vector<stk::mesh::Entity> localElmts;
232  mesh.getMyElements(localElmts);
233 
234  // loop over elements adding edges to sidesets
235  std::vector<stk::mesh::Entity>::const_iterator itr;
236  for(itr=localElmts.begin();itr!=localElmts.end();++itr) {
237  stk::mesh::Entity element = (*itr);
238  stk::mesh::EntityId gid = mesh.elementGlobalId(element);
239 
240  std::size_t nx,ny;
241  ny = (gid-1) / totalXElems;
242  nx = gid-ny*totalXElems-1;
243 
244  // vertical boundaries
246 
247  if(nx+1==totalXElems) {
248  stk::mesh::Entity edge = mesh.findConnectivityById(element, stk::topology::EDGE_RANK, 1);
249 
250  // on the right
251  if(mesh.entityOwnerRank(edge)==machRank_)
252  mesh.addEntityToSideset(edge,right);
253  }
254 
255  if(nx==0) {
256  stk::mesh::Entity edge = mesh.findConnectivityById(element, stk::topology::EDGE_RANK, 3);
257 
258  // on the left
259  if(mesh.entityOwnerRank(edge)==machRank_)
260  mesh.addEntityToSideset(edge,left);
261  }
262 
263  // horizontal boundaries
265 
266  if(ny==0) {
267  stk::mesh::Entity edge = mesh.findConnectivityById(element, stk::topology::EDGE_RANK, 0);
268 
269  // on the bottom
270  if(mesh.entityOwnerRank(edge)==machRank_)
271  mesh.addEntityToSideset(edge,bottom);
272  }
273 
274  if(ny+1==totalYElems) {
275  stk::mesh::Entity edge = mesh.findConnectivityById(element, stk::topology::EDGE_RANK, 2);
276 
277  // on the top
278  if(mesh.entityOwnerRank(edge)==machRank_)
279  mesh.addEntityToSideset(edge,top);
280  }
281  }
282 
283  mesh.endModification();
284 }
285 
286 } // end panzer_stk
virtual Teuchos::RCP< STK_Interface > buildUncommitedMesh(stk::ParallelMachine parallelMach) const
unsigned entityOwnerRank(stk::mesh::Entity entity) const
stk::mesh::Part * getElementBlockPart(const std::string &name) const
get the block part
void addEntityToSideset(stk::mesh::Entity entity, stk::mesh::Part *sideset)
ParameterList & set(std::string const &name, T &&value, std::string const &docString="", RCP< const ParameterEntryValidator > const &validator=null)
stk::mesh::Part * getSideset(const std::string &name) const
void buildBlock(stk::ParallelMachine parallelMach, int xBlock, int yBlock, 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 initialize(stk::ParallelMachine parallelMach, bool setupIO=true, const bool buildRefinementSupport=false)
void rebalance(STK_Interface &mesh) const
stk::mesh::EntityId elementGlobalId(std::size_t lid) const
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const
From ParameterListAcceptor.
Teuchos::RCP< STK_Interface > buildMesh(stk::ParallelMachine parallelMach) const
Build the mesh object.
void addNode(stk::mesh::EntityId gid, const std::vector< double > &coord)
void buildSubcells()
force the mesh to build subcells: edges and faces
bool isInitialized() const
Has initialize been called on this mesh object?
void buildElements(stk::ParallelMachine parallelMach, STK_Interface &mesh) const
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)
std::pair< int, int > determineYElemSizeAndStart(int yBlock, unsigned int size, unsigned int rank) const
std::pair< int, int > determineXElemSizeAndStart(int xBlock, unsigned int size, unsigned int rank) const
void setParameterList(const Teuchos::RCP< Teuchos::ParameterList > &paramList)
From ParameterListAcceptor.
virtual void completeMeshConstruction(STK_Interface &mesh, stk::ParallelMachine parallelMach) const
void getMyElements(std::vector< stk::mesh::Entity > &elements) const
#define TEUCHOS_ASSERT(assertion_test)
void buildMetaData(stk::ParallelMachine parallelMach, STK_Interface &mesh) const
void addElementBlock(const std::string &name, const CellTopologyData *ctData)