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

The FiniteElement object is used to represent a mesh element \( \Omega^e \) corresponding to a mesh \( \mathcal{M} \) . More...

#include </home/docs/checkouts/readthedocs.org/user_builds/axom/checkouts/v0.5.0/src/axom/mint/fem/FiniteElement.hpp>

Public Member Functions

 FiniteElement ()=delete
 Default constructor. Disabled. More...
 
 FiniteElement (numerics::Matrix< double > &M, CellType cellType, bool useExternal=false)
 Constructs a FiniteElement instance from a matrix consisting of the element coordinates and a cell type. More...
 
 ~FiniteElement ()
 Destructor. More...
 
bool usesExternalBuffer () const
 Checks if this instance is pointing to an external, user-supplied, buffer for the element coordinate information. More...
 
void setMaxSolverIterations (int numIters)
 Overrides the max number of iterations for the Newton-Raphson, used for the isoparametric inverse mapping from physical coordinates to reference coordinates. More...
 
int getMaxSolverIterations () const
 Returns the max number of iterations used for the Newton-Raphson. More...
 
CellType getCellType () const
 Returns the cell type associated with this FiniteElement instance. More...
 
int getBasisType () const
 Returns the Basis associated with this FiniteElement instance. More...
 
int getPhysicalDimension () const
 Returns the physical dimension for this FiniteElement instance. More...
 
int getReferenceDimension () const
 Returns the dimension of the element in reference space. More...
 
int getNumNodes () const
 Returns the number of nodes, i.e., degrees-of-freedom (dofs) of the FiniteElement. More...
 
int getNumDofs () const
 Returns the number of degrees of freedom associated with the Finite Element basis bound to this FiniteElement instance. More...
 
int computeReferenceCoords (const double *xp, double *xr, double TOL=1.e-12)
 Given a point in physical space, \( \bar{x_p} \), this method computes the corresponding reference coordinates, \( \bar{xi} \) in the reference space of the element. More...
 
void computePhysicalCoords (const double *xr, double *xp)
 Given reference coordinates, \( \xi \in \bar{\Omega}^e \), this method computes the corresponding physical coordinates, given by \( \mathbf{x}(\xi)=\sum\limits_{i=0}^n{N}_i^e({\xi}){x_i} \). More...
 
void jacobian (const double *xr, numerics::Matrix< double > &J)
 Given reference coordinates, \( \xi \in \bar{\Omega}^e \), this method computes the jacobian, \( \mathcal{J}(\xi) \). More...
 
void evaluateShapeFunctions (const double *xr, double *phi)
 Given reference coordinates, \( \xi \in \bar{\Omega}^e \), this method evaluates the shape functions at each degree of freedom, \( \phi_i(\xi)=N_i^e(\xi) \). More...
 
void evaluateDerivatives (const double *xr, double *phidot)
 Given reference coordinates, \( \xi \in \bar{\Omega}^e \), this method evaluates the first derivatives of the shape function at each degree of freedom, \( \partial\phi_i(\xi)=N_i^e(\xi) \). More...
 
double * getPhysicalNodes ()
 Returns a pointer to the nodes of the element in physical space, \( \forall\hat{x_i} \in \Omega^e \). More...
 
const double * getPhysicalNodes () const
 
double * getReferenceNodes ()
 Returns a pointer to the nodes of the element in reference space, \( \forall\hat{\xi_i} \in \bar{\Omega}^e \). More...
 
const double * getReferenceNodes () const
 
double * getReferenceCenter ()
 Returns a pointer to the centroid of the reference element \( \xi_c \in \bar{\Omega}^e \). More...
 
const double * getReferenceCenter () const
 

Friends

Friend Functions
template<int BasisType, CellType CellType>
void bind_basis (FiniteElement &fe)
 Binds the given FiniteElement instance to a FiniteElement basis. More...
 

Detailed Description

The FiniteElement object is used to represent a mesh element \( \Omega^e \) corresponding to a mesh \( \mathcal{M} \) .

The FiniteElement is associated with a Finite Element Basis (FEBasis), which defines a set of shape functions \( \mathcal{N}_i( \hat{\xi} ) \), on a corresponding element, \( \bar{\Omega}^e \) in reference space. This class provides the following key operations:

  • Evaluate the shape function and its derivatives
  • Forward/Inverse mapping between physical/reference coordinates
  • Compute the jacobian \( \mathcal{J}(\hat{\xi} ) \)

The FiniteElement class is the primary object that application codes will use to perform these operations. Once a FiniteElement object is instantiated it can then be bound to a Finite Element Basis by invoking the bind_basis( ) method on the target FiniteElement instance, which is templated on two enum values: (a) the BasisType, defined in (FEBasisTypes.hpp) and (b) the CellType defined in(CellType.hpp). The rationale behind this is to insulate the user from the low-level classes and provide a simple and unified interface to Finite Element operations.

A simple example illustrating how to use the FiniteElement class is given below:

using axom;
double xp[2]; // test query point
double xr[2]; // query point reference coordinates
double wgts[4]; // buffer to store interpolation weights
...
mint::Mesh* m = get_mesh();
const int N = m->getMeshNumberOfCells();
const bool zero_copy = true;
for ( int idx=0; idx < N; ++idx ) {
...
// gather cell coordinates in a buffer
double coords[ ] = {
x1[ idx ], y1[ idx ],
x2[ idx ], y2[ idx ],
x3[ idx ], y3[ idx ],
x4[ idx ], y4[ idx ],
};
// put cell coordinates in matrix form
numeric::Matrix< double > m( 2, 4, coords, zero_copy );
// construct finite element and bind to basis
mint::FiniteElement fe( m, mint::QUAD, zero_copy );
mint::bind_basis< MINT_LAGRANGE_BASIS, mint::QUAD >( fe );
// compute reference coordinates
int status = fe.computeReferenceCoords( xp, xr );
if ( status == mint::INSIDE_ELEMENT ) {
fe.evaluateShapeFunction( xr, wgts );
....
} // END if point is inside fe
} // END for all cells
...
See also
FEBasis
CellType.hpp
FEBasisTypes.hpp
ShapeFunction

Constructor & Destructor Documentation

◆ FiniteElement() [1/2]

axom::mint::FiniteElement::FiniteElement ( )
delete

Default constructor. Disabled.

◆ FiniteElement() [2/2]

axom::mint::FiniteElement::FiniteElement ( numerics::Matrix< double > &  M,
CellType  cellType,
bool  useExternal = false 
)

Constructs a FiniteElement instance from a matrix consisting of the element coordinates and a cell type.

Parameters
[in]M(ndims x nnodes) matrix consisting of the element coordinates
[in]cellTypethe celltype, e.g., MINT_QUAD, etc.
[in]useExternaloptional argument that indicates whether the coordinate data will be shallow copied to the internal data-structures. Defaults to false if argument is not specified.
Note
The element coordinates are arranged in the supplied matrix, M, such that the number of rows corresponds to the dimension of the element and number of columns corresponds to the number of nodes.

◆ ~FiniteElement()

axom::mint::FiniteElement::~FiniteElement ( )

Destructor.

Member Function Documentation

◆ usesExternalBuffer()

bool axom::mint::FiniteElement::usesExternalBuffer ( ) const
inline

Checks if this instance is pointing to an external, user-supplied, buffer for the element coordinate information.

Returns
status true if an external buffer is used, otherwise, false.

◆ setMaxSolverIterations()

void axom::mint::FiniteElement::setMaxSolverIterations ( int  numIters)
inline

Overrides the max number of iterations for the Newton-Raphson, used for the isoparametric inverse mapping from physical coordinates to reference coordinates.

Note
Each Basis/Element pair prescribes a default value for the maximum number of Newton-Raphson iterations. This method provides the flexibility of overriding this value.
Warning
If bind_basis() is called after invoking this method, the max newton iterations would be set to the nominal value. Typically, this method should be called after bind_basis() is invoked.
Parameters
[in]Nuser-supplied number for

◆ getMaxSolverIterations()

int axom::mint::FiniteElement::getMaxSolverIterations ( ) const
inline

Returns the max number of iterations used for the Newton-Raphson.

Returns
N max number of iterations for the Newton-Raphson

◆ getCellType()

CellType axom::mint::FiniteElement::getCellType ( ) const
inline

Returns the cell type associated with this FiniteElement instance.

Returns
ctype the corresponding cell type.
Postcondition
(ctype != MINT_UNDEFINED_CELL) && (ctype < MINT_NUM_CELL_TYPES)
See also
CellType.hpp for a list of possible return values.

Referenced by axom::mint::bind_basis().

◆ getBasisType()

int axom::mint::FiniteElement::getBasisType ( ) const
inline

Returns the Basis associated with this FiniteElement instance.

Returns
basis the basis type, may be MINT_UNDEFINED_BASIS if no basis is bound to this FiniteElement.
See also
FEBasis

◆ getPhysicalDimension()

int axom::mint::FiniteElement::getPhysicalDimension ( ) const
inline

Returns the physical dimension for this FiniteElement instance.

Returns
dim the dimension of this FiniteElement.

◆ getReferenceDimension()

int axom::mint::FiniteElement::getReferenceDimension ( ) const
inline

Returns the dimension of the element in reference space.

Returns
ref_dim the reference dimension of the element.
Postcondition
ref_dim \( \in [1,3] \)

◆ getNumNodes()

int axom::mint::FiniteElement::getNumNodes ( ) const
inline

Returns the number of nodes, i.e., degrees-of-freedom (dofs) of the FiniteElement.

Returns
N the number of nodes of the element

◆ getNumDofs()

int axom::mint::FiniteElement::getNumDofs ( ) const
inline

Returns the number of degrees of freedom associated with the Finite Element basis bound to this FiniteElement instance.

Returns
N the number of degrees of freedom, or -1 if no basis is set.

◆ getPhysicalNodes() [1/2]

double* axom::mint::FiniteElement::getPhysicalNodes ( )
inline

Returns a pointer to the nodes of the element in physical space, \( \forall\hat{x_i} \in \Omega^e \).

Returns
ptr pointer to the nodes of the elements in physical space.
Note
The nodes of the element are arranged in a flat array using a column-major layout. It is convenient to access the nodes of the element by wrapping the return pointer in a matrix object with getNumNodes() columns and getPhysicalDimension() rows. This can be achieved as follows:
...
double* nodesptr = fe->getPhysicalNodes( );
numerics::Matrix< double > pnodes( ndims, nnodes, nodesptr, true );
for ( int inode=0; inode < nnodes; ++inode ) {
double* node = pnodes.getColumn( inode );
std::cout << "Node [" << inode << "]: ";
for ( int idim=0; idim < ndims << ++idim ) {
std::cout << nodes[ idim ] << " ";
} // END
std::cout << std::endl;
} // END for all nodes
...
See also
numerics::Matrix
Postcondition
ptr != nullptr

◆ getPhysicalNodes() [2/2]

const double* axom::mint::FiniteElement::getPhysicalNodes ( ) const
inline

◆ getReferenceNodes() [1/2]

double* axom::mint::FiniteElement::getReferenceNodes ( )
inline

Returns a pointer to the nodes of the element in reference space, \( \forall\hat{\xi_i} \in \bar{\Omega}^e \).

Returns
ptr pointer to the nodes of the element in reference space.
Note
The nodes of the element are arranged in a flat array using a column-major layout. It is convenient to access the nodes of the element by wrapping the return pointer in a matrix object as follows:
...
double* nodesptr = fe->getReferenceNodes( );
numerics::Matrix< double > pnodes( ndims, nnodes, nodesptr, true );
for ( int inode=0; inode < nnodes; ++inode ) {
double* node = pnodes.getColumn( inode );
std::cout << "Node [" << inode << "]: ";
for ( int idim=0; idim < ndims << ++idim ) {
std::cout << nodes[ idim ] << " ";
} // END
std::cout << std::endl;
} // END for all nodes
...
Postcondition
ptr == nullptr iff getBasisType()==MINT_UNDEFINED_BASIS

◆ getReferenceNodes() [2/2]

const double* axom::mint::FiniteElement::getReferenceNodes ( ) const
inline

◆ getReferenceCenter() [1/2]

double* axom::mint::FiniteElement::getReferenceCenter ( )
inline

Returns a pointer to the centroid of the reference element \( \xi_c \in \bar{\Omega}^e \).

Returns
ptr pointer to the reference center.
Postcondition
ptr == nullptr iff getBasisType()==MINT_UNDEFINED_BASIS
Note
ptr points to a buffer that is ndims long, where ndims is the dimension of the reference element.

◆ getReferenceCenter() [2/2]

const double* axom::mint::FiniteElement::getReferenceCenter ( ) const
inline

◆ computeReferenceCoords()

int axom::mint::FiniteElement::computeReferenceCoords ( const double *  xp,
double *  xr,
double  TOL = 1.e-12 
)

Given a point in physical space, \( \bar{x_p} \), this method computes the corresponding reference coordinates, \( \bar{xi} \) in the reference space of the element.

Parameters
[in]xpphysical coordinates of the point in query.
[out]xrcomputed reference coordinates \( \bar{\xi} \)
[in]TOLoptional tolerance for Newton-Raphson. Default is 1.e-12.
Returns
rc return code
Precondition
xp != nullptr
xr != nullptr
this->getBasisType() != MINT_UNDEFINED_BASIS

Referenced by getReferenceCenter().

◆ computePhysicalCoords()

void axom::mint::FiniteElement::computePhysicalCoords ( const double *  xr,
double *  xp 
)

Given reference coordinates, \( \xi \in \bar{\Omega}^e \), this method computes the corresponding physical coordinates, given by \( \mathbf{x}(\xi)=\sum\limits_{i=0}^n{N}_i^e({\xi}){x_i} \).

Parameters
[in]xrbuffer consisting of the reference coordinates
[out]ptbuffer to store the computed physical coordinates
Precondition
xr != nullptr
xp != nullptr
this->getBasisType() != MINT_UNDEFINED_BASIS

Referenced by getReferenceCenter().

◆ jacobian()

void axom::mint::FiniteElement::jacobian ( const double *  xr,
numerics::Matrix< double > &  J 
)

Given reference coordinates, \( \xi \in \bar{\Omega}^e \), this method computes the jacobian, \( \mathcal{J}(\xi) \).

Parameters
[in]xrbuffer consisting of the reference coordinates
[out]Jthe jacobian matrix evaluated at the given location
Precondition
xr != nullptr
this->getBasisType() != MINT_UNDEFINED_BASIS

Referenced by getReferenceCenter().

◆ evaluateShapeFunctions()

void axom::mint::FiniteElement::evaluateShapeFunctions ( const double *  xr,
double *  phi 
)

Given reference coordinates, \( \xi \in \bar{\Omega}^e \), this method evaluates the shape functions at each degree of freedom, \( \phi_i(\xi)=N_i^e(\xi) \).

Parameters
[in]xrreference coordinates where the shape functions is evaluated
[out]phibuffer to store the shape functions \( \phi_i \)
Note
The user-supplied output buffer has to be at least ndofs long
Precondition
xr != nullptr
phi != nullptr

Referenced by getReferenceCenter().

◆ evaluateDerivatives()

void axom::mint::FiniteElement::evaluateDerivatives ( const double *  xr,
double *  phidot 
)

Given reference coordinates, \( \xi \in \bar{\Omega}^e \), this method evaluates the first derivatives of the shape function at each degree of freedom, \( \partial\phi_i(\xi)=N_i^e(\xi) \).

Parameters
[in]xrreference coordinates where derivatives are computed
[out]phidotbuffer to store the shape function derivatives
Note
The output buffer
The derivatives in the output buffer, phidot, are arranged using a column-major flat array layout. It is therefore convenient to think of the data as a Matrix, where each row represents a degree of freedom and each column a dimension.
The matrix class may be used as a way to access the derivative information more conveniently as illustrated below:
...
fe->evaluateDerivatives( xr, phidot);
numerics::Matrix< double > shape_derivs( ndofs, ndims, phidot, true );
for ( int idim=0; idim < ndims; ++idim ) {
double* derivs = pnodes.getColumn( idim );
std::cout << "dimension => " << idim << ":\n";
for ( int idof=0; idof < ndofs << ++idof ) {
std::cout << "\t dof[" << idof << "]=" << derivs[ idof ] << "\n";
} // END for each dof
} // END for each dimension
...
Precondition
xr != nullptr
phidot != nullptr

Referenced by getReferenceCenter().

Friends And Related Function Documentation

◆ bind_basis

template<int BasisType, CellType CellType>
void bind_basis ( FiniteElement fe)
friend

Binds the given FiniteElement instance to a FiniteElement basis.

Parameters
[in]fereference to a finite element object.
Template Parameters
BasisTypethe basis type, e.g., MINT_LAGRANGE_BASIS
CellTypethe cell type, e.g., MINT_QUAD, etc.
Postcondition
fe.getBasisType() != MINT_UNDEFINED_BASIS

Referenced by getReferenceCenter().


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