Panzer
Version of the Day
|
#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 (const bool find_and_set_shared_nodes_in_stk=true) |
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 VectorFieldType & | getCoordinatesField () const |
const VectorFieldType & | getEdgesField () const |
const VectorFieldType & | getFacesField () 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 | setupExodusFile (const std::string &filename, const std::vector< Ioss::Property > &ioss_properties, const bool append=false, const bool append_after_restart_time=false, const double restart_time=0.0) |
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 |
ProcIdFieldType * | getProcessorIdField () |
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 ¶ms) |
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_ |
VectorFieldType * | coordinatesField_ |
VectorFieldType * | edgesField_ |
VectorFieldType * | facesField_ |
ProcIdFieldType * | processorIdField_ |
SolutionFieldType * | loadBalField_ |
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_ |
Definition at line 70 of file Panzer_STK_Interface.hpp.
typedef double panzer_stk::STK_Interface::ProcIdData |
Definition at line 72 of file Panzer_STK_Interface.hpp.
typedef stk::mesh::Field<double> panzer_stk::STK_Interface::SolutionFieldType |
Definition at line 73 of file Panzer_STK_Interface.hpp.
typedef stk::mesh::Field<double> panzer_stk::STK_Interface::VectorFieldType |
Definition at line 74 of file Panzer_STK_Interface.hpp.
typedef stk::mesh::Field<ProcIdData> panzer_stk::STK_Interface::ProcIdFieldType |
Definition at line 75 of file Panzer_STK_Interface.hpp.
panzer_stk::STK_Interface::STK_Interface | ( | ) |
Definition at line 78 of file Panzer_STK_Interface.cpp.
panzer_stk::STK_Interface::STK_Interface | ( | unsigned | dim | ) |
Default constructor
Definition at line 92 of file Panzer_STK_Interface.cpp.
panzer_stk::STK_Interface::STK_Interface | ( | Teuchos::RCP< stk::mesh::MetaData > | metaData | ) |
Definition at line 85 of file Panzer_STK_Interface.cpp.
void panzer_stk::STK_Interface::addElementBlock | ( | const std::string & | name, |
const CellTopologyData * | ctData | ||
) |
Add an element block with a string name
Definition at line 1568 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 1599 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 1620 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 1644 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 1665 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 104 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 116 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 130 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 149 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 169 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 190 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.
[in] | blockId | Element block to use |
[in] | coordPrefix | Prefix for solution coordinate field (assumes using "X","Y","Z" as postfix) |
[in] | dispPrefix | Prefix for displacment (output) field |
Definition at line 211 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 255 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.
[in] | parallelMach | Communicator |
[in] | setupIO | If set to true and IOSS is enabled, the output mesh will be initialized. |
[in] | buildRefinementSupport | If true, build percept uniform refinement objects. |
Definition at line 260 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 392 of file Panzer_STK_Interface.cpp.
void panzer_stk::STK_Interface::beginModification | ( | ) |
Put the bulk data manager in modification mode.
Definition at line 402 of file Panzer_STK_Interface.cpp.
void panzer_stk::STK_Interface::endModification | ( | const bool | find_and_set_shared_nodes_in_stk = true | ) |
Take the bulk data manager out of modification mode.
find_and_set_shared_nodes_in_stk | (bool) If set to true, this performs a geometric search to figure out which nodes are shared entities and marks them as shared in stk. This is very inefficient and communication intensive. Set this to false for exodus as the sharing is already known. Must be true for inline meshes. We should fix the inline mesh factories to mark sharing at node construction time to eliminate the need for searching for shared nodes. Defaults to true to maintain backwards compatibility. |
Definition at line 410 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.
coord.size()==getDimension()
isModifiable()==true
Definition at line 433 of file Panzer_STK_Interface.cpp.
void panzer_stk::STK_Interface::addElement | ( | const Teuchos::RCP< ElementDescriptor > & | ed, |
stk::mesh::Part * | block | ||
) |
Definition at line 501 of file Panzer_STK_Interface.cpp.
void panzer_stk::STK_Interface::addEdges | ( | ) |
Definition at line 524 of file Panzer_STK_Interface.cpp.
void panzer_stk::STK_Interface::addFaces | ( | ) |
Definition at line 552 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 451 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 459 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 467 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 474 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 484 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 491 of file Panzer_STK_Interface.cpp.
|
inline |
Grab the coordinates field
Definition at line 250 of file Panzer_STK_Interface.hpp.
|
inline |
Grab the edges field
Definition at line 255 of file Panzer_STK_Interface.hpp.
|
inline |
Definition at line 258 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 1072 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 1078 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 1083 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 1104 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 1114 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 1129 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 1139 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 1153 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()
.
[in] | edgeBlockName | Name of edge block |
[in,out] | edges | Vector of entities containing the requested edges. |
Definition at line 1162 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()
.
[in] | edgeBlockName | Name of edge block |
[in] | blockName | Name of block |
[in,out] | edges | Vector of entities containing the requested edges. |
Definition at line 1175 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()
.
[in] | edgeBlockName | Name of edge block |
[in,out] | edges | Vector of entities containing the requested edges. |
Definition at line 1192 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()
.
[in] | edgeBlockName | Name of edge block |
[in] | blockName | Name of block |
[in,out] | edges | Vector of entities containing the requested edges. |
Definition at line 1204 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 1221 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()
.
[in] | faceBlockName | Name of face block |
[in,out] | faces | Vector of entities containing the requested faces. |
Definition at line 1231 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()
.
[in] | faceBlockName | Name of face block |
[in] | blockName | Name of block |
[in,out] | faces | Vector of entities containing the requested faces. |
Definition at line 1244 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()
.
[in] | faceBlockName | Name of face block |
[in,out] | faces | Vector of entities containing the requested faces. |
Definition at line 1261 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()
.
[in] | faceBlockName | Name of face block |
[in] | blockName | Name of block |
[in,out] | faces | Vector of entities containing the requested faces. |
Definition at line 1273 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()
.
[in] | sideName | Name of side set |
[in,out] | sides | Vector of entities containing the requested sides. |
Definition at line 1290 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()
.
[in] | sideName | Name of side set |
[in] | blockName | Name of block |
[in,out] | sides | Vector of entities containing the requested sides. |
Definition at line 1303 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()
.
[in] | sideName | Name of side set |
[in,out] | sides | Vector of entities containing the requested sides. |
Definition at line 1320 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()
.
[in] | sideName | Name of side set |
[in] | blockName | Name of block |
[in,out] | sides | Vector of entities containing the requested sides. |
Definition at line 1332 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
.
[in] | sideName | Name of side set |
[in] | blockName | Name of block |
[in,out] | nodes | Vector of entities containing the requested nodes. |
Definition at line 1349 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.
[in] | src | The handle to the source entity (the 'from' part of the relation) |
[in] | tgt_rank | The entity rank of relations to search |
[in] | rel_id | The id of the relation to search for |
Definition at line 582 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.
[in] | filename | The name of the output Exodus file. |
[in] | append | If 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 603 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.
writeToExodus(double timestep)
.[in] | filename | The name of the output Exodus file. |
[in] | append | If 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_time | If 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_time | If append_after_restart_time is true, this is the time value to append after. Otherwise this value is ignored. |
<tt>std::logic_error</tt> | If the STK_Interface does not yet have a MPI communicator. |
Definition at line 617 of file Panzer_STK_Interface.cpp.
void panzer_stk::STK_Interface::setupExodusFile | ( | const std::string & | filename, |
const std::vector< Ioss::Property > & | ioss_properties, | ||
const bool | append = false , |
||
const bool | append_after_restart_time = false , |
||
const double | restart_time = 0.0 |
||
) |
Definition at line 632 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()
.
setupExodusFile(const std::string& filename)
must be called before any invocations of this routine.[in] | timestep | The simulation time at which you're writing the results. |
Definition at line 696 of file Panzer_STK_Interface.cpp.
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.
writeToExodus()
is called.[in] | key | The 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] | value | The value of the global variable you'll be adding. |
Definition at line 847 of file Panzer_STK_Interface.cpp.
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.
writeToExodus()
is called.[in] | key | The 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] | value | The value of the global variable you'll be adding. |
Definition at line 861 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.
writeToExodus()
is called.[in] | key | The 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] | value | The value of the global variable you'll be adding. |
Definition at line 875 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.
writeToExodus()
is called.[in] | key | The 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] | value | The value of the global variable you'll be adding. |
Definition at line 890 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 2019 of file Panzer_STK_Interface.cpp.
|
inline |
Definition at line 597 of file Panzer_STK_Interface.hpp.
|
inline |
Definition at line 598 of file Panzer_STK_Interface.hpp.
bool panzer_stk::STK_Interface::isWritable | ( | ) | const |
Definition at line 898 of file Panzer_STK_Interface.cpp.
|
inline |
Definition at line 608 of file Panzer_STK_Interface.hpp.
|
inline |
get the dimension
Definition at line 613 of file Panzer_STK_Interface.hpp.
|
inline |
get the block count
Definition at line 617 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
.
[in,out] | names | Vector of names of the element blocks. |
Definition at line 1366 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
.
[in,out] | names | Vector of names of the side sets |
Definition at line 1378 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
.
[in,out] | names | Vector of names of the node sets. |
Definition at line 1390 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
.
[in,out] | names | Vector of names of the edge blocks. |
Definition at line 1400 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
.
[in,out] | names | Vector of names of the face blocks. |
Definition at line 1410 of file Panzer_STK_Interface.cpp.
|
inline |
Get a pointer to the locally owned part.
Definition at line 666 of file Panzer_STK_Interface.hpp.
|
inline |
get the block part
Definition at line 670 of file Panzer_STK_Interface.hpp.
|
inline |
get the block part
Definition at line 678 of file Panzer_STK_Interface.hpp.
|
inline |
get the block part
Definition at line 686 of file Panzer_STK_Interface.hpp.
|
inline |
get the side set count
Definition at line 694 of file Panzer_STK_Interface.hpp.
|
inline |
Definition at line 697 of file Panzer_STK_Interface.hpp.
|
inline |
get the side set count
Definition at line 704 of file Panzer_STK_Interface.hpp.
|
inline |
Definition at line 707 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 1043 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 1051 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 907 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 990 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 920 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 940 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 958 of file Panzer_STK_Interface.cpp.
void panzer_stk::STK_Interface::buildSubcells | ( | ) |
force the mesh to build subcells: edges and faces
Definition at line 1059 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 1420 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 1427 of file Panzer_STK_Interface.cpp.
|
inline |
Get an elements global index
Definition at line 755 of file Panzer_STK_Interface.hpp.
|
inline |
Get an elements global index
Definition at line 760 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 1438 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 1443 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 1452 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 1457 of file Panzer_STK_Interface.cpp.
|
inline |
Get an edge's global index
Definition at line 781 of file Panzer_STK_Interface.hpp.
|
inline |
Get an edge's global index
Definition at line 786 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 1464 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 1469 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 1478 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 1483 of file Panzer_STK_Interface.cpp.
|
inline |
Get a face's global index
Definition at line 807 of file Panzer_STK_Interface.hpp.
|
inline |
Get a face's global index
Definition at line 812 of file Panzer_STK_Interface.hpp.
|
inline |
Get an Entity's parallel owner (process rank)
Definition at line 817 of file Panzer_STK_Interface.hpp.
|
inline |
Check if entity handle is valid
Definition at line 822 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 1491 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 1499 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 1513 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 1527 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 1541 of file Panzer_STK_Interface.cpp.
|
inline |
Definition at line 857 of file Panzer_STK_Interface.hpp.
|
inline |
Has initialize
been called on this mesh object?
Definition at line 860 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 1555 of file Panzer_STK_Interface.cpp.
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.
[in] | fieldName | Name of field to be filled |
[in] | blockId | Name of block this set of elements belongs to |
[in] | localElementIds | Local element IDs for this set of solution values |
[in] | solutionValues | An two dimensional array object sized by (Cells,Basis Count) |
Definition at line 1567 of file Panzer_STK_Interface.hpp.
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.
[in] | fieldName | Name of field to be filled |
[in] | blockId | Name of block this set of elements belongs to |
[in] | localElementIds | Local element IDs for this set of solution values |
[in] | solutionValues | An two dimensional array object sized by (Cells,Basis Count) |
Definition at line 1629 of file Panzer_STK_Interface.hpp.
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.
[in] | fieldName | Name of field to be filled |
[in] | blockId | Name of block this set of elements belongs to |
[in] | localElementIds | Local element IDs for this set of solution values |
[in] | solutionValues | A one dimensional array object sized by (Cells) |
Definition at line 1659 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 1586 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 1631 of file Panzer_STK_Interface.cpp.
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.
[in] | fieldName | Name of field to be filled |
[in] | blockId | Name of block this set of elements belongs to |
[in] | localEdgeIds | Local edge IDs for this set of solution values |
[in] | edgeValues | A one dimensional array object sized by (edges) |
Definition at line 1680 of file Panzer_STK_Interface.hpp.
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.
[in] | fieldName | Name of field to be filled |
[in] | blockId | Name of block this set of elements belongs to |
[in] | localFaceIds | Local face IDs for this set of solution values |
[in] | faceValues | A one dimensional array object sized by (faces) |
Definition at line 1701 of file Panzer_STK_Interface.hpp.
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.
[in] | localIds | Element local IDs to construct vertices |
[out] | vertices | Output array that will be sized (localIds.size() ,#Vertices,#Dim) |
localIds
is 0, the function will silently return Definition at line 1724 of file Panzer_STK_Interface.hpp.
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.
[in] | elements | Element entities to construct vertices |
[out] | vertices | Output array that will be sized (localIds.size() ,#Vertices,#Dim) |
localIds
is 0, the function will silently return Definition at line 1748 of file Panzer_STK_Interface.hpp.
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.
[in] | localIds | Element local IDs to construct vertices |
[in] | eBlock | Element block the elements are in |
[out] | vertices | Output array that will be sized (localIds.size() ,#Vertices,#Dim) |
localIds
is 0, the function will silently return Definition at line 1772 of file Panzer_STK_Interface.hpp.
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.
[in] | elements | Element entities to construct vertices |
[in] | eBlock | Element block the elements are in |
[out] | vertices | Output array that will be sized (localIds.size() ,#Vertices,#Dim) |
localIds
is 0, the function will silently return Definition at line 1761 of file Panzer_STK_Interface.hpp.
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.
[in] | localIds | Element local IDs to construct vertices |
[out] | vertices | Output array that will be sized (localIds.size() ,#Vertices,#Dim) |
localIds
is 0, the function will silently return Definition at line 1790 of file Panzer_STK_Interface.hpp.
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.
[in] | elements | Element entities to construct vertices |
[out] | vertices | Output array that will be sized (localIds.size() ,#Vertices,#Dim) |
localIds
is 0, the function will silently return Definition at line 1814 of file Panzer_STK_Interface.hpp.
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.
[in] | localIds | Element local IDs to construct vertices |
[in] | eBlock | Element block the elements are in |
[out] | vertices | Output array that will be sized (localIds.size() ,#Vertices,#Dim) |
localIds
is 0, the function will silently return Definition at line 1838 of file Panzer_STK_Interface.hpp.
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.
[in] | elements | Element entities to construct vertices |
[in] | eBlock | Element block the elements are in |
[out] | vertices | Output array that will be sized (localIds.size() ,#Vertices,#Dim) |
localIds
is 0, the function will silently return Definition at line 1827 of file Panzer_STK_Interface.hpp.
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.
[in] | localIds | Element local IDs to construct nodes |
[out] | nodes | Output array that will be sized (localIds.size() ,#Nodes,#Dim) |
localIds
is 0, the function will silently return Definition at line 2053 of file Panzer_STK_Interface.hpp.
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.
[in] | elements | Element entities to construct nodes |
[out] | nodes | Output array that will be sized (localIds.size() ,#Nodes,#Dim) |
localIds
is 0, the function will silently return Definition at line 2077 of file Panzer_STK_Interface.hpp.
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.
[in] | localIds | Element local IDs to construct nodes |
[in] | eBlock | Element block the elements are in |
[out] | nodes | Output array that will be sized (localIds.size() ,#Nodes,#Dim) |
localIds
is 0, the function will silently return Definition at line 2101 of file Panzer_STK_Interface.hpp.
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.
[in] | elements | Element entities to construct nodes |
[in] | eBlock | Element block the elements are in |
[out] | nodes | Output array that will be sized (localIds.size() ,#Nodes,#Dim) |
localIds
is 0, the function will silently return Definition at line 2090 of file Panzer_STK_Interface.hpp.
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.
[in] | localIds | Element local IDs to construct nodes |
[out] | nodes | Output array that will be sized (localIds.size() ,#Nodes,#Dim) |
localIds
is 0, the function will silently return Definition at line 2119 of file Panzer_STK_Interface.hpp.
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.
[in] | elements | Element entities to construct nodes |
[out] | nodes | Output array that will be sized (localIds.size() ,#Nodes,#Dim) |
localIds
is 0, the function will silently return Definition at line 2143 of file Panzer_STK_Interface.hpp.
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.
[in] | localIds | Element local IDs to construct nodes |
[in] | eBlock | Element block the elements are in |
[out] | nodes | Output array that will be sized (localIds.size() ,#Nodes,#Dim) |
localIds
is 0, the function will silently return Definition at line 2167 of file Panzer_STK_Interface.hpp.
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.
[in] | elements | Element entities to construct nodes |
[in] | eBlock | Element block the elements are in |
[out] | nodes | Output array that will be sized (localIds.size() ,#Nodes,#Dim) |
localIds
is 0, the function will silently return Definition at line 2156 of file Panzer_STK_Interface.hpp.
|
inline |
Definition at line 1159 of file Panzer_STK_Interface.hpp.
|
inline |
Definition at line 1160 of file Panzer_STK_Interface.hpp.
|
inline |
Definition at line 1161 of file Panzer_STK_Interface.hpp.
|
inline |
Definition at line 1162 of file Panzer_STK_Interface.hpp.
|
inline |
Definition at line 1163 of file Panzer_STK_Interface.hpp.
void panzer_stk::STK_Interface::initializeFromMetaData | ( | ) |
Build fields and parts from the meta data
Definition at line 1676 of file Panzer_STK_Interface.cpp.
void panzer_stk::STK_Interface::buildLocalElementIDs | ( | ) |
Setup local element IDs
Definition at line 1700 of file Panzer_STK_Interface.cpp.
void panzer_stk::STK_Interface::buildLocalEdgeIDs | ( | ) |
Setup local edge IDs
Definition at line 1764 of file Panzer_STK_Interface.cpp.
void panzer_stk::STK_Interface::buildLocalFaceIDs | ( | ) |
Setup local face IDs
Definition at line 1784 of file Panzer_STK_Interface.cpp.
|
inline |
Return a vector containing all the periodic boundary conditions.
Definition at line 1184 of file Panzer_STK_Interface.hpp.
|
inline |
Return a vector containing all the periodic boundary conditions.
Definition at line 1190 of file Panzer_STK_Interface.hpp.
|
inline |
Return a flag indicating if the bounding box search is used when matching periodic Ids.
Definition at line 1195 of file Panzer_STK_Interface.hpp.
|
inline |
Set the periodic search flag. Indicates if the bounding box search is used
Definition at line 1199 of file Panzer_STK_Interface.hpp.
|
inline |
Add a periodic boundary condition.
Definition at line 1207 of file Panzer_STK_Interface.hpp.
|
inline |
Add a set of periodic boundary conditions.
Definition at line 1215 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 1827 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 1866 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 1873 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 1923 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 1963 of file Panzer_STK_Interface.cpp.
|
inline |
Get the value of the time-stamp the last time this object was written to Exodus (default 0.0)
Definition at line 1243 of file Panzer_STK_Interface.hpp.
|
inline |
Get the value of the time-stamp when this object was created (default 0.0) or when setInitialStateTime was called.
Definition at line 1250 of file Panzer_STK_Interface.hpp.
|
inline |
Set the value of the initial time-stamp
Definition at line 1256 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 1990 of file Panzer_STK_Interface.cpp.
|
inline |
Set the weight for a particular element block. Larger means more costly to assemble and evaluate.
Definition at line 1265 of file Panzer_STK_Interface.hpp.
|
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 1274 of file Panzer_STK_Interface.hpp.
|
inline |
Return the use field coordinates flag
Definition at line 1278 of file Panzer_STK_Interface.hpp.
|
inline |
Use lower case (or not) for I/O
Definition at line 1282 of file Panzer_STK_Interface.hpp.
|
inline |
Use lower case (or not) for I/O
Definition at line 1286 of file Panzer_STK_Interface.hpp.
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).
[in] | elements | Element entities to construct vertices |
[in] | eBlock | Element block the elements are in |
[out] | vertices | Output array that will be sized (localIds.size() ,#Vertices,#Dim) |
localIds
is 0, the function will silently return Definition at line 1953 of file Panzer_STK_Interface.hpp.
void panzer_stk::STK_Interface::getElementVertices_FromFieldNoResize | ( | const std::vector< stk::mesh::Entity > & | elements, |
const std::string & | eBlock, | ||
ArrayT & | vertices | ||
) | const |
Definition at line 2006 of file Panzer_STK_Interface.hpp.
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.
[in] | elements | Element entities to construct vertices |
[out] | vertices | Output array that will be sized (localIds.size() ,#Vertices,#Dim) |
localIds
is 0, the function will silently return Definition at line 1856 of file Panzer_STK_Interface.hpp.
void panzer_stk::STK_Interface::getElementVertices_FromCoordsNoResize | ( | const std::vector< stk::mesh::Entity > & | elements, |
ArrayT & | vertices | ||
) | const |
Definition at line 1907 of file Panzer_STK_Interface.hpp.
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).
[in] | elements | Element entities to construct nodes |
[in] | eBlock | Element block the elements are in |
[out] | nodes | Output array that will be sized (localIds.size() ,#Nodes,#Dim) |
localIds
is 0, the function will silently return Definition at line 2282 of file Panzer_STK_Interface.hpp.
void panzer_stk::STK_Interface::getElementNodes_FromFieldNoResize | ( | const std::vector< stk::mesh::Entity > & | elements, |
const std::string & | eBlock, | ||
ArrayT & | nodes | ||
) | const |
Definition at line 2336 of file Panzer_STK_Interface.hpp.
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.
[in] | elements | Element entities to construct nodes |
[out] | nodes | Output array that will be sized (localIds.size() ,#Nodes,#Dim) |
localIds
is 0, the function will silently return Definition at line 2185 of file Panzer_STK_Interface.hpp.
void panzer_stk::STK_Interface::getElementNodes_FromCoordsNoResize | ( | const std::vector< stk::mesh::Entity > & | elements, |
ArrayT & | nodes | ||
) | const |
Definition at line 2236 of file Panzer_STK_Interface.hpp.
void panzer_stk::STK_Interface::refineMesh | ( | const int | numberOfLevels, |
const bool | deleteParentElements | ||
) |
Uniformly refine the mesh using Percept
[in] | numberOfLevels | Number of uniform refinement levels to apply. Must be >=1. |
[in] | deleteParentElements | If true, deletes the parent elements from the mesh to save memory. |
Definition at line 2025 of file Panzer_STK_Interface.cpp.
|
protected |
Compute global entity counts.
Definition at line 1001 of file Panzer_STK_Interface.cpp.
|
protected |
Compute global entity counts.
Definition at line 1007 of file Panzer_STK_Interface.cpp.
|
protected |
Initialize STK fields using a map (allocate space for storage and writing) to a specific entity rank.
Definition at line 368 of file Panzer_STK_Interface.cpp.
|
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 1979 of file Panzer_STK_Interface.cpp.
|
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 1743 of file Panzer_STK_Interface.cpp.
|
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 1805 of file Panzer_STK_Interface.cpp.
|
protected |
Writes a particular field to an array as a coordinate diplacement. Notice this is setup to work with the worksets associated with Panzer.
[in] | fieldName | Name of field to be filled |
[in] | blockId | Name of block this set of elements belongs to |
[in] | dimension | What coordinate axis to write to |
[in] | localElementIds | Local element IDs for this set of solution values |
[in] | solutionValues | An two dimensional array object sized by (Cells,Basis Count) |
Definition at line 1600 of file Panzer_STK_Interface.hpp.
|
static |
Definition at line 1365 of file Panzer_STK_Interface.hpp.
|
static |
Definition at line 1366 of file Panzer_STK_Interface.hpp.
|
static |
Definition at line 1367 of file Panzer_STK_Interface.hpp.
|
static |
Definition at line 1368 of file Panzer_STK_Interface.hpp.
|
static |
Definition at line 1369 of file Panzer_STK_Interface.hpp.
|
static |
Definition at line 1370 of file Panzer_STK_Interface.hpp.
|
protected |
Definition at line 1425 of file Panzer_STK_Interface.hpp.
|
protected |
Definition at line 1426 of file Panzer_STK_Interface.hpp.
|
protected |
Definition at line 1428 of file Panzer_STK_Interface.hpp.
|
protected |
Definition at line 1429 of file Panzer_STK_Interface.hpp.
|
protected |
Definition at line 1435 of file Panzer_STK_Interface.hpp.
|
protected |
Definition at line 1436 of file Panzer_STK_Interface.hpp.
|
protected |
Definition at line 1437 of file Panzer_STK_Interface.hpp.
|
protected |
Definition at line 1438 of file Panzer_STK_Interface.hpp.
|
protected |
Definition at line 1439 of file Panzer_STK_Interface.hpp.
|
protected |
Definition at line 1441 of file Panzer_STK_Interface.hpp.
|
protected |
Definition at line 1444 of file Panzer_STK_Interface.hpp.
|
protected |
Definition at line 1445 of file Panzer_STK_Interface.hpp.
|
protected |
Definition at line 1446 of file Panzer_STK_Interface.hpp.
|
protected |
Definition at line 1447 of file Panzer_STK_Interface.hpp.
|
protected |
Definition at line 1448 of file Panzer_STK_Interface.hpp.
|
protected |
Definition at line 1449 of file Panzer_STK_Interface.hpp.
|
protected |
Definition at line 1451 of file Panzer_STK_Interface.hpp.
|
protected |
Definition at line 1452 of file Panzer_STK_Interface.hpp.
|
protected |
Definition at line 1453 of file Panzer_STK_Interface.hpp.
|
protected |
Definition at line 1454 of file Panzer_STK_Interface.hpp.
|
protected |
Definition at line 1455 of file Panzer_STK_Interface.hpp.
|
protected |
Definition at line 1458 of file Panzer_STK_Interface.hpp.
|
protected |
Definition at line 1459 of file Panzer_STK_Interface.hpp.
|
protected |
Definition at line 1460 of file Panzer_STK_Interface.hpp.
|
protected |
Definition at line 1461 of file Panzer_STK_Interface.hpp.
|
protected |
Definition at line 1464 of file Panzer_STK_Interface.hpp.
|
protected |
Definition at line 1466 of file Panzer_STK_Interface.hpp.
|
protected |
Definition at line 1468 of file Panzer_STK_Interface.hpp.
|
protected |
Definition at line 1471 of file Panzer_STK_Interface.hpp.
|
protected |
Definition at line 1474 of file Panzer_STK_Interface.hpp.
|
protected |
Definition at line 1476 of file Panzer_STK_Interface.hpp.
|
protected |
Definition at line 1477 of file Panzer_STK_Interface.hpp.
|
protected |
Definition at line 1479 of file Panzer_STK_Interface.hpp.
|
protected |
Definition at line 1481 of file Panzer_STK_Interface.hpp.
|
protected |
Definition at line 1482 of file Panzer_STK_Interface.hpp.
|
mutableprotected |
Definition at line 1526 of file Panzer_STK_Interface.hpp.
|
mutableprotected |
Definition at line 1529 of file Panzer_STK_Interface.hpp.
|
mutableprotected |
Definition at line 1532 of file Panzer_STK_Interface.hpp.
|
protected |
Definition at line 1535 of file Panzer_STK_Interface.hpp.
|
protected |
Definition at line 1537 of file Panzer_STK_Interface.hpp.
|
protected |
Definition at line 1538 of file Panzer_STK_Interface.hpp.
|
protected |
Definition at line 1539 of file Panzer_STK_Interface.hpp.
|
protected |
Definition at line 1544 of file Panzer_STK_Interface.hpp.
|
protected |
Definition at line 1545 of file Panzer_STK_Interface.hpp.
|
protected |
Definition at line 1547 of file Panzer_STK_Interface.hpp.
|
protected |
Definition at line 1549 of file Panzer_STK_Interface.hpp.