Panzer Version of the Day
Loading...
Searching...
No Matches
Panzer_STK_QuadraticToLinearMeshFactory.hpp
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#ifndef __Panzer_STK_QuadraticToLinearMeshFactory_hpp__
12#define __Panzer_STK_QuadraticToLinearMeshFactory_hpp__
13
16#include <vector>
17#include <string>
18
19namespace panzer_stk {
20
21class STK_Interface;
22
28public:
29
30 QuadraticToLinearMeshFactory(const std::string& quadMeshFileName,
31 stk::ParallelMachine mpi_comm = MPI_COMM_WORLD, // CHECK: ALLOW MPI_COMM_WORLD
32 const bool print_debug = false);
33
34 QuadraticToLinearMeshFactory(const Teuchos::RCP<panzer_stk::STK_Interface>& quadMesh,
35 const bool print_debug = false);
36
37 Teuchos::RCP<STK_Interface> buildMesh(stk::ParallelMachine parallelMach) const;
38 virtual Teuchos::RCP<STK_Interface> buildUncommitedMesh(stk::ParallelMachine parallelMach) const;
39 virtual void completeMeshConstruction(STK_Interface & mesh,stk::ParallelMachine parallelMach) const;
40
42 void setParameterList(const Teuchos::RCP<Teuchos::ParameterList> & paramList);
43
45 Teuchos::RCP<const Teuchos::ParameterList> getValidParameters() const;
46
47protected:
48
49 void buildMetaData(stk::ParallelMachine parallelMach,STK_Interface & mesh) const;
50 void buildElements(stk::ParallelMachine parallelMach,STK_Interface & mesh) const;
51
52 void addSideSets(STK_Interface & mesh) const;
53 void addNodeSets(STK_Interface & mesh) const;
54 void addEdgeBlocks(STK_Interface & mesh) const;
55 void copyCellFieldData(STK_Interface & mesh) const;
56
62 void getOutputTopology();
63
64 Teuchos::RCP<panzer_stk::STK_Interface> quadMesh_;
65
66 mutable unsigned int machRank_, machSize_;
67
69
72 mutable stk::mesh::EntityId offset_;
73
75
76 std::string edgeBlockName_;
77
78 const CellTopologyData * outputTopoData_;
79
80 unsigned int nDim_;
81 unsigned int nNodes_;
82
84 std::vector<shards::CellTopology> supportedInputTopos_ = {
85 shards::CellTopology(shards::getCellTopologyData<shards::Quadrilateral<8>>()),
86 shards::CellTopology(shards::getCellTopologyData<shards::Triangle<6>>()),
87 shards::CellTopology(shards::getCellTopologyData<shards::Tetrahedron<10>>()),
88 shards::CellTopology(shards::getCellTopologyData<shards::Hexahedron<20>>())
89 };
90
93 std::map<const std::string,const CellTopologyData *> outputTopoMap_ = {
94 {shards::getCellTopologyData<shards::Quadrilateral<8>>()->name,
95 shards::getCellTopologyData<shards::Quadrilateral<4>>()},
96 {shards::getCellTopologyData<shards::Triangle<6>>()->name,
97 shards::getCellTopologyData<shards::Triangle<3>>()},
98 {shards::getCellTopologyData<shards::Tetrahedron<10>>()->name,
99 shards::getCellTopologyData<shards::Tetrahedron<4>>()},
100 {shards::getCellTopologyData<shards::Hexahedron<20>>()->name,
101 shards::getCellTopologyData<shards::Hexahedron<8>>()}
102 };
103};
104
105}
106
107#endif
bool offsetGIDs_
If true, offset mesh GIDs to exercise 32-bit limits.
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_
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.