Serac  0.1
Serac is an implicit thermal strucural mechanics simulation code.
Public Member Functions | Protected Attributes | List of all members
serac::FiniteElementState Class Reference

Class for encapsulating the critical MFEM components of a primal finite element field. More...

#include <finite_element_state.hpp>

Collaboration diagram for serac::FiniteElementState:
Collaboration graph
[legend]

Public Member Functions

 FiniteElementState (const FiniteElementState &rhs)
 Copy constructor. More...
 
 FiniteElementState (FiniteElementState &&rhs)
 Move construct a new Finite Element State object. More...
 
FiniteElementStateoperator= (const FiniteElementState &rhs)
 Copy assignment. More...
 
FiniteElementStateoperator= (FiniteElementState &&rhs)
 Move assignment. More...
 
FiniteElementStateoperator= (const mfem::HypreParVector &rhs)
 Copy assignment with HypreParVector. More...
 
FiniteElementStateoperator= (const mfem::Vector &rhs)
 Copy assignment with mfem::Vector. More...
 
FiniteElementStateoperator= (double rhs)
 Copy assignment with double. More...
 
void fillGridFunction (mfem::ParGridFunction &grid_function) const
 Fill a user-provided grid function based on the underlying true vector. More...
 
void setFromGridFunction (const mfem::ParGridFunction &grid_function)
 Initialize the true vector in the FiniteElementState based on an input grid function. More...
 
void project (mfem::VectorCoefficient &coef, mfem::Array< int > &dof_list)
 Project a vector coefficient onto a set of dofs. More...
 
void project (mfem::Coefficient &coef, mfem::Array< int > &dof_list, std::optional< int > component={})
 Project a scalar coefficient onto a set of dofs. More...
 
void project (const GeneralCoefficient &coef)
 
void project (mfem::Coefficient &coef)
 This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.
 
void project (mfem::VectorCoefficient &coef)
 This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.
 
void projectOnBoundary (mfem::Coefficient &coef, const mfem::Array< int > &markers)
 Project a coefficient on a specific set of marked boundaries. More...
 
void projectOnBoundary (mfem::VectorCoefficient &coef, const mfem::Array< int > &markers)
 This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.
 
mfem::ParGridFunction & gridFunction () const
 Construct a grid function from the finite element state true vector. More...
 
 FiniteElementVector (const mfem::ParFiniteElementSpace &space, const std::string &name="")
 Minimal constructor for a FiniteElementVector given a finite element space. More...
 
template<typename FunctionSpace >
 FiniteElementVector (mfem::ParMesh &mesh, FunctionSpace, const std::string &name="")
 Construct a new Finite Element Vector object given a templated function space. More...
 
 FiniteElementVector (const FiniteElementVector &rhs)
 Copy constructor. More...
 
 FiniteElementVector (FiniteElementVector &&rhs)
 Move construct a new Finite Element Vector object. More...
 
- Public Member Functions inherited from serac::FiniteElementVector
 FiniteElementVector (const mfem::ParFiniteElementSpace &space, const std::string &name="")
 Minimal constructor for a FiniteElementVector given a finite element space. More...
 
template<typename FunctionSpace >
 FiniteElementVector (mfem::ParMesh &mesh, FunctionSpace, const std::string &name="")
 Construct a new Finite Element Vector object given a templated function space. More...
 
 FiniteElementVector (const FiniteElementVector &rhs)
 Copy constructor. More...
 
 FiniteElementVector (FiniteElementVector &&rhs)
 Move construct a new Finite Element Vector object. More...
 
FiniteElementVectoroperator= (const FiniteElementVector &rhs)
 Copy assignment. More...
 
FiniteElementVectoroperator= (FiniteElementVector &&rhs)
 Move assignment. More...
 
FiniteElementVectoroperator= (const mfem::HypreParVector &rhs)
 Copy assignment from a hypre par vector. More...
 
FiniteElementVectoroperator= (const mfem::Vector &rhs)
 Copy assignment from a hypre par vector. More...
 
MPI_Comm comm () const
 Returns the MPI communicator for the state. More...
 
mfem::ParMesh & mesh ()
 Returns a non-owning reference to the internal mesh object. More...
 
const mfem::ParMesh & mesh () const
 This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.
 
mfem::ParFiniteElementSpace & space ()
 Returns a non-owning reference to the internal FESpace. More...
 
const mfem::ParFiniteElementSpace & space () const
 This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.
 
std::string name () const
 Returns the name of the FEState (field) More...
 
FiniteElementVectoroperator= (const double value)
 Set a finite element state to a constant value. More...
 
virtual ~FiniteElementVector ()
 Destroy the Finite Element Vector object.
 

Protected Attributes

std::unique_ptr< mfem::ParGridFunction > grid_func_
 An optional container for a grid function (L-vector) view of the finite element state. More...
 
- Protected Attributes inherited from serac::FiniteElementVector
std::reference_wrapper< mfem::ParMesh > mesh_
 A reference to the mesh object on which the field is defined.
 
std::unique_ptr< mfem::FiniteElementCollection > coll_
 Handle to the FiniteElementCollection, which is owned by MFEMSidreDataCollection. More...
 
std::unique_ptr< mfem::ParFiniteElementSpace > space_
 Handle to the mfem::ParFiniteElementSpace, which is owned by MFEMSidreDataCollection.
 
std::string name_ = ""
 The name of the finite element vector.
 

Detailed Description

Class for encapsulating the critical MFEM components of a primal finite element field.

Namely: Mesh, FiniteElementCollection, FiniteElementState, and the true vector of the solution

Definition at line 48 of file finite_element_state.hpp.

Constructor & Destructor Documentation

◆ FiniteElementState() [1/2]

serac::FiniteElementState::FiniteElementState ( const FiniteElementState rhs)
inline

Copy constructor.

Parameters
[in]rhsThe input state used for construction

Definition at line 58 of file finite_element_state.hpp.

◆ FiniteElementState() [2/2]

serac::FiniteElementState::FiniteElementState ( FiniteElementState &&  rhs)
inline

Move construct a new Finite Element State object.

Parameters
[in]rhsThe input vector used for construction

Definition at line 65 of file finite_element_state.hpp.

Member Function Documentation

◆ fillGridFunction()

void serac::FiniteElementState::fillGridFunction ( mfem::ParGridFunction &  grid_function) const
inline

Fill a user-provided grid function based on the underlying true vector.

This distributes true vector dofs to the finite element (local) dofs by multiplying the true dofs by the prolongation operator.

See also
MFEM documentation for details

Definition at line 136 of file finite_element_state.hpp.

◆ FiniteElementVector() [1/4]

serac::FiniteElementVector::FiniteElementVector
inline

Copy constructor.

Parameters
[in]rhsThe input vector used for construction

Definition at line 109 of file finite_element_vector.hpp.

◆ FiniteElementVector() [2/4]

serac::FiniteElementVector::FiniteElementVector

Minimal constructor for a FiniteElementVector given a finite element space.

Parameters
[in]spaceThe space to use for the finite element state. This space is deep copied into the new FE state
[in]nameThe name of the field

Definition at line 56 of file finite_element_vector.cpp.

◆ FiniteElementVector() [3/4]

serac::FiniteElementVector::FiniteElementVector

Move construct a new Finite Element Vector object.

Parameters
[in]rhsThe input vector used for construction

Definition at line 119 of file finite_element_vector.cpp.

◆ FiniteElementVector() [4/4]

template<typename FunctionSpace >
serac::FiniteElementVector::FiniteElementVector ( typename FunctionSpace  )
inline

Construct a new Finite Element Vector object given a templated function space.

Template Parameters
FunctionSpacewhat kind of interpolating functions to use
Parameters
meshThe mesh used to construct the finite element state
nameThe name of the new finite element state field

Definition at line 66 of file finite_element_vector.hpp.

◆ gridFunction()

mfem::ParGridFunction & serac::FiniteElementState::gridFunction ( ) const

Construct a grid function from the finite element state true vector.

Returns
The constructed grid function

Definition at line 76 of file finite_element_state.cpp.

◆ operator=() [1/5]

FiniteElementState& serac::FiniteElementState::operator= ( const FiniteElementState rhs)
inline

Copy assignment.

Parameters
rhsThe right hand side input state
Returns
The assigned FiniteElementState

Definition at line 73 of file finite_element_state.hpp.

◆ operator=() [2/5]

FiniteElementState& serac::FiniteElementState::operator= ( const mfem::HypreParVector &  rhs)
inline

Copy assignment with HypreParVector.

Parameters
rhsThe right hand side input HypreParVector
Returns
The assigned FiniteElementState

Definition at line 97 of file finite_element_state.hpp.

◆ operator=() [3/5]

FiniteElementState& serac::FiniteElementState::operator= ( const mfem::Vector &  rhs)
inline

Copy assignment with mfem::Vector.

Parameters
rhsThe right hand side input State
Returns
The assigned FiniteElementState

Definition at line 109 of file finite_element_state.hpp.

◆ operator=() [4/5]

FiniteElementState& serac::FiniteElementState::operator= ( double  rhs)
inline

Copy assignment with double.

Parameters
rhsThe right hand side input double
Returns
The assigned FiniteElementState

Definition at line 121 of file finite_element_state.hpp.

◆ operator=() [5/5]

FiniteElementState& serac::FiniteElementState::operator= ( FiniteElementState &&  rhs)
inline

Move assignment.

Parameters
rhsThe right hand side input State
Returns
The assigned FiniteElementState

Definition at line 85 of file finite_element_state.hpp.

◆ project() [1/3]

void serac::FiniteElementState::project ( const GeneralCoefficient coef)

Projects a coefficient (vector or scalar) onto the field

Parameters
[in]coefThe coefficient to project
Note
This only sets nodal values based on the coefficient at that point. It does not perform a full least squares projection.

Definition at line 32 of file finite_element_state.cpp.

◆ project() [2/3]

void serac::FiniteElementState::project ( mfem::Coefficient &  coef,
mfem::Array< int > &  dof_list,
std::optional< int >  component = {} 
)

Project a scalar coefficient onto a set of dofs.

Parameters
coefThe vector coefficient to project
dof_listA list of true degrees of freedom to set. Note this is the scalar dof (not vdof) numbering.
componentThe component to set
Note
This only sets nodal values based on the coefficient at that point. It does not perform a full least squares projection.

Definition at line 19 of file finite_element_state.cpp.

◆ project() [3/3]

void serac::FiniteElementState::project ( mfem::VectorCoefficient &  coef,
mfem::Array< int > &  dof_list 
)

Project a vector coefficient onto a set of dofs.

Parameters
coefThe vector coefficient to project
dof_listA list of true degrees of freedom to set. Note this is the scalar dof (not vdof) numbering.
Note
This only sets nodal values based on the coefficient at that point. It does not perform a full least squares projection.

Definition at line 12 of file finite_element_state.cpp.

◆ projectOnBoundary()

void serac::FiniteElementState::projectOnBoundary ( mfem::Coefficient &  coef,
const mfem::Array< int > &  markers 
)

Project a coefficient on a specific set of marked boundaries.

Parameters
coefThe coefficient to project
markersA marker array of the boundaries to set
Note
This only sets nodal values based on the coefficient at that point. It does not perform a full least squares projection.

Definition at line 60 of file finite_element_state.cpp.

◆ setFromGridFunction()

void serac::FiniteElementState::setFromGridFunction ( const mfem::ParGridFunction &  grid_function)
inline

Initialize the true vector in the FiniteElementState based on an input grid function.

This distributes the grid function dofs to the true vector dofs by multiplying by the restriction operator.

See also
MFEM documentation for details
Parameters
grid_functionThe grid function used to initialize the underlying true vector.

Definition at line 148 of file finite_element_state.hpp.

Member Data Documentation

◆ grid_func_

std::unique_ptr<mfem::ParGridFunction> serac::FiniteElementState::grid_func_
mutableprotected

An optional container for a grid function (L-vector) view of the finite element state.

If a user requests it, it is constructed and potentially reused during subsequent calls. It is not updated unless specifically requested via the gridFunction method.

Definition at line 216 of file finite_element_state.hpp.


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