Panzer Version of the Day
Loading...
Searching...
No Matches
Panzer_STK_QuadraticToLinearMeshFactory.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
52using Teuchos::RCP;
53using Teuchos::rcp;
54
55namespace panzer_stk {
56
58 stk::ParallelMachine mpi_comm,
59 const bool print_debug)
60 : createEdgeBlocks_(false),
61 print_debug_(print_debug)
62{
63 panzer_stk::STK_ExodusReaderFactory factory;
64 RCP<Teuchos::ParameterList> pl = rcp(new Teuchos::ParameterList);
65 pl->set("File Name",quadMeshFileName);
66 factory.setParameterList(pl);
67 quadMesh_ = factory.buildMesh(mpi_comm);
68
69 // TODO for currently supported input topologies, this should always be true
70 // but may need to be generalized in the future
72
73 // check the conversion is supported
74 // and get output topology
75 this->getOutputTopology();
76}
77
78QuadraticToLinearMeshFactory::QuadraticToLinearMeshFactory(const Teuchos::RCP<panzer_stk::STK_Interface>& quadMesh,
79 const bool print_debug)
80 : quadMesh_(quadMesh),
81 createEdgeBlocks_(false),
82 print_debug_(print_debug)
83{
84 // TODO for currently supported input topologies, this should always be true
85 // but may need to be generalized in the future
87
88 // check the conversion is supported
89 // and get output topology
90 this->getOutputTopology();
91}
92
94Teuchos::RCP<STK_Interface> QuadraticToLinearMeshFactory::buildMesh(stk::ParallelMachine parallelMach) const
95{
96 PANZER_FUNC_TIME_MONITOR("panzer::QuadraticToLinearMeshFactory::buildMesh()");
97
98 // Make sure the Comms match between the meshes
99 {
100 int result = MPI_UNEQUAL;
101 MPI_Comm_compare(parallelMach, quadMesh_->getBulkData()->parallel(), &result);
102 TEUCHOS_ASSERT(result != MPI_UNEQUAL);
103 }
104
105 // build all meta data
106 RCP<STK_Interface> mesh = this->buildUncommitedMesh(parallelMach);
107
108 // commit meta data
109 mesh->initialize(parallelMach);
110
111 // build bulk data
112 this->completeMeshConstruction(*mesh,parallelMach);
113
114 return mesh;
115}
116
118{
119 bool errFlag = false;
120
121 std::vector<std::string> eblock_names;
122 quadMesh_->getElementBlockNames(eblock_names);
123
124 // check that we have a supported topology
125 auto inputTopo = quadMesh_->getCellTopology(eblock_names[0]);
126 if (std::find(supportedInputTopos_.begin(),
127 supportedInputTopos_.end(),*inputTopo) == supportedInputTopos_.end()) errFlag = true;
128 TEUCHOS_TEST_FOR_EXCEPTION(errFlag,std::logic_error,
129 "ERROR :: Input topology " << inputTopo->getName() << " currently unsupported by QuadraticToLinearMeshFactory!");
130
131 // check that the topology is the same over blocks
132 // not sure this is 100% foolproof
133 for (auto & eblock : eblock_names) {
134 auto cellTopo = quadMesh_->getCellTopology(eblock);
135 if (*cellTopo != *inputTopo) errFlag = true;
136 }
137
138 TEUCHOS_TEST_FOR_EXCEPTION(errFlag, std::logic_error,
139 "ERROR :: The mesh has different topologies on different blocks!");
140
141 outputTopoData_ = outputTopoMap_[inputTopo->getName()];
142
143 nDim_ = outputTopoData_->dimension;
144 nNodes_ = outputTopoData_->node_count;
145
146 return;
147}
148
149Teuchos::RCP<STK_Interface> QuadraticToLinearMeshFactory::buildUncommitedMesh(stk::ParallelMachine parallelMach) const
150{
151 PANZER_FUNC_TIME_MONITOR("panzer::QuadraticToLinearMeshFactory::buildUncomittedMesh()");
152
153 RCP<STK_Interface> mesh = rcp(new STK_Interface(nDim_));
154
155 machRank_ = stk::parallel_machine_rank(parallelMach);
156 machSize_ = stk::parallel_machine_size(parallelMach);
157
158 // build meta information: blocks and side set setups
159 this->buildMetaData(parallelMach,*mesh);
160
161 mesh->addPeriodicBCs(periodicBCVec_);
162 mesh->setBoundingBoxSearchFlag(useBBoxSearch_);
163
164 return mesh;
165}
166
167void QuadraticToLinearMeshFactory::completeMeshConstruction(STK_Interface & mesh,stk::ParallelMachine parallelMach) const
168{
169 PANZER_FUNC_TIME_MONITOR("panzer::QuadraticToLinearMeshFactory::completeMeshConstruction()");
170
171 if(not mesh.isInitialized())
172 mesh.initialize(parallelMach);
173
174 // add node and element information
175 this->buildElements(parallelMach,mesh);
176
177 // finish up the edges
178#ifndef ENABLE_UNIFORM
179 mesh.buildSubcells();
180#endif
183 mesh.buildLocalEdgeIDs();
184 }
185
186 // now that edges are built, sidesets can be added
187#ifndef ENABLE_UNIFORM
188 this->addSideSets(mesh);
189#endif
190
191 // add nodesets
192 this->addNodeSets(mesh);
193
194 // TODO this functionality may be untested
196 this->addEdgeBlocks(mesh);
197 }
198
199 // Copy field data
200 this->copyCellFieldData(mesh);
201
202 // calls Stk_MeshFactory::rebalance
203 // this->rebalance(mesh);
204}
205
207void QuadraticToLinearMeshFactory::setParameterList(const Teuchos::RCP<Teuchos::ParameterList> & paramList)
208{
209 paramList->validateParametersAndSetDefaults(*getValidParameters(),0);
210
211 setMyParamList(paramList);
212
213 // offsetGIDs_ = (paramList->get<std::string>("Offset mesh GIDs above 32-bit int limit") == "ON") ? true : false;
214
215 createEdgeBlocks_ = paramList->get<bool>("Create Edge Blocks");
216
217 // read in periodic boundary conditions
218 parsePeriodicBCList(Teuchos::rcpFromRef(paramList->sublist("Periodic BCs")),periodicBCVec_,useBBoxSearch_);
219}
220
222Teuchos::RCP<const Teuchos::ParameterList> QuadraticToLinearMeshFactory::getValidParameters() const
223{
224 static RCP<Teuchos::ParameterList> defaultParams;
225
226 // fill with default values
227 if(defaultParams == Teuchos::null) {
228 defaultParams = rcp(new Teuchos::ParameterList);
229
230 Teuchos::setStringToIntegralParameter<int>(
231 "Offset mesh GIDs above 32-bit int limit",
232 "OFF",
233 "If 64-bit GIDs are supported, the mesh element and node global indices will start at a value greater than 32-bit limit.",
234 Teuchos::tuple<std::string>("OFF", "ON"),
235 defaultParams.get());
236
237 // default to false for backward compatibility
238 defaultParams->set<bool>("Create Edge Blocks",false,"Create edge blocks in the mesh");
239
240 Teuchos::ParameterList & bcs = defaultParams->sublist("Periodic BCs");
241 bcs.set<int>("Count",0); // no default periodic boundary conditions
242 }
243
244 return defaultParams;
245}
246
247void QuadraticToLinearMeshFactory::buildMetaData(stk::ParallelMachine /* parallelMach */, STK_Interface & mesh) const
248{
249 if (print_debug_) {
250 std::cout << "\n\n**** DEBUG: begin printing source quad mesh exodus file metadata ****\n";
251 stk::mesh::impl::dump_all_meta_info(*(quadMesh_->getMetaData()), std::cout);
252 std::cout << "\n\n**** DEBUG: end printing source quad mesh exodus file metadata ****\n";
253 }
254
255 // Required topologies
256 auto ctd = outputTopoData_;
257 // will only work for 2D and 3D meshes
258 const CellTopologyData * side_ctd = shards::CellTopology(ctd).getBaseCellTopologyData(nDim_-1,0);
259
260 // Add in element blocks
261 std::vector<std::string> element_block_names;
262 quadMesh_->getElementBlockNames(element_block_names);
263 {
264 for (const auto& n : element_block_names)
265 mesh.addElementBlock(n,ctd);
266 }
267
268 // Add in sidesets
269 {
270 std::vector<std::string> sideset_names;
271 quadMesh_->getSidesetNames(sideset_names);
272 for (const auto& n : sideset_names)
273 mesh.addSideset(n,side_ctd);
274 }
275
276 // Add in nodesets
277 {
278 std::vector<std::string> nodeset_names;
279 quadMesh_->getNodesetNames(nodeset_names);
280 for (const auto& n : nodeset_names)
281 mesh.addNodeset(n);
282 }
283
285 const CellTopologyData * edge_ctd = shards::CellTopology(ctd).getBaseCellTopologyData(1,0);
286 std::vector<std::string> element_block_names;
287 quadMesh_->getElementBlockNames(element_block_names);
288 for (const auto& block_name : element_block_names)
289 mesh.addEdgeBlock(block_name,edgeBlockName_,edge_ctd);
290 }
291
292 // Add in nodal fields
293 {
294 const auto& fields = quadMesh_->getMetaData()->get_fields(mesh.getNodeRank());
295 for (const auto& f : fields) {
296 if (print_debug_)
297 std::cout << "Field=" << f->name() << ", rank=" << f->entity_rank() << std::endl;
298
299 // Cull out the coordinate fields. That is a default field in
300 // stk that is automatically created by stk.
301 if (f->name() != "coordinates") {
302 for (const auto& n : element_block_names)
303 mesh.addSolutionField(f->name(),n);
304 }
305 }
306 }
307
308 // Add in element fields
309 {
310 const auto& fields = quadMesh_->getMetaData()->get_fields(mesh.getElementRank());
311 for (const auto& f : fields) {
312 if (print_debug_)
313 std::cout << "Add Cell Field=" << f->name() << ", rank=" << f->entity_rank() << std::endl;
314
315 for (const auto& n : element_block_names)
316 mesh.addCellField(f->name(),n);
317 }
318 }
319
320 // NOTE: skipping edge and face fields for now. Can't error out since sidesets count as edge fields.
321 // TEUCHOS_TEST_FOR_EXCEPTION(quadMesh_->getMetaData()->get_fields(mesh.getEdgeRank()).size() != 0,std::runtime_error,
322 // "ERROR: the Quad8 mesh contains edge fields\""
323 // << quadMesh_->getMetaData()->get_fields(mesh.getEdgeRank())[0]->name()
324 // << "\". Edge fields are not supported in Quad8ToQuad4!");
325
326 if (print_debug_) {
327 std::cout << "\n\n**** DEBUG: begin printing source linear mesh exodus file metadata ****\n";
328 stk::mesh::impl::dump_all_meta_info(*(mesh.getMetaData()), std::cout);
329 std::cout << "\n\n**** DEBUG: end printing source linear mesh exodus file metadata ****\n";
330 }
331}
332
333void QuadraticToLinearMeshFactory::buildElements(stk::ParallelMachine parallelMach,STK_Interface & mesh) const
334{
335 mesh.beginModification();
336
337 auto metadata = mesh.getMetaData();
338 auto bulkdata = mesh.getBulkData();
339
340 // Loop over element blocks
341 std::vector<std::string> block_names;
342 quadMesh_->getElementBlockNames(block_names);
343 for (const auto& block_name : block_names) {
344
345 // Get the elements on this process
346 std::vector<stk::mesh::Entity> elements;
347 quadMesh_->getMyElements(block_name,elements);
348
349 if (print_debug_) {
350 std::cout << "*************************************************" << std::endl;
351 std::cout << "block_name=" << block_name << ", num_my_elements=" << elements.size() << std::endl;
352 std::cout << "*************************************************" << std::endl;
353 }
354
355 for (const auto& element : elements) {
356
357 const auto element_gid = quadMesh_->getBulkData()->identifier(element);
358
359 if (print_debug_) {
360 std::cout << "rank=" << machRank_
361 << ", block=" << block_name
362 << ", element entity_id=" << element
363 << ", gid=" << element_gid << std::endl;
364 }
365
366 // Register nodes with the mesh
367 std::vector<stk::mesh::EntityId> nodes(nNodes_);
368 for (size_t i=0; i < nNodes_; ++i) {
369 // NOTE: this assumes that the linear topology is nested in
370 // the quadratic topology as an extended topology - i.e. the first n
371 // nodes of the quadratic topology are the corresponding linear
372 // nodes. Shards topologies enforce this through the concept
373 // of extended topologies.
374 stk::mesh::Entity node_entity = quadMesh_->findConnectivityById(element,mesh.getNodeRank(),i);
375 TEUCHOS_ASSERT(node_entity.is_local_offset_valid());
376 const auto node_gid = quadMesh_->getBulkData()->identifier(node_entity);
377 const double* node_coords = quadMesh_->getNodeCoordinates(node_entity);
378 std::vector<double> coords_vec(nDim_);
379 for (size_t j=0; j < nDim_; ++j) coords_vec[j] = node_coords[j];
380 mesh.addNode(node_gid,coords_vec);
381 nodes[i]=node_gid;
382 if (print_debug_) {
383 if (nDim_==2) {
384 std::cout << "elem gid=" << quadMesh_->getBulkData()->identifier(element)
385 << ", node_gid=" << node_gid << ", ("
386 << coords_vec[0] << "," << coords_vec[1] << ")\n";
387 } else {
388 std::cout << "elem gid=" << quadMesh_->getBulkData()->identifier(element)
389 << ", node_gid=" << node_gid << ", ("
390 << coords_vec[0] << "," << coords_vec[1] << "," << coords_vec[2] << ")\n";
391 }
392 }
393 }
394
395 // Register element with the element block
396 auto element_descriptor = panzer_stk::buildElementDescriptor(element_gid,nodes);
397 auto element_block_part = mesh.getMetaData()->get_part(block_name);
398 mesh.addElement(element_descriptor,element_block_part);
399 }
400 }
401 mesh.endModification();
402}
403
405{
406 mesh.beginModification();
407
408 // Loop over sidesets
409 std::vector<std::string> sideset_names;
410 quadMesh_->getSidesetNames(sideset_names);
411 for (const auto& sideset_name : sideset_names) {
412
413 stk::mesh::Part* sideset_part = mesh.getSideset(sideset_name);
414
415 std::vector<stk::mesh::Entity> q_sides;
416 quadMesh_->getMySides(sideset_name,q_sides);
417
418 // Loop over edges
419 for (const auto q_ent : q_sides) {
420 // The edge numbering scheme uses the element/node gids, so it
421 // should be consistent between the quadratic and linear meshes
422 // since we used the same gids. We use this fact to populate
423 // the quadratic sidesets.
424 stk::mesh::EntityId ent_gid = quadMesh_->getBulkData()->identifier(q_ent);
425 stk::mesh::Entity lin_ent = mesh.getBulkData()->get_entity(mesh.getSideRank(),ent_gid);
426 mesh.addEntityToSideset(lin_ent,sideset_part);
427 }
428 }
429
430 mesh.endModification();
431}
432
435
437{
438 mesh.beginModification();
439
440 Teuchos::RCP<stk::mesh::BulkData> bulkData = mesh.getBulkData();
441 Teuchos::RCP<stk::mesh::MetaData> metaData = mesh.getMetaData();
442
443 stk::mesh::Part * edge_block = mesh.getEdgeBlock(edgeBlockName_);
444
445 stk::mesh::Selector owned_block = metaData->locally_owned_part();
446
447 std::vector<stk::mesh::Entity> edges;
448 bulkData->get_entities(mesh.getEdgeRank(), owned_block, edges);
449 mesh.addEntitiesToEdgeBlock(edges, edge_block);
450
451 mesh.endModification();
452}
453
455{
456 // Vector of pointers to field data
457 const auto& fields = lin_mesh.getMetaData()->get_fields(lin_mesh.getElementRank());
458
459 // loop over fields and add the data to the new mesh.
460 for (const auto& field : fields) {
461
462 if (print_debug_)
463 std::cout << "Adding field values for \"" << *field << "\" to the linear mesh!\n";
464
465 auto field_name = field->name();
466
467 // Divide into scalar and vector fields, ignoring any other types
468 // for now.
469 if (field->type_is<double>() &&
470 field_name != "coordinates" &&
471 field_name != "PROC_ID" &&
472 field_name != "LOAD_BAL") {
473
474 // Loop over element blocks
475 std::vector<std::string> block_names;
476 quadMesh_->getElementBlockNames(block_names);
477 for (const auto& block : block_names) {
478
479 auto* lin_field = lin_mesh.getCellField(field_name,block);
480 TEUCHOS_ASSERT(lin_field != nullptr);
481 // The q mesh doesn't have the field names set up, so a query
482 // fails. Go to stk directly in this case.
483 auto* q_field = quadMesh_->getMetaData()->get_field(quadMesh_->getElementRank(),field_name);
484#ifdef PANZER_DEBUG
485 TEUCHOS_ASSERT(q_field != nullptr);
486#endif
487
488 // Get the elements for this block.
489 std::vector<stk::mesh::Entity> lin_elements;
490 lin_mesh.getMyElements(block,lin_elements);
491 std::vector<stk::mesh::Entity> q_elements;
492 quadMesh_->getMyElements(block,q_elements);
493 TEUCHOS_ASSERT(lin_elements.size() == q_elements.size());
494
495 for (size_t i=0; i < lin_elements.size(); ++i) {
496#ifdef PANZER_DEBUG
497 TEUCHOS_ASSERT(lin_mesh.getBulkData()->identifier(lin_elements[i]) ==
498 quadMesh_->getBulkData()->identifier(q_elements[i]));
499#endif
500
501 double* const lin_val = static_cast<double*>(stk::mesh::field_data(*lin_field,lin_elements[i]));
502 const double* const q_val = static_cast<double*>(stk::mesh::field_data(*q_field,q_elements[i]));
503 *lin_val = *q_val;
504
505 if (print_debug_) {
506 std::cout << "field=" << field_name << ", block=" << block
507 << ", lin_e(" << lin_mesh.getBulkData()->identifier(lin_elements[i]) << ") = " << *lin_val
508 << ", q_e(" << quadMesh_->getBulkData()->identifier(q_elements[i]) << ") = " << *q_val
509 << std::endl;
510 }
511
512 }
513 }
514 }
515 }
516}
517
518} // end panzer_stk
virtual void completeMeshConstruction(STK_Interface &mesh, stk::ParallelMachine parallelMach) const
void buildMetaData(stk::ParallelMachine parallelMach, STK_Interface &mesh) const
std::map< const std::string, const CellTopologyData * > outputTopoMap_
QuadraticToLinearMeshFactory(const std::string &quadMeshFileName, stk::ParallelMachine mpi_comm=MPI_COMM_WORLD, const bool print_debug=false)
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const
Derived from ParameterListAcceptor.
void setParameterList(const Teuchos::RCP< Teuchos::ParameterList > &paramList)
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
virtual Teuchos::RCP< STK_Interface > buildUncommitedMesh(stk::ParallelMachine parallelMach) const
std::vector< shards::CellTopology > supportedInputTopos_
Nodes in one element of the linear basis.
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)
stk::mesh::EntityRank getSideRank() const
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)