Panzer  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Panzer_STK_LineMeshFactory.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> LineMeshFactory::buildMesh(stk::ParallelMachine parallelMach) const
64 {
65  PANZER_FUNC_TIME_MONITOR("panzer::LineMeshFactory::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> LineMeshFactory::buildUncommitedMesh(stk::ParallelMachine parallelMach) const
80 {
81  PANZER_FUNC_TIME_MONITOR("panzer::LineMeshFactory::buildUncomittedMesh()");
82 
83  RCP<STK_Interface> mesh = rcp(new STK_Interface(1));
84 
85  machRank_ = stk::parallel_machine_rank(parallelMach);
86  machSize_ = stk::parallel_machine_size(parallelMach);
87 
89 
90  // build meta information: blocks and side set setups
91  buildMetaData(parallelMach,*mesh);
92 
93  mesh->addPeriodicBCs(periodicBCVec_);
94 
95  return mesh;
96 }
97 
98 void LineMeshFactory::completeMeshConstruction(STK_Interface & mesh,stk::ParallelMachine parallelMach) const
99 {
100  PANZER_FUNC_TIME_MONITOR("panzer::LineMeshFactory::completeMeshConstruction()");
101 
102  if(not mesh.isInitialized())
103  mesh.initialize(parallelMach);
104 
105  // add node and element information
106  buildElements(parallelMach,mesh);
107 
108  // finish up the edges
109  mesh.buildSubcells();
110  mesh.buildLocalElementIDs();
111 
112  // now that edges are built, sidets can be added
113  addSideSets(mesh);
114 
115  // calls Stk_MeshFactory::rebalance
116  this->rebalance(mesh);
117 }
118 
121 {
123 
124  setMyParamList(paramList);
125 
126  x0_ = paramList->get<double>("X0");
127  xf_ = paramList->get<double>("Xf");
128  xBlocks_ = paramList->get<int>("X Blocks");
129  nXElems_ = paramList->get<int>("X Elements");
130 
131  // read in periodic boundary conditions
132  parsePeriodicBCList(Teuchos::rcpFromRef(paramList->sublist("Periodic BCs")),periodicBCVec_);
133 }
134 
137 {
138  static RCP<Teuchos::ParameterList> defaultParams;
139 
140  // fill with default values
141  if(defaultParams == Teuchos::null) {
142  defaultParams = rcp(new Teuchos::ParameterList);
143 
144  defaultParams->set<double>("X0",0.0);
145  defaultParams->set<double>("Xf",1.0);
146  defaultParams->set<int>("X Blocks",1);
147  defaultParams->set<int>("X Elements",5);
148 
149  Teuchos::ParameterList & bcs = defaultParams->sublist("Periodic BCs");
150  bcs.set<int>("Count",0); // no default periodic boundary conditions
151  }
152 
153  return defaultParams;
154 }
155 
157 {
158  // get valid parameters
160 
161  // set that parameter list
162  setParameterList(validParams);
163 }
164 
165 void LineMeshFactory::buildMetaData(stk::ParallelMachine /* parallelMach */, STK_Interface & mesh) const
166 {
167  typedef shards::Line<2> LineTopo;
168  const CellTopologyData * ctd = shards::getCellTopologyData<LineTopo>();
169  const CellTopologyData * side_ctd = shards::CellTopology(ctd).getBaseCellTopologyData(0,0);
170 
171  // build meta data
172  //mesh.setDimension(2);
173  for(int bx=0;bx<xBlocks_;bx++) {
174  // add this element block
175  {
176  std::stringstream ebPostfix;
177  ebPostfix << "-" << bx;
178 
179  // add element blocks
180  mesh.addElementBlock("eblock"+ebPostfix.str(),ctd);
181  }
182  }
183 
184  // add sidesets
185  mesh.addSideset("left",side_ctd);
186  mesh.addSideset("right",side_ctd);
187 }
188 
189 void LineMeshFactory::buildElements(stk::ParallelMachine parallelMach,STK_Interface & mesh) const
190 {
191  mesh.beginModification();
192  // build each block
193  for(int xBlock=0;xBlock<xBlocks_;xBlock++) {
194  buildBlock(parallelMach,xBlock,mesh);
195  }
196  mesh.endModification();
197 }
198 
199 void LineMeshFactory::buildBlock(stk::ParallelMachine /* parallelMach */, int xBlock, STK_Interface& mesh) const
200 {
201  // grab this processors rank and machine size
202  std::pair<int,int> sizeAndStartX = determineXElemSizeAndStart(xBlock,machSize_,machRank_);
203 
204  int myXElems_start = sizeAndStartX.first;
205  int myXElems_end = myXElems_start+sizeAndStartX.second;
206  int totalXElems = nXElems_*xBlocks_;
207 
208  double deltaX = (xf_-x0_)/double(totalXElems);
209 
210  // build the nodes
211  std::vector<double> coord(1,0.0);
212  for(int nx=myXElems_start;nx<myXElems_end+1;++nx) {
213  coord[0] = double(nx)*deltaX+x0_;
214  mesh.addNode(nx+1,coord);
215  }
216 
217  std::stringstream blockName;
218  blockName << "eblock-" << xBlock;
219  stk::mesh::Part * block = mesh.getElementBlockPart(blockName.str());
220 
221  // build the elements
222  for(int nx=myXElems_start;nx<myXElems_end;++nx) {
223  stk::mesh::EntityId gid = nx+1;
224  std::vector<stk::mesh::EntityId> nodes(2);
225  nodes[0] = nx+1;
226  nodes[1] = nodes[0]+1;
227 
228  RCP<ElementDescriptor> ed = rcp(new ElementDescriptor(gid,nodes));
229  mesh.addElement(ed,block);
230  }
231 }
232 
233 std::pair<int,int> LineMeshFactory::determineXElemSizeAndStart(int xBlock,unsigned int size,unsigned int /* rank */) const
234 {
235  std::size_t xProcLoc = procTuple_[0];
236  unsigned int minElements = nXElems_/size;
237  unsigned int extra = nXElems_ - minElements*size;
238 
239  TEUCHOS_ASSERT(minElements>0);
240 
241  // first "extra" elements get an extra column of elements
242  // this determines the starting X index and number of elements
243  int nume=0, start=0;
244  if(xProcLoc<extra) {
245  nume = minElements+1;
246  start = xProcLoc*(minElements+1);
247  }
248  else {
249  nume = minElements;
250  start = extra*(minElements+1)+(xProcLoc-extra)*minElements;
251  }
252 
253  return std::make_pair(start+nXElems_*xBlock,nume);
254 }
255 
257 {
258  mesh.beginModification();
259  const stk::mesh::EntityRank sideRank = mesh.getSideRank();
260 
261  std::size_t totalXElems = nXElems_*xBlocks_;
262 
263  // get all part vectors
264  stk::mesh::Part * left = mesh.getSideset("left");
265  stk::mesh::Part * right = mesh.getSideset("right");
266 
267  std::vector<stk::mesh::Entity> localElmts;
268  mesh.getMyElements(localElmts);
269 
270  // loop over elements adding edges to sidesets
271  std::vector<stk::mesh::Entity>::const_iterator itr;
272  for(itr=localElmts.begin();itr!=localElmts.end();++itr) {
273  stk::mesh::Entity element = (*itr);
274  stk::mesh::EntityId gid = mesh.elementGlobalId(element);
275 
276  std::size_t nx = gid-1;
277 
278  // vertical boundaries
280 
281  if(nx+1==totalXElems) {
282  stk::mesh::Entity edge = mesh.findConnectivityById(element, sideRank, 1);
283 
284  // on the right
285  if(mesh.entityOwnerRank(edge)==machRank_)
286  mesh.addEntityToSideset(edge,right);
287  }
288 
289  if(nx==0) {
290  stk::mesh::Entity edge = mesh.findConnectivityById(element, sideRank, 0);
291 
292  // on the left
293  if(mesh.entityOwnerRank(edge)==machRank_)
294  mesh.addEntityToSideset(edge,left);
295  }
296  }
297 
298  mesh.endModification();
299 }
300 
303 {
304  std::size_t i=0,j=0;
305 
306  j = procRank/machSize_;
307  procRank = procRank % machSize_;
308  i = procRank;
309 
310  return Teuchos::tuple(i,j);
311 }
312 
313 } // end panzer_stk
void buildMetaData(stk::ParallelMachine parallelMach, STK_Interface &mesh) const
Teuchos::Tuple< std::size_t, 2 > procRankToProcTuple(std::size_t procRank) const
what is the 2D tuple describe this processor distribution
unsigned entityOwnerRank(stk::mesh::Entity entity) const
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 addEntityToSideset(stk::mesh::Entity entity, stk::mesh::Part *sideset)
Teuchos::Tuple< std::size_t, 2 > procTuple_
virtual Teuchos::RCP< STK_Interface > buildUncommitedMesh(stk::ParallelMachine parallelMach) const
void initialize(stk::ParallelMachine parallelMach, bool setupIO=true)
void buildBlock(stk::ParallelMachine machRank, int xBlock, STK_Interface &mesh) const
stk::mesh::Part * getSideset(const std::string &name) 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)
stk::mesh::EntityRank getSideRank() const
Teuchos::RCP< STK_Interface > buildMesh(stk::ParallelMachine parallelMach) const
Build the mesh object.
void rebalance(STK_Interface &mesh) const
std::pair< int, int > determineXElemSizeAndStart(int xBlock, unsigned int size, unsigned int rank) const
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const
From ParameterListAcceptor.
virtual void completeMeshConstruction(STK_Interface &mesh, stk::ParallelMachine parallelMach) 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
bool isInitialized() const
Has initialize been called on this mesh object?
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::vector< Teuchos::RCP< const PeriodicBC_MatcherBase > > periodicBCVec_
void buildElements(stk::ParallelMachine parallelMach, STK_Interface &mesh) 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)
void addElementBlock(const std::string &name, const CellTopologyData *ctData)
void addSideSets(STK_Interface &mesh) const
void setParameterList(const Teuchos::RCP< Teuchos::ParameterList > &paramList)
From ParameterListAcceptor.
static void parsePeriodicBCList(const Teuchos::RCP< Teuchos::ParameterList > &pl, std::vector< Teuchos::RCP< const PeriodicBC_MatcherBase > > &periodicBC)