11#include "PanzerAdaptersSTK_config.hpp"
15#include "Teuchos_TimeMonitor.hpp"
16#include "Teuchos_StandardParameterEntryValidators.hpp"
17#include <stk_mesh/base/FieldBase.hpp>
18#include <stk_mesh/base/DumpMeshInfo.hpp>
28 stk::ParallelMachine mpi_comm,
29 const bool print_debug)
30 : createEdgeBlocks_(false),
31 print_debug_(print_debug)
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);
43 const bool print_debug)
44 : quad8Mesh_(quad8Mesh),
45 createEdgeBlocks_(false),
46 print_debug_(print_debug)
54 PANZER_FUNC_TIME_MONITOR(
"panzer::Quad8ToQuad4MeshFactory::buildMesh()");
59 MPI_Comm_compare(parallelMach,
quad8Mesh_->getBulkData()->parallel(), &
result);
60 TEUCHOS_ASSERT(
result != MPI_UNEQUAL);
67 mesh->initialize(parallelMach);
77 PANZER_FUNC_TIME_MONITOR(
"panzer::Quad8ToQuad4MeshFactory::buildUncomittedMesh()");
81 machRank_ = stk::parallel_machine_rank(parallelMach);
82 machSize_ = stk::parallel_machine_size(parallelMach);
95 PANZER_FUNC_TIME_MONITOR(
"panzer::Quad8ToQuad4MeshFactory::completeMeshConstruction()");
104#ifndef ENABLE_UNIFORM
113#ifndef ENABLE_UNIFORM
136 setMyParamList(paramList);
149 static RCP<Teuchos::ParameterList> defaultParams;
152 if(defaultParams == Teuchos::null) {
153 defaultParams = rcp(
new Teuchos::ParameterList);
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"))));
160 defaultParams->set<
bool>(
"Create Edge Blocks",
false,
"Create edge blocks in the mesh");
162 Teuchos::ParameterList & bcs = defaultParams->sublist(
"Periodic BCs");
163 bcs.set<
int>(
"Count",0);
166 return defaultParams;
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";
178 using QuadTopo = shards::Quadrilateral<4>;
179 const CellTopologyData * ctd = shards::getCellTopologyData<QuadTopo>();
180 const CellTopologyData * side_ctd = shards::CellTopology(ctd).getBaseCellTopologyData(1,0);
183 std::vector<std::string> element_block_names;
184 quad8Mesh_->getElementBlockNames(element_block_names);
186 for (
const auto& n : element_block_names)
192 std::vector<std::string> sideset_names;
194 for (
const auto& n : sideset_names)
200 std::vector<std::string> nodeset_names;
202 for (
const auto& n : nodeset_names)
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)
217 for (
const auto& f : fields) {
219 std::cout <<
"Field=" << f->name() <<
", rank=" << f->entity_rank() << std::endl;
223 if (f->name() !=
"coordinates") {
224 for (
const auto& n : element_block_names)
233 for (
const auto& f : fields) {
235 std::cout <<
"Add Cell Field=" << f->name() <<
", rank=" << f->entity_rank() << std::endl;
237 for (
const auto& n : element_block_names)
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";
263 std::vector<std::string> block_names;
264 quad8Mesh_->getElementBlockNames(block_names);
265 for (
const auto& block_name : block_names) {
268 std::vector<stk::mesh::Entity> elements;
269 quad8Mesh_->getMyElements(block_name,elements);
272 std::cout <<
"*************************************************" << std::endl;
273 std::cout <<
"block_name=" << block_name <<
", num_my_elements=" << elements.size() << std::endl;
274 std::cout <<
"*************************************************" << std::endl;
277 for (
const auto& element : elements) {
279 const auto element_gid =
quad8Mesh_->getBulkData()->identifier(element);
283 <<
", block=" << block_name
284 <<
", element entity_id=" << element
285 <<
", gid=" << element_gid << std::endl;
289 std::vector<stk::mesh::EntityId> nodes(4);
290 for (
int i=0; i < 4; ++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);
306 std::cout <<
"elem gid=" <<
quad8Mesh_->getBulkData()->identifier(element)
307 <<
", node_gid=" << node_gid <<
", (" << coords_vec[0] <<
"," << coords_vec[1] <<
")\n";
313 auto element_block_part = mesh.
getMetaData()->get_part(block_name);
314 mesh.
addElement(element_descriptor,element_block_part);
325 std::vector<std::string> sideset_names;
327 for (
const auto& sideset_name : sideset_names) {
329 stk::mesh::Part* sideset_part = mesh.
getSideset(sideset_name);
331 std::vector<stk::mesh::Entity> q8_sides;
332 quad8Mesh_->getMySides(sideset_name,q8_sides);
335 for (
const auto q8_edge : q8_sides) {
340 stk::mesh::EntityId edge_gid =
quad8Mesh_->getBulkData()->identifier(q8_edge);
356 Teuchos::RCP<stk::mesh::BulkData> bulkData = mesh.
getBulkData();
357 Teuchos::RCP<stk::mesh::MetaData> metaData = mesh.
getMetaData();
361 stk::mesh::Selector owned_block = metaData->locally_owned_part();
363 std::vector<stk::mesh::Entity> edges;
364 bulkData->get_entities(mesh.
getEdgeRank(), owned_block, edges);
376 for (
const auto&
field : fields) {
379 std::cout <<
"Adding field values for \"" << *
field <<
"\" to the Quad4 mesh!\n";
381 auto field_name =
field->name();
385 if (
field->type_is<
double>() &&
386 field_name !=
"coordinates" &&
387 field_name !=
"PROC_ID" &&
388 field_name !=
"LOAD_BAL") {
391 std::vector<std::string> block_names;
392 quad8Mesh_->getElementBlockNames(block_names);
393 for (
const auto& block : block_names) {
395 auto* q4_field = q4_mesh.
getCellField(field_name,block);
396 TEUCHOS_ASSERT(q4_field !=
nullptr);
401 TEUCHOS_ASSERT(q8_field !=
nullptr);
405 std::vector<stk::mesh::Entity> q4_elements;
407 std::vector<stk::mesh::Entity> q8_elements;
409 TEUCHOS_ASSERT(q4_elements.size() == q8_elements.size());
411 for (
size_t i=0; i < q4_elements.size(); ++i) {
413 TEUCHOS_ASSERT(q4_mesh.
getBulkData()->identifier(q4_elements[i]) ==
414 quad8Mesh_->getBulkData()->identifier(q8_elements[i]));
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]));
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
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
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 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)