Panzer  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Namespaces | Classes | Functions
panzer_stk Namespace Reference

Namespaces

 workset_utils
 
 periodic_helpers
 

Classes

class  GatherFields
 
class  GatherRefCoords
 
class  ProjectField
 Given a field, perform a local L2 projection onto the desired basis. More...
 
class  ScatterCellAvgQuantity
 
class  ScatterCellAvgVector
 
class  ScatterCellQuantity
 
class  ScatterFields
 
class  ScatterVectorFields
 
class  IOClosureModelFactory
 
class  IOClosureModelFactory_TemplateBuilder
 
class  ModelEvaluatorFactory
 
class  NOXObserverFactory
 
class  TempusObserverFactory
 
class  WorksetFactory
 
class  STKConnManager
 
class  ResponseEvaluatorFactory_SolutionWriter
 
struct  RespFactorySolnWriter_Builder
 
class  CubeHexMeshFactory
 
class  CubeTetMeshFactory
 
class  CustomMeshFactory
 
class  ElementDescriptor
 
class  STK_Interface
 
class  LineMeshFactory
 
class  STK_MeshFactory
 
class  MultiBlockMeshFactory
 
class  CoordMatcher
 
class  PlaneMatcher
 
class  QuarterPlaneMatcher
 
class  WedgeMatcher
 
class  PeriodicBC_MatcherBase
 
class  PeriodicBC_Matcher
 
class  PeriodicBC_Parser
 
class  Quad8ToQuad4MeshFactory
 
class  QuadraticToLinearMeshFactory
 
class  SculptMeshFactory
 
class  SquareQuadMeshFactory
 
class  SquareTriMeshFactory
 

Functions

bool checkSidesetOverlap (const std::string &side_a_name, const std::string &side_b_name, const panzer_stk::STK_Interface &mesh)
 Returns true if the sidesets overlap. More...
 
Teuchos::RCP
< panzer::LocalMeshInfo
generateLocalMeshInfo (const panzer_stk::STK_Interface &mesh)
 Create a structure containing information about the local portion of a given element block. More...
 
Teuchos::RCP
< Thyra::LinearOpWithSolveFactoryBase
< double > > 
buildLOWSFactory (bool blockedAssembly, const Teuchos::RCP< const panzer::GlobalIndexer > &globalIndexer, const Teuchos::RCP< panzer_stk::STKConnManager > &stkConn_manager, int spatialDim, const Teuchos::RCP< const Teuchos::MpiComm< int > > &mpi_comm, const Teuchos::RCP< Teuchos::ParameterList > &strat_params, bool writeCoordinates, bool writeTopo, const Teuchos::RCP< const panzer::GlobalIndexer > &auxGlobalIndexer, bool useCoordinates)
 
Teuchos::RCP
< Thyra::LinearOpWithSolveFactoryBase
< double > > 
buildLOWSFactory (bool blockedAssembly, const Teuchos::RCP< const panzer::GlobalIndexer > &globalIndexer, const Teuchos::RCP< panzer::ConnManager > &conn_manager, int spatialDim, const Teuchos::RCP< const Teuchos::MpiComm< int > > &mpi_comm, const Teuchos::RCP< Teuchos::ParameterList > &strat_params, bool writeCoordinates, bool writeTopo, const Teuchos::RCP< const panzer::GlobalIndexer > &auxGlobalIndexer, bool useCoordinates)
 
Teuchos::RCP< std::vector
< panzer::Workset > > 
buildWorksets (const panzer_stk::STK_Interface &mesh, const std::string &eBlock, const panzer::WorksetNeeds &needs)
 
Teuchos::RCP< std::vector
< panzer::Workset > > 
buildWorksets (const panzer_stk::STK_Interface &mesh, const panzer::WorksetNeeds &needs, const std::string &sideset, const std::string &eBlock, bool useCascade)
 
Teuchos::RCP< std::map
< unsigned, panzer::Workset > > 
buildBCWorksets (const panzer_stk::STK_Interface &mesh, const panzer::WorksetNeeds &needs_a, const std::string &blockid_a, const panzer::WorksetNeeds &needs_b, const std::string &blockid_b, const std::string &sideset)
 
Teuchos::RCP< std::map
< unsigned, panzer::Workset > > 
buildBCWorksets (const panzer_stk::STK_Interface &mesh, const panzer::WorksetNeeds &needs, const std::string &eblockID, const std::string &sidesetID)
 
void computeSidesetNodeNormals (std::unordered_map< unsigned, std::vector< double > > &normals, const Teuchos::RCP< const panzer_stk::STK_Interface > &mesh, const std::string &sidesetName, const std::string &elementBlockName, std::ostream *out=NULL, std::ostream *pout=NULL)
 Computes the normals for all nodes associated with a sideset surface. More...
 
void computeSidesetNodeNormals (std::unordered_map< std::size_t, Kokkos::DynRankView< double, PHX::Device > > &elementToNormalMap, const Teuchos::RCP< const panzer_stk::STK_Interface > &mesh, const std::string &sidesetName, const std::string &elementBlockName, std::ostream *out=NULL, std::ostream *pout=NULL)
 Computes the normals for all nodes associated with a sideset surface. More...
 
std::string transformBCNameForIOSS (std::string &bc_name)
 
std::string version ()
 
std::size_t getElementIdx (const std::vector< stk::mesh::Entity > &elements, stk::mesh::Entity const e)
 
Teuchos::RCP< ElementDescriptorbuildElementDescriptor (stk::mesh::EntityId elmtId, std::vector< stk::mesh::EntityId > &nodes)
 
template<typename Matcher >
Teuchos::RCP
< PeriodicBC_MatcherBase
buildPeriodicBC_Matcher (const std::string &left, const std::string &right, const Matcher &matcher, const std::string type="coord")
 
static std::string trim_left (const std::string &s)
 
static std::string trim_right (const std::string &s)
 
static std::string trim (const std::string &s)
 

Function Documentation

bool panzer_stk::checkSidesetOverlap ( const std::string &  side_a_name,
const std::string &  side_b_name,
const panzer_stk::STK_Interface mesh 
)

Returns true if the sidesets overlap.

Definition at line 50 of file Panzer_STK_CheckSidesetOverlap.cpp.

Teuchos::RCP< panzer::LocalMeshInfo > panzer_stk::generateLocalMeshInfo ( const panzer_stk::STK_Interface mesh)

Create a structure containing information about the local portion of a given element block.

Parameters
[in]meshReference to STK mesh interface
Returns
Structure containing local mesh information

Definition at line 455 of file Panzer_STK_LocalMeshUtilities.cpp.

Teuchos::RCP< Thyra::LinearOpWithSolveFactoryBase< double > > panzer_stk::buildLOWSFactory ( bool  blockedAssembly,
const Teuchos::RCP< const panzer::GlobalIndexer > &  globalIndexer,
const Teuchos::RCP< panzer_stk::STKConnManager > &  stkConn_manager,
int  spatialDim,
const Teuchos::RCP< const Teuchos::MpiComm< int > > &  mpi_comm,
const Teuchos::RCP< Teuchos::ParameterList > &  strat_params,
bool  writeCoordinates = false,
bool  writeTopo = false,
const Teuchos::RCP< const panzer::GlobalIndexer > &  auxGlobalIndexer = Teuchos::null,
bool  useCoordinates = true 
)

Build LOWS factory.

Definition at line 148 of file Panzer_STK_SetupLOWSFactory.cpp.

Teuchos::RCP< Thyra::LinearOpWithSolveFactoryBase< double > > panzer_stk::buildLOWSFactory ( bool  blockedAssembly,
const Teuchos::RCP< const panzer::GlobalIndexer > &  globalIndexer,
const Teuchos::RCP< panzer::ConnManager > &  conn_manager,
int  spatialDim,
const Teuchos::RCP< const Teuchos::MpiComm< int > > &  mpi_comm,
const Teuchos::RCP< Teuchos::ParameterList > &  strat_params,
bool  writeCoordinates = false,
bool  writeTopo = false,
const Teuchos::RCP< const panzer::GlobalIndexer > &  auxGlobalIndexer = Teuchos::null,
bool  useCoordinates = true 
)

Build LOWS factory.

Definition at line 485 of file Panzer_STK_SetupLOWSFactory.cpp.

Teuchos::RCP< std::vector< panzer::Workset > > panzer_stk::buildWorksets ( const panzer_stk::STK_Interface mesh,
const std::string &  eBlock,
const panzer::WorksetNeeds needs 
)

Build volumetric worksets for a STK mesh

Parameters
[in]meshA pointer to the STK_Interface used to construct the worksets
[in]needsPhysics block associated with a particular element block
[in]eBlockElement block to build worksets over (the descriptor information)
Returns
vector of worksets for the corresponding element block.

Definition at line 53 of file Panzer_STK_SetupUtilities.cpp.

Teuchos::RCP< std::vector< panzer::Workset > > panzer_stk::buildWorksets ( const panzer_stk::STK_Interface mesh,
const panzer::WorksetNeeds needs,
const std::string &  sideset,
const std::string &  eBlock,
bool  useCascade = false 
)

Build volumetric worksets for a STK mesh with elements that touch a particular sideset.

Parameters
[in]meshA pointer to the STK_Interface used to construct the worksets
[in]needsNeeds associated with the element block
[in]workset_sizeThe size of each workset measured in the number of elements
[in]sidesetThe sideset id used to locate volume elements associated with the sideset
[in]eBlockElement block to build worksets over (the descriptor information)
[in]useCascadeIf true, worksets will be built for every local node, edge and face that touches the side set. Note that this implies that the workset will have repeated elements. This is useful for higher-order surface flux calculations.
Returns
vector of worksets for the corresponding element block.

Definition at line 74 of file Panzer_STK_SetupUtilities.cpp.

Teuchos::RCP< std::map< unsigned, panzer::Workset > > panzer_stk::buildBCWorksets ( const panzer_stk::STK_Interface mesh,
const panzer::WorksetNeeds needs_a,
const std::string &  blockid_a,
const panzer::WorksetNeeds needs_b,
const std::string &  blockid_b,
const std::string &  sideset 
)

Build side worksets with elements on both sides (this is for DG)

Parameters
[in]meshA pointer to the STK_Interface used to construct the worksets
[in]needs_aNeeds of these worksets
[in]blockid_aElement block id
[in]needs_bNeeds of these worksets
[in]blockid_bElement block id
[in]workset_sizeThe size of each workset measured in the number of elements
[in]sidesetThe sideset id used to locate volume elements associated with the sideset
Returns
vector of worksets for the corresponding edge

Definition at line 193 of file Panzer_STK_SetupUtilities.cpp.

Teuchos::RCP< std::map< unsigned, panzer::Workset > > panzer_stk::buildBCWorksets ( const panzer_stk::STK_Interface mesh,
const panzer::WorksetNeeds needs,
const std::string &  eblockID,
const std::string &  sidesetID 
)

Build boundary condition worksets for a STK mesh

Parameters
[in]meshA pointer to the STK_Interface used to construct the worksets
[in]needsPhysics block associated with the element block
[in]eblockIDName of sideset
[in]sidesetIDName of sideset
Returns
Map relating local element side IDs to a workset.
Note
All elements for a bc that are associated with a particular side ID are grouped into a single workset

Definition at line 282 of file Panzer_STK_SetupUtilities.cpp.

void panzer_stk::computeSidesetNodeNormals ( std::unordered_map< unsigned, std::vector< double > > &  normals,
const Teuchos::RCP< const panzer_stk::STK_Interface > &  mesh,
const std::string &  sidesetName,
const std::string &  elementBlockName,
std::ostream *  out = NULL,
std::ostream *  pout = NULL 
)

Computes the normals for all nodes associated with a sideset surface.

Computes the node normals for a given side set. This computation will use ghosted element contributions for all surfaces associated with the node, but will only use elements associated with the sideset. So if two sidesets adjoin at the seam, the nodes on the seam will NOT have the same node normal. A simple future addition would be to allow for a list of sidesets to be provided to allow for the union computation.

Parameters
[out]normalsMap of the node normals (including ghosted nodes). Key is the global id from the stk mesh node entity id and value is vector of the size of the parent element dimension containing the normal vector components. Vector components will be resized on calling this method.
[in]mesh(Required) Panzer stk mesh
[in]sidesetName(Required) Name of the sideset that the normals will be computed on
[in]elementBlockName(Required) Name of the element block that the outward facing normals will be computed on
[in]out(Optional) The ostream used for serial debug output on print process only. If non-null this will print debug info.
[in]pout(Optional) The ostream used for parallel debug output by all processes. If non-null this will print debug info.

Definition at line 67 of file Panzer_STK_SurfaceNodeNormals.cpp.

void panzer_stk::computeSidesetNodeNormals ( std::unordered_map< std::size_t, Kokkos::DynRankView< double, PHX::Device > > &  elementToNormalMap,
const Teuchos::RCP< const panzer_stk::STK_Interface > &  mesh,
const std::string &  sidesetName,
const std::string &  elementBlockName,
std::ostream *  out = NULL,
std::ostream *  pout = NULL 
)

Computes the normals for all nodes associated with a sideset surface.

Computes the node normals for a given side set. This computation will use ghosted element contributions for all surfaces associated with the node, but will only use elements associated with the sideset. So if two sidesets adjoin at the seam, the nodes on the seam will NOT have the same node normal. A simple future addition would be to allow for a list of sidesets to be provided to allow for the union computation.

Parameters
[out]normalsMap of the node normals (including ghosted nodes). Key is the panzer_stk::STK_Interface local element id and value is a multidimensional array of the size of the number of nodes times the parent element dimension containing the normal vector components. Vector components will be resized on calling this method.
[in]mesh(Required) Panzer stk mesh
[in]sidesetName(Required) Name of the sideset that the normals will be computed on
[in]elementBlockName(Required) Name of the element block that the outward facing normals will be computed on
[in]out(Optional) The ostream used for serial debug output on print process only. If non-null this will print debug info.
[in]pout(Optional) The ostream used for parallel debug output by all processes. If non-null this will print debug info.

Definition at line 241 of file Panzer_STK_SurfaceNodeNormals.cpp.

std::string panzer_stk::transformBCNameForIOSS ( std::string &  bc_name)

Takes a sideset or nodeset name and transforms it for IOSS restrictions. Cubit replaces whitespace in sideset and nodeset names with underscores when writing to exodus. When IOSS reads exodus files into STK, it replaces all capital letters with lowercase letters. This function allows you to use the name labels used in Cubit and tie them to the correct nodeset or sideset in a STK mesh database.

Definition at line 6 of file Panzer_STK_TransformBCNameForIOSS.cpp.

std::string panzer_stk::version ( )

Definition at line 48 of file Panzer_STK_Version.cpp.

std::size_t panzer_stk::getElementIdx ( const std::vector< stk::mesh::Entity > &  elements,
stk::mesh::Entity const  e 
)
inline

Definition at line 415 of file Panzer_STKConnManager.cpp.

Teuchos::RCP< ElementDescriptor > panzer_stk::buildElementDescriptor ( stk::mesh::EntityId  elmtId,
std::vector< stk::mesh::EntityId > &  nodes 
)

Constructor function for building the element descriptors.

Definition at line 96 of file Panzer_STK_Interface.cpp.

template<typename Matcher >
Teuchos::RCP<PeriodicBC_MatcherBase> panzer_stk::buildPeriodicBC_Matcher ( const std::string &  left,
const std::string &  right,
const Matcher &  matcher,
const std::string  type = "coord" 
)

A simple constructor function for building a matcher object. This prevents the need to directly instantiate the templated derived class. It is a convenience.

Definition at line 406 of file Panzer_STK_PeriodicBC_Matcher.hpp.

static std::string panzer_stk::trim_left ( const std::string &  s)
static

Definition at line 176 of file Panzer_STK_PeriodicBC_Parser.cpp.

static std::string panzer_stk::trim_right ( const std::string &  s)
static

Definition at line 183 of file Panzer_STK_PeriodicBC_Parser.cpp.

static std::string panzer_stk::trim ( const std::string &  s)
static

Definition at line 190 of file Panzer_STK_PeriodicBC_Parser.cpp.