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

#include <Panzer_STK_Interface.hpp>

Classes

struct  EdgeBlockException
 
struct  ElementBlockException
 
struct  FaceBlockException
 
class  LocalIdCompare
 
struct  SidesetException
 

Public Types

typedef double ProcIdData
 
typedef stk::mesh::Field< double > SolutionFieldType
 
typedef stk::mesh::Field< double > VectorFieldType
 
typedef stk::mesh::Field
< ProcIdData
ProcIdFieldType
 

Public Member Functions

 STK_Interface ()
 
 STK_Interface (unsigned dim)
 
 STK_Interface (Teuchos::RCP< stk::mesh::MetaData > metaData)
 
void addElementBlock (const std::string &name, const CellTopologyData *ctData)
 
void addEdgeBlock (const std::string &elemBlockName, const std::string &edgeBlockName, const stk::topology &topology)
 
void addEdgeBlock (const std::string &elemBlockName, const std::string &edgeBlockName, const CellTopologyData *ctData)
 
void addFaceBlock (const std::string &elemBlockName, const std::string &faceBlockName, const stk::topology &topology)
 
void addFaceBlock (const std::string &elemBlockName, const std::string &faceBlockName, const CellTopologyData *ctData)
 
void addSideset (const std::string &name, const CellTopologyData *ctData)
 
void addNodeset (const std::string &name)
 
void addSolutionField (const std::string &fieldName, const std::string &blockId)
 
void addCellField (const std::string &fieldName, const std::string &blockId)
 
void addEdgeField (const std::string &fieldName, const std::string &blockId)
 
void addFaceField (const std::string &fieldName, const std::string &blockId)
 
void addMeshCoordFields (const std::string &blockId, const std::vector< std::string > &coordField, const std::string &dispPrefix)
 
void addInformationRecords (const std::vector< std::string > &info_records)
 
void initialize (stk::ParallelMachine parallelMach, bool setupIO=true, const bool buildRefinementSupport=false)
 
void instantiateBulkData (stk::ParallelMachine parallelMach)
 
void beginModification ()
 
void endModification ()
 
void addNode (stk::mesh::EntityId gid, const std::vector< double > &coord)
 
void addElement (const Teuchos::RCP< ElementDescriptor > &ed, stk::mesh::Part *block)
 
void addEdges ()
 
void addFaces ()
 
void addEntityToSideset (stk::mesh::Entity entity, stk::mesh::Part *sideset)
 
void addEntityToNodeset (stk::mesh::Entity entity, stk::mesh::Part *nodeset)
 
void addEntityToEdgeBlock (stk::mesh::Entity entity, stk::mesh::Part *edgeblock)
 
void addEntitiesToEdgeBlock (std::vector< stk::mesh::Entity > entities, stk::mesh::Part *edgeblock)
 
void addEntityToFaceBlock (stk::mesh::Entity entity, stk::mesh::Part *faceblock)
 
void addEntitiesToFaceBlock (std::vector< stk::mesh::Entity > entities, stk::mesh::Part *faceblock)
 
const VectorFieldTypegetCoordinatesField () const
 
const VectorFieldTypegetEdgesField () const
 
const VectorFieldTypegetFacesField () const
 
const double * getNodeCoordinates (stk::mesh::EntityId nodeId) const
 
const double * getNodeCoordinates (stk::mesh::Entity node) const
 
void getSubcellIndices (unsigned entityRank, stk::mesh::EntityId elementId, std::vector< stk::mesh::EntityId > &subcellIds) const
 
void getMyElements (std::vector< stk::mesh::Entity > &elements) const
 
void getMyElements (const std::string &blockID, std::vector< stk::mesh::Entity > &elements) const
 
void getNeighborElements (std::vector< stk::mesh::Entity > &elements) const
 
void getNeighborElements (const std::string &blockID, std::vector< stk::mesh::Entity > &elements) const
 
void getMyEdges (std::vector< stk::mesh::Entity > &edges) const
 
void getMyEdges (const std::string &edgeBlockName, std::vector< stk::mesh::Entity > &edges) const
 
void getMyEdges (const std::string &edgeBlockName, const std::string &blockName, std::vector< stk::mesh::Entity > &edges) const
 
void getAllEdges (const std::string &edgeBlockName, std::vector< stk::mesh::Entity > &edges) const
 
void getAllEdges (const std::string &edgeBlockName, const std::string &blockName, std::vector< stk::mesh::Entity > &edges) const
 
void getMyFaces (std::vector< stk::mesh::Entity > &faces) const
 
void getMyFaces (const std::string &faceBlockName, std::vector< stk::mesh::Entity > &faces) const
 
void getMyFaces (const std::string &faceBlockName, const std::string &blockName, std::vector< stk::mesh::Entity > &faces) const
 
void getAllFaces (const std::string &faceBlockName, std::vector< stk::mesh::Entity > &faces) const
 
void getAllFaces (const std::string &faceBlockName, const std::string &blockName, std::vector< stk::mesh::Entity > &faces) const
 
void getMySides (const std::string &sideName, std::vector< stk::mesh::Entity > &sides) const
 
void getMySides (const std::string &sideName, const std::string &blockName, std::vector< stk::mesh::Entity > &sides) const
 
void getAllSides (const std::string &sideName, std::vector< stk::mesh::Entity > &sides) const
 
void getAllSides (const std::string &sideName, const std::string &blockName, std::vector< stk::mesh::Entity > &sides) const
 
void getMyNodes (const std::string &sideName, const std::string &blockName, std::vector< stk::mesh::Entity > &nodes) const
 
stk::mesh::Entity findConnectivityById (stk::mesh::Entity src, stk::mesh::EntityRank tgt_rank, unsigned rel_id) const
 
void writeToExodus (const std::string &filename, const bool append=false)
 Write this mesh and associated fields to the given output file. More...
 
void setupExodusFile (const std::string &filename, const bool append=false, const bool append_after_restart_time=false, const double restart_time=0.0)
 Set up an output Exodus file for writing results. More...
 
void writeToExodus (double timestep)
 Write this mesh and associated fields at the given timestep. More...
 
void addGlobalToExodus (const std::string &key, const int &value)
 Add an int global variable to the information to be written to the Exodus output file. More...
 
void addGlobalToExodus (const std::string &key, const double &value)
 Add a double global variable to the information to be written to the Exodus output file. More...
 
void addGlobalToExodus (const std::string &key, const std::vector< int > &value)
 Add a std::vector<int> global variable to the information to be written to the Exodus output file. More...
 
void addGlobalToExodus (const std::string &key, const std::vector< double > &value)
 Add a std::vector<double> global variable to the information to be written to the Exodus output file. More...
 
Teuchos::RCP< const
Teuchos::Comm< int > > 
getComm () const
 get the comm associated with this mesh More...
 
Teuchos::RCP< stk::mesh::BulkData > getBulkData () const
 
Teuchos::RCP< stk::mesh::MetaData > getMetaData () const
 
bool isWritable () const
 
bool isModifiable () const
 
unsigned getDimension () const
 get the dimension More...
 
std::size_t getNumElementBlocks () const
 get the block count More...
 
void getElementBlockNames (std::vector< std::string > &names) const
 
void getSidesetNames (std::vector< std::string > &name) const
 
void getNodesetNames (std::vector< std::string > &name) const
 
void getEdgeBlockNames (std::vector< std::string > &names) const
 
void getFaceBlockNames (std::vector< std::string > &names) const
 
stk::mesh::Part * getOwnedPart () const
 Get a pointer to the locally owned part. More...
 
stk::mesh::Part * getElementBlockPart (const std::string &name) const
 get the block part More...
 
stk::mesh::Part * getEdgeBlock (const std::string &name) const
 get the block part More...
 
stk::mesh::Part * getFaceBlock (const std::string &name) const
 get the block part More...
 
std::size_t getNumSidesets () const
 get the side set count More...
 
stk::mesh::Part * getSideset (const std::string &name) const
 
std::size_t getNumNodesets () const
 get the side set count More...
 
stk::mesh::Part * getNodeset (const std::string &name) const
 
std::size_t getEntityCounts (unsigned entityRank) const
 get the global counts for the entity of specified rank More...
 
stk::mesh::EntityId getMaxEntityId (unsigned entityRank) const
 get max entity ID of type entityRank More...
 
void getElementsSharingNode (stk::mesh::EntityId nodeId, std::vector< stk::mesh::Entity > &elements) const
 get a set of elements sharing a single node More...
 
void getNodeIdsForElement (stk::mesh::Entity element, std::vector< stk::mesh::EntityId > &nodeIds) const
 get a list of node ids for nodes connected to an element More...
 
void getOwnedElementsSharingNode (stk::mesh::Entity node, std::vector< stk::mesh::Entity > &elements, std::vector< int > &relIds) const
 
void getOwnedElementsSharingNode (stk::mesh::EntityId nodeId, std::vector< stk::mesh::Entity > &elements, std::vector< int > &relIds, unsigned int matchType) const
 
void getElementsSharingNodes (const std::vector< stk::mesh::EntityId > nodeId, std::vector< stk::mesh::Entity > &elements) const
 get a set of elements sharing multiple nodes More...
 
void buildSubcells ()
 force the mesh to build subcells: edges and faces More...
 
std::size_t elementLocalId (stk::mesh::Entity elmt) const
 
std::size_t elementLocalId (stk::mesh::EntityId gid) const
 
stk::mesh::EntityId elementGlobalId (std::size_t lid) const
 
stk::mesh::EntityId elementGlobalId (stk::mesh::Entity elmt) const
 
bool isEdgeLocal (stk::mesh::Entity edge) const
 
bool isEdgeLocal (stk::mesh::EntityId gid) const
 
std::size_t edgeLocalId (stk::mesh::Entity elmt) const
 
std::size_t edgeLocalId (stk::mesh::EntityId gid) const
 
stk::mesh::EntityId edgeGlobalId (std::size_t lid) const
 
stk::mesh::EntityId edgeGlobalId (stk::mesh::Entity edge) const
 
bool isFaceLocal (stk::mesh::Entity face) const
 
bool isFaceLocal (stk::mesh::EntityId gid) const
 
std::size_t faceLocalId (stk::mesh::Entity elmt) const
 
std::size_t faceLocalId (stk::mesh::EntityId gid) const
 
stk::mesh::EntityId faceGlobalId (std::size_t lid) const
 
stk::mesh::EntityId faceGlobalId (stk::mesh::Entity face) const
 
unsigned entityOwnerRank (stk::mesh::Entity entity) const
 
bool isValid (stk::mesh::Entity entity) const
 
std::string containingBlockId (stk::mesh::Entity elmt) const
 
stk::mesh::Field< double > * getSolutionField (const std::string &fieldName, const std::string &blockId) const
 
stk::mesh::Field< double > * getCellField (const std::string &fieldName, const std::string &blockId) const
 
stk::mesh::Field< double > * getEdgeField (const std::string &fieldName, const std::string &blockId) const
 
stk::mesh::Field< double > * getFaceField (const std::string &fieldName, const std::string &blockId) const
 
ProcIdFieldTypegetProcessorIdField ()
 
bool isInitialized () const
 Has initialize been called on this mesh object? More...
 
Teuchos::RCP< const
std::vector< stk::mesh::Entity > > 
getElementsOrderedByLID () const
 
template<typename ArrayT >
void setSolutionFieldData (const std::string &fieldName, const std::string &blockId, const std::vector< std::size_t > &localElementIds, const ArrayT &solutionValues, double scaleValue=1.0)
 
template<typename ArrayT >
void getSolutionFieldData (const std::string &fieldName, const std::string &blockId, const std::vector< std::size_t > &localElementIds, ArrayT &solutionValues) const
 
template<typename ArrayT >
void setCellFieldData (const std::string &fieldName, const std::string &blockId, const std::vector< std::size_t > &localElementIds, const ArrayT &solutionValues, double scaleValue=1.0)
 
Teuchos::RCP< const
std::vector< stk::mesh::Entity > > 
getEdgesOrderedByLID () const
 
Teuchos::RCP< const
std::vector< stk::mesh::Entity > > 
getFacesOrderedByLID () const
 
template<typename ArrayT >
void setEdgeFieldData (const std::string &fieldName, const std::string &blockId, const std::vector< std::size_t > &localEdgeIds, const ArrayT &edgeValues, double scaleValue=1.0)
 
template<typename ArrayT >
void setFaceFieldData (const std::string &fieldName, const std::string &blockId, const std::vector< std::size_t > &localFaceIds, const ArrayT &faceValues, double scaleValue=1.0)
 
template<typename ArrayT >
void getElementVertices (const std::vector< std::size_t > &localIds, ArrayT &vertices) const
 
template<typename ArrayT >
void getElementVertices (const std::vector< stk::mesh::Entity > &elements, ArrayT &vertices) const
 
template<typename ArrayT >
void getElementVertices (const std::vector< std::size_t > &localIds, const std::string &eBlock, ArrayT &vertices) const
 
template<typename ArrayT >
void getElementVertices (const std::vector< stk::mesh::Entity > &elements, const std::string &eBlock, ArrayT &vertices) const
 
template<typename ArrayT >
void getElementVerticesNoResize (const std::vector< std::size_t > &localIds, ArrayT &vertices) const
 
template<typename ArrayT >
void getElementVerticesNoResize (const std::vector< stk::mesh::Entity > &elements, ArrayT &vertices) const
 
template<typename ArrayT >
void getElementVerticesNoResize (const std::vector< std::size_t > &localIds, const std::string &eBlock, ArrayT &vertices) const
 
template<typename ArrayT >
void getElementVerticesNoResize (const std::vector< stk::mesh::Entity > &elements, const std::string &eBlock, ArrayT &vertices) const
 
template<typename ArrayT >
void getElementNodes (const std::vector< std::size_t > &localIds, ArrayT &nodes) const
 
template<typename ArrayT >
void getElementNodes (const std::vector< stk::mesh::Entity > &elements, ArrayT &nodes) const
 
template<typename ArrayT >
void getElementNodes (const std::vector< std::size_t > &localIds, const std::string &eBlock, ArrayT &nodes) const
 
template<typename ArrayT >
void getElementNodes (const std::vector< stk::mesh::Entity > &elements, const std::string &eBlock, ArrayT &nodes) const
 
template<typename ArrayT >
void getElementNodesNoResize (const std::vector< std::size_t > &localIds, ArrayT &nodes) const
 
template<typename ArrayT >
void getElementNodesNoResize (const std::vector< stk::mesh::Entity > &elements, ArrayT &nodes) const
 
template<typename ArrayT >
void getElementNodesNoResize (const std::vector< std::size_t > &localIds, const std::string &eBlock, ArrayT &nodes) const
 
template<typename ArrayT >
void getElementNodesNoResize (const std::vector< stk::mesh::Entity > &elements, const std::string &eBlock, ArrayT &nodes) const
 
stk::mesh::EntityRank getElementRank () const
 
stk::mesh::EntityRank getSideRank () const
 
stk::mesh::EntityRank getFaceRank () const
 
stk::mesh::EntityRank getEdgeRank () const
 
stk::mesh::EntityRank getNodeRank () const
 
void initializeFromMetaData ()
 
void buildLocalElementIDs ()
 
void buildLocalEdgeIDs ()
 
void buildLocalFaceIDs ()
 
const std::vector
< Teuchos::RCP< const
PeriodicBC_MatcherBase > > & 
getPeriodicBCVector () const
 
std::vector< Teuchos::RCP
< const PeriodicBC_MatcherBase > > & 
getPeriodicBCVector ()
 
const bool & useBoundingBoxSearch () const
 
void setBoundingBoxSearchFlag (const bool &searchFlag)
 
void addPeriodicBC (const Teuchos::RCP< const PeriodicBC_MatcherBase > &bc)
 
void addPeriodicBCs (const std::vector< Teuchos::RCP< const PeriodicBC_MatcherBase > > &bc_vec)
 
std::pair< Teuchos::RCP
< std::vector< std::pair
< std::size_t, std::size_t >
> >, Teuchos::RCP
< std::vector< unsigned int > > > 
getPeriodicNodePairing () const
 
bool validBlockId (const std::string &blockId) const
 
void print (std::ostream &os) const
 
void printMetaData (std::ostream &os) const
 
Teuchos::RCP< const
shards::CellTopology > 
getCellTopology (const std::string &eBlock) const
 
double getCurrentStateTime () const
 
double getInitialStateTime () const
 
void setInitialStateTime (double value)
 
void rebalance (const Teuchos::ParameterList &params)
 
void setBlockWeight (const std::string &blockId, double weight)
 
void setUseFieldCoordinates (bool useFieldCoordinates)
 
bool getUseFieldCoordinates () const
 
void setUseLowerCaseForIO (bool useLowerCase)
 
bool getUseLowerCaseForIO () const
 
template<typename ArrayT >
void getElementVertices_FromField (const std::vector< stk::mesh::Entity > &elements, const std::string &eBlock, ArrayT &vertices) const
 
template<typename ArrayT >
void getElementVertices_FromFieldNoResize (const std::vector< stk::mesh::Entity > &elements, const std::string &eBlock, ArrayT &vertices) const
 
template<typename ArrayT >
void getElementVertices_FromCoords (const std::vector< stk::mesh::Entity > &elements, ArrayT &vertices) const
 
template<typename ArrayT >
void getElementVertices_FromCoordsNoResize (const std::vector< stk::mesh::Entity > &elements, ArrayT &vertices) const
 
template<typename ArrayT >
void getElementNodes_FromField (const std::vector< stk::mesh::Entity > &elements, const std::string &eBlock, ArrayT &nodes) const
 
template<typename ArrayT >
void getElementNodes_FromFieldNoResize (const std::vector< stk::mesh::Entity > &elements, const std::string &eBlock, ArrayT &nodes) const
 
template<typename ArrayT >
void getElementNodes_FromCoords (const std::vector< stk::mesh::Entity > &elements, ArrayT &nodes) const
 
template<typename ArrayT >
void getElementNodes_FromCoordsNoResize (const std::vector< stk::mesh::Entity > &elements, ArrayT &nodes) const
 
void refineMesh (const int numberOfLevels, const bool deleteParentElements)
 

Static Public Attributes

static const std::string coordsString = "coordinates"
 
static const std::string nodesString = "nodes"
 
static const std::string edgesString = "edges"
 
static const std::string edgeBlockString = "edge_block"
 
static const std::string faceBlockString = "face_block"
 
static const std::string facesString = "faces"
 

Protected Member Functions

void buildEntityCounts ()
 
void buildMaxEntityIds ()
 
void initializeFieldsInSTK (const std::map< std::pair< std::string, std::string >, SolutionFieldType * > &nameToField, bool setupIO)
 
Teuchos::RCP< Teuchos::MpiComm
< int > > 
getSafeCommunicator (stk::ParallelMachine parallelMach) const
 
void applyElementLoadBalanceWeights ()
 
bool isMeshCoordField (const std::string &eBlock, const std::string &fieldName, int &axis) const
 
template<typename ArrayT >
void setDispFieldData (const std::string &fieldName, const std::string &blockId, int axis, const std::vector< std::size_t > &localElementIds, const ArrayT &solutionValues)
 

Protected Attributes

std::vector< Teuchos::RCP
< const PeriodicBC_MatcherBase > > 
periodicBCs_
 
bool useBBoxSearch_ = false
 
Teuchos::RCP< stk::mesh::MetaData > metaData_
 
Teuchos::RCP< stk::mesh::BulkData > bulkData_
 
std::map< std::string,
stk::mesh::Part * > 
elementBlocks_
 
std::map< std::string,
stk::mesh::Part * > 
sidesets_
 
std::map< std::string,
stk::mesh::Part * > 
nodesets_
 
std::map< std::string,
stk::mesh::Part * > 
edgeBlocks_
 
std::map< std::string,
stk::mesh::Part * > 
faceBlocks_
 
std::map< std::string,
Teuchos::RCP
< shards::CellTopology > > 
elementBlockCT_
 
stk::mesh::Part * nodesPart_
 
std::vector< stk::mesh::Part * > nodesPartVec_
 
stk::mesh::Part * edgesPart_
 
std::vector< stk::mesh::Part * > edgesPartVec_
 
stk::mesh::Part * facesPart_
 
std::vector< stk::mesh::Part * > facesPartVec_
 
VectorFieldTypecoordinatesField_
 
VectorFieldTypeedgesField_
 
VectorFieldTypefacesField_
 
ProcIdFieldTypeprocessorIdField_
 
SolutionFieldTypeloadBalField_
 
std::map< std::pair
< std::string, std::string >
, SolutionFieldType * > 
fieldNameToSolution_
 
std::map< std::pair
< std::string, std::string >
, SolutionFieldType * > 
fieldNameToCellField_
 
std::map< std::pair
< std::string, std::string >
, SolutionFieldType * > 
fieldNameToEdgeField_
 
std::map< std::pair
< std::string, std::string >
, SolutionFieldType * > 
fieldNameToFaceField_
 
std::set< std::string > informationRecords_
 
unsigned dimension_
 
bool initialized_
 
std::vector< std::size_t > entityCounts_
 
std::vector< stk::mesh::EntityId > maxEntityId_
 
unsigned procRank_
 
std::size_t currentLocalId_
 
Teuchos::RCP< Teuchos::MpiComm
< int > > 
mpiComm_
 
double initialStateTime_
 
double currentStateTime_
 
Teuchos::RCP< std::vector
< stk::mesh::Entity > > 
orderedElementVector_
 
Teuchos::RCP< std::vector
< stk::mesh::Entity > > 
orderedEdgeVector_
 
Teuchos::RCP< std::vector
< stk::mesh::Entity > > 
orderedFaceVector_
 
std::map< std::string, double > blockWeights_
 
std::unordered_map
< stk::mesh::EntityId,
std::size_t > 
localIDHash_
 
std::unordered_map
< stk::mesh::EntityId,
std::size_t > 
localEdgeIDHash_
 
std::unordered_map
< stk::mesh::EntityId,
std::size_t > 
localFaceIDHash_
 
std::map< std::string,
std::vector< std::string > > 
meshCoordFields_
 
std::map< std::string,
std::vector< std::string > > 
meshDispFields_
 
bool useFieldCoordinates_
 
bool useLowerCase_
 

Detailed Description

Definition at line 102 of file Panzer_STK_Interface.hpp.

Member Typedef Documentation

Definition at line 104 of file Panzer_STK_Interface.hpp.

typedef stk::mesh::Field<double> panzer_stk::STK_Interface::SolutionFieldType

Definition at line 105 of file Panzer_STK_Interface.hpp.

typedef stk::mesh::Field<double> panzer_stk::STK_Interface::VectorFieldType

Definition at line 106 of file Panzer_STK_Interface.hpp.

Definition at line 107 of file Panzer_STK_Interface.hpp.

Constructor & Destructor Documentation

panzer_stk::STK_Interface::STK_Interface ( )

Definition at line 108 of file Panzer_STK_Interface.cpp.

panzer_stk::STK_Interface::STK_Interface ( unsigned  dim)

Default constructor

Definition at line 122 of file Panzer_STK_Interface.cpp.

panzer_stk::STK_Interface::STK_Interface ( Teuchos::RCP< stk::mesh::MetaData >  metaData)

Definition at line 115 of file Panzer_STK_Interface.cpp.

Member Function Documentation

void panzer_stk::STK_Interface::addElementBlock ( const std::string &  name,
const CellTopologyData *  ctData 
)

Add an element block with a string name

Definition at line 1612 of file Panzer_STK_Interface.cpp.

void panzer_stk::STK_Interface::addEdgeBlock ( const std::string &  elemBlockName,
const std::string &  edgeBlockName,
const stk::topology &  topology 
)

Add an edge block with a string name

Definition at line 1643 of file Panzer_STK_Interface.cpp.

void panzer_stk::STK_Interface::addEdgeBlock ( const std::string &  elemBlockName,
const std::string &  edgeBlockName,
const CellTopologyData *  ctData 
)

Definition at line 1664 of file Panzer_STK_Interface.cpp.

void panzer_stk::STK_Interface::addFaceBlock ( const std::string &  elemBlockName,
const std::string &  faceBlockName,
const stk::topology &  topology 
)

Add a face block with a string name

Definition at line 1688 of file Panzer_STK_Interface.cpp.

void panzer_stk::STK_Interface::addFaceBlock ( const std::string &  elemBlockName,
const std::string &  faceBlockName,
const CellTopologyData *  ctData 
)

Definition at line 1709 of file Panzer_STK_Interface.cpp.

void panzer_stk::STK_Interface::addSideset ( const std::string &  name,
const CellTopologyData *  ctData 
)

Add a side set with a string name

Definition at line 134 of file Panzer_STK_Interface.cpp.

void panzer_stk::STK_Interface::addNodeset ( const std::string &  name)

Add a node set with a string name

Definition at line 146 of file Panzer_STK_Interface.cpp.

void panzer_stk::STK_Interface::addSolutionField ( const std::string &  fieldName,
const std::string &  blockId 
)

Add a solution field

Definition at line 160 of file Panzer_STK_Interface.cpp.

void panzer_stk::STK_Interface::addCellField ( const std::string &  fieldName,
const std::string &  blockId 
)

Add a solution field

Definition at line 179 of file Panzer_STK_Interface.cpp.

void panzer_stk::STK_Interface::addEdgeField ( const std::string &  fieldName,
const std::string &  blockId 
)

Add an edge field

Definition at line 199 of file Panzer_STK_Interface.cpp.

void panzer_stk::STK_Interface::addFaceField ( const std::string &  fieldName,
const std::string &  blockId 
)

Add a face field

Definition at line 220 of file Panzer_STK_Interface.cpp.

void panzer_stk::STK_Interface::addMeshCoordFields ( const std::string &  blockId,
const std::vector< std::string > &  coordField,
const std::string &  dispPrefix 
)

Add a solution field for coordinates with a particular prefix, force it to be outputed as a to be mesh displacement field. This is really only relevant for I/O and how the field is stored internally in the mesh.

Parameters
[in]blockIdElement block to use
[in]coordPrefixPrefix for solution coordinate field (assumes using "X","Y","Z" as postfix)
[in]dispPrefixPrefix for displacment (output) field

Definition at line 241 of file Panzer_STK_Interface.cpp.

void panzer_stk::STK_Interface::addInformationRecords ( const std::vector< std::string > &  info_records)

Add a vector of strings to the Information Records block. Each string will be it's own record. The info records will be deduped before they are added to IOSS.

Definition at line 285 of file Panzer_STK_Interface.cpp.

void panzer_stk::STK_Interface::initialize ( stk::ParallelMachine  parallelMach,
bool  setupIO = true,
const bool  buildRefinementSupport = false 
)

Initialize the mesh with the current dimension This also calls commit on the meta data causing it to be frozen. Information about elements blocks has to be commited before this. If parallel machine has already been specified through instantiateBulkData that communicator is used. Otherwise a new copy is constructed and will be used through out this mesh object's lifetime.

Parameters
[in]parallelMachCommunicator
[in]setupIOIf set to true and IOSS is enabled, the output mesh will be initialized.
[in]buildRefinementSupportIf true, build percept uniform refinement objects.

Definition at line 290 of file Panzer_STK_Interface.cpp.

void panzer_stk::STK_Interface::instantiateBulkData ( stk::ParallelMachine  parallelMach)

Build a bulk data object but don't do anything with it. If parallel machine has already been specified through initialize that communicator is used. Otherwise a new copy is constructed and will be used by default throughout the lifetime of this object.

Definition at line 422 of file Panzer_STK_Interface.cpp.

void panzer_stk::STK_Interface::beginModification ( )

Put the bulk data manager in modification mode.

Definition at line 432 of file Panzer_STK_Interface.cpp.

void panzer_stk::STK_Interface::endModification ( )

Take the bulk data manager out of modification mode.

Definition at line 440 of file Panzer_STK_Interface.cpp.

void panzer_stk::STK_Interface::addNode ( stk::mesh::EntityId  gid,
const std::vector< double > &  coord 
)

Add a node to the mesh with a specific set of coordinates to the mesh.

Precondition
coord.size()==getDimension()
isModifiable()==true

Definition at line 496 of file Panzer_STK_Interface.cpp.

void panzer_stk::STK_Interface::addElement ( const Teuchos::RCP< ElementDescriptor > &  ed,
stk::mesh::Part *  block 
)

Definition at line 564 of file Panzer_STK_Interface.cpp.

void panzer_stk::STK_Interface::addEdges ( )

Definition at line 587 of file Panzer_STK_Interface.cpp.

void panzer_stk::STK_Interface::addFaces ( )

Definition at line 615 of file Panzer_STK_Interface.cpp.

void panzer_stk::STK_Interface::addEntityToSideset ( stk::mesh::Entity  entity,
stk::mesh::Part *  sideset 
)

Adds an entity to a specified side set.

Definition at line 514 of file Panzer_STK_Interface.cpp.

void panzer_stk::STK_Interface::addEntityToNodeset ( stk::mesh::Entity  entity,
stk::mesh::Part *  nodeset 
)

Adds an entity to a specified node set.

Definition at line 522 of file Panzer_STK_Interface.cpp.

void panzer_stk::STK_Interface::addEntityToEdgeBlock ( stk::mesh::Entity  entity,
stk::mesh::Part *  edgeblock 
)

Adds an entity to a specified edge block.

Definition at line 530 of file Panzer_STK_Interface.cpp.

void panzer_stk::STK_Interface::addEntitiesToEdgeBlock ( std::vector< stk::mesh::Entity >  entities,
stk::mesh::Part *  edgeblock 
)

Adds a vector of entities to a specified edge block.

Definition at line 537 of file Panzer_STK_Interface.cpp.

void panzer_stk::STK_Interface::addEntityToFaceBlock ( stk::mesh::Entity  entity,
stk::mesh::Part *  faceblock 
)

Adds an entity to a specified face block.

Definition at line 547 of file Panzer_STK_Interface.cpp.

void panzer_stk::STK_Interface::addEntitiesToFaceBlock ( std::vector< stk::mesh::Entity >  entities,
stk::mesh::Part *  faceblock 
)

Adds a vector of entities to a specified face block.

Definition at line 554 of file Panzer_STK_Interface.cpp.

const VectorFieldType& panzer_stk::STK_Interface::getCoordinatesField ( ) const
inline

Grab the coordinates field

Definition at line 272 of file Panzer_STK_Interface.hpp.

const VectorFieldType& panzer_stk::STK_Interface::getEdgesField ( ) const
inline

Grab the edges field

Definition at line 277 of file Panzer_STK_Interface.hpp.

const VectorFieldType& panzer_stk::STK_Interface::getFacesField ( ) const
inline

Definition at line 280 of file Panzer_STK_Interface.hpp.

const double * panzer_stk::STK_Interface::getNodeCoordinates ( stk::mesh::EntityId  nodeId) const

Look up a global node and get the coordinate.

Definition at line 1116 of file Panzer_STK_Interface.cpp.

const double * panzer_stk::STK_Interface::getNodeCoordinates ( stk::mesh::Entity  node) const

Look up a global node and get the coordinate.

Definition at line 1122 of file Panzer_STK_Interface.cpp.

void panzer_stk::STK_Interface::getSubcellIndices ( unsigned  entityRank,
stk::mesh::EntityId  elementId,
std::vector< stk::mesh::EntityId > &  subcellIds 
) const

Get subcell global IDs

Definition at line 1127 of file Panzer_STK_Interface.cpp.

void panzer_stk::STK_Interface::getMyElements ( std::vector< stk::mesh::Entity > &  elements) const

Get a vector of elements owned by this processor

Definition at line 1148 of file Panzer_STK_Interface.cpp.

void panzer_stk::STK_Interface::getMyElements ( const std::string &  blockID,
std::vector< stk::mesh::Entity > &  elements 
) const

Get a vector of elements owned by this processor on a particular block ID

Definition at line 1158 of file Panzer_STK_Interface.cpp.

void panzer_stk::STK_Interface::getNeighborElements ( std::vector< stk::mesh::Entity > &  elements) const

Get a vector of elements that share an edge/face with an owned element. Note that these elements are not owned.

Definition at line 1173 of file Panzer_STK_Interface.cpp.

void panzer_stk::STK_Interface::getNeighborElements ( const std::string &  blockID,
std::vector< stk::mesh::Entity > &  elements 
) const

Get a vector of elements not owned by this processor but in a particular block

Definition at line 1183 of file Panzer_STK_Interface.cpp.

void panzer_stk::STK_Interface::getMyEdges ( std::vector< stk::mesh::Entity > &  edges) const

Get a vector of edges owned by this processor

Definition at line 1197 of file Panzer_STK_Interface.cpp.

void panzer_stk::STK_Interface::getMyEdges ( const std::string &  edgeBlockName,
std::vector< stk::mesh::Entity > &  edges 
) const

Get Entities corresponding to the edge block requested. The Entites in the vector should be a dimension lower than getDimension().

Parameters
[in]edgeBlockNameName of edge block
[in,out]edgesVector of entities containing the requested edges.

Definition at line 1206 of file Panzer_STK_Interface.cpp.

void panzer_stk::STK_Interface::getMyEdges ( const std::string &  edgeBlockName,
const std::string &  blockName,
std::vector< stk::mesh::Entity > &  edges 
) const

Get Entities corresponding to the locally owned part of the edge block requested. This also limits the entities to be in a particular element block. The Entites in the vector should be a dimension lower than getDimension().

Parameters
[in]edgeBlockNameName of edge block
[in]blockNameName of block
[in,out]edgesVector of entities containing the requested edges.

Definition at line 1219 of file Panzer_STK_Interface.cpp.

void panzer_stk::STK_Interface::getAllEdges ( const std::string &  edgeBlockName,
std::vector< stk::mesh::Entity > &  edges 
) const

Get Entities corresponding to the locally owned part of the edge block requested. The Entites in the vector should be a dimension lower than getDimension().

Parameters
[in]edgeBlockNameName of edge block
[in,out]edgesVector of entities containing the requested edges.

Definition at line 1236 of file Panzer_STK_Interface.cpp.

void panzer_stk::STK_Interface::getAllEdges ( const std::string &  edgeBlockName,
const std::string &  blockName,
std::vector< stk::mesh::Entity > &  edges 
) const

Get Entities corresponding to the edge block requested. This also limits the entities to be in a particular element block. The Entites in the vector should be a dimension lower than getDimension().

Parameters
[in]edgeBlockNameName of edge block
[in]blockNameName of block
[in,out]edgesVector of entities containing the requested edges.

Definition at line 1248 of file Panzer_STK_Interface.cpp.

void panzer_stk::STK_Interface::getMyFaces ( std::vector< stk::mesh::Entity > &  faces) const

Get a vector of faces owned by this processor

Definition at line 1265 of file Panzer_STK_Interface.cpp.

void panzer_stk::STK_Interface::getMyFaces ( const std::string &  faceBlockName,
std::vector< stk::mesh::Entity > &  faces 
) const

Get Entities corresponding to the face block requested. The Entites in the vector should be a dimension lower than getDimension().

Parameters
[in]faceBlockNameName of face block
[in,out]facesVector of entities containing the requested faces.

Definition at line 1275 of file Panzer_STK_Interface.cpp.

void panzer_stk::STK_Interface::getMyFaces ( const std::string &  faceBlockName,
const std::string &  blockName,
std::vector< stk::mesh::Entity > &  faces 
) const

Get Entities corresponding to the locally owned part of the face block requested. This also limits the entities to be in a particular element block. The Entites in the vector should be a dimension lower than getDimension().

Parameters
[in]faceBlockNameName of face block
[in]blockNameName of block
[in,out]facesVector of entities containing the requested faces.

Definition at line 1288 of file Panzer_STK_Interface.cpp.

void panzer_stk::STK_Interface::getAllFaces ( const std::string &  faceBlockName,
std::vector< stk::mesh::Entity > &  faces 
) const

Get Entities corresponding to the locally owned part of the face block requested. The Entites in the vector should be a dimension lower than getDimension().

Parameters
[in]faceBlockNameName of face block
[in,out]facesVector of entities containing the requested faces.

Definition at line 1305 of file Panzer_STK_Interface.cpp.

void panzer_stk::STK_Interface::getAllFaces ( const std::string &  faceBlockName,
const std::string &  blockName,
std::vector< stk::mesh::Entity > &  faces 
) const

Get Entities corresponding to the face block requested. This also limits the entities to be in a particular element block. The Entites in the vector should be a dimension lower than getDimension().

Parameters
[in]faceBlockNameName of face block
[in]blockNameName of block
[in,out]facesVector of entities containing the requested faces.

Definition at line 1317 of file Panzer_STK_Interface.cpp.

void panzer_stk::STK_Interface::getMySides ( const std::string &  sideName,
std::vector< stk::mesh::Entity > &  sides 
) const

Get Entities corresponding to the side set requested. The Entites in the vector should be a dimension lower then getDimension().

Parameters
[in]sideNameName of side set
[in,out]sidesVector of entities containing the requested sides.

Definition at line 1334 of file Panzer_STK_Interface.cpp.

void panzer_stk::STK_Interface::getMySides ( const std::string &  sideName,
const std::string &  blockName,
std::vector< stk::mesh::Entity > &  sides 
) const

Get Entities corresponding to the locally owned part of the side set requested. This also limits the entities to be in a particular element block. The Entites in the vector should be a dimension lower then getDimension().

Parameters
[in]sideNameName of side set
[in]blockNameName of block
[in,out]sidesVector of entities containing the requested sides.

Definition at line 1347 of file Panzer_STK_Interface.cpp.

void panzer_stk::STK_Interface::getAllSides ( const std::string &  sideName,
std::vector< stk::mesh::Entity > &  sides 
) const

Get Entities corresponding to the locally owned part of the side set requested. The Entites in the vector should be a dimension lower then getDimension().

Parameters
[in]sideNameName of side set
[in,out]sidesVector of entities containing the requested sides.

Definition at line 1364 of file Panzer_STK_Interface.cpp.

void panzer_stk::STK_Interface::getAllSides ( const std::string &  sideName,
const std::string &  blockName,
std::vector< stk::mesh::Entity > &  sides 
) const

Get Entities corresponding to the side set requested. This also limits the entities to be in a particular element block. The Entites in the vector should be a dimension lower then getDimension().

Parameters
[in]sideNameName of side set
[in]blockNameName of block
[in,out]sidesVector of entities containing the requested sides.

Definition at line 1376 of file Panzer_STK_Interface.cpp.

void panzer_stk::STK_Interface::getMyNodes ( const std::string &  sideName,
const std::string &  blockName,
std::vector< stk::mesh::Entity > &  nodes 
) const

Get Entities corresponding to the node set requested. This also limits the entities to be in a particular element block. The Entites in the vector should be of dimension 0.

Parameters
[in]sideNameName of side set
[in]blockNameName of block
[in,out]nodesVector of entities containing the requested nodes.

Definition at line 1393 of file Panzer_STK_Interface.cpp.

stk::mesh::Entity panzer_stk::STK_Interface::findConnectivityById ( stk::mesh::Entity  src,
stk::mesh::EntityRank  tgt_rank,
unsigned  rel_id 
) const

Searches for connected entity by rank and relation id. Returns invalid entity on failure.

Parameters
[in]srcThe handle to the source entity (the 'from' part of the relation)
[in]tgt_rankThe entity rank of relations to search
[in]rel_idThe id of the relation to search for

Definition at line 645 of file Panzer_STK_Interface.cpp.

void panzer_stk::STK_Interface::writeToExodus ( const std::string &  filename,
const bool  append = false 
)

Write this mesh and associated fields to the given output file.

Note
This will only write a single timestep at time = 0.
Parameters
[in]filenameThe name of the output Exodus file.
[in]appendIf set to true, the output will be appended to the output Exodus file. If set to false, output file will be overwritten. Default is false.

Definition at line 666 of file Panzer_STK_Interface.cpp.

void panzer_stk::STK_Interface::setupExodusFile ( const std::string &  filename,
const bool  append = false,
const bool  append_after_restart_time = false,
const double  restart_time = 0.0 
)

Set up an output Exodus file for writing results.

Create an output mesh associated with the given filename and add the fields to it.

Note
No writing is actually done at this stage. That happens with a call to writeToExodus(double timestep).
Parameters
[in]filenameThe name of the output Exodus file.
[in]appendIf set to true, the output will be appended to the output Exodus file. If set to false, output file will be overwritten. Default is false.
[in]append_after_restart_timeIf set to true, instead of appending to the end of the Exodus file, this option will append new writes after a specified time and overwrite subsequent time steps that were in the file. Allows users to restart anywhere in the exodus file and write consistently. If set to false, the next write will append to the end of the file. Default is false.
[in]restart_timeIf append_after_restart_time is true, this is the time value to append after. Otherwise this value is ignored.
Exceptions
<tt>std::logic_error</tt>If the STK_Interface does not yet have a MPI communicator.

Definition at line 680 of file Panzer_STK_Interface.cpp.

void panzer_stk::STK_Interface::writeToExodus ( double  timestep)

Write this mesh and associated fields at the given timestep.

Write this mesh and the current state of the associated fields for time = timestep to the output Exodus file created via setupExodusFile().

Note
setupExodusFile(const std::string& filename) must be called before any invocations of this routine.
Parameters
[in]timestepThe simulation time at which you're writing the results.

Definition at line 740 of file Panzer_STK_Interface.cpp.

void panzer_stk::STK_Interface::addGlobalToExodus ( const std::string &  key,
const int &  value 
)

Add an int global variable to the information to be written to the Exodus output file.

This allows you to add any arbitrary global data (that is, data not associated with nodes, elements, etc.) to the Exodus output file tagged with a string key.

Note
Data is not actually written to the output Exodus file until writeToExodus() is called.
Parameters
[in]keyThe name of the global variable you'll be adding. Note that keys should be unique if you're adding multiple global variables with repeated calls to addGlobalToExodus(). Repeated calls with identical keys will result in the value associated with the key simply being overwritten.
[in]valueThe value of the global variable you'll be adding.

Definition at line 891 of file Panzer_STK_Interface.cpp.

void panzer_stk::STK_Interface::addGlobalToExodus ( const std::string &  key,
const double &  value 
)

Add a double global variable to the information to be written to the Exodus output file.

This allows you to add any arbitrary global data (that is, data not associated with nodes, elements, etc.) to the Exodus output file tagged with a string key.

Note
Data is not actually written to the output Exodus file until writeToExodus() is called.
Parameters
[in]keyThe name of the global variable you'll be adding. Note that keys should be unique if you're adding multiple global variables with repeated calls to addGlobalToExodus(). Repeated calls with identical keys will result in the value associated with the key simply being overwritten.
[in]valueThe value of the global variable you'll be adding.

Definition at line 905 of file Panzer_STK_Interface.cpp.

void panzer_stk::STK_Interface::addGlobalToExodus ( const std::string &  key,
const std::vector< int > &  value 
)

Add a std::vector<int> global variable to the information to be written to the Exodus output file.

This allows you to add any arbitrary global data (that is, data not associated with nodes, elements, etc.) to the Exodus output file tagged with a string key.

Note
Data is not actually written to the output Exodus file until writeToExodus() is called.
Parameters
[in]keyThe name of the global variable you'll be adding. Note that keys should be unique if you're adding multiple global variables with repeated calls to addGlobalToExodus(). Repeated calls with identical keys will result in the value associated with the key simply being overwritten.
[in]valueThe value of the global variable you'll be adding.

Definition at line 919 of file Panzer_STK_Interface.cpp.

void panzer_stk::STK_Interface::addGlobalToExodus ( const std::string &  key,
const std::vector< double > &  value 
)

Add a std::vector<double> global variable to the information to be written to the Exodus output file.

This allows you to add any arbitrary global data (that is, data not associated with nodes, elements, etc.) to the Exodus output file tagged with a string key.

Note
Data is not actually written to the output Exodus file until writeToExodus() is called.
Parameters
[in]keyThe name of the global variable you'll be adding. Note that keys should be unique if you're adding multiple global variables with repeated calls to addGlobalToExodus(). Repeated calls with identical keys will result in the value associated with the key simply being overwritten.
[in]valueThe value of the global variable you'll be adding.

Definition at line 934 of file Panzer_STK_Interface.cpp.

Teuchos::RCP< const Teuchos::Comm< int > > panzer_stk::STK_Interface::getComm ( ) const

get the comm associated with this mesh

Definition at line 2063 of file Panzer_STK_Interface.cpp.

Teuchos::RCP<stk::mesh::BulkData> panzer_stk::STK_Interface::getBulkData ( ) const
inline

Definition at line 612 of file Panzer_STK_Interface.hpp.

Teuchos::RCP<stk::mesh::MetaData> panzer_stk::STK_Interface::getMetaData ( ) const
inline

Definition at line 613 of file Panzer_STK_Interface.hpp.

bool panzer_stk::STK_Interface::isWritable ( ) const

Definition at line 942 of file Panzer_STK_Interface.cpp.

bool panzer_stk::STK_Interface::isModifiable ( ) const
inline

Definition at line 623 of file Panzer_STK_Interface.hpp.

unsigned panzer_stk::STK_Interface::getDimension ( ) const
inline

get the dimension

Definition at line 628 of file Panzer_STK_Interface.hpp.

std::size_t panzer_stk::STK_Interface::getNumElementBlocks ( ) const
inline

get the block count

Definition at line 632 of file Panzer_STK_Interface.hpp.

void panzer_stk::STK_Interface::getElementBlockNames ( std::vector< std::string > &  names) const

Get a vector containing the names of the element blocks. This function always returns the current set of element blocks in lexiographic order (uses the sorting built into the std::map). This method can only be called after initialize.

Parameters
[in,out]namesVector of names of the element blocks.

Definition at line 1410 of file Panzer_STK_Interface.cpp.

void panzer_stk::STK_Interface::getSidesetNames ( std::vector< std::string > &  name) const

Get a vector containing the names of the side sets. This function always returns the current set of side sets in lexiographic order (uses the sorting built into the std::map). This method can only be called after initialize.

Parameters
[in,out]namesVector of names of the side sets

Definition at line 1422 of file Panzer_STK_Interface.cpp.

void panzer_stk::STK_Interface::getNodesetNames ( std::vector< std::string > &  name) const

Get a vector containing the names of the node sets. This function always returns the current set of node sets in lexiographic order (uses the sorting built into the std::map). This method can only be called after initialize.

Parameters
[in,out]namesVector of names of the node sets.

Definition at line 1434 of file Panzer_STK_Interface.cpp.

void panzer_stk::STK_Interface::getEdgeBlockNames ( std::vector< std::string > &  names) const

Get a vector containing the names of the edge blocks. This function always returns the current set of edge blocks in lexiographic order (uses the sorting built into the std::map). This method can only be called after initialize.

Parameters
[in,out]namesVector of names of the edge blocks.

Definition at line 1444 of file Panzer_STK_Interface.cpp.

void panzer_stk::STK_Interface::getFaceBlockNames ( std::vector< std::string > &  names) const

Get a vector containing the names of the face blocks. This function always returns the current set of face blocks in lexiographic order (uses the sorting built into the std::map). This method can only be called after initialize.

Parameters
[in,out]namesVector of names of the face blocks.

Definition at line 1454 of file Panzer_STK_Interface.cpp.

stk::mesh::Part* panzer_stk::STK_Interface::getOwnedPart ( ) const
inline

Get a pointer to the locally owned part.

Definition at line 681 of file Panzer_STK_Interface.hpp.

stk::mesh::Part* panzer_stk::STK_Interface::getElementBlockPart ( const std::string &  name) const
inline

get the block part

Definition at line 685 of file Panzer_STK_Interface.hpp.

stk::mesh::Part* panzer_stk::STK_Interface::getEdgeBlock ( const std::string &  name) const
inline

get the block part

Definition at line 693 of file Panzer_STK_Interface.hpp.

stk::mesh::Part* panzer_stk::STK_Interface::getFaceBlock ( const std::string &  name) const
inline

get the block part

Definition at line 701 of file Panzer_STK_Interface.hpp.

std::size_t panzer_stk::STK_Interface::getNumSidesets ( ) const
inline

get the side set count

Definition at line 709 of file Panzer_STK_Interface.hpp.

stk::mesh::Part* panzer_stk::STK_Interface::getSideset ( const std::string &  name) const
inline

Definition at line 712 of file Panzer_STK_Interface.hpp.

std::size_t panzer_stk::STK_Interface::getNumNodesets ( ) const
inline

get the side set count

Definition at line 719 of file Panzer_STK_Interface.hpp.

stk::mesh::Part* panzer_stk::STK_Interface::getNodeset ( const std::string &  name) const
inline

Definition at line 722 of file Panzer_STK_Interface.hpp.

std::size_t panzer_stk::STK_Interface::getEntityCounts ( unsigned  entityRank) const

get the global counts for the entity of specified rank

Definition at line 1087 of file Panzer_STK_Interface.cpp.

stk::mesh::EntityId panzer_stk::STK_Interface::getMaxEntityId ( unsigned  entityRank) const

get max entity ID of type entityRank

Definition at line 1095 of file Panzer_STK_Interface.cpp.

void panzer_stk::STK_Interface::getElementsSharingNode ( stk::mesh::EntityId  nodeId,
std::vector< stk::mesh::Entity > &  elements 
) const

get a set of elements sharing a single node

Definition at line 951 of file Panzer_STK_Interface.cpp.

void panzer_stk::STK_Interface::getNodeIdsForElement ( stk::mesh::Entity  element,
std::vector< stk::mesh::EntityId > &  nodeIds 
) const

get a list of node ids for nodes connected to an element

Definition at line 1034 of file Panzer_STK_Interface.cpp.

void panzer_stk::STK_Interface::getOwnedElementsSharingNode ( stk::mesh::Entity  node,
std::vector< stk::mesh::Entity > &  elements,
std::vector< int > &  relIds 
) const

Get set of element sharing a single node and its local node id.

Definition at line 964 of file Panzer_STK_Interface.cpp.

void panzer_stk::STK_Interface::getOwnedElementsSharingNode ( stk::mesh::EntityId  nodeId,
std::vector< stk::mesh::Entity > &  elements,
std::vector< int > &  relIds,
unsigned int  matchType 
) const

Get set of element sharing a single node and its local node id.

Definition at line 984 of file Panzer_STK_Interface.cpp.

void panzer_stk::STK_Interface::getElementsSharingNodes ( const std::vector< stk::mesh::EntityId >  nodeId,
std::vector< stk::mesh::Entity > &  elements 
) const

get a set of elements sharing multiple nodes

Definition at line 1002 of file Panzer_STK_Interface.cpp.

void panzer_stk::STK_Interface::buildSubcells ( )

force the mesh to build subcells: edges and faces

Definition at line 1103 of file Panzer_STK_Interface.cpp.

std::size_t panzer_stk::STK_Interface::elementLocalId ( stk::mesh::Entity  elmt) const

Get an elements local index

Definition at line 1464 of file Panzer_STK_Interface.cpp.

std::size_t panzer_stk::STK_Interface::elementLocalId ( stk::mesh::EntityId  gid) const

Get an elements local index

Definition at line 1471 of file Panzer_STK_Interface.cpp.

stk::mesh::EntityId panzer_stk::STK_Interface::elementGlobalId ( std::size_t  lid) const
inline

Get an elements global index

Definition at line 770 of file Panzer_STK_Interface.hpp.

stk::mesh::EntityId panzer_stk::STK_Interface::elementGlobalId ( stk::mesh::Entity  elmt) const
inline

Get an elements global index

Definition at line 775 of file Panzer_STK_Interface.hpp.

bool panzer_stk::STK_Interface::isEdgeLocal ( stk::mesh::Entity  edge) const

Is an edge local to this processor?

Definition at line 1482 of file Panzer_STK_Interface.cpp.

bool panzer_stk::STK_Interface::isEdgeLocal ( stk::mesh::EntityId  gid) const

Is an edge local to this processor?

Definition at line 1487 of file Panzer_STK_Interface.cpp.

std::size_t panzer_stk::STK_Interface::edgeLocalId ( stk::mesh::Entity  elmt) const

Get an edge's local index

Definition at line 1496 of file Panzer_STK_Interface.cpp.

std::size_t panzer_stk::STK_Interface::edgeLocalId ( stk::mesh::EntityId  gid) const

Get an edge's local index

Definition at line 1501 of file Panzer_STK_Interface.cpp.

stk::mesh::EntityId panzer_stk::STK_Interface::edgeGlobalId ( std::size_t  lid) const
inline

Get an edge's global index

Definition at line 796 of file Panzer_STK_Interface.hpp.

stk::mesh::EntityId panzer_stk::STK_Interface::edgeGlobalId ( stk::mesh::Entity  edge) const
inline

Get an edge's global index

Definition at line 801 of file Panzer_STK_Interface.hpp.

bool panzer_stk::STK_Interface::isFaceLocal ( stk::mesh::Entity  face) const

Is a face local to this processor?

Definition at line 1508 of file Panzer_STK_Interface.cpp.

bool panzer_stk::STK_Interface::isFaceLocal ( stk::mesh::EntityId  gid) const

Is a face local to this processor?

Definition at line 1513 of file Panzer_STK_Interface.cpp.

std::size_t panzer_stk::STK_Interface::faceLocalId ( stk::mesh::Entity  elmt) const

Get a face's local index

Definition at line 1522 of file Panzer_STK_Interface.cpp.

std::size_t panzer_stk::STK_Interface::faceLocalId ( stk::mesh::EntityId  gid) const

Get a face's local index

Definition at line 1527 of file Panzer_STK_Interface.cpp.

stk::mesh::EntityId panzer_stk::STK_Interface::faceGlobalId ( std::size_t  lid) const
inline

Get a face's global index

Definition at line 822 of file Panzer_STK_Interface.hpp.

stk::mesh::EntityId panzer_stk::STK_Interface::faceGlobalId ( stk::mesh::Entity  face) const
inline

Get a face's global index

Definition at line 827 of file Panzer_STK_Interface.hpp.

unsigned panzer_stk::STK_Interface::entityOwnerRank ( stk::mesh::Entity  entity) const
inline

Get an Entity's parallel owner (process rank)

Definition at line 832 of file Panzer_STK_Interface.hpp.

bool panzer_stk::STK_Interface::isValid ( stk::mesh::Entity  entity) const
inline

Check if entity handle is valid

Definition at line 837 of file Panzer_STK_Interface.hpp.

std::string panzer_stk::STK_Interface::containingBlockId ( stk::mesh::Entity  elmt) const

Get the containing block ID of this element.

Definition at line 1535 of file Panzer_STK_Interface.cpp.

stk::mesh::Field< double > * panzer_stk::STK_Interface::getSolutionField ( const std::string &  fieldName,
const std::string &  blockId 
) const

Get the stk mesh field pointer associated with a particular solution value Assumes there is a field associated with "fieldName,blockId" pair. If none is found an exception (std::runtime_error) is raised.

Definition at line 1543 of file Panzer_STK_Interface.cpp.

stk::mesh::Field< double > * panzer_stk::STK_Interface::getCellField ( const std::string &  fieldName,
const std::string &  blockId 
) const

Get the stk mesh field pointer associated with a particular value Assumes there is a field associated with "fieldName,blockId" pair. If none is found an exception (std::runtime_error) is raised.

Definition at line 1557 of file Panzer_STK_Interface.cpp.

stk::mesh::Field< double > * panzer_stk::STK_Interface::getEdgeField ( const std::string &  fieldName,
const std::string &  blockId 
) const

Get the stk mesh field pointer associated with a particular value Assumes there is a field associated with "fieldName,blockId" pair. If none is found an exception (std::runtime_error) is raised.

Definition at line 1571 of file Panzer_STK_Interface.cpp.

stk::mesh::Field< double > * panzer_stk::STK_Interface::getFaceField ( const std::string &  fieldName,
const std::string &  blockId 
) const

Get the stk mesh field pointer associated with a particular value Assumes there is a field associated with "fieldName,blockId" pair. If none is found an exception (std::runtime_error) is raised.

Definition at line 1585 of file Panzer_STK_Interface.cpp.

ProcIdFieldType* panzer_stk::STK_Interface::getProcessorIdField ( )
inline

Definition at line 872 of file Panzer_STK_Interface.hpp.

bool panzer_stk::STK_Interface::isInitialized ( ) const
inline

Has initialize been called on this mesh object?

Definition at line 875 of file Panzer_STK_Interface.hpp.

Teuchos::RCP< const std::vector< stk::mesh::Entity > > panzer_stk::STK_Interface::getElementsOrderedByLID ( ) const

Get Vector of element entities ordered by their LID, returns an RCP so that it is easily stored by the caller.

Definition at line 1599 of file Panzer_STK_Interface.cpp.

template<typename ArrayT >
void panzer_stk::STK_Interface::setSolutionFieldData ( const std::string &  fieldName,
const std::string &  blockId,
const std::vector< std::size_t > &  localElementIds,
const ArrayT &  solutionValues,
double  scaleValue = 1.0 
)

Writes a particular field to an array. Notice this is setup to work with the worksets associated with Panzer.

Parameters
[in]fieldNameName of field to be filled
[in]blockIdName of block this set of elements belongs to
[in]localElementIdsLocal element IDs for this set of solution values
[in]solutionValuesAn two dimensional array object sized by (Cells,Basis Count)
Note
The block ID is not strictly needed in this context. However forcing the user to provide it does permit an additional level of safety. The implicit assumption is that the elements being "set" are part of the specified block. This prevents the need to perform a null pointer check on the field data, because the STK_Interface construction of the fields should force it to be nonnull...

Definition at line 1582 of file Panzer_STK_Interface.hpp.

template<typename ArrayT >
void panzer_stk::STK_Interface::getSolutionFieldData ( const std::string &  fieldName,
const std::string &  blockId,
const std::vector< std::size_t > &  localElementIds,
ArrayT &  solutionValues 
) const

Reads a particular field into an array. Notice this is setup to work with the worksets associated with Panzer.

Parameters
[in]fieldNameName of field to be filled
[in]blockIdName of block this set of elements belongs to
[in]localElementIdsLocal element IDs for this set of solution values
[in]solutionValuesAn two dimensional array object sized by (Cells,Basis Count)
Note
The block ID is not strictly needed in this context. However forcing the user to provide it does permit an additional level of safety. The implicit assumption is that the elements being retrieved are part of the specified block. This prevents the need to perform a null pointer check on the field data, because the STK_Interface construction of the fields should force it to be nonnull...

Definition at line 1644 of file Panzer_STK_Interface.hpp.

template<typename ArrayT >
void panzer_stk::STK_Interface::setCellFieldData ( const std::string &  fieldName,
const std::string &  blockId,
const std::vector< std::size_t > &  localElementIds,
const ArrayT &  solutionValues,
double  scaleValue = 1.0 
)

Writes a particular field to a cell array. Notice this is setup to work with the worksets associated with Panzer.

Parameters
[in]fieldNameName of field to be filled
[in]blockIdName of block this set of elements belongs to
[in]localElementIdsLocal element IDs for this set of solution values
[in]solutionValuesA one dimensional array object sized by (Cells)
Note
The block ID is not strictly needed in this context. However forcing the user to provide it does permit an additional level of safety. The implicit assumption is that the elements being "set" are part of the specified block. This prevents the need to perform a null pointer check on the field data, because the STK_Interface construction of the fields should force it to be nonnull...

Definition at line 1674 of file Panzer_STK_Interface.hpp.

Teuchos::RCP< const std::vector< stk::mesh::Entity > > panzer_stk::STK_Interface::getEdgesOrderedByLID ( ) const

Get Vector of edge entities ordered by their LID, returns an RCP so that it is easily stored by the caller.

Definition at line 1630 of file Panzer_STK_Interface.cpp.

Teuchos::RCP< const std::vector< stk::mesh::Entity > > panzer_stk::STK_Interface::getFacesOrderedByLID ( ) const

Get Vector of face entities ordered by their LID, returns an RCP so that it is easily stored by the caller.

Definition at line 1675 of file Panzer_STK_Interface.cpp.

template<typename ArrayT >
void panzer_stk::STK_Interface::setEdgeFieldData ( const std::string &  fieldName,
const std::string &  blockId,
const std::vector< std::size_t > &  localEdgeIds,
const ArrayT &  edgeValues,
double  scaleValue = 1.0 
)

Writes a particular field to an edge array. Notice this is setup to work with the worksets associated with Panzer.

Parameters
[in]fieldNameName of field to be filled
[in]blockIdName of block this set of elements belongs to
[in]localEdgeIdsLocal edge IDs for this set of solution values
[in]edgeValuesA one dimensional array object sized by (edges)
Note
The block ID is not strictly needed in this context. However forcing the user to provide it does permit an additional level of safety. The implicit assumption is that the elements being "set" are part of the specified block. This prevents the need to perform a null pointer check on the field data, because the STK_Interface construction of the fields should force it to be nonnull...

Definition at line 1695 of file Panzer_STK_Interface.hpp.

template<typename ArrayT >
void panzer_stk::STK_Interface::setFaceFieldData ( const std::string &  fieldName,
const std::string &  blockId,
const std::vector< std::size_t > &  localFaceIds,
const ArrayT &  faceValues,
double  scaleValue = 1.0 
)

Writes a particular field to a face array. Notice this is setup to work with the worksets associated with Panzer.

Parameters
[in]fieldNameName of field to be filled
[in]blockIdName of block this set of elements belongs to
[in]localFaceIdsLocal face IDs for this set of solution values
[in]faceValuesA one dimensional array object sized by (faces)
Note
The block ID is not strictly needed in this context. However forcing the user to provide it does permit an additional level of safety. The implicit assumption is that the elements being "set" are part of the specified block. This prevents the need to perform a null pointer check on the field data, because the STK_Interface construction of the fields should force it to be nonnull...

Definition at line 1716 of file Panzer_STK_Interface.hpp.

template<typename ArrayT >
void panzer_stk::STK_Interface::getElementVertices ( const std::vector< std::size_t > &  localIds,
ArrayT &  vertices 
) const

Get vertices associated with a number of elements of the same geometry.

Parameters
[in]localIdsElement local IDs to construct vertices
[out]verticesOutput array that will be sized (localIds.size(),#Vertices,#Dim)
Note
If not all elements have the same number of vertices an exception is thrown. If the size of localIds is 0, the function will silently return

Definition at line 1739 of file Panzer_STK_Interface.hpp.

template<typename ArrayT >
void panzer_stk::STK_Interface::getElementVertices ( const std::vector< stk::mesh::Entity > &  elements,
ArrayT &  vertices 
) const

Get vertices associated with a number of elements of the same geometry.

Parameters
[in]elementsElement entities to construct vertices
[out]verticesOutput array that will be sized (localIds.size(),#Vertices,#Dim)
Note
If not all elements have the same number of vertices an exception is thrown. If the size of localIds is 0, the function will silently return

Definition at line 1763 of file Panzer_STK_Interface.hpp.

template<typename ArrayT >
void panzer_stk::STK_Interface::getElementVertices ( const std::vector< std::size_t > &  localIds,
const std::string &  eBlock,
ArrayT &  vertices 
) const

Get vertices associated with a number of elements of the same geometry.

Parameters
[in]localIdsElement local IDs to construct vertices
[in]eBlockElement block the elements are in
[out]verticesOutput array that will be sized (localIds.size(),#Vertices,#Dim)
Note
If not all elements have the same number of vertices an exception is thrown. If the size of localIds is 0, the function will silently return

Definition at line 1787 of file Panzer_STK_Interface.hpp.

template<typename ArrayT >
void panzer_stk::STK_Interface::getElementVertices ( const std::vector< stk::mesh::Entity > &  elements,
const std::string &  eBlock,
ArrayT &  vertices 
) const

Get vertices associated with a number of elements of the same geometry.

Parameters
[in]elementsElement entities to construct vertices
[in]eBlockElement block the elements are in
[out]verticesOutput array that will be sized (localIds.size(),#Vertices,#Dim)
Note
If not all elements have the same number of vertices an exception is thrown. If the size of localIds is 0, the function will silently return

Definition at line 1776 of file Panzer_STK_Interface.hpp.

template<typename ArrayT >
void panzer_stk::STK_Interface::getElementVerticesNoResize ( const std::vector< std::size_t > &  localIds,
ArrayT &  vertices 
) const

Get vertices associated with a number of elements of the same geometry. The vertex array will not be resized.

Parameters
[in]localIdsElement local IDs to construct vertices
[out]verticesOutput array that will be sized (localIds.size(),#Vertices,#Dim)
Note
If not all elements have the same number of vertices an exception is thrown. If the size of localIds is 0, the function will silently return

Definition at line 1805 of file Panzer_STK_Interface.hpp.

template<typename ArrayT >
void panzer_stk::STK_Interface::getElementVerticesNoResize ( const std::vector< stk::mesh::Entity > &  elements,
ArrayT &  vertices 
) const

Get vertices associated with a number of elements of the same geometry. The vertex array will not be resized.

Parameters
[in]elementsElement entities to construct vertices
[out]verticesOutput array that will be sized (localIds.size(),#Vertices,#Dim)
Note
If not all elements have the same number of vertices an exception is thrown. If the size of localIds is 0, the function will silently return

Definition at line 1829 of file Panzer_STK_Interface.hpp.

template<typename ArrayT >
void panzer_stk::STK_Interface::getElementVerticesNoResize ( const std::vector< std::size_t > &  localIds,
const std::string &  eBlock,
ArrayT &  vertices 
) const

Get vertices associated with a number of elements of the same geometry. The vertex array will not be resized.

Parameters
[in]localIdsElement local IDs to construct vertices
[in]eBlockElement block the elements are in
[out]verticesOutput array that will be sized (localIds.size(),#Vertices,#Dim)
Note
If not all elements have the same number of vertices an exception is thrown. If the size of localIds is 0, the function will silently return

Definition at line 1853 of file Panzer_STK_Interface.hpp.

template<typename ArrayT >
void panzer_stk::STK_Interface::getElementVerticesNoResize ( const std::vector< stk::mesh::Entity > &  elements,
const std::string &  eBlock,
ArrayT &  vertices 
) const

Get vertices associated with a number of elements of the same geometry. The vertex array will not be resized.

Parameters
[in]elementsElement entities to construct vertices
[in]eBlockElement block the elements are in
[out]verticesOutput array that will be sized (localIds.size(),#Vertices,#Dim)
Note
If not all elements have the same number of vertices an exception is thrown. If the size of localIds is 0, the function will silently return

Definition at line 1842 of file Panzer_STK_Interface.hpp.

template<typename ArrayT >
void panzer_stk::STK_Interface::getElementNodes ( const std::vector< std::size_t > &  localIds,
ArrayT &  nodes 
) const

Get nodes associated with a number of elements of the same geometry.

Parameters
[in]localIdsElement local IDs to construct nodes
[out]nodesOutput array that will be sized (localIds.size(),#Nodes,#Dim)
Note
If not all elements have the same number of nodes an exception is thrown. If the size of localIds is 0, the function will silently return

Definition at line 2068 of file Panzer_STK_Interface.hpp.

template<typename ArrayT >
void panzer_stk::STK_Interface::getElementNodes ( const std::vector< stk::mesh::Entity > &  elements,
ArrayT &  nodes 
) const

Get nodes associated with a number of elements of the same geometry.

Parameters
[in]elementsElement entities to construct nodes
[out]nodesOutput array that will be sized (localIds.size(),#Nodes,#Dim)
Note
If not all elements have the same number of nodes an exception is thrown. If the size of localIds is 0, the function will silently return

Definition at line 2092 of file Panzer_STK_Interface.hpp.

template<typename ArrayT >
void panzer_stk::STK_Interface::getElementNodes ( const std::vector< std::size_t > &  localIds,
const std::string &  eBlock,
ArrayT &  nodes 
) const

Get nodes associated with a number of elements of the same geometry.

Parameters
[in]localIdsElement local IDs to construct nodes
[in]eBlockElement block the elements are in
[out]nodesOutput array that will be sized (localIds.size(),#Nodes,#Dim)
Note
If not all elements have the same number of nodes an exception is thrown. If the size of localIds is 0, the function will silently return

Definition at line 2116 of file Panzer_STK_Interface.hpp.

template<typename ArrayT >
void panzer_stk::STK_Interface::getElementNodes ( const std::vector< stk::mesh::Entity > &  elements,
const std::string &  eBlock,
ArrayT &  nodes 
) const

Get nodes associated with a number of elements of the same geometry.

Parameters
[in]elementsElement entities to construct nodes
[in]eBlockElement block the elements are in
[out]nodesOutput array that will be sized (localIds.size(),#Nodes,#Dim)
Note
If not all elements have the same number of nodes an exception is thrown. If the size of localIds is 0, the function will silently return

Definition at line 2105 of file Panzer_STK_Interface.hpp.

template<typename ArrayT >
void panzer_stk::STK_Interface::getElementNodesNoResize ( const std::vector< std::size_t > &  localIds,
ArrayT &  nodes 
) const

Get nodes associated with a number of elements of the same geometry. The node array will not be resized.

Parameters
[in]localIdsElement local IDs to construct nodes
[out]nodesOutput array that will be sized (localIds.size(),#Nodes,#Dim)
Note
If not all elements have the same number of nodes an exception is thrown. If the size of localIds is 0, the function will silently return

Definition at line 2134 of file Panzer_STK_Interface.hpp.

template<typename ArrayT >
void panzer_stk::STK_Interface::getElementNodesNoResize ( const std::vector< stk::mesh::Entity > &  elements,
ArrayT &  nodes 
) const

Get nodes associated with a number of elements of the same geometry. The node array will not be resized.

Parameters
[in]elementsElement entities to construct nodes
[out]nodesOutput array that will be sized (localIds.size(),#Nodes,#Dim)
Note
If not all elements have the same number of nodes an exception is thrown. If the size of localIds is 0, the function will silently return

Definition at line 2158 of file Panzer_STK_Interface.hpp.

template<typename ArrayT >
void panzer_stk::STK_Interface::getElementNodesNoResize ( const std::vector< std::size_t > &  localIds,
const std::string &  eBlock,
ArrayT &  nodes 
) const

Get nodes associated with a number of elements of the same geometry. The node array will not be resized.

Parameters
[in]localIdsElement local IDs to construct nodes
[in]eBlockElement block the elements are in
[out]nodesOutput array that will be sized (localIds.size(),#Nodes,#Dim)
Note
If not all elements have the same number of nodes an exception is thrown. If the size of localIds is 0, the function will silently return

Definition at line 2182 of file Panzer_STK_Interface.hpp.

template<typename ArrayT >
void panzer_stk::STK_Interface::getElementNodesNoResize ( const std::vector< stk::mesh::Entity > &  elements,
const std::string &  eBlock,
ArrayT &  nodes 
) const

Get nodes associated with a number of elements of the same geometry. The node array will not be resized.

Parameters
[in]elementsElement entities to construct nodes
[in]eBlockElement block the elements are in
[out]nodesOutput array that will be sized (localIds.size(),#Nodes,#Dim)
Note
If not all elements have the same number of nodes an exception is thrown. If the size of localIds is 0, the function will silently return

Definition at line 2171 of file Panzer_STK_Interface.hpp.

stk::mesh::EntityRank panzer_stk::STK_Interface::getElementRank ( ) const
inline

Definition at line 1174 of file Panzer_STK_Interface.hpp.

stk::mesh::EntityRank panzer_stk::STK_Interface::getSideRank ( ) const
inline

Definition at line 1175 of file Panzer_STK_Interface.hpp.

stk::mesh::EntityRank panzer_stk::STK_Interface::getFaceRank ( ) const
inline

Definition at line 1176 of file Panzer_STK_Interface.hpp.

stk::mesh::EntityRank panzer_stk::STK_Interface::getEdgeRank ( ) const
inline

Definition at line 1177 of file Panzer_STK_Interface.hpp.

stk::mesh::EntityRank panzer_stk::STK_Interface::getNodeRank ( ) const
inline

Definition at line 1178 of file Panzer_STK_Interface.hpp.

void panzer_stk::STK_Interface::initializeFromMetaData ( )

Build fields and parts from the meta data

Definition at line 1720 of file Panzer_STK_Interface.cpp.

void panzer_stk::STK_Interface::buildLocalElementIDs ( )

Setup local element IDs

Definition at line 1744 of file Panzer_STK_Interface.cpp.

void panzer_stk::STK_Interface::buildLocalEdgeIDs ( )

Setup local edge IDs

Definition at line 1808 of file Panzer_STK_Interface.cpp.

void panzer_stk::STK_Interface::buildLocalFaceIDs ( )

Setup local face IDs

Definition at line 1828 of file Panzer_STK_Interface.cpp.

const std::vector<Teuchos::RCP<const PeriodicBC_MatcherBase> >& panzer_stk::STK_Interface::getPeriodicBCVector ( ) const
inline

Return a vector containing all the periodic boundary conditions.

Definition at line 1199 of file Panzer_STK_Interface.hpp.

std::vector<Teuchos::RCP<const PeriodicBC_MatcherBase> >& panzer_stk::STK_Interface::getPeriodicBCVector ( )
inline

Return a vector containing all the periodic boundary conditions.

Definition at line 1205 of file Panzer_STK_Interface.hpp.

const bool& panzer_stk::STK_Interface::useBoundingBoxSearch ( ) const
inline

Return a flag indicating if the bounding box search is used when matching periodic Ids.

Definition at line 1210 of file Panzer_STK_Interface.hpp.

void panzer_stk::STK_Interface::setBoundingBoxSearchFlag ( const bool &  searchFlag)
inline

Set the periodic search flag. Indicates if the bounding box search is used

Definition at line 1214 of file Panzer_STK_Interface.hpp.

void panzer_stk::STK_Interface::addPeriodicBC ( const Teuchos::RCP< const PeriodicBC_MatcherBase > &  bc)
inline

Add a periodic boundary condition.

Note
This does not actually change the underlying mesh. The object itself only communciates the matched IDs (currently nodes)

Definition at line 1222 of file Panzer_STK_Interface.hpp.

void panzer_stk::STK_Interface::addPeriodicBCs ( const std::vector< Teuchos::RCP< const PeriodicBC_MatcherBase > > &  bc_vec)
inline

Add a set of periodic boundary conditions.

Note
This does not actually change the underlying mesh. The object itself only communciates the matched IDs (currently nodes)

Definition at line 1230 of file Panzer_STK_Interface.hpp.

std::pair< Teuchos::RCP< std::vector< std::pair< std::size_t, std::size_t > > >, Teuchos::RCP< std::vector< unsigned int > > > panzer_stk::STK_Interface::getPeriodicNodePairing ( ) const

Pairs DOFs on periodic entities

Definition at line 1871 of file Panzer_STK_Interface.cpp.

bool panzer_stk::STK_Interface::validBlockId ( const std::string &  blockId) const

check for a valid block id

Definition at line 1910 of file Panzer_STK_Interface.cpp.

void panzer_stk::STK_Interface::print ( std::ostream &  os) const

Print a brief description of the STK mesh to a stream.

Definition at line 1917 of file Panzer_STK_Interface.cpp.

void panzer_stk::STK_Interface::printMetaData ( std::ostream &  os) const

Print a brief description of the STK mesh to a stream.

Definition at line 1967 of file Panzer_STK_Interface.cpp.

Teuchos::RCP< const shards::CellTopology > panzer_stk::STK_Interface::getCellTopology ( const std::string &  eBlock) const

Get the cell topology from the element block.

Definition at line 2007 of file Panzer_STK_Interface.cpp.

double panzer_stk::STK_Interface::getCurrentStateTime ( ) const
inline

Get the value of the time-stamp the last time this object was written to Exodus (default 0.0)

Note
The initial state time is completely disconnected from the current state time

Definition at line 1258 of file Panzer_STK_Interface.hpp.

double panzer_stk::STK_Interface::getInitialStateTime ( ) const
inline

Get the value of the time-stamp when this object was created (default 0.0) or when setInitialStateTime was called.

Note
The initial state time is completely disconnected from the current state time

Definition at line 1265 of file Panzer_STK_Interface.hpp.

void panzer_stk::STK_Interface::setInitialStateTime ( double  value)
inline

Set the value of the initial time-stamp

Note
The initial state time is completely disconnected from the current state time

Definition at line 1271 of file Panzer_STK_Interface.hpp.

void panzer_stk::STK_Interface::rebalance ( const Teuchos::ParameterList params)

Rebalance using zoltan. Note that this will void the local element ids.

Definition at line 2034 of file Panzer_STK_Interface.cpp.

void panzer_stk::STK_Interface::setBlockWeight ( const std::string &  blockId,
double  weight 
)
inline

Set the weight for a particular element block. Larger means more costly to assemble and evaluate.

Definition at line 1280 of file Panzer_STK_Interface.hpp.

void panzer_stk::STK_Interface::setUseFieldCoordinates ( bool  useFieldCoordinates)
inline

When coordinates are returned in the getElementNodes method, extract coordinates using a specified field (not the intrinsic coordinates) where available (where unavailable the intrinsic coordinates are used. Note that this does not change the behavior of getNodeCoordinates. This is set to false by default.

Definition at line 1289 of file Panzer_STK_Interface.hpp.

bool panzer_stk::STK_Interface::getUseFieldCoordinates ( ) const
inline

Return the use field coordinates flag

Definition at line 1293 of file Panzer_STK_Interface.hpp.

void panzer_stk::STK_Interface::setUseLowerCaseForIO ( bool  useLowerCase)
inline

Use lower case (or not) for I/O

Definition at line 1297 of file Panzer_STK_Interface.hpp.

bool panzer_stk::STK_Interface::getUseLowerCaseForIO ( ) const
inline

Use lower case (or not) for I/O

Definition at line 1301 of file Panzer_STK_Interface.hpp.

template<typename ArrayT >
void panzer_stk::STK_Interface::getElementVertices_FromField ( const std::vector< stk::mesh::Entity > &  elements,
const std::string &  eBlock,
ArrayT &  vertices 
) const

Get vertices associated with a number of elements of the same geometry, note that a coordinate field will be used (if not is available an exception will be thrown).

Parameters
[in]elementsElement entities to construct vertices
[in]eBlockElement block the elements are in
[out]verticesOutput array that will be sized (localIds.size(),#Vertices,#Dim)
Note
If not all elements have the same number of vertices an exception is thrown. If the size of localIds is 0, the function will silently return

Definition at line 1968 of file Panzer_STK_Interface.hpp.

template<typename ArrayT >
void panzer_stk::STK_Interface::getElementVertices_FromFieldNoResize ( const std::vector< stk::mesh::Entity > &  elements,
const std::string &  eBlock,
ArrayT &  vertices 
) const

Definition at line 2021 of file Panzer_STK_Interface.hpp.

template<typename ArrayT >
void panzer_stk::STK_Interface::getElementVertices_FromCoords ( const std::vector< stk::mesh::Entity > &  elements,
ArrayT &  vertices 
) const

Get vertices associated with a number of elements of the same geometry. This access the true mesh coordinates array.

Parameters
[in]elementsElement entities to construct vertices
[out]verticesOutput array that will be sized (localIds.size(),#Vertices,#Dim)
Note
If not all elements have the same number of vertices an exception is thrown. If the size of localIds is 0, the function will silently return

Definition at line 1871 of file Panzer_STK_Interface.hpp.

template<typename ArrayT >
void panzer_stk::STK_Interface::getElementVertices_FromCoordsNoResize ( const std::vector< stk::mesh::Entity > &  elements,
ArrayT &  vertices 
) const

Definition at line 1922 of file Panzer_STK_Interface.hpp.

template<typename ArrayT >
void panzer_stk::STK_Interface::getElementNodes_FromField ( const std::vector< stk::mesh::Entity > &  elements,
const std::string &  eBlock,
ArrayT &  nodes 
) const

Get nodes associated with a number of elements of the same geometry, note that a coordinate field will be used (if not is available an exception will be thrown).

Parameters
[in]elementsElement entities to construct nodes
[in]eBlockElement block the elements are in
[out]nodesOutput array that will be sized (localIds.size(),#Nodes,#Dim)
Note
If not all elements have the same number of nodes an exception is thrown. If the size of localIds is 0, the function will silently return

Definition at line 2297 of file Panzer_STK_Interface.hpp.

template<typename ArrayT >
void panzer_stk::STK_Interface::getElementNodes_FromFieldNoResize ( const std::vector< stk::mesh::Entity > &  elements,
const std::string &  eBlock,
ArrayT &  nodes 
) const

Definition at line 2351 of file Panzer_STK_Interface.hpp.

template<typename ArrayT >
void panzer_stk::STK_Interface::getElementNodes_FromCoords ( const std::vector< stk::mesh::Entity > &  elements,
ArrayT &  nodes 
) const

Get nodes associated with a number of elements of the same geometry. This access the true mesh coordinates array.

Parameters
[in]elementsElement entities to construct nodes
[out]nodesOutput array that will be sized (localIds.size(),#Nodes,#Dim)
Note
If not all elements have the same number of nodes an exception is thrown. If the size of localIds is 0, the function will silently return

Definition at line 2200 of file Panzer_STK_Interface.hpp.

template<typename ArrayT >
void panzer_stk::STK_Interface::getElementNodes_FromCoordsNoResize ( const std::vector< stk::mesh::Entity > &  elements,
ArrayT &  nodes 
) const

Definition at line 2251 of file Panzer_STK_Interface.hpp.

void panzer_stk::STK_Interface::refineMesh ( const int  numberOfLevels,
const bool  deleteParentElements 
)

Uniformly refine the mesh using Percept

Parameters
[in]numberOfLevelsNumber of uniform refinement levels to apply. Must be >=1.
[in]deleteParentElementsIf true, deletes the parent elements from the mesh to save memory.

Definition at line 2069 of file Panzer_STK_Interface.cpp.

void panzer_stk::STK_Interface::buildEntityCounts ( )
protected

Compute global entity counts.

Definition at line 1045 of file Panzer_STK_Interface.cpp.

void panzer_stk::STK_Interface::buildMaxEntityIds ( )
protected

Compute global entity counts.

Definition at line 1051 of file Panzer_STK_Interface.cpp.

void panzer_stk::STK_Interface::initializeFieldsInSTK ( const std::map< std::pair< std::string, std::string >, SolutionFieldType * > &  nameToField,
bool  setupIO 
)
protected

Initialize STK fields using a map (allocate space for storage and writing) to a specific entity rank.

Definition at line 398 of file Panzer_STK_Interface.cpp.

Teuchos::RCP< Teuchos::MpiComm< int > > panzer_stk::STK_Interface::getSafeCommunicator ( stk::ParallelMachine  parallelMach) const
protected

Build a safely handled Teuchos MPI communicator from a parallel machine. This object asserts ownership of the communicator so that we can gurantee existence so the outer user can do whatever they'd like with the original.

Definition at line 2023 of file Panzer_STK_Interface.cpp.

void panzer_stk::STK_Interface::applyElementLoadBalanceWeights ( )
protected

In a pure local operation apply the user specified block weights for each element block to the field that defines the load balance weighting. This uses the blockWeights_ member to determine the user value that has been set for a particular element block.

Definition at line 1787 of file Panzer_STK_Interface.cpp.

bool panzer_stk::STK_Interface::isMeshCoordField ( const std::string &  eBlock,
const std::string &  fieldName,
int &  axis 
) const
protected

Determine if a particular field in an element block is a displacement field. Return the displacement field name in the requested argument if so.

Definition at line 1849 of file Panzer_STK_Interface.cpp.

template<typename ArrayT >
void panzer_stk::STK_Interface::setDispFieldData ( const std::string &  fieldName,
const std::string &  blockId,
int  axis,
const std::vector< std::size_t > &  localElementIds,
const ArrayT &  solutionValues 
)
protected

Writes a particular field to an array as a coordinate diplacement. Notice this is setup to work with the worksets associated with Panzer.

Parameters
[in]fieldNameName of field to be filled
[in]blockIdName of block this set of elements belongs to
[in]dimensionWhat coordinate axis to write to
[in]localElementIdsLocal element IDs for this set of solution values
[in]solutionValuesAn two dimensional array object sized by (Cells,Basis Count)
Note
The block ID is not strictly needed in this context. However forcing the user to provide it does permit an additional level of safety. The implicit assumption is that the elements being "set" are part of the specified block. This prevents the need to perform a null pointer check on the field data, because the STK_Interface construction of the fields should force it to be nonnull...

Definition at line 1615 of file Panzer_STK_Interface.hpp.

Member Data Documentation

const std::string panzer_stk::STK_Interface::coordsString = "coordinates"
static

Definition at line 1380 of file Panzer_STK_Interface.hpp.

const std::string panzer_stk::STK_Interface::nodesString = "nodes"
static

Definition at line 1381 of file Panzer_STK_Interface.hpp.

const std::string panzer_stk::STK_Interface::edgesString = "edges"
static

Definition at line 1382 of file Panzer_STK_Interface.hpp.

const std::string panzer_stk::STK_Interface::edgeBlockString = "edge_block"
static

Definition at line 1383 of file Panzer_STK_Interface.hpp.

const std::string panzer_stk::STK_Interface::faceBlockString = "face_block"
static

Definition at line 1384 of file Panzer_STK_Interface.hpp.

const std::string panzer_stk::STK_Interface::facesString = "faces"
static

Definition at line 1385 of file Panzer_STK_Interface.hpp.

std::vector<Teuchos::RCP<const PeriodicBC_MatcherBase> > panzer_stk::STK_Interface::periodicBCs_
protected

Definition at line 1440 of file Panzer_STK_Interface.hpp.

bool panzer_stk::STK_Interface::useBBoxSearch_ = false
protected

Definition at line 1441 of file Panzer_STK_Interface.hpp.

Teuchos::RCP<stk::mesh::MetaData> panzer_stk::STK_Interface::metaData_
protected

Definition at line 1443 of file Panzer_STK_Interface.hpp.

Teuchos::RCP<stk::mesh::BulkData> panzer_stk::STK_Interface::bulkData_
protected

Definition at line 1444 of file Panzer_STK_Interface.hpp.

std::map<std::string, stk::mesh::Part*> panzer_stk::STK_Interface::elementBlocks_
protected

Definition at line 1450 of file Panzer_STK_Interface.hpp.

std::map<std::string, stk::mesh::Part*> panzer_stk::STK_Interface::sidesets_
protected

Definition at line 1451 of file Panzer_STK_Interface.hpp.

std::map<std::string, stk::mesh::Part*> panzer_stk::STK_Interface::nodesets_
protected

Definition at line 1452 of file Panzer_STK_Interface.hpp.

std::map<std::string, stk::mesh::Part*> panzer_stk::STK_Interface::edgeBlocks_
protected

Definition at line 1453 of file Panzer_STK_Interface.hpp.

std::map<std::string, stk::mesh::Part*> panzer_stk::STK_Interface::faceBlocks_
protected

Definition at line 1454 of file Panzer_STK_Interface.hpp.

std::map<std::string, Teuchos::RCP<shards::CellTopology> > panzer_stk::STK_Interface::elementBlockCT_
protected

Definition at line 1456 of file Panzer_STK_Interface.hpp.

stk::mesh::Part* panzer_stk::STK_Interface::nodesPart_
protected

Definition at line 1459 of file Panzer_STK_Interface.hpp.

std::vector<stk::mesh::Part*> panzer_stk::STK_Interface::nodesPartVec_
protected

Definition at line 1460 of file Panzer_STK_Interface.hpp.

stk::mesh::Part* panzer_stk::STK_Interface::edgesPart_
protected

Definition at line 1461 of file Panzer_STK_Interface.hpp.

std::vector<stk::mesh::Part*> panzer_stk::STK_Interface::edgesPartVec_
protected

Definition at line 1462 of file Panzer_STK_Interface.hpp.

stk::mesh::Part* panzer_stk::STK_Interface::facesPart_
protected

Definition at line 1463 of file Panzer_STK_Interface.hpp.

std::vector<stk::mesh::Part*> panzer_stk::STK_Interface::facesPartVec_
protected

Definition at line 1464 of file Panzer_STK_Interface.hpp.

VectorFieldType* panzer_stk::STK_Interface::coordinatesField_
protected

Definition at line 1466 of file Panzer_STK_Interface.hpp.

VectorFieldType* panzer_stk::STK_Interface::edgesField_
protected

Definition at line 1467 of file Panzer_STK_Interface.hpp.

VectorFieldType* panzer_stk::STK_Interface::facesField_
protected

Definition at line 1468 of file Panzer_STK_Interface.hpp.

ProcIdFieldType* panzer_stk::STK_Interface::processorIdField_
protected

Definition at line 1469 of file Panzer_STK_Interface.hpp.

SolutionFieldType* panzer_stk::STK_Interface::loadBalField_
protected

Definition at line 1470 of file Panzer_STK_Interface.hpp.

std::map<std::pair<std::string,std::string>,SolutionFieldType*> panzer_stk::STK_Interface::fieldNameToSolution_
protected

Definition at line 1473 of file Panzer_STK_Interface.hpp.

std::map<std::pair<std::string,std::string>,SolutionFieldType*> panzer_stk::STK_Interface::fieldNameToCellField_
protected

Definition at line 1474 of file Panzer_STK_Interface.hpp.

std::map<std::pair<std::string,std::string>,SolutionFieldType*> panzer_stk::STK_Interface::fieldNameToEdgeField_
protected

Definition at line 1475 of file Panzer_STK_Interface.hpp.

std::map<std::pair<std::string,std::string>,SolutionFieldType*> panzer_stk::STK_Interface::fieldNameToFaceField_
protected

Definition at line 1476 of file Panzer_STK_Interface.hpp.

std::set<std::string> panzer_stk::STK_Interface::informationRecords_
protected

Definition at line 1479 of file Panzer_STK_Interface.hpp.

unsigned panzer_stk::STK_Interface::dimension_
protected

Definition at line 1481 of file Panzer_STK_Interface.hpp.

bool panzer_stk::STK_Interface::initialized_
protected

Definition at line 1483 of file Panzer_STK_Interface.hpp.

std::vector<std::size_t> panzer_stk::STK_Interface::entityCounts_
protected

Definition at line 1486 of file Panzer_STK_Interface.hpp.

std::vector<stk::mesh::EntityId> panzer_stk::STK_Interface::maxEntityId_
protected

Definition at line 1489 of file Panzer_STK_Interface.hpp.

unsigned panzer_stk::STK_Interface::procRank_
protected

Definition at line 1491 of file Panzer_STK_Interface.hpp.

std::size_t panzer_stk::STK_Interface::currentLocalId_
protected

Definition at line 1492 of file Panzer_STK_Interface.hpp.

Teuchos::RCP<Teuchos::MpiComm<int> > panzer_stk::STK_Interface::mpiComm_
protected

Definition at line 1494 of file Panzer_STK_Interface.hpp.

double panzer_stk::STK_Interface::initialStateTime_
protected

Definition at line 1496 of file Panzer_STK_Interface.hpp.

double panzer_stk::STK_Interface::currentStateTime_
protected

Definition at line 1497 of file Panzer_STK_Interface.hpp.

Teuchos::RCP<std::vector<stk::mesh::Entity> > panzer_stk::STK_Interface::orderedElementVector_
mutableprotected

Definition at line 1541 of file Panzer_STK_Interface.hpp.

Teuchos::RCP<std::vector<stk::mesh::Entity> > panzer_stk::STK_Interface::orderedEdgeVector_
mutableprotected

Definition at line 1544 of file Panzer_STK_Interface.hpp.

Teuchos::RCP<std::vector<stk::mesh::Entity> > panzer_stk::STK_Interface::orderedFaceVector_
mutableprotected

Definition at line 1547 of file Panzer_STK_Interface.hpp.

std::map<std::string,double> panzer_stk::STK_Interface::blockWeights_
protected

Definition at line 1550 of file Panzer_STK_Interface.hpp.

std::unordered_map<stk::mesh::EntityId,std::size_t> panzer_stk::STK_Interface::localIDHash_
protected

Definition at line 1552 of file Panzer_STK_Interface.hpp.

std::unordered_map<stk::mesh::EntityId,std::size_t> panzer_stk::STK_Interface::localEdgeIDHash_
protected

Definition at line 1553 of file Panzer_STK_Interface.hpp.

std::unordered_map<stk::mesh::EntityId,std::size_t> panzer_stk::STK_Interface::localFaceIDHash_
protected

Definition at line 1554 of file Panzer_STK_Interface.hpp.

std::map<std::string,std::vector<std::string> > panzer_stk::STK_Interface::meshCoordFields_
protected

Definition at line 1559 of file Panzer_STK_Interface.hpp.

std::map<std::string,std::vector<std::string> > panzer_stk::STK_Interface::meshDispFields_
protected

Definition at line 1560 of file Panzer_STK_Interface.hpp.

bool panzer_stk::STK_Interface::useFieldCoordinates_
protected

Definition at line 1562 of file Panzer_STK_Interface.hpp.

bool panzer_stk::STK_Interface::useLowerCase_
protected

Definition at line 1564 of file Panzer_STK_Interface.hpp.


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