Panzer Version of the Day
Loading...
Searching...
No Matches
Panzer_STK_LineMeshFactory.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
15using Teuchos::RCP;
16using Teuchos::rcp;
17
18namespace panzer_stk {
19
24
29
31Teuchos::RCP<STK_Interface> LineMeshFactory::buildMesh(stk::ParallelMachine parallelMach) const
32{
33 PANZER_FUNC_TIME_MONITOR("panzer::LineMeshFactory::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
47Teuchos::RCP<STK_Interface> LineMeshFactory::buildUncommitedMesh(stk::ParallelMachine parallelMach) const
48{
49 PANZER_FUNC_TIME_MONITOR("panzer::LineMeshFactory::buildUncomittedMesh()");
50
51 RCP<STK_Interface> mesh = rcp(new STK_Interface(1));
52
53 machRank_ = stk::parallel_machine_rank(parallelMach);
54 machSize_ = stk::parallel_machine_size(parallelMach);
55
57
58 // build meta information: blocks and side set setups
59 buildMetaData(parallelMach,*mesh);
60
61 mesh->addPeriodicBCs(periodicBCVec_);
62 mesh->setBoundingBoxSearchFlag(useBBoxSearch_);
63
64 return mesh;
65}
66
67void LineMeshFactory::completeMeshConstruction(STK_Interface & mesh,stk::ParallelMachine parallelMach) const
68{
69 PANZER_FUNC_TIME_MONITOR("panzer::LineMeshFactory::completeMeshConstruction()");
70
71 if(not mesh.isInitialized())
72 mesh.initialize(parallelMach);
73
74 // add node and element information
75 buildElements(parallelMach,mesh);
76
77 // finish up the edges
78 mesh.buildSubcells();
80
81 // now that edges are built, sidets can be added
82 addSideSets(mesh);
83
84 // calls Stk_MeshFactory::rebalance
85 this->rebalance(mesh);
86}
87
89void LineMeshFactory::setParameterList(const Teuchos::RCP<Teuchos::ParameterList> & paramList)
90{
91 paramList->validateParametersAndSetDefaults(*getValidParameters(),0);
92
93 setMyParamList(paramList);
94
95 x0_ = paramList->get<double>("X0");
96 xf_ = paramList->get<double>("Xf");
97 xBlocks_ = paramList->get<int>("X Blocks");
98 nXElems_ = paramList->get<int>("X Elements");
99
100 // read in periodic boundary conditions
101 parsePeriodicBCList(Teuchos::rcpFromRef(paramList->sublist("Periodic BCs")),periodicBCVec_,useBBoxSearch_);
102}
103
105Teuchos::RCP<const Teuchos::ParameterList> LineMeshFactory::getValidParameters() const
106{
107 static RCP<Teuchos::ParameterList> defaultParams;
108
109 // fill with default values
110 if(defaultParams == Teuchos::null) {
111 defaultParams = rcp(new Teuchos::ParameterList);
112
113 defaultParams->set<double>("X0",0.0);
114 defaultParams->set<double>("Xf",1.0);
115 defaultParams->set<int>("X Blocks",1);
116 defaultParams->set<int>("X Elements",5);
117
118 Teuchos::ParameterList & bcs = defaultParams->sublist("Periodic BCs");
119 bcs.set<int>("Count",0); // no default periodic boundary conditions
120 }
121
122 return defaultParams;
123}
124
126{
127 // get valid parameters
128 RCP<Teuchos::ParameterList> validParams = rcp(new Teuchos::ParameterList(*getValidParameters()));
129
130 // set that parameter list
131 setParameterList(validParams);
132}
133
134void LineMeshFactory::buildMetaData(stk::ParallelMachine /* parallelMach */, STK_Interface & mesh) const
135{
136 typedef shards::Line<2> LineTopo;
137 const CellTopologyData * ctd = shards::getCellTopologyData<LineTopo>();
138 const CellTopologyData * side_ctd = shards::CellTopology(ctd).getBaseCellTopologyData(0,0);
139
140 // build meta data
141 //mesh.setDimension(2);
142 for(int bx=0;bx<xBlocks_;bx++) {
143 // add this element block
144 {
145 std::stringstream ebPostfix;
146 ebPostfix << "-" << bx;
147
148 // add element blocks
149 mesh.addElementBlock("eblock"+ebPostfix.str(),ctd);
150 }
151 }
152
153 // add sidesets
154 mesh.addSideset("left",side_ctd);
155 mesh.addSideset("right",side_ctd);
156}
157
158void LineMeshFactory::buildElements(stk::ParallelMachine parallelMach,STK_Interface & mesh) const
159{
160 mesh.beginModification();
161 // build each block
162 for(int xBlock=0;xBlock<xBlocks_;xBlock++) {
163 buildBlock(parallelMach,xBlock,mesh);
164 }
165 mesh.endModification();
166}
167
168void LineMeshFactory::buildBlock(stk::ParallelMachine /* parallelMach */, int xBlock, STK_Interface& mesh) const
169{
170 // grab this processors rank and machine size
171 std::pair<int,int> sizeAndStartX = determineXElemSizeAndStart(xBlock,machSize_,machRank_);
172
173 int myXElems_start = sizeAndStartX.first;
174 int myXElems_end = myXElems_start+sizeAndStartX.second;
175 int totalXElems = nXElems_*xBlocks_;
176
177 double deltaX = (xf_-x0_)/static_cast<double>(totalXElems);
178
179 // build the nodes
180 std::vector<double> coord(1,0.0);
181 for(int nx=myXElems_start;nx<myXElems_end+1;++nx) {
182 coord[0] = this->getMeshCoord(nx, deltaX, x0_);
183 mesh.addNode(nx+1,coord);
184 }
185
186 std::stringstream blockName;
187 blockName << "eblock-" << xBlock;
188 stk::mesh::Part * block = mesh.getElementBlockPart(blockName.str());
189
190 // build the elements
191 for(int nx=myXElems_start;nx<myXElems_end;++nx) {
192 stk::mesh::EntityId gid = nx+1;
193 std::vector<stk::mesh::EntityId> nodes(2);
194 nodes[0] = nx+1;
195 nodes[1] = nodes[0]+1;
196
197 RCP<ElementDescriptor> ed = rcp(new ElementDescriptor(gid,nodes));
198 mesh.addElement(ed,block);
199 }
200}
201
202std::pair<int,int> LineMeshFactory::determineXElemSizeAndStart(int xBlock,unsigned int size,unsigned int /* rank */) const
203{
204 std::size_t xProcLoc = procTuple_[0];
205 unsigned int minElements = nXElems_/size;
206 unsigned int extra = nXElems_ - minElements*size;
207
208 TEUCHOS_ASSERT(minElements>0);
209
210 // first "extra" elements get an extra column of elements
211 // this determines the starting X index and number of elements
212 int nume=0, start=0;
213 if(xProcLoc<extra) {
214 nume = minElements+1;
215 start = xProcLoc*(minElements+1);
216 }
217 else {
218 nume = minElements;
219 start = extra*(minElements+1)+(xProcLoc-extra)*minElements;
220 }
221
222 return std::make_pair(start+nXElems_*xBlock,nume);
223}
224
226{
227 mesh.beginModification();
228 const stk::mesh::EntityRank sideRank = mesh.getSideRank();
229
230 std::size_t totalXElems = nXElems_*xBlocks_;
231
232 // get all part vectors
233 stk::mesh::Part * left = mesh.getSideset("left");
234 stk::mesh::Part * right = mesh.getSideset("right");
235
236 std::vector<stk::mesh::Entity> localElmts;
237 mesh.getMyElements(localElmts);
238
239 // loop over elements adding edges to sidesets
240 std::vector<stk::mesh::Entity>::const_iterator itr;
241 for(itr=localElmts.begin();itr!=localElmts.end();++itr) {
242 stk::mesh::Entity element = (*itr);
243 stk::mesh::EntityId gid = mesh.elementGlobalId(element);
244
245 std::size_t nx = gid-1;
246
247 // vertical boundaries
249
250 if(nx+1==totalXElems) {
251 stk::mesh::Entity edge = mesh.findConnectivityById(element, sideRank, 1);
252
253 // on the right
254 if(mesh.entityOwnerRank(edge)==machRank_)
255 mesh.addEntityToSideset(edge,right);
256 }
257
258 if(nx==0) {
259 stk::mesh::Entity edge = mesh.findConnectivityById(element, sideRank, 0);
260
261 // on the left
262 if(mesh.entityOwnerRank(edge)==machRank_)
263 mesh.addEntityToSideset(edge,left);
264 }
265 }
266
267 mesh.endModification();
268}
269
271Teuchos::Tuple<std::size_t,2> LineMeshFactory::procRankToProcTuple(std::size_t procRank) const
272{
273 std::size_t i=0,j=0;
274
275 j = procRank/machSize_;
276 procRank = procRank % machSize_;
277 i = procRank;
278
279 return Teuchos::tuple(i,j);
280}
281
282} // end panzer_stk
Teuchos::RCP< STK_Interface > buildMesh(stk::ParallelMachine parallelMach) const
Build the mesh object.
void buildMetaData(stk::ParallelMachine parallelMach, STK_Interface &mesh) const
std::pair< int, int > determineXElemSizeAndStart(int xBlock, unsigned int size, unsigned int rank) const
void buildElements(stk::ParallelMachine parallelMach, STK_Interface &mesh) const
void buildBlock(stk::ParallelMachine machRank, int xBlock, STK_Interface &mesh) const
void addSideSets(STK_Interface &mesh) const
virtual Teuchos::RCP< STK_Interface > buildUncommitedMesh(stk::ParallelMachine parallelMach) const
Teuchos::Tuple< std::size_t, 2 > procRankToProcTuple(std::size_t procRank) const
what is the 2D tuple describe this processor distribution
virtual void completeMeshConstruction(STK_Interface &mesh, stk::ParallelMachine parallelMach) const
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const
From ParameterListAcceptor.
Teuchos::Tuple< std::size_t, 2 > procTuple_
void setParameterList(const Teuchos::RCP< Teuchos::ParameterList > &paramList)
From ParameterListAcceptor.
stk::mesh::Entity findConnectivityById(stk::mesh::Entity src, stk::mesh::EntityRank tgt_rank, unsigned rel_id) const
void initialize(stk::ParallelMachine parallelMach, bool setupIO=true, const bool buildRefinementSupport=false)
stk::mesh::Part * getElementBlockPart(const std::string &name) const
get the block part
stk::mesh::EntityId elementGlobalId(std::size_t lid) const
bool isInitialized() const
Has initialize been called on this mesh object?
void buildSubcells()
force the mesh to build subcells: edges and faces
void addElement(const Teuchos::RCP< ElementDescriptor > &ed, stk::mesh::Part *block)
void addNode(stk::mesh::EntityId gid, const std::vector< double > &coord)
void addSideset(const std::string &name, const CellTopologyData *ctData)
void endModification(const bool find_and_set_shared_nodes_in_stk=true)
unsigned entityOwnerRank(stk::mesh::Entity entity) const
void addEntityToSideset(stk::mesh::Entity entity, stk::mesh::Part *sideset)
stk::mesh::EntityRank getSideRank() const
void getMyElements(std::vector< stk::mesh::Entity > &elements) const
void addElementBlock(const std::string &name, const CellTopologyData *ctData)
stk::mesh::Part * getSideset(const std::string &name) const
void rebalance(STK_Interface &mesh) const
static void parsePeriodicBCList(const Teuchos::RCP< Teuchos::ParameterList > &pl, std::vector< Teuchos::RCP< const PeriodicBC_MatcherBase > > &periodicBC, bool &useBBoxSearch)
double getMeshCoord(const int nx, const double deltaX, const double x0) const
std::vector< Teuchos::RCP< const PeriodicBC_MatcherBase > > periodicBCVec_