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

Namespaces

 workset_utils
 
 periodic_helpers
 

Classes

class  GatherFields
 
class  GatherRefCoords
 
class  ScatterFields
 
class  IOClosureModelFactory
 
class  IOClosureModelFactory_TemplateBuilder
 
class  STKConnManager
 
class  ModelEvaluatorFactory
 
class  NOXObserverFactory
 
class  RythmosObserverFactory
 
struct  PermFunctor
 
class  WorksetFactory
 
class  LocalIdCompare
 
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  PeriodicBC_MatcherBase
 
class  PeriodicBC_Matcher
 
class  PeriodicBC_Parser
 
class  SculptMeshFactory
 
class  SquareQuadMeshFactory
 
class  SquareTriMeshFactory
 

Typedefs

typedef
panzer_stk_classic::STK_Interface::SolutionFieldType 
VariableField
 

Functions

 PHX_EVALUATOR_CTOR (ScatterCellAvgQuantity, p)
 
 PHX_POST_REGISTRATION_SETUP (ScatterCellAvgQuantity, d, fm)
 
 PHX_EVALUATE_FIELDS (ScatterCellAvgQuantity, workset)
 
 PHX_EVALUATOR_CTOR (ScatterCellAvgVector, p)
 
 PHX_POST_REGISTRATION_SETUP (ScatterCellAvgVector, d, fm)
 
 PHX_EVALUATE_FIELDS (ScatterCellAvgVector, workset)
 
 PHX_EVALUATOR_CTOR (ScatterCellQuantity, p)
 
 PHX_POST_REGISTRATION_SETUP (ScatterCellQuantity, d, fm)
 
 PHX_EVALUATE_FIELDS (ScatterCellQuantity, workset)
 
 ScatterVectorFields (const std::string &scatterName, const Teuchos::RCP< STK_Interface > mesh, const Teuchos::RCP< const panzer::PointRule > &pointRule, const std::vector< std::string > &names)
 
 PHX_EVALUATOR_CTOR (ScatterVectorFields, p)
 
 PHX_POST_REGISTRATION_SETUP (ScatterVectorFields, d, fm)
 
 PHX_EVALUATE_FIELDS (ScatterVectorFields, workset)
 
template Teuchos::RCP
< Thyra::LinearOpWithSolveFactoryBase
< double > > 
buildLOWSFactory< int > (bool blockedAssembly, const Teuchos::RCP< const panzer::UniqueGlobalIndexerBase > &globalIndexer, const Teuchos::RCP< panzer_stk_classic::STKConnManager< int > > &stkConn_manager, int spatialDim, const Teuchos::RCP< const Teuchos::MpiComm< int > > &mpi_comm, const Teuchos::RCP< Teuchos::ParameterList > &strat_params, bool writeCoordinates, bool writeTopo)
 
template Teuchos::RCP
< Thyra::LinearOpWithSolveFactoryBase
< double > > 
buildLOWSFactory< panzer::Ordinal64 > (bool blockedAssembly, const Teuchos::RCP< const panzer::UniqueGlobalIndexerBase > &globalIndexer, const Teuchos::RCP< panzer_stk_classic::STKConnManager< panzer::Ordinal64 > > &stkConn_manager, int spatialDim, const Teuchos::RCP< const Teuchos::MpiComm< int > > &mpi_comm, const Teuchos::RCP< Teuchos::ParameterList > &strat_params, bool writeCoordinates, bool writeTopo)
 
Teuchos::RCP
< Thyra::LinearOpWithSolveFactoryBase
< double > > 
buildLOWSFactory (bool blockedAssembly, const Teuchos::RCP< const panzer::UniqueGlobalIndexerBase > &globalIndexer, const Teuchos::RCP< panzer::ConnManagerBase< int > > &conn_manager, int spatialDim, const Teuchos::RCP< const Teuchos::MpiComm< int > > &mpi_comm, const Teuchos::RCP< Teuchos::ParameterList > &strat_params, bool writeCoordinates, bool writeTopo)
 
template<typename GO >
Teuchos::RCP
< Thyra::LinearOpWithSolveFactoryBase
< double > > 
buildLOWSFactory (bool blockedAssembly, const Teuchos::RCP< const panzer::UniqueGlobalIndexerBase > &globalIndexer, const Teuchos::RCP< panzer_stk_classic::STKConnManager< GO > > &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)
 
template<typename GO >
void writeTopology (const panzer::BlockedDOFManager< int, GO > &blkDofs)
 
Teuchos::RCP< std::vector
< panzer::Workset > > 
buildWorksets (const panzer_stk_classic::STK_Interface &mesh, const panzer::PhysicsBlock &pb)
 
Teuchos::RCP< std::vector
< panzer::Workset > > 
buildWorksets (const panzer_stk_classic::STK_Interface &mesh, const std::string &eBlock, const panzer::WorksetNeeds &needs)
 
Teuchos::RCP< std::vector
< panzer::Workset > > 
buildWorksets (const panzer_stk_classic::STK_Interface &mesh, const panzer::PhysicsBlock &pb, const std::string &sideset, bool useCascade)
 
Teuchos::RCP< std::vector
< panzer::Workset > > 
buildWorksets (const panzer_stk_classic::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_classic::STK_Interface &mesh, const panzer::PhysicsBlock &pb_a, const panzer::PhysicsBlock &pb_b, const std::string &sideset)
 
Teuchos::RCP< std::map
< unsigned, panzer::Workset > > 
buildBCWorksets (const panzer_stk_classic::STK_Interface &mesh, const panzer::PhysicsBlock &pb, const std::string &sidesetID)
 
Teuchos::RCP< std::map
< unsigned, panzer::Workset > > 
buildBCWorksets (const panzer_stk_classic::STK_Interface &mesh, const panzer::WorksetNeeds &needs, const std::string &eblockID, const std::string &sidesetID)
 
void computeSidesetNodeNormals (boost::unordered_map< unsigned, std::vector< double > > &normals, const Teuchos::RCP< const panzer_stk_classic::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 (boost::unordered_map< std::size_t, Intrepid::FieldContainer< double > > &elementToNormalMap, const Teuchos::RCP< const panzer_stk_classic::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...
 
template<typename GlobalOrdinal >
static void gather_in_block (const std::string &blockId, const panzer::UniqueGlobalIndexer< int, GlobalOrdinal > &dofMngr, const Epetra_Vector &x, const std::vector< std::size_t > &localCellIds, std::map< std::string, Intrepid::FieldContainer< double > > &fc)
 
static void build_local_ids (const panzer_stk_classic::STK_Interface &mesh, std::map< std::string, Teuchos::RCP< std::vector< std::size_t > > > &localIds)
 
void write_cell_data (panzer_stk_classic::STK_Interface &mesh, const std::vector< double > &data, const std::string &fieldName)
 
template<typename GlobalOrdinal >
void write_solution_data (const panzer::UniqueGlobalIndexer< int, GlobalOrdinal > &dofMngr, panzer_stk_classic::STK_Interface &mesh, const Epetra_MultiVector &x, const std::string &prefix, const std::string &postfix)
 
template void write_solution_data< int > (const panzer::UniqueGlobalIndexer< int, int > &dofMngr, panzer_stk_classic::STK_Interface &mesh, const Epetra_MultiVector &x, const std::string &prefix, const std::string &postfix)
 
template<typename GlobalOrdinal >
void write_solution_data (const panzer::UniqueGlobalIndexer< int, GlobalOrdinal > &dofMngr, panzer_stk_classic::STK_Interface &mesh, const Epetra_Vector &x, const std::string &prefix, const std::string &postfix)
 
template void write_solution_data< int > (const panzer::UniqueGlobalIndexer< int, int > &dofMngr, panzer_stk_classic::STK_Interface &mesh, const Epetra_Vector &x, const std::string &prefix, const std::string &postfix)
 
template<typename GlobalOrdinal >
void gather_in_block (const std::string &blockId, const panzer::UniqueGlobalIndexer< int, GlobalOrdinal > &dofMngr, const Epetra_Vector &x, const std::vector< std::size_t > &localCellIds, std::map< std::string, Intrepid::FieldContainer< double > > &fc)
 
template<typename RAContainer , class Compare >
void sorted_permutation (const RAContainer &cont, std::vector< std::size_t > &permutation, const Compare &comp)
 
template<typename RAContainer >
void sorted_permutation (const RAContainer &cont, std::vector< std::size_t > &permutation)
 
std::string version ()
 
std::size_t getElementIdx (const std::vector< stk_classic::mesh::Entity * > &elements, const stk_classic::mesh::Entity *const e)
 
Teuchos::RCP< ElementDescriptorbuildElementDescriptor (stk_classic::mesh::EntityId elmtId, std::vector< stk_classic::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)
 

Variables

std::size_t numValues_
 
std::vector< PHX::MDField
< const ScalarT, panzer::Cell,
panzer::Point > > 
scatterFields_
 
Teuchos::RCP< STK_Interfacemesh_
 
std::vector< VariableField * > stkFields_
 
std::vector< std::string > names_
 
PHX::MDField< const ScalarT,
panzer::Cell, panzer::IP,
panzer::Dim > 
pointField_
 
int spatialDimension_
 

Typedef Documentation

This class is a scatter operation to the mesh. It takes a set of field names and an integration rule those quantities are then averaged over the cell and written to the mesh.

The constructor takes a STK_Interface RCP and parameter list that is required to contain the following two fields "Scatter Name" string specifying the name of this evaulator "Field Names" of type this is a comma seperated list of strings, "IR" of type Teuchos::RCP<panzer::IntegrationRule> and "Mesh" of type Teuchos::RCP<const panzer_stk_classic::STK_Interface>.

This class is a scatter operation to the mesh. It takes a set of field names and an integration rule. Those quantities are components of vector fields. They are averaged over the cell and written to the mesh.

The constructor takes a STK_Interface RCP and parameter list that is required to contain the following fields: "Scatter Name" string specifying the name of this evaulator "Field Names" of type this is a comma seperated list of strings, "IR" of type Teuchos::RCP<panzer::IntegrationRule> and "Mesh" of type Teuchos::RCP<const panzer_stk_classic::STK_Interface>.

This class is a scatter operation to the mesh. It takes a set of field names and basis objects and then writes them to the mesh object.

The constructor takes a STK_Interface RCP and parameter list that is required to contain the following two fields "Field Names" of type Teuchos::RCP<std::vector<std::string> >, "Basis" of type Teuchos::RCP<panzer::BasisIRLayout> and "Mesh" of type Teuchos::RCP<const panzer_stk_classic::STK_Interface>.

Definition at line 74 of file Panzer_STK_ScatterCellAvgQuantity_decl.hpp.

Function Documentation

panzer_stk_classic::PHX_EVALUATOR_CTOR ( ScatterCellAvgQuantity  ,
 
)

Definition at line 62 of file Panzer_STK_ScatterCellAvgQuantity_impl.hpp.

panzer_stk_classic::PHX_POST_REGISTRATION_SETUP ( ScatterCellAvgQuantity  ,
,
fm   
)

Definition at line 92 of file Panzer_STK_ScatterCellAvgQuantity_impl.hpp.

panzer_stk_classic::PHX_EVALUATE_FIELDS ( ScatterCellAvgQuantity  ,
workset   
)

Definition at line 104 of file Panzer_STK_ScatterCellAvgQuantity_impl.hpp.

panzer_stk_classic::PHX_EVALUATOR_CTOR ( ScatterCellAvgVector  ,
 
)

Definition at line 62 of file Panzer_STK_ScatterCellAvgVector_impl.hpp.

panzer_stk_classic::PHX_POST_REGISTRATION_SETUP ( ScatterCellAvgVector  ,
,
fm   
)

Definition at line 94 of file Panzer_STK_ScatterCellAvgVector_impl.hpp.

panzer_stk_classic::PHX_EVALUATE_FIELDS ( ScatterCellAvgVector  ,
workset   
)

Definition at line 108 of file Panzer_STK_ScatterCellAvgVector_impl.hpp.

panzer_stk_classic::PHX_EVALUATOR_CTOR ( ScatterCellQuantity  ,
 
)

Definition at line 62 of file Panzer_STK_ScatterCellQuantity_impl.hpp.

panzer_stk_classic::PHX_POST_REGISTRATION_SETUP ( ScatterCellQuantity  ,
,
fm   
)

Definition at line 90 of file Panzer_STK_ScatterCellQuantity_impl.hpp.

panzer_stk_classic::PHX_EVALUATE_FIELDS ( ScatterCellQuantity  ,
workset   
)

Definition at line 100 of file Panzer_STK_ScatterCellQuantity_impl.hpp.

panzer_stk_classic::ScatterVectorFields::ScatterVectorFields ( const std::string &  scatterName,
const Teuchos::RCP< STK_Interface >  mesh,
const Teuchos::RCP< const panzer::PointRule > &  pointRule,
const std::vector< std::string > &  names 
)

Definition at line 69 of file Panzer_STK_ScatterVectorFields_impl.hpp.

panzer_stk_classic::PHX_EVALUATOR_CTOR ( ScatterVectorFields  ,
 
)

Definition at line 61 of file Panzer_STK_ScatterVectorFields_impl.hpp.

panzer_stk_classic::PHX_POST_REGISTRATION_SETUP ( ScatterVectorFields  ,
,
fm   
)

Definition at line 99 of file Panzer_STK_ScatterVectorFields_impl.hpp.

panzer_stk_classic::PHX_EVALUATE_FIELDS ( ScatterVectorFields  ,
workset   
)

Definition at line 109 of file Panzer_STK_ScatterVectorFields_impl.hpp.

template Teuchos::RCP<Thyra::LinearOpWithSolveFactoryBase<double> > panzer_stk_classic::buildLOWSFactory< int > ( bool  blockedAssembly,
const Teuchos::RCP< const panzer::UniqueGlobalIndexerBase > &  globalIndexer,
const Teuchos::RCP< panzer_stk_classic::STKConnManager< int > > &  stkConn_manager,
int  spatialDim,
const Teuchos::RCP< const Teuchos::MpiComm< int > > &  mpi_comm,
const Teuchos::RCP< Teuchos::ParameterList > &  strat_params,
bool  writeCoordinates,
bool  writeTopo 
)
template Teuchos::RCP<Thyra::LinearOpWithSolveFactoryBase<double> > panzer_stk_classic::buildLOWSFactory< panzer::Ordinal64 > ( bool  blockedAssembly,
const Teuchos::RCP< const panzer::UniqueGlobalIndexerBase > &  globalIndexer,
const Teuchos::RCP< panzer_stk_classic::STKConnManager< panzer::Ordinal64 > > &  stkConn_manager,
int  spatialDim,
const Teuchos::RCP< const Teuchos::MpiComm< int > > &  mpi_comm,
const Teuchos::RCP< Teuchos::ParameterList > &  strat_params,
bool  writeCoordinates,
bool  writeTopo 
)
Teuchos::RCP< Thyra::LinearOpWithSolveFactoryBase< double > > panzer_stk_classic::buildLOWSFactory ( bool  blockedAssembly,
const Teuchos::RCP< const panzer::UniqueGlobalIndexerBase > &  globalIndexer,
const Teuchos::RCP< panzer::ConnManagerBase< int > > &  conn_manager,
int  spatialDim,
const Teuchos::RCP< const Teuchos::MpiComm< int > > &  mpi_comm,
const Teuchos::RCP< Teuchos::ParameterList > &  strat_params,
bool  writeCoordinates,
bool  writeTopo 
)

Definition at line 83 of file Panzer_STK_SetupLOWSFactory.cpp.

template<typename GO >
Teuchos::RCP< Thyra::LinearOpWithSolveFactoryBase< double > > panzer_stk_classic::buildLOWSFactory ( bool  blockedAssembly,
const Teuchos::RCP< const panzer::UniqueGlobalIndexerBase > &  globalIndexer,
const Teuchos::RCP< panzer_stk_classic::STKConnManager< GO > > &  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 
)

Build LOWS factory.

Definition at line 203 of file Panzer_STK_SetupLOWSFactory_impl.hpp.

template<typename GO >
void panzer_stk_classic::writeTopology ( const panzer::BlockedDOFManager< int, GO > &  blkDofs)

Definition at line 420 of file Panzer_STK_SetupLOWSFactory_impl.hpp.

Teuchos::RCP< std::vector< panzer::Workset > > panzer_stk_classic::buildWorksets ( const panzer_stk_classic::STK_Interface mesh,
const panzer::PhysicsBlock pb 
)

Build volumetric worksets for a STK mesh

Parameters
[in]meshA pointer to the STK_Interface used to construct the worksets
[in]pbPhysics block associated with a particular element block
Returns
vector of worksets for the corresponding element block.

Definition at line 52 of file Panzer_STK_SetupUtilities.cpp.

Teuchos::RCP< std::vector< panzer::Workset > > panzer_stk_classic::buildWorksets ( const panzer_stk_classic::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 72 of file Panzer_STK_SetupUtilities.cpp.

Teuchos::RCP< std::vector< panzer::Workset > > panzer_stk_classic::buildWorksets ( const panzer_stk_classic::STK_Interface mesh,
const panzer::PhysicsBlock pb,
const std::string &  sideset,
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]pbPhysics block 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]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 93 of file Panzer_STK_SetupUtilities.cpp.

Teuchos::RCP< std::vector< panzer::Workset > > panzer_stk_classic::buildWorksets ( const panzer_stk_classic::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 211 of file Panzer_STK_SetupUtilities.cpp.

Teuchos::RCP< std::map< unsigned, panzer::Workset > > panzer_stk_classic::buildBCWorksets ( const panzer_stk_classic::STK_Interface mesh,
const panzer::PhysicsBlock pb_a,
const panzer::PhysicsBlock pb_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]pb_aPhysics block associated with the element block
[in]pb_bPhysics block 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
Returns
vector of worksets for the corresponding edge

Definition at line 330 of file Panzer_STK_SetupUtilities.cpp.

Teuchos::RCP< std::map< unsigned, panzer::Workset > > panzer_stk_classic::buildBCWorksets ( const panzer_stk_classic::STK_Interface mesh,
const panzer::PhysicsBlock pb,
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]pbPhysics block associated with the element block
[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 417 of file Panzer_STK_SetupUtilities.cpp.

Teuchos::RCP< std::map< unsigned, panzer::Workset > > panzer_stk_classic::buildBCWorksets ( const panzer_stk_classic::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 493 of file Panzer_STK_SetupUtilities.cpp.

void panzer_stk_classic::computeSidesetNodeNormals ( boost::unordered_map< unsigned, std::vector< double > > &  normals,
const Teuchos::RCP< const panzer_stk_classic::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 66 of file Panzer_STK_SurfaceNodeNormals.cpp.

void panzer_stk_classic::computeSidesetNodeNormals ( boost::unordered_map< std::size_t, Intrepid::FieldContainer< double > > &  elementToNormalMap,
const Teuchos::RCP< const panzer_stk_classic::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_classic::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 211 of file Panzer_STK_SurfaceNodeNormals.cpp.

template<typename GlobalOrdinal >
static void panzer_stk_classic::gather_in_block ( const std::string &  blockId,
const panzer::UniqueGlobalIndexer< int, GlobalOrdinal > &  dofMngr,
const Epetra_Vector x,
const std::vector< std::size_t > &  localCellIds,
std::map< std::string, Intrepid::FieldContainer< double > > &  fc 
)
static

Definition at line 178 of file Panzer_STK_Utilities.cpp.

void panzer_stk_classic::build_local_ids ( const panzer_stk_classic::STK_Interface mesh,
std::map< std::string, Teuchos::RCP< std::vector< std::size_t > > > &  localIds 
)
static

Definition at line 253 of file Panzer_STK_Utilities.cpp.

void panzer_stk_classic::write_cell_data ( panzer_stk_classic::STK_Interface mesh,
const std::vector< double > &  data,
const std::string &  fieldName 
)

Write a vector to the cell data of a STK mesh. This will look up the cell field prefix+fieldName+postfix, which is assumed to be in the STK mesh object. If not an assertion exeption will be thrown.

Parameters
[in]meshSTK mesh object
[in]dataVector of doubles equatl to the total number of elements on this processor
[in]fieldNameName of field to be written (must be a STK field)

Definition at line 71 of file Panzer_STK_Utilities.cpp.

template<typename GlobalOrdinal >
void panzer_stk_classic::write_solution_data ( const panzer::UniqueGlobalIndexer< int, GlobalOrdinal > &  dofMngr,
panzer_stk_classic::STK_Interface mesh,
const Epetra_MultiVector x,
const std::string &  prefix,
const std::string &  postfix 
)

Definition at line 95 of file Panzer_STK_Utilities.cpp.

template void panzer_stk_classic::write_solution_data< int > ( const panzer::UniqueGlobalIndexer< int, int > &  dofMngr,
panzer_stk_classic::STK_Interface mesh,
const Epetra_MultiVector x,
const std::string &  prefix,
const std::string &  postfix 
)
template<typename GlobalOrdinal >
void panzer_stk_classic::write_solution_data ( const panzer::UniqueGlobalIndexer< int, GlobalOrdinal > &  dofMngr,
panzer_stk_classic::STK_Interface mesh,
const Epetra_Vector x,
const std::string &  prefix,
const std::string &  postfix 
)

Definition at line 108 of file Panzer_STK_Utilities.cpp.

template void panzer_stk_classic::write_solution_data< int > ( const panzer::UniqueGlobalIndexer< int, int > &  dofMngr,
panzer_stk_classic::STK_Interface mesh,
const Epetra_Vector x,
const std::string &  prefix,
const std::string &  postfix 
)
template<typename GlobalOrdinal >
void panzer_stk_classic::gather_in_block ( const std::string &  blockId,
const panzer::UniqueGlobalIndexer< int, GlobalOrdinal > &  dofMngr,
const Epetra_Vector x,
const std::vector< std::size_t > &  localCellIds,
std::map< std::string, Intrepid::FieldContainer< double > > &  fc 
)

Definition at line 178 of file Panzer_STK_Utilities.cpp.

template<typename RAContainer , class Compare >
void panzer_stk_classic::sorted_permutation ( const RAContainer &  cont,
std::vector< std::size_t > &  permutation,
const Compare &  comp 
)

Using a container, compute the sorted permutation vector do not modifiy the original container.

Motivated by this board on StackOverflow: http://stackoverflow.com/questions/4523220/sorting-a-vector-of-double-precision-reals-and-obtain-their-order

Definition at line 127 of file Panzer_STK_Utilities.hpp.

template<typename RAContainer >
void panzer_stk_classic::sorted_permutation ( const RAContainer &  cont,
std::vector< std::size_t > &  permutation 
)

Using a container, compute the sorted permutation vector do not modifiy the original container.

Motivated by this board on StackOverflow: http://stackoverflow.com/questions/4523220/sorting-a-vector-of-double-precision-reals-and-obtain-their-order

Definition at line 120 of file Panzer_STK_Utilities.hpp.

std::string panzer_stk_classic::version ( )

Definition at line 48 of file Panzer_STK_Version.cpp.

std::size_t panzer_stk_classic::getElementIdx ( const std::vector< stk_classic::mesh::Entity * > &  elements,
const stk_classic::mesh::Entity *const  e 
)
inline

Definition at line 383 of file Panzer_STKConnManager_impl.hpp.

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

Constructor function for building the element descriptors.

Definition at line 86 of file Panzer_STK_Interface.cpp.

template<typename Matcher >
Teuchos::RCP<PeriodicBC_MatcherBase> panzer_stk_classic::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 233 of file Panzer_STK_PeriodicBC_Matcher.hpp.

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

Definition at line 145 of file Panzer_STK_PeriodicBC_Parser.cpp.

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

Definition at line 152 of file Panzer_STK_PeriodicBC_Parser.cpp.

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

Definition at line 159 of file Panzer_STK_PeriodicBC_Parser.cpp.

Variable Documentation

std::size_t panzer_stk_classic::numValues_

Definition at line 76 of file Panzer_STK_ScatterCellAvgQuantity_decl.hpp.

std::vector< PHX::MDField< const ScalarT, panzer::Cell, panzer::IP, panzer::Dim > > panzer_stk_classic::scatterFields_

This class is a scatter operation to the mesh. It takes a set of field names on cells and writes them to the mesh.

The constructor takes a STK_Interface RCP and parameter list that is required to contain the following fields "Scatter Name" string specifying the name of this evaulator "Field Names" of type this is a comma seperated list of strings, "Workset Size" of type int "Mesh" of type Teuchos::RCP<const panzer_stk_classic::STK_Interface>.

Definition at line 78 of file Panzer_STK_ScatterCellAvgQuantity_decl.hpp.

Teuchos::RCP< STK_Interface > panzer_stk_classic::mesh_

Definition at line 79 of file Panzer_STK_ScatterCellAvgQuantity_decl.hpp.

std::vector< VariableField * > panzer_stk_classic::stkFields_

Definition at line 81 of file Panzer_STK_ScatterCellAvgQuantity_decl.hpp.

std::vector<std::string> panzer_stk_classic::names_

Definition at line 74 of file Panzer_STK_ScatterVectorFields_decl.hpp.

PHX::MDField<const ScalarT,panzer::Cell,panzer::IP,panzer::Dim> panzer_stk_classic::pointField_

Definition at line 76 of file Panzer_STK_ScatterVectorFields_decl.hpp.

int panzer_stk_classic::spatialDimension_

Definition at line 79 of file Panzer_STK_ScatterVectorFields_decl.hpp.