Panzer Version of the Day
Loading...
Searching...
No Matches
Panzer_STK_Quad8ToQuad4MeshFactory.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
43#include "PanzerAdaptersSTK_config.hpp"
47#include "Teuchos_TimeMonitor.hpp"
48#include "Teuchos_StandardParameterEntryValidators.hpp" // for plist validation
49#include <stk_mesh/base/FieldBase.hpp>
50#include <stk_mesh/base/DumpMeshInfo.hpp>
51
52// #define ENABLE_UNIFORM
53
54using Teuchos::RCP;
55using Teuchos::rcp;
56
57namespace panzer_stk {
58
59Quad8ToQuad4MeshFactory::Quad8ToQuad4MeshFactory(const std::string& quad8MeshFileName,
60 stk::ParallelMachine mpi_comm,
61 const bool print_debug)
62 : createEdgeBlocks_(false),
63 print_debug_(print_debug)
64{
65 panzer_stk::STK_ExodusReaderFactory factory;
66 RCP<Teuchos::ParameterList> pl = rcp(new Teuchos::ParameterList);
67 pl->set("File Name",quad8MeshFileName);
68 factory.setParameterList(pl);
69 quad8Mesh_ = factory.buildMesh(mpi_comm);
70
72}
73
74Quad8ToQuad4MeshFactory::Quad8ToQuad4MeshFactory(const Teuchos::RCP<panzer_stk::STK_Interface>& quad8Mesh,
75 const bool print_debug)
76 : quad8Mesh_(quad8Mesh),
77 createEdgeBlocks_(false),
78 print_debug_(print_debug)
79{
81}
82
84Teuchos::RCP<STK_Interface> Quad8ToQuad4MeshFactory::buildMesh(stk::ParallelMachine parallelMach) const
85{
86 PANZER_FUNC_TIME_MONITOR("panzer::Quad8ToQuad4MeshFactory::buildMesh()");
87
88 // Make sure the Quad8 and Quad4 Comms match
89 {
90 int result = MPI_UNEQUAL;
91 MPI_Comm_compare(parallelMach, quad8Mesh_->getBulkData()->parallel(), &result);
92 TEUCHOS_ASSERT(result != MPI_UNEQUAL);
93 }
94
95 // build all meta data
96 RCP<STK_Interface> mesh = this->buildUncommitedMesh(parallelMach);
97
98 // commit meta data
99 mesh->initialize(parallelMach);
100
101 // build bulk data
102 this->completeMeshConstruction(*mesh,parallelMach);
103
104 return mesh;
105}
106
107Teuchos::RCP<STK_Interface> Quad8ToQuad4MeshFactory::buildUncommitedMesh(stk::ParallelMachine parallelMach) const
108{
109 PANZER_FUNC_TIME_MONITOR("panzer::Quad8ToQuad4MeshFactory::buildUncomittedMesh()");
110
111 RCP<STK_Interface> mesh = rcp(new STK_Interface(2));
112
113 machRank_ = stk::parallel_machine_rank(parallelMach);
114 machSize_ = stk::parallel_machine_size(parallelMach);
115
116 // build meta information: blocks and side set setups
117 this->buildMetaData(parallelMach,*mesh);
118
119 mesh->addPeriodicBCs(periodicBCVec_);
120 mesh->setBoundingBoxSearchFlag(useBBoxSearch_);
121
122 return mesh;
123}
124
125void Quad8ToQuad4MeshFactory::completeMeshConstruction(STK_Interface & mesh,stk::ParallelMachine parallelMach) const
126{
127 PANZER_FUNC_TIME_MONITOR("panzer::Quad8ToQuad4MeshFactory::completeMeshConstruction()");
128
129 if(not mesh.isInitialized())
130 mesh.initialize(parallelMach);
131
132 // add node and element information
133 this->buildElements(parallelMach,mesh);
134
135 // finish up the edges
136#ifndef ENABLE_UNIFORM
137 mesh.buildSubcells();
138#endif
141 mesh.buildLocalEdgeIDs();
142 }
143
144 // now that edges are built, sidsets can be added
145#ifndef ENABLE_UNIFORM
146 this->addSideSets(mesh);
147#endif
148
149 // add nodesets
150 this->addNodeSets(mesh);
151
153 this->addEdgeBlocks(mesh);
154 }
155
156 // Copy field data
157 this->copyCellFieldData(mesh);
158
159 // calls Stk_MeshFactory::rebalance
160 // this->rebalance(mesh);
161}
162
164void Quad8ToQuad4MeshFactory::setParameterList(const Teuchos::RCP<Teuchos::ParameterList> & paramList)
165{
166 paramList->validateParametersAndSetDefaults(*getValidParameters(),0);
167
168 setMyParamList(paramList);
169
170 // offsetGIDs_ = (paramList->get<std::string>("Offset mesh GIDs above 32-bit int limit") == "ON") ? true : false;
171
172 createEdgeBlocks_ = paramList->get<bool>("Create Edge Blocks");
173
174 // read in periodic boundary conditions
175 parsePeriodicBCList(Teuchos::rcpFromRef(paramList->sublist("Periodic BCs")),periodicBCVec_,useBBoxSearch_);
176}
177
179Teuchos::RCP<const Teuchos::ParameterList> Quad8ToQuad4MeshFactory::getValidParameters() const
180{
181 static RCP<Teuchos::ParameterList> defaultParams;
182
183 // fill with default values
184 if(defaultParams == Teuchos::null) {
185 defaultParams = rcp(new Teuchos::ParameterList);
186
187 Teuchos::setStringToIntegralParameter<int>(
188 "Offset mesh GIDs above 32-bit int limit",
189 "OFF",
190 "If 64-bit GIDs are supported, the mesh element and node global indices will start at a value greater than 32-bit limit.",
191 Teuchos::tuple<std::string>("OFF", "ON"),
192 defaultParams.get());
193
194 // default to false for backward compatibility
195 defaultParams->set<bool>("Create Edge Blocks",false,"Create edge blocks in the mesh");
196
197 Teuchos::ParameterList & bcs = defaultParams->sublist("Periodic BCs");
198 bcs.set<int>("Count",0); // no default periodic boundary conditions
199 }
200
201 return defaultParams;
202}
203
204void Quad8ToQuad4MeshFactory::buildMetaData(stk::ParallelMachine /* parallelMach */, STK_Interface & mesh) const
205{
206 if (print_debug_) {
207 std::cout << "\n\n**** DEBUG: begin printing source Quad8 exodus file metadata ****\n";
208 stk::mesh::impl::dump_all_meta_info(*(quad8Mesh_->getMetaData()), std::cout);
209 std::cout << "\n\n**** DEBUG: end printing source Quad8 exodus file metadata ****\n";
210 }
211
212 // Required topologies
213 using QuadTopo = shards::Quadrilateral<4>;
214 const CellTopologyData * ctd = shards::getCellTopologyData<QuadTopo>();
215 const CellTopologyData * side_ctd = shards::CellTopology(ctd).getBaseCellTopologyData(1,0);
216
217 // Add in element blocks
218 std::vector<std::string> element_block_names;
219 quad8Mesh_->getElementBlockNames(element_block_names);
220 {
221 for (const auto& n : element_block_names)
222 mesh.addElementBlock(n,ctd);
223 }
224
225 // Add in sidesets
226 {
227 std::vector<std::string> sideset_names;
228 quad8Mesh_->getSidesetNames(sideset_names);
229 for (const auto& n : sideset_names)
230 mesh.addSideset(n,side_ctd);
231 }
232
233 // Add in nodesets
234 {
235 std::vector<std::string> nodeset_names;
236 quad8Mesh_->getNodesetNames(nodeset_names);
237 for (const auto& n : nodeset_names)
238 mesh.addNodeset(n);
239 }
240
242 const CellTopologyData * edge_ctd = shards::CellTopology(ctd).getBaseCellTopologyData(1,0);
243 std::vector<std::string> element_block_names;
244 quad8Mesh_->getElementBlockNames(element_block_names);
245 for (const auto& block_name : element_block_names)
246 mesh.addEdgeBlock(block_name,edgeBlockName_,edge_ctd);
247 }
248
249 // Add in nodal fields
250 {
251 const auto& fields = quad8Mesh_->getMetaData()->get_fields(mesh.getNodeRank());
252 for (const auto& f : fields) {
253 if (print_debug_)
254 std::cout << "Field=" << f->name() << ", rank=" << f->entity_rank() << std::endl;
255
256 // Cull out the coordinate fields. That is a default field in
257 // stk that is automatically created by stk.
258 if (f->name() != "coordinates") {
259 for (const auto& n : element_block_names)
260 mesh.addSolutionField(f->name(),n);
261 }
262 }
263 }
264
265 // Add in element fields
266 {
267 const auto& fields = quad8Mesh_->getMetaData()->get_fields(mesh.getElementRank());
268 for (const auto& f : fields) {
269 if (print_debug_)
270 std::cout << "Add Cell Field=" << f->name() << ", rank=" << f->entity_rank() << std::endl;
271
272 for (const auto& n : element_block_names)
273 mesh.addCellField(f->name(),n);
274 }
275 }
276
277 // NOTE: skipping edge and face fields for now. Can't error out since sidesets count as edge fields.
278 // TEUCHOS_TEST_FOR_EXCEPTION(quad8Mesh_->getMetaData()->get_fields(mesh.getEdgeRank()).size() != 0,std::runtime_error,
279 // "ERROR: the Quad8 mesh contains edge fields\""
280 // << quad8Mesh_->getMetaData()->get_fields(mesh.getEdgeRank())[0]->name()
281 // << "\". Edge fields are not supported in Quad8ToQuad4!");
282
283 if (print_debug_) {
284 std::cout << "\n\n**** DEBUG: begin printing source Quad4 exodus file metadata ****\n";
285 stk::mesh::impl::dump_all_meta_info(*(mesh.getMetaData()), std::cout);
286 std::cout << "\n\n**** DEBUG: end printing source Quad4 exodus file metadata ****\n";
287 }
288}
289
290void Quad8ToQuad4MeshFactory::buildElements(stk::ParallelMachine parallelMach,STK_Interface & mesh) const
291{
292 mesh.beginModification();
293
294 auto metadata = mesh.getMetaData();
295 auto bulkdata = mesh.getBulkData();
296
297 // Loop over element blocks
298 std::vector<std::string> block_names;
299 quad8Mesh_->getElementBlockNames(block_names);
300 for (const auto& block_name : block_names) {
301
302 // Get the elements on this process
303 std::vector<stk::mesh::Entity> elements;
304 quad8Mesh_->getMyElements(block_name,elements);
305
306 if (print_debug_) {
307 std::cout << "*************************************************" << std::endl;
308 std::cout << "block_name=" << block_name << ", num_my_elements=" << elements.size() << std::endl;
309 std::cout << "*************************************************" << std::endl;
310 }
311
312 for (const auto& element : elements) {
313
314 const auto element_gid = quad8Mesh_->getBulkData()->identifier(element);
315
316 if (print_debug_) {
317 std::cout << "rank=" << machRank_
318 << ", block=" << block_name
319 << ", element entity_id=" << element
320 << ", gid=" << element_gid << std::endl;
321 }
322
323 // Register nodes with the mesh
324 std::vector<stk::mesh::EntityId> nodes(4);
325 for (int i=0; i < 4; ++i) {
326 // NOTE: this assumes that the Quad4 topology is nested in
327 // the Quad8 as an extended topology - i.e. the first four
328 // nodes of the quad8 topology are the corresponding quad4
329 // nodes. Shards topologies enforce this throught the concept
330 // of extended topologies.
331 stk::mesh::Entity node_entity = quad8Mesh_->findConnectivityById(element,mesh.getNodeRank(),i);
332 TEUCHOS_ASSERT(node_entity.is_local_offset_valid());
333 const auto node_gid = quad8Mesh_->getBulkData()->identifier(node_entity);
334 const double* node_coords = quad8Mesh_->getNodeCoordinates(node_entity);
335 std::vector<double> coords_vec(2);
336 coords_vec[0] = node_coords[0];
337 coords_vec[1] = node_coords[1];
338 mesh.addNode(node_gid,coords_vec);
339 nodes[i]=node_gid;
340 if (print_debug_) {
341 std::cout << "elem gid=" << quad8Mesh_->getBulkData()->identifier(element)
342 << ", node_gid=" << node_gid << ", (" << coords_vec[0] << "," << coords_vec[1] << ")\n";
343 }
344 }
345
346 // Register element with the element block
347 auto element_descriptor = panzer_stk::buildElementDescriptor(element_gid,nodes);
348 auto element_block_part = mesh.getMetaData()->get_part(block_name);
349 mesh.addElement(element_descriptor,element_block_part);
350 }
351 }
352 mesh.endModification();
353}
354
356{
357 mesh.beginModification();
358
359 // Loop over sidesets
360 std::vector<std::string> sideset_names;
361 quad8Mesh_->getSidesetNames(sideset_names);
362 for (const auto& sideset_name : sideset_names) {
363
364 stk::mesh::Part* sideset_part = mesh.getSideset(sideset_name);
365
366 std::vector<stk::mesh::Entity> q8_sides;
367 quad8Mesh_->getMySides(sideset_name,q8_sides);
368
369 // Loop over edges
370 for (const auto q8_edge : q8_sides) {
371 // The edge numbering scheme uses the element/node gids, so it
372 // should be consistent between the quad8 and quad4 meshes
373 // since we used the same gids. We use this fact to populate
374 // the quad4 sidesets.
375 stk::mesh::EntityId edge_gid = quad8Mesh_->getBulkData()->identifier(q8_edge);
376 stk::mesh::Entity q4_edge = mesh.getBulkData()->get_entity(mesh.getEdgeRank(),edge_gid);
377 mesh.addEntityToSideset(q4_edge,sideset_part);
378 }
379 }
380
381 mesh.endModification();
382}
383
386
388{
389 mesh.beginModification();
390
391 Teuchos::RCP<stk::mesh::BulkData> bulkData = mesh.getBulkData();
392 Teuchos::RCP<stk::mesh::MetaData> metaData = mesh.getMetaData();
393
394 stk::mesh::Part * edge_block = mesh.getEdgeBlock(edgeBlockName_);
395
396 stk::mesh::Selector owned_block = metaData->locally_owned_part();
397
398 std::vector<stk::mesh::Entity> edges;
399 bulkData->get_entities(mesh.getEdgeRank(), owned_block, edges);
400 mesh.addEntitiesToEdgeBlock(edges, edge_block);
401
402 mesh.endModification();
403}
404
406{
407 // Vector of pointers to field data
408 const auto& fields = q4_mesh.getMetaData()->get_fields(q4_mesh.getElementRank());
409
410 // loop over fields and add the data to the new mesh.
411 for (const auto& field : fields) {
412
413 if (print_debug_)
414 std::cout << "Adding field values for \"" << *field << "\" to the Quad4 mesh!\n";
415
416 auto field_name = field->name();
417
418 // Divide into scalar and vector fields, ignoring any other types
419 // for now.
420 if (field->type_is<double>() &&
421 field_name != "coordinates" &&
422 field_name != "PROC_ID" &&
423 field_name != "LOAD_BAL") {
424
425 // Loop over element blocks
426 std::vector<std::string> block_names;
427 quad8Mesh_->getElementBlockNames(block_names);
428 for (const auto& block : block_names) {
429
430 auto* q4_field = q4_mesh.getCellField(field_name,block);
431 TEUCHOS_ASSERT(q4_field != nullptr);
432 // The q8 mesh doesn't have the field names set up, so a query
433 // fails. Go to stk directly in this case.
434 auto* q8_field = quad8Mesh_->getMetaData()->get_field(quad8Mesh_->getElementRank(),field_name);
435#ifdef PANZER_DEBUG
436 TEUCHOS_ASSERT(q8_field != nullptr);
437#endif
438
439 // Get the elements for this block.
440 std::vector<stk::mesh::Entity> q4_elements;
441 q4_mesh.getMyElements(block,q4_elements);
442 std::vector<stk::mesh::Entity> q8_elements;
443 quad8Mesh_->getMyElements(block,q8_elements);
444 TEUCHOS_ASSERT(q4_elements.size() == q8_elements.size());
445
446 for (size_t i=0; i < q4_elements.size(); ++i) {
447#ifdef PANZER_DEBUG
448 TEUCHOS_ASSERT(q4_mesh.getBulkData()->identifier(q4_elements[i]) ==
449 quad8Mesh_->getBulkData()->identifier(q8_elements[i]));
450#endif
451
452 double* const q4_val = static_cast<double*>(stk::mesh::field_data(*q4_field,q4_elements[i]));
453 const double* const q8_val = static_cast<double*>(stk::mesh::field_data(*q8_field,q8_elements[i]));
454 *q4_val = *q8_val;
455
456 if (print_debug_) {
457 std::cout << "field=" << field_name << ", block=" << block
458 << ", q4e(" << q4_mesh.getBulkData()->identifier(q4_elements[i]) << ") = " << *q4_val
459 << ", q8e(" << quad8Mesh_->getBulkData()->identifier(q8_elements[i]) << ") = " << *q8_val
460 << std::endl;
461 }
462
463 }
464 }
465 }
466 }
467}
468
469} // end panzer_stk
virtual void completeMeshConstruction(STK_Interface &mesh, stk::ParallelMachine parallelMach) const
Quad8ToQuad4MeshFactory(const std::string &quad8MeshFileName, stk::ParallelMachine mpi_comm=MPI_COMM_WORLD, const bool print_debug=false)
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const
Derived from ParameterListAcceptor.
Teuchos::RCP< STK_Interface > buildMesh(stk::ParallelMachine parallelMach) const
Build the mesh object.
void buildElements(stk::ParallelMachine parallelMach, STK_Interface &mesh) const
void buildMetaData(stk::ParallelMachine parallelMach, STK_Interface &mesh) const
Teuchos::RCP< panzer_stk::STK_Interface > quad8Mesh_
virtual Teuchos::RCP< STK_Interface > buildUncommitedMesh(stk::ParallelMachine parallelMach) const
void setParameterList(const Teuchos::RCP< Teuchos::ParameterList > &paramList)
Derived from ParameterListAcceptor.
void initialize(stk::ParallelMachine parallelMach, bool setupIO=true, const bool buildRefinementSupport=false)
static const std::string edgeBlockString
stk::mesh::Part * getEdgeBlock(const std::string &name) const
get the block part
bool isInitialized() const
Has initialize been called on this mesh object?
stk::mesh::EntityRank getNodeRank() const
void addEdgeBlock(const std::string &elemBlockName, const std::string &edgeBlockName, const stk::topology &topology)
void buildSubcells()
force the mesh to build subcells: edges and faces
stk::mesh::EntityRank getElementRank() const
void addElement(const Teuchos::RCP< ElementDescriptor > &ed, stk::mesh::Part *block)
void addNode(stk::mesh::EntityId gid, const std::vector< double > &coord)
void addCellField(const std::string &fieldName, const std::string &blockId)
void addNodeset(const std::string &name)
void addSideset(const std::string &name, const CellTopologyData *ctData)
Teuchos::RCP< stk::mesh::MetaData > getMetaData() const
stk::mesh::Field< double > * getCellField(const std::string &fieldName, const std::string &blockId) const
void addEntitiesToEdgeBlock(std::vector< stk::mesh::Entity > entities, stk::mesh::Part *edgeblock)
void addSolutionField(const std::string &fieldName, const std::string &blockId)
void addEntityToSideset(stk::mesh::Entity entity, stk::mesh::Part *sideset)
void getMyElements(std::vector< stk::mesh::Entity > &elements) const
void addElementBlock(const std::string &name, const CellTopologyData *ctData)
Teuchos::RCP< stk::mesh::BulkData > getBulkData() const
stk::mesh::EntityRank getEdgeRank() const
stk::mesh::Part * getSideset(const std::string &name) const
static void parsePeriodicBCList(const Teuchos::RCP< Teuchos::ParameterList > &pl, std::vector< Teuchos::RCP< const PeriodicBC_MatcherBase > > &periodicBC, bool &useBBoxSearch)
std::vector< Teuchos::RCP< const PeriodicBC_MatcherBase > > periodicBCVec_
Teuchos::RCP< ElementDescriptor > buildElementDescriptor(stk::mesh::EntityId elmtId, std::vector< stk::mesh::EntityId > &nodes)