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