AXOM
Axom provides a robust, flexible software infrastructure for the development of multi-physics applications and computational tools.
axom::mint::Mesh Class Referenceabstract

Base class that defines the core API common to all Mesh types. More...

#include </home/docs/checkouts/readthedocs.org/user_builds/axom/checkouts/latest/src/axom/mint/mesh/Mesh.hpp>

Inheritance diagram for axom::mint::Mesh:

Public Member Functions

 Mesh ()=delete
 Default constructor. Disabled. More...
 
Virtual methods
virtual ~Mesh ()
 Destructor. More...
 
Cells
virtual IndexType getNumberOfCells () const =0
 Returns the number of cells in this mesh instance. More...
 
virtual IndexType getCellCapacity () const
 Returns the capacity for number of cell in this mesh instance. More...
 
virtual CellType getCellType (IndexType cellID=0) const =0
 Return the type of the given cell. More...
 
virtual IndexType getNumberOfCellNodes (IndexType cellID=0) const =0
 Return the number of nodes associated with the given cell. More...
 
virtual IndexType getCellNodeIDs (IndexType AXOM_UNUSED_PARAM(cellID), IndexType *AXOM_UNUSED_PARAM(nodes)) const =0
 Copy the connectivity of the given cell into the provided buffer. The buffer must be of length at least getNumberOfCellNodes( cellID ). More...
 
virtual IndexType getNumberOfCellFaces (IndexType AXOM_UNUSED_PARAM(cellID)=0) const =0
 Return the number of faces associated with the given cell. More...
 
virtual IndexType getCellFaceIDs (IndexType AXOM_UNUSED_PARAM(cellID), IndexType *AXOM_UNUSED_PARAM(faces)) const =0
 Populates the given buffer with the IDs of the faces of the given cell and returns the number of faces. More...
 
Nodes
virtual IndexType getNumberOfNodes () const =0
 Returns the number of nodes in this mesh instance. More...
 
virtual IndexType getNodeCapacity () const
 Returns the capacity for number of nodes in this mesh instance. More...
 
virtual void getNode (IndexType nodeID, double *node) const =0
 Copy the coordinates of the given node into the provided buffer. More...
 
virtual double * getCoordinateArray (int dim)=0
 Returns pointer to the requested mesh coordinate buffer. More...
 
virtual const double * getCoordinateArray (int dim) const =0
 
Faces
virtual IndexType getNumberOfFaces () const =0
 Returns the number of faces in this mesh instance. More...
 
virtual IndexType getFaceCapacity () const
 Returns the capacity for number of faces in this mesh instance. More...
 
virtual CellType getFaceType (IndexType AXOM_UNUSED_PARAM(faceID)) const =0
 Return the type of the given face. More...
 
virtual IndexType getNumberOfFaceNodes (IndexType AXOM_UNUSED_PARAM(faceID)) const =0
 Return the number of nodes associated with the given face. More...
 
virtual IndexType getFaceNodeIDs (IndexType AXOM_UNUSED_PARAM(faceID), IndexType *AXOM_UNUSED_PARAM(nodes)) const =0
 Copy the IDs of the nodes that compose the given face into the provided buffer. More...
 
virtual void getFaceCellIDs (IndexType AXOM_UNUSED_PARAM(faceID), IndexType &AXOM_UNUSED_PARAM(cellIDOne), IndexType &AXOM_UNUSED_PARAM(cellIDTwo)) const =0
 Copy the IDs of the cells adjacent to the given face into the provided indices. More...
 
Edges
virtual IndexType getNumberOfEdges () const =0
 Returns the number of edges in this mesh instance. More...
 
virtual IndexType getEdgeCapacity () const
 Returns the capacity for number of edges in this mesh instance. More...
 
virtual bool isExternal () const =0
 Returns true iff the mesh was constructed with external arrays. More...
 
Mesh Attribute get/set Methods
int getDimension () const
 Returns the dimension for this mesh instance. More...
 
int getBlockId () const
 Returns the ID of this mesh instance. More...
 
void setBlockId (int ID)
 set the block ID of this mesh instance. More...
 
int getPartitionId () const
 Returns the partition ID of this mesh instance. More...
 
void setPartitionId (int ID)
 set the partition ID of this mesh instance. More...
 
int getMeshType () const
 Returns the mesh type of this mesh instance. More...
 
bool hasExplicitCoordinates () const
 Checks if this mesh instance has explicit coordinates. More...
 
bool hasExplicitConnectivity () const
 Checks if this mesh instance has explicit connectivity. More...
 
bool hasMixedCellTypes () const
 Checks if the mesh has mixed cell types, e.g., consisting of both triangle and quad elements or hex,pyramid,prisms and tets in 3-D. More...
 
bool isStructured () const
 Returns true if the mesh type is structured. More...
 
bool isUnstructured () const
 Returns true if the mesh type is unstructured. More...
 
bool hasSidreGroup () const
 Checks if this Mesh instance is associated with a Sidre Group. More...
 
Methods to Create, Access & Remove Fields from a Mesh
const FieldDatagetFieldData (int association) const
 Returns const pointer to the FieldData instance with the specified mesh field association, e.g., NODE_CENTERED, CELL_CENTERED, etc. More...
 
bool hasField (const std::string &name, int association=ANY_CENTERING) const
 Check if a field with the given name and association exists. More...
 
template<typename T >
T * createField (const std::string &name, int association, IndexType num_components=1, bool storeInSidre=true)
 Creates a new field with the given name and specified mesh field association, e.g., NODE_CENTERED, CELL_CENTERED, etc. More...
 
template<typename T >
T * createField (const std::string &name, int association, T *data, IndexType num_components=1, IndexType capacity=USE_DEFAULT)
 Creates a new field from an external buffer that has the given name and specified mesh field association, e.g., NODE_CENTERED, CELL_CENTERED, etc. More...
 
bool removeField (const std::string &name, int association)
 Removes the field with the given name and specified association. More...
 
template<typename T >
T * getFieldPtr (const std::string &name, int association, IndexType &num_components)
 Returns pointer to buffer of the field with the given ane and specified mesh field association. More...
 
template<typename T >
T * getFieldPtr (const std::string &name, int association)
 
template<typename T >
const T * getFieldPtr (const std::string &name, int association, IndexType &num_components) const
 
template<typename T >
const T * getFieldPtr (const std::string &name, int association) const
 

Protected Member Functions

Protected Constructors (used in derived classes )
 Mesh (int ndims, int type)
 Mesh Constructor. More...
 

Protected Attributes

Protected Members
int m_ndims
 
int m_type
 
int m_block_idx
 
int m_part_idx
 
bool m_explicit_coords
 
bool m_explicit_connectivity
 
bool m_has_mixed_topology
 
FieldDatam_mesh_fields [NUM_FIELD_ASSOCIATIONS]
 

Detailed Description

Base class that defines the core API common to all Mesh types.

A Mesh, \( \mathcal{M}(\Omega) \), provides an approximation of a physical domain, \( \Omega \in \mathcal{R}^d \), where \( d \in [1,3] \) . The Mesh is essentially a discrete representation of a problem and is used to facilitate the analysis, e.g., FEA, CFD, etc. The solution domain is approximated by dividing it into a finite number of nodes and cells at which the variables of the underlying mathematical model (i.e., a PDE) are then computed via a numerical method, such as, Finite Difference (FD), Finite Volume (FV) or Finite Element (FE), chief among them.

There are a variety of mesh types. Mint supports the following mesh types:

  • Structured (Curvilinear) Mesh

    A structured mesh divides the solution domain according to a logical grid where each node/cell of the mesh can be uniquely identified by a corresponding logical ijk-index. The nodes of the mesh are found at the intersection of the grid lines, but, are explicitly defined via a mapping onto the physical Cartesian coordinate system of the domain. For this reason, these types of meshes are also called mapped or body-fitted meshes.

    However, the mesh topology (e.g., connectivity, neighbor information) is implicitly defined by the logical indexing scheme. For example, a structured mesh is composed of quadrilateral cells in 2-D and hexahedron cells in 3-D. Given the logical index of a cell one can compute the node indices of the cell and neighbor information by performing simple shift operations on the associated index.

  • Rectilinear Mesh

    A rectilinear mesh , also known as a product mesh, is similar to the structured mesh in that it is also defined according to a logical indexing scheme and has implicit topology.

    However, the nodes and cells on a rectilinear mesh are arranged on a regular lattice that is axis-aligned with the Cartesian coordinate system. In contrast to the general structured mesh , the rectilinear mesh does not explicitly define all the nodes of the mesh. Instead, the nodes are only defined along each coordinate axis and may have variable spacing between nodes. Given a logical index, the corresponding physical position of a node can be evaluated by taking the cartesian product of the corresponding coordinate along each coordinate axis.

  • Uniform Mesh

    A uniform mesh , also called a regular mesh, subdivides the domain in cells that have uniform spacing across each coordinate axis. Similar to the structured mesh , a uniform mesh adheres to the same logical indexing scheme and implicit topology representation. Moreover, The nodes and cells of a uniform mesh are arranged on a regular lattice as with the rectilinear mesh . However, both topology and geometry is implicitly defined on a uniform mesh . The geometry is solely defined by an origin, \( \hat{x_0} \) and spacing, \( \hat{h} \), along each axis. The coordinates of a node can be evaluated algebraically by the following: \( \hat{p} = \hat{x_0} + \hat{i} \times \hat{h} \), where \(\hat{i}\) is the logical ijk-index of the corresponding node.

  • Unstructured Mesh

    An unstructured mesh stores both node and topology information explicitly. This allows the flexibility of discretizing the solution domain using a variety of cell types not just quadrilateral (in 2D) or hexahedral (in 3D) cells. Due to this added flexibility, the use of unstructured meshes is more common when dealing with complex geometries. However, unstructured meshes require additional storage and generally incur some performance penalty to store, create and access mesh topology information respectively.

    Mint classifies unstructured meshes in two basic types based on the underlying mesh topology:

    • Single Cell Topology

      In this case, the unstructured mesh consists of a single cell type, e.g., a quad or triangle mesh in 2D, or, a hex or tet mesh in 3D. In this case the underlying implementation is optimized for the specified cell type (specified in the constructor).

    • Mixed Cell Topology

      When mixed cell topology is specified, the unstructured mesh can be composed of any of the supported cell types, e.g., a mesh consisting of both quads and triangles. This mode incurs additional overhead for storage and access to mesh topology information, since it requires indirection.

    The list of supported cell types for an unstructured mesh is available in CellTypes.hpp

  • Particle Mesh

    A particle mesh discretizes the solution domain using a set of discrete particle elements which correspond to the the nodes of the mesh. There is no ordering imposed on the particles and the coordinates of the particles are explicitly defined. A particle mesh has no connectivity information connecting the particles, which is why in literature methods using a particle discretization are also referred to as meshless or meshfree methods.

The Mesh class provides the means to create, access and remove fields on a mesh given the field name and its association. The field association designates the corresponding mesh entity at which the field is stored, e.g. whether the field is stored at the nodes or cell centers. A Field may be a scalar quantity, e.g., pressure, or a vector field, such as, velocity.

Warning
When using Sidre, field names have to be unique. For example, if there exists a "pressure" node-centered field, there cannot be a corresponding cell-centered field.
Note
Mesh is a base class and cannot be instantiated directly
Typically, the computational mesh can be defined across one or more blocks , e.g., for multi-block problems, where each block is subsequently decomposed into several partitions for parallel computation. A Mesh instance represents a single partition for a given block.
See also
mint::UnstructuredMesh
mint::StructuredMesh
mint::CurvilinearMesh
mint::RectilinearMesh
mint::UniformMesh
mint::Field
mint::FieldData
mint::MeshTypes

Constructor & Destructor Documentation

◆ Mesh() [1/2]

axom::mint::Mesh::Mesh ( )
delete

Default constructor. Disabled.

◆ ~Mesh()

virtual axom::mint::Mesh::~Mesh ( )
virtual

Destructor.

◆ Mesh() [2/2]

axom::mint::Mesh::Mesh ( int  ndims,
int  type 
)
protected

Mesh Constructor.

Parameters
[in]ndimsthe number of dimensions
[in]typethe mesh type.
[in]blockIdthe block ID for this mesh instance.
[in]partIdthe partition ID for this mesh instance.

Member Function Documentation

◆ getNumberOfCells()

virtual IndexType axom::mint::Mesh::getNumberOfCells ( ) const
pure virtual

Returns the number of cells in this mesh instance.

Returns
N the number of cells
Postcondition
N >= 0

Implemented in axom::mint::UnstructuredMesh< TOPO >, axom::mint::StructuredMesh, and axom::mint::ParticleMesh.

◆ getCellCapacity()

virtual IndexType axom::mint::Mesh::getCellCapacity ( ) const
inlinevirtual

Returns the capacity for number of cell in this mesh instance.

Returns
N the cell capacity
Postcondition
N >= 0

Reimplemented in axom::mint::UnstructuredMesh< TOPO >, and axom::mint::ParticleMesh.

References getNumberOfCells().

◆ getCellType()

virtual CellType axom::mint::Mesh::getCellType ( IndexType  cellID = 0) const
pure virtual

Return the type of the given cell.

Parameters
[in]cellIDthe ID of the cell in question, this parameter is ignored if hasMixedCellTypes() == false.
Precondition
0 <= cellID < getNumberOfCells()

Implemented in axom::mint::UnstructuredMesh< TOPO >.

◆ getNumberOfCellNodes()

virtual IndexType axom::mint::Mesh::getNumberOfCellNodes ( IndexType  cellID = 0) const
pure virtual

Return the number of nodes associated with the given cell.

Parameters
[in]cellIDthe ID of the cell in question, this parameter is ignored unless hasMixedCellTypes() == true.
Precondition
0 <= cellID < getNumberOfCells()

Implemented in axom::mint::UnstructuredMesh< TOPO >.

◆ getCellNodeIDs()

virtual IndexType axom::mint::Mesh::getCellNodeIDs ( IndexType   AXOM_UNUSED_PARAMcellID,
IndexType AXOM_UNUSED_PARAMnodes 
) const
pure virtual

Copy the connectivity of the given cell into the provided buffer. The buffer must be of length at least getNumberOfCellNodes( cellID ).

Parameters
[in]cellIDthe ID of the cell in question.
[out]nodesthe buffer into which the connectivity is copied, must be of length at least getNumberOfCellNodes( cellID ).
Returns
The number of nodes for the given cell.
Precondition
nodes != nullptr
0 <= cellID < getNumberOfCells()

◆ getNumberOfCellFaces()

virtual IndexType axom::mint::Mesh::getNumberOfCellFaces ( IndexType   AXOM_UNUSED_PARAMcellID = 0) const
pure virtual

Return the number of faces associated with the given cell.

Parameters
[in]cellIDthe ID of the cell in question.

Implemented in axom::mint::StructuredMesh, and axom::mint::ParticleMesh.

◆ getCellFaceIDs()

virtual IndexType axom::mint::Mesh::getCellFaceIDs ( IndexType   AXOM_UNUSED_PARAMcellID,
IndexType AXOM_UNUSED_PARAMfaces 
) const
pure virtual

Populates the given buffer with the IDs of the faces of the given cell and returns the number of faces.

Parameters
[in]cellIDthe ID of the cellID in question.
[out]facesbuffer to populate with the face IDs. Must be of length at least getNumberOfCellFaces( cellID ).
Precondition
faces != nullptr
0 <= cellID < getNumberOfCells()

Implemented in axom::mint::ParticleMesh.

◆ getNumberOfNodes()

virtual IndexType axom::mint::Mesh::getNumberOfNodes ( ) const
pure virtual

Returns the number of nodes in this mesh instance.

Returns
N the number of nodes
Postcondition
N >= 0

Implemented in axom::mint::UnstructuredMesh< TOPO >, axom::mint::StructuredMesh, and axom::mint::ParticleMesh.

◆ getNodeCapacity()

virtual IndexType axom::mint::Mesh::getNodeCapacity ( ) const
inlinevirtual

Returns the capacity for number of nodes in this mesh instance.

Returns
N the node capacity
Postcondition
N >= 0

Reimplemented in axom::mint::UnstructuredMesh< TOPO >, and axom::mint::ParticleMesh.

References getNumberOfNodes().

◆ getNode()

virtual void axom::mint::Mesh::getNode ( IndexType  nodeID,
double *  node 
) const
pure virtual

Copy the coordinates of the given node into the provided buffer.

Parameters
[in]nodeIDthe ID of the node in question.
[in]coordsthe buffer to copy the coordinates into, of length at least getDimension().
Precondition
0 <= nodeID < getNumberOfNodes()
coords != nullptr

Implemented in axom::mint::StructuredMesh, axom::mint::UniformMesh, axom::mint::RectilinearMesh, axom::mint::ParticleMesh, axom::mint::UnstructuredMesh< TOPO >, and axom::mint::CurvilinearMesh.

◆ getCoordinateArray() [1/2]

virtual double* axom::mint::Mesh::getCoordinateArray ( int  dim)
pure virtual

Returns pointer to the requested mesh coordinate buffer.

Parameters
[in]dimthe dimension of the requested coordinate buffer
Returns
ptr pointer to the coordinate buffer.
Note
if hasExplicitCoordinates() == true then the length of the returned buffer is getNumberOfNodes(). Otherwise the UniformMesh returns nullptr and the RectilinearMesh returns a pointer to the associated dimension scale which is of length static_cast< RectilinearMesh* >( this )->getNodeResolution().
Precondition
dim >= 0 && dim < dimension()
dim == X_COORDINATE || dim == Y_COORDINATE || dim == Z_COORDINATE

Implemented in axom::mint::StructuredMesh, axom::mint::UnstructuredMesh< TOPO >, axom::mint::RectilinearMesh, axom::mint::ParticleMesh, and axom::mint::CurvilinearMesh.

◆ getCoordinateArray() [2/2]

virtual const double* axom::mint::Mesh::getCoordinateArray ( int  dim) const
pure virtual

◆ getNumberOfFaces()

virtual IndexType axom::mint::Mesh::getNumberOfFaces ( ) const
pure virtual

Returns the number of faces in this mesh instance.

Returns
N the number of faces
Postcondition
N >= 0

Implemented in axom::mint::UnstructuredMesh< TOPO >, axom::mint::StructuredMesh, and axom::mint::ParticleMesh.

◆ getFaceCapacity()

virtual IndexType axom::mint::Mesh::getFaceCapacity ( ) const
inlinevirtual

Returns the capacity for number of faces in this mesh instance.

Returns
N the face capacity
Postcondition
N >= 0

Reimplemented in axom::mint::UnstructuredMesh< TOPO >.

References getNumberOfFaces().

◆ getFaceType()

virtual CellType axom::mint::Mesh::getFaceType ( IndexType   AXOM_UNUSED_PARAMfaceID) const
pure virtual

Return the type of the given face.

Parameters
[in]faceIDthe ID of the face in question.

Implemented in axom::mint::StructuredMesh, and axom::mint::ParticleMesh.

◆ getNumberOfFaceNodes()

virtual IndexType axom::mint::Mesh::getNumberOfFaceNodes ( IndexType   AXOM_UNUSED_PARAMfaceID) const
pure virtual

Return the number of nodes associated with the given face.

Parameters
[in]faceIDthe ID of the face in question.

Implemented in axom::mint::StructuredMesh, and axom::mint::ParticleMesh.

◆ getFaceNodeIDs()

virtual IndexType axom::mint::Mesh::getFaceNodeIDs ( IndexType   AXOM_UNUSED_PARAMfaceID,
IndexType AXOM_UNUSED_PARAMnodes 
) const
pure virtual

Copy the IDs of the nodes that compose the given face into the provided buffer.

Parameters
[in]faceIDthe ID of the face in question.
[out]nodesthe buffer into which the node IDs are copied, must be of length at least getNumberOfFaceNodes().
Returns
The number of nodes for the given face.
Precondition
nodes != nullptr
0 <= faceID < getNumberOfFaces()

Implemented in axom::mint::ParticleMesh.

◆ getFaceCellIDs()

virtual void axom::mint::Mesh::getFaceCellIDs ( IndexType   AXOM_UNUSED_PARAMfaceID,
IndexType AXOM_UNUSED_PARAMcellIDOne,
IndexType AXOM_UNUSED_PARAMcellIDTwo 
) const
pure virtual

Copy the IDs of the cells adjacent to the given face into the provided indices.

Parameters
[in]faceIDthe ID of the face in question.
[out]cellIDOnethe ID of the first cell.
[out]cellIDTwothe ID of the second cell.
Note
If no cell exists (the face is external) then the ID will be set to -1.
Precondition
0 <= faceID < getNumberOfFaces()

Implemented in axom::mint::ParticleMesh.

◆ getNumberOfEdges()

virtual IndexType axom::mint::Mesh::getNumberOfEdges ( ) const
pure virtual

Returns the number of edges in this mesh instance.

Returns
N the number of edges
Postcondition
N >= 0

Implemented in axom::mint::UnstructuredMesh< TOPO >, axom::mint::StructuredMesh, and axom::mint::ParticleMesh.

◆ getEdgeCapacity()

virtual IndexType axom::mint::Mesh::getEdgeCapacity ( ) const
inlinevirtual

Returns the capacity for number of edges in this mesh instance.

Returns
N the edge capacity
Postcondition
N >= 0

Reimplemented in axom::mint::UnstructuredMesh< TOPO >, and axom::mint::ParticleMesh.

References getNumberOfEdges().

◆ isExternal()

virtual bool axom::mint::Mesh::isExternal ( ) const
pure virtual

Returns true iff the mesh was constructed with external arrays.

Returns
status true if the mesh points to external buffers, else, false.

Implemented in axom::mint::StructuredMesh, axom::mint::UnstructuredMesh< TOPO >, axom::mint::RectilinearMesh, axom::mint::ParticleMesh, and axom::mint::CurvilinearMesh.

◆ getDimension()

int axom::mint::Mesh::getDimension ( ) const
inline

Returns the dimension for this mesh instance.

Returns
ndims the dimension of this mesh instance.
Postcondition
ndims >= 1 && ndims <= 3

References m_ndims.

◆ getBlockId()

int axom::mint::Mesh::getBlockId ( ) const
inline

Returns the ID of this mesh instance.

Returns
Id the ID of the mesh.

References m_block_idx.

◆ setBlockId()

void axom::mint::Mesh::setBlockId ( int  ID)

set the block ID of this mesh instance.

Parameters
[in]IDthe new block ID.
Postcondition
getBlockId() == ID

◆ getPartitionId()

int axom::mint::Mesh::getPartitionId ( ) const
inline

Returns the partition ID of this mesh instance.

Returns
partitionId the partition ID of the mesh.

References m_part_idx.

◆ setPartitionId()

void axom::mint::Mesh::setPartitionId ( int  ID)

set the partition ID of this mesh instance.

Parameters
[in]IDthe new partition ID.
Postcondition
getPartitionId() == ID

◆ getMeshType()

int axom::mint::Mesh::getMeshType ( ) const
inline

Returns the mesh type of this mesh instance.

Returns
meshType the mesh type
See also
MeshType

References m_type.

◆ hasExplicitCoordinates()

bool axom::mint::Mesh::hasExplicitCoordinates ( ) const
inline

Checks if this mesh instance has explicit coordinates.

Returns
status true iff the mesh defines coordinates explicitly.

References m_explicit_coords.

◆ hasExplicitConnectivity()

bool axom::mint::Mesh::hasExplicitConnectivity ( ) const
inline

Checks if this mesh instance has explicit connectivity.

Returns
status true iff the mesh defines cell connectivity explicitly.

References m_explicit_connectivity.

◆ hasMixedCellTypes()

bool axom::mint::Mesh::hasMixedCellTypes ( ) const
inline

Checks if the mesh has mixed cell types, e.g., consisting of both triangle and quad elements or hex,pyramid,prisms and tets in 3-D.

Returns
status true iff the mesh has mixed cell types.

References m_has_mixed_topology.

◆ isStructured()

bool axom::mint::Mesh::isStructured ( ) const
inline

Returns true if the mesh type is structured.

Returns
status true if the mesh type is structured, else, false.

References m_type, axom::mint::STRUCTURED_CURVILINEAR_MESH, axom::mint::STRUCTURED_RECTILINEAR_MESH, and axom::mint::STRUCTURED_UNIFORM_MESH.

◆ isUnstructured()

bool axom::mint::Mesh::isUnstructured ( ) const
inline

Returns true if the mesh type is unstructured.

Returns
status true if the mesh type is unstructured, else, false.

References m_type, and axom::mint::UNSTRUCTURED_MESH.

◆ hasSidreGroup()

bool axom::mint::Mesh::hasSidreGroup ( ) const
inline

Checks if this Mesh instance is associated with a Sidre Group.

Returns
status true if the Mesh is associated with a group in a Sidre hierarchy, else, false.

◆ getFieldData()

const FieldData * axom::mint::Mesh::getFieldData ( int  association) const
inline

Returns const pointer to the FieldData instance with the specified mesh field association, e.g., NODE_CENTERED, CELL_CENTERED, etc.

Parameters
[in]associationthe specified mesh field association
Returns
fd pointer to the requested FieldData instance
Precondition
association >= 0 && association < NUM_FIELD_ASSOCIATION
Postcondition
fd != nullptr
See also
FieldAssociation
FieldData

References m_mesh_fields, m_type, axom::mint::NODE_CENTERED, axom::mint::NUM_FIELD_ASSOCIATIONS, axom::mint::PARTICLE_MESH, and SLIC_ERROR_IF.

◆ hasField()

bool axom::mint::Mesh::hasField ( const std::string &  name,
int  association = ANY_CENTERING 
) const
inline

Check if a field with the given name and association exists.

Parameters
[in]namethe name of the field in query.
[in]associationthe field association (optional)
Returns
status true if the field exists, else, false.
Note
If an association is not explicitly specified, the code will check if a field by the given name exists in any available centeering.
Precondition
name.empty()==false
association >= 0 && association < NUM_FIELD_ASSOCIATION
See also
FieldAssociation

References axom::mint::ANY_CENTERING, getFieldData(), axom::mint::FieldData::hasField(), m_type, axom::mint::NUM_FIELD_ASSOCIATIONS, axom::mint::PARTICLE_MESH, and SLIC_ASSERT.

◆ createField() [1/2]

template<typename T >
T * axom::mint::Mesh::createField ( const std::string &  name,
int  association,
IndexType  num_components = 1,
bool  storeInSidre = true 
)
inline

Creates a new field with the given name and specified mesh field association, e.g., NODE_CENTERED, CELL_CENTERED, etc.

Parameters
[in]namethe name of the new field.
[in]associationthe mesh field association.
[in]num_componentsnumber of components of the field (optional).
[in]storeInSidreindicates whether to store the field in the corresponding Sidre group (optional).
[in]capacity
Returns
ptr raw pointer to the data buffer of the new field.
Note
This method throws an error and aborts if any of the pre-conditions is not satisfied.
Precondition
name.empty() == false
hasField( name ) == false
association >= 0 && association < NUM_FIELD_ASSOCIATION
Postcondition
ptr != nullptr
hasField( name ) == true
See also
FieldAssociation

References axom::mint::FieldData::createField(), getFieldData(), hasField(), SLIC_ASSERT, and SLIC_ERROR_IF.

◆ createField() [2/2]

template<typename T >
T * axom::mint::Mesh::createField ( const std::string &  name,
int  association,
T *  data,
IndexType  num_components = 1,
IndexType  capacity = USE_DEFAULT 
)
inline

Creates a new field from an external buffer that has the given name and specified mesh field association, e.g., NODE_CENTERED, CELL_CENTERED, etc.

Parameters
[in]namethe name of the new field.
[in]associationthe mesh field association.
[in]datapointer to the external data buffer.
[in]num_componentsnumber of components of the field (optional).
Returns
ptr raw pointer to the data buffer of the new field.
Note
This method throws an error and aborts if any of the pre-conditions is not satisfied.
Precondition
name.empty() == false
hasField( name ) == false
data != nullptr
association >= 0 && association < NUM_FIELD_ASSOCIATION
Postcondition
ptr != nullptr
ptr == data
hasField( name ) == true
See also
FieldAssociation

References axom::mint::FieldData::createField(), getFieldData(), hasField(), SLIC_ASSERT, and SLIC_ERROR_IF.

◆ removeField()

bool axom::mint::Mesh::removeField ( const std::string &  name,
int  association 
)
inline

Removes the field with the given name and specified association.

Parameters
[in]namethe name of the field to remove.
[in]associationthe mesh field association.
Returns
status true if the field is removed successfully, else, false.
Precondition
name.emtpy() == false
association >= 0 && association < NUM_FIELD_ASSOCIATION
See also
FieldAssociation

References getFieldData(), axom::mint::FieldData::hasField(), hasField(), axom::mint::FieldData::removeField(), and SLIC_WARNING_IF.

◆ getFieldPtr() [1/4]

template<typename T >
T * axom::mint::Mesh::getFieldPtr ( const std::string &  name,
int  association,
IndexType num_components 
)
inline

Returns pointer to buffer of the field with the given ane and specified mesh field association.

Parameters
[in]namethe name of the requested field.
[in]associationthe mesh field association.
[out]num_componentsthe number of components per tuple (optional).
Returns
ptr raw pointer to the data buffer of the requested field.
Precondition
name.empty() == false
hasField( name )
association >= 0 && association < NUM_FIELD_ASSOCIATION
See also
FieldAssociation

References getFieldPtr().

◆ getFieldPtr() [2/4]

template<typename T >
T * axom::mint::Mesh::getFieldPtr ( const std::string &  name,
int  association 
)
inline

◆ getFieldPtr() [3/4]

template<typename T >
const T * axom::mint::Mesh::getFieldPtr ( const std::string &  name,
int  association,
IndexType num_components 
) const
inline

◆ getFieldPtr() [4/4]

template<typename T >
const T * axom::mint::Mesh::getFieldPtr ( const std::string &  name,
int  association 
) const
inline

Member Data Documentation

◆ m_ndims

int axom::mint::Mesh::m_ndims
protected

◆ m_type

int axom::mint::Mesh::m_type
protected

mesh dimension

◆ m_block_idx

int axom::mint::Mesh::m_block_idx
protected

the type of the mesh

◆ m_part_idx

int axom::mint::Mesh::m_part_idx
protected

the Block ID of the mesh

◆ m_explicit_coords

bool axom::mint::Mesh::m_explicit_coords
protected

the partition ID of the mesh

◆ m_explicit_connectivity

bool axom::mint::Mesh::m_explicit_connectivity
protected

◆ m_has_mixed_topology

bool axom::mint::Mesh::m_has_mixed_topology
protected

◆ m_mesh_fields

FieldData* axom::mint::Mesh::m_mesh_fields[NUM_FIELD_ASSOCIATIONS]
protected

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