Panzer  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
panzer_stk::QuadraticToLinearMeshFactory Class Reference

#include <Panzer_STK_QuadraticToLinearMeshFactory.hpp>

Inheritance diagram for panzer_stk::QuadraticToLinearMeshFactory:
Inheritance graph
[legend]

Public Member Functions

 QuadraticToLinearMeshFactory (const std::string &quadMeshFileName, stk::ParallelMachine mpi_comm=MPI_COMM_WORLD, const bool print_debug=false)
 
 QuadraticToLinearMeshFactory (const Teuchos::RCP< panzer_stk::STK_Interface > &quadMesh, const bool print_debug=false)
 
Teuchos::RCP< STK_InterfacebuildMesh (stk::ParallelMachine parallelMach) const
 Build the mesh object. More...
 
virtual Teuchos::RCP
< STK_Interface
buildUncommitedMesh (stk::ParallelMachine parallelMach) const
 
virtual void completeMeshConstruction (STK_Interface &mesh, stk::ParallelMachine parallelMach) const
 
void setParameterList (const Teuchos::RCP< Teuchos::ParameterList > &paramList)
 Derived from ParameterListAcceptor. More...
 
Teuchos::RCP< const
Teuchos::ParameterList
getValidParameters () const
 Derived from ParameterListAcceptor. More...
 
- Public Member Functions inherited from panzer_stk::STK_MeshFactory
 STK_MeshFactory ()
 
void enableRebalance (bool enable, const Teuchos::RCP< const Teuchos::ParameterList > &rebalanceList=Teuchos::null)
 
void rebalance (STK_Interface &mesh) const
 
double getMeshCoord (const int nx, const double deltaX, const double x0) const
 
- Public Member Functions inherited from Teuchos::ParameterListAcceptorDefaultBase
RCP< ParameterListgetNonconstParameterList ()
 
RCP< ParameterListunsetParameterList ()
 
RCP< const ParameterListgetParameterList () const
 

Protected Member Functions

void buildMetaData (stk::ParallelMachine parallelMach, STK_Interface &mesh) const
 
void buildElements (stk::ParallelMachine parallelMach, STK_Interface &mesh) const
 
void addSideSets (STK_Interface &mesh) const
 
void addNodeSets (STK_Interface &mesh) const
 
void addEdgeBlocks (STK_Interface &mesh) const
 
void copyCellFieldData (STK_Interface &mesh) const
 
void getOutputTopology ()
 
- Protected Member Functions inherited from Teuchos::ParameterListAcceptorDefaultBase
void setMyParamList (const RCP< ParameterList > &paramList)
 
RCP< ParameterListgetMyNonconstParamList ()
 
RCP< const ParameterListgetMyParamList () const
 

Protected Attributes

Teuchos::RCP
< panzer_stk::STK_Interface
quadMesh_
 
unsigned int machRank_
 Second order mesh. More...
 
unsigned int machSize_
 
bool createEdgeBlocks_
 
bool offsetGIDs_
 If true, offset mesh GIDs to exercise 32-bit limits. More...
 
stk::mesh::EntityId offset_
 
bool print_debug_
 
std::string edgeBlockName_
 
const CellTopologyData * outputTopoData_
 
unsigned int nDim_
 Output mesh topology data. More...
 
unsigned int nNodes_
 Dimension of the mesh. More...
 
std::vector< shards::CellTopology > supportedInputTopos_
 Nodes in one element of the linear basis. More...
 
std::map< const std::string,
const CellTopologyData * > 
outputTopoMap_
 
- Protected Attributes inherited from panzer_stk::STK_MeshFactory
std::vector< Teuchos::RCP
< const PeriodicBC_MatcherBase > > 
periodicBCVec_
 
bool useBBoxSearch_
 
bool enableRebalance_
 
Teuchos::RCP< const
Teuchos::ParameterList
rebalanceList_
 

Additional Inherited Members

- Static Public Member Functions inherited from panzer_stk::STK_MeshFactory
static void parsePeriodicBCList (const Teuchos::RCP< Teuchos::ParameterList > &pl, std::vector< Teuchos::RCP< const PeriodicBC_MatcherBase > > &periodicBC, bool &useBBoxSearch)
 

Detailed Description

This class reads in a second-order (quadratic) mesh and converts it to a first-order (linear) mesh. It will create the cells/nodes/edges, copy the sideset and nodeset data, and copy a list of user provided cell (not nodal) field data into the mesh.

Definition at line 59 of file Panzer_STK_QuadraticToLinearMeshFactory.hpp.

Constructor & Destructor Documentation

panzer_stk::QuadraticToLinearMeshFactory::QuadraticToLinearMeshFactory ( const std::string &  quadMeshFileName,
stk::ParallelMachine  mpi_comm = MPI_COMM_WORLD,
const bool  print_debug = false 
)

Definition at line 57 of file Panzer_STK_QuadraticToLinearMeshFactory.cpp.

panzer_stk::QuadraticToLinearMeshFactory::QuadraticToLinearMeshFactory ( const Teuchos::RCP< panzer_stk::STK_Interface > &  quadMesh,
const bool  print_debug = false 
)

Definition at line 78 of file Panzer_STK_QuadraticToLinearMeshFactory.cpp.

Member Function Documentation

Teuchos::RCP< STK_Interface > panzer_stk::QuadraticToLinearMeshFactory::buildMesh ( stk::ParallelMachine  parallelMach) const
virtual

Build the mesh object.

Implements panzer_stk::STK_MeshFactory.

Definition at line 94 of file Panzer_STK_QuadraticToLinearMeshFactory.cpp.

Teuchos::RCP< STK_Interface > panzer_stk::QuadraticToLinearMeshFactory::buildUncommitedMesh ( stk::ParallelMachine  parallelMach) const
virtual

This builds all the meta data of the mesh. Does not call metaData->commit. Allows user to add solution fields and other pieces. The mesh can be "completed" by calling completeMeshConstruction.

Implements panzer_stk::STK_MeshFactory.

Definition at line 149 of file Panzer_STK_QuadraticToLinearMeshFactory.cpp.

void panzer_stk::QuadraticToLinearMeshFactory::completeMeshConstruction ( STK_Interface mesh,
stk::ParallelMachine  parallelMach 
) const
virtual

Finishes building a mesh object started by buildUncommitedMesh.

Implements panzer_stk::STK_MeshFactory.

Definition at line 167 of file Panzer_STK_QuadraticToLinearMeshFactory.cpp.

void panzer_stk::QuadraticToLinearMeshFactory::setParameterList ( const Teuchos::RCP< Teuchos::ParameterList > &  paramList)
virtual

Derived from ParameterListAcceptor.

From ParameterListAcceptor.

Implements Teuchos::ParameterListAcceptor.

Definition at line 207 of file Panzer_STK_QuadraticToLinearMeshFactory.cpp.

Teuchos::RCP< const Teuchos::ParameterList > panzer_stk::QuadraticToLinearMeshFactory::getValidParameters ( ) const
virtual

Derived from ParameterListAcceptor.

From ParameterListAcceptor.

Reimplemented from Teuchos::ParameterListAcceptor.

Definition at line 222 of file Panzer_STK_QuadraticToLinearMeshFactory.cpp.

void panzer_stk::QuadraticToLinearMeshFactory::buildMetaData ( stk::ParallelMachine  parallelMach,
STK_Interface mesh 
) const
protected

Definition at line 244 of file Panzer_STK_QuadraticToLinearMeshFactory.cpp.

void panzer_stk::QuadraticToLinearMeshFactory::buildElements ( stk::ParallelMachine  parallelMach,
STK_Interface mesh 
) const
protected

Definition at line 330 of file Panzer_STK_QuadraticToLinearMeshFactory.cpp.

void panzer_stk::QuadraticToLinearMeshFactory::addSideSets ( STK_Interface mesh) const
protected

Definition at line 401 of file Panzer_STK_QuadraticToLinearMeshFactory.cpp.

void panzer_stk::QuadraticToLinearMeshFactory::addNodeSets ( STK_Interface mesh) const
protected

Definition at line 430 of file Panzer_STK_QuadraticToLinearMeshFactory.cpp.

void panzer_stk::QuadraticToLinearMeshFactory::addEdgeBlocks ( STK_Interface mesh) const
protected

Definition at line 433 of file Panzer_STK_QuadraticToLinearMeshFactory.cpp.

void panzer_stk::QuadraticToLinearMeshFactory::copyCellFieldData ( STK_Interface mesh) const
protected

Definition at line 451 of file Panzer_STK_QuadraticToLinearMeshFactory.cpp.

void panzer_stk::QuadraticToLinearMeshFactory::getOutputTopology ( )
protected

Infer the output topology from the given input mesh. Also, ensure the input mesh topology matches the parameter list. Currently, each element block must have the same topology.

Definition at line 117 of file Panzer_STK_QuadraticToLinearMeshFactory.cpp.

Member Data Documentation

Teuchos::RCP<panzer_stk::STK_Interface> panzer_stk::QuadraticToLinearMeshFactory::quadMesh_
protected

Definition at line 96 of file Panzer_STK_QuadraticToLinearMeshFactory.hpp.

unsigned int panzer_stk::QuadraticToLinearMeshFactory::machRank_
mutableprotected

Second order mesh.

Definition at line 98 of file Panzer_STK_QuadraticToLinearMeshFactory.hpp.

unsigned int panzer_stk::QuadraticToLinearMeshFactory::machSize_
mutableprotected

Definition at line 98 of file Panzer_STK_QuadraticToLinearMeshFactory.hpp.

bool panzer_stk::QuadraticToLinearMeshFactory::createEdgeBlocks_
protected

Definition at line 100 of file Panzer_STK_QuadraticToLinearMeshFactory.hpp.

bool panzer_stk::QuadraticToLinearMeshFactory::offsetGIDs_
protected

If true, offset mesh GIDs to exercise 32-bit limits.

Definition at line 103 of file Panzer_STK_QuadraticToLinearMeshFactory.hpp.

stk::mesh::EntityId panzer_stk::QuadraticToLinearMeshFactory::offset_
mutableprotected

Definition at line 104 of file Panzer_STK_QuadraticToLinearMeshFactory.hpp.

bool panzer_stk::QuadraticToLinearMeshFactory::print_debug_
protected

Definition at line 106 of file Panzer_STK_QuadraticToLinearMeshFactory.hpp.

std::string panzer_stk::QuadraticToLinearMeshFactory::edgeBlockName_
protected

Definition at line 108 of file Panzer_STK_QuadraticToLinearMeshFactory.hpp.

const CellTopologyData* panzer_stk::QuadraticToLinearMeshFactory::outputTopoData_
protected

Definition at line 110 of file Panzer_STK_QuadraticToLinearMeshFactory.hpp.

unsigned int panzer_stk::QuadraticToLinearMeshFactory::nDim_
protected

Output mesh topology data.

Definition at line 112 of file Panzer_STK_QuadraticToLinearMeshFactory.hpp.

unsigned int panzer_stk::QuadraticToLinearMeshFactory::nNodes_
protected

Dimension of the mesh.

Definition at line 113 of file Panzer_STK_QuadraticToLinearMeshFactory.hpp.

std::vector<shards::CellTopology> panzer_stk::QuadraticToLinearMeshFactory::supportedInputTopos_
protected
Initial value:
= {
shards::CellTopology(shards::getCellTopologyData<shards::Quadrilateral<8>>()),
shards::CellTopology(shards::getCellTopologyData<shards::Triangle<6>>()),
shards::CellTopology(shards::getCellTopologyData<shards::Tetrahedron<10>>()),
shards::CellTopology(shards::getCellTopologyData<shards::Hexahedron<20>>())
}

Nodes in one element of the linear basis.

List of currently supported input topologies.

Definition at line 116 of file Panzer_STK_QuadraticToLinearMeshFactory.hpp.

std::map<const std::string,const CellTopologyData *> panzer_stk::QuadraticToLinearMeshFactory::outputTopoMap_
protected
Initial value:
= {
{shards::getCellTopologyData<shards::Quadrilateral<8>>()->name,
shards::getCellTopologyData<shards::Quadrilateral<4>>()},
{shards::getCellTopologyData<shards::Triangle<6>>()->name,
shards::getCellTopologyData<shards::Triangle<3>>()},
{shards::getCellTopologyData<shards::Tetrahedron<10>>()->name,
shards::getCellTopologyData<shards::Tetrahedron<4>>()},
{shards::getCellTopologyData<shards::Hexahedron<20>>()->name,
shards::getCellTopologyData<shards::Hexahedron<8>>()}
}

Map from input topology to the output shards topology data. The list here is currently supported. Right now, this is one-to-one and may need to be expanded.

Definition at line 125 of file Panzer_STK_QuadraticToLinearMeshFactory.hpp.


The documentation for this class was generated from the following files: