43#include "PanzerAdaptersSTK_config.hpp"
47#include "Teuchos_TimeMonitor.hpp"
48#include "Teuchos_StandardParameterEntryValidators.hpp"
49#include <stk_mesh/base/FieldBase.hpp>
50#include <stk_mesh/base/DumpMeshInfo.hpp>
60 stk::ParallelMachine mpi_comm,
61 const bool print_debug)
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);
75 const bool print_debug)
86 PANZER_FUNC_TIME_MONITOR(
"panzer::Quad8ToQuad4MeshFactory::buildMesh()");
90 int result = MPI_UNEQUAL;
91 MPI_Comm_compare(parallelMach,
quad8Mesh_->getBulkData()->parallel(), &result);
92 TEUCHOS_ASSERT(result != MPI_UNEQUAL);
99 mesh->initialize(parallelMach);
109 PANZER_FUNC_TIME_MONITOR(
"panzer::Quad8ToQuad4MeshFactory::buildUncomittedMesh()");
113 machRank_ = stk::parallel_machine_rank(parallelMach);
114 machSize_ = stk::parallel_machine_size(parallelMach);
127 PANZER_FUNC_TIME_MONITOR(
"panzer::Quad8ToQuad4MeshFactory::completeMeshConstruction()");
136#ifndef ENABLE_UNIFORM
145#ifndef ENABLE_UNIFORM
168 setMyParamList(paramList);
181 static RCP<Teuchos::ParameterList> defaultParams;
184 if(defaultParams == Teuchos::null) {
185 defaultParams = rcp(
new Teuchos::ParameterList);
187 Teuchos::setStringToIntegralParameter<int>(
188 "Offset mesh GIDs above 32-bit int limit",
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());
195 defaultParams->set<
bool>(
"Create Edge Blocks",
false,
"Create edge blocks in the mesh");
197 Teuchos::ParameterList & bcs = defaultParams->sublist(
"Periodic BCs");
198 bcs.set<
int>(
"Count",0);
201 return defaultParams;
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";
213 using QuadTopo = shards::Quadrilateral<4>;
214 const CellTopologyData * ctd = shards::getCellTopologyData<QuadTopo>();
215 const CellTopologyData * side_ctd = shards::CellTopology(ctd).getBaseCellTopologyData(1,0);
218 std::vector<std::string> element_block_names;
219 quad8Mesh_->getElementBlockNames(element_block_names);
221 for (
const auto& n : element_block_names)
227 std::vector<std::string> sideset_names;
229 for (
const auto& n : sideset_names)
235 std::vector<std::string> nodeset_names;
237 for (
const auto& n : nodeset_names)
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)
252 for (
const auto& f : fields) {
254 std::cout <<
"Field=" << f->name() <<
", rank=" << f->entity_rank() << std::endl;
258 if (f->name() !=
"coordinates") {
259 for (
const auto& n : element_block_names)
268 for (
const auto& f : fields) {
270 std::cout <<
"Add Cell Field=" << f->name() <<
", rank=" << f->entity_rank() << std::endl;
272 for (
const auto& n : element_block_names)
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";
298 std::vector<std::string> block_names;
299 quad8Mesh_->getElementBlockNames(block_names);
300 for (
const auto& block_name : block_names) {
303 std::vector<stk::mesh::Entity> elements;
304 quad8Mesh_->getMyElements(block_name,elements);
307 std::cout <<
"*************************************************" << std::endl;
308 std::cout <<
"block_name=" << block_name <<
", num_my_elements=" << elements.size() << std::endl;
309 std::cout <<
"*************************************************" << std::endl;
312 for (
const auto& element : elements) {
314 const auto element_gid =
quad8Mesh_->getBulkData()->identifier(element);
318 <<
", block=" << block_name
319 <<
", element entity_id=" << element
320 <<
", gid=" << element_gid << std::endl;
324 std::vector<stk::mesh::EntityId> nodes(4);
325 for (
int i=0; i < 4; ++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);
341 std::cout <<
"elem gid=" <<
quad8Mesh_->getBulkData()->identifier(element)
342 <<
", node_gid=" << node_gid <<
", (" << coords_vec[0] <<
"," << coords_vec[1] <<
")\n";
348 auto element_block_part = mesh.
getMetaData()->get_part(block_name);
349 mesh.
addElement(element_descriptor,element_block_part);
360 std::vector<std::string> sideset_names;
362 for (
const auto& sideset_name : sideset_names) {
364 stk::mesh::Part* sideset_part = mesh.
getSideset(sideset_name);
366 std::vector<stk::mesh::Entity> q8_sides;
367 quad8Mesh_->getMySides(sideset_name,q8_sides);
370 for (
const auto q8_edge : q8_sides) {
375 stk::mesh::EntityId edge_gid =
quad8Mesh_->getBulkData()->identifier(q8_edge);
391 Teuchos::RCP<stk::mesh::BulkData> bulkData = mesh.
getBulkData();
392 Teuchos::RCP<stk::mesh::MetaData> metaData = mesh.
getMetaData();
396 stk::mesh::Selector owned_block = metaData->locally_owned_part();
398 std::vector<stk::mesh::Entity> edges;
399 bulkData->get_entities(mesh.
getEdgeRank(), owned_block, edges);
411 for (
const auto& field : fields) {
414 std::cout <<
"Adding field values for \"" << *field <<
"\" to the Quad4 mesh!\n";
416 auto field_name = field->name();
420 if (field->type_is<
double>() &&
421 field_name !=
"coordinates" &&
422 field_name !=
"PROC_ID" &&
423 field_name !=
"LOAD_BAL") {
426 std::vector<std::string> block_names;
427 quad8Mesh_->getElementBlockNames(block_names);
428 for (
const auto& block : block_names) {
430 auto* q4_field = q4_mesh.
getCellField(field_name,block);
431 TEUCHOS_ASSERT(q4_field !=
nullptr);
436 TEUCHOS_ASSERT(q8_field !=
nullptr);
440 std::vector<stk::mesh::Entity> q4_elements;
442 std::vector<stk::mesh::Entity> q8_elements;
444 TEUCHOS_ASSERT(q4_elements.size() == q8_elements.size());
446 for (
size_t i=0; i < q4_elements.size(); ++i) {
448 TEUCHOS_ASSERT(q4_mesh.
getBulkData()->identifier(q4_elements[i]) ==
449 quad8Mesh_->getBulkData()->identifier(q8_elements[i]));
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]));
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
virtual void completeMeshConstruction(STK_Interface &mesh, stk::ParallelMachine parallelMach) const
void addSideSets(STK_Interface &mesh) const
void addNodeSets(STK_Interface &mesh) 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_
void addEdgeBlocks(STK_Interface &mesh) const
virtual Teuchos::RCP< STK_Interface > buildUncommitedMesh(stk::ParallelMachine parallelMach) const
void copyCellFieldData(STK_Interface &mesh) const
std::string edgeBlockName_
void setParameterList(const Teuchos::RCP< Teuchos::ParameterList > ¶mList)
Derived from ParameterListAcceptor.
void buildLocalElementIDs()
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)