v0.15.0
Loading...
Searching...
No Matches
Public Member Functions | Public Attributes | List of all members
MoFEM::FEMethod Struct Reference

Structure for user loop methods on finite elements. More...

#include "src/interfaces/LoopMethods.hpp"

Inheritance diagram for MoFEM::FEMethod:
[legend]
Collaboration diagram for MoFEM::FEMethod:
[legend]

Public Member Functions

MoFEMErrorCode query_interface (boost::typeindex::type_index type_index, UnknownInterface **iface) const
 Query interface for type casting.
 
 FEMethod ()=default
 Default constructor.
 
auto getFEName () const
 Get the name of the current finite element.
 
auto getDataDofsPtr () const
 Get pointer to DOF data for the current finite element.
 
auto getDataVectorDofsPtr () const
 Get pointer to vector DOF data for the current finite element.
 
const FieldEntity_vector_viewgetDataFieldEnts () const
 Get reference to data field entities for the current finite element.
 
boost::shared_ptr< FieldEntity_vector_view > & getDataFieldEntsPtr () const
 Get shared pointer to data field entities for the current finite element.
 
auto & getRowFieldEnts () const
 Get reference to row field entities for the current finite element.
 
auto & getRowFieldEntsPtr () const
 Get shared pointer to row field entities for the current finite element.
 
auto & getColFieldEnts () const
 Get reference to column field entities for the current finite element.
 
auto & getColFieldEntsPtr () const
 Get shared pointer to column field entities for the current finite element.
 
auto getRowDofsPtr () const
 Get pointer to row DOFs for the current finite element.
 
auto getColDofsPtr () const
 Get pointer to column DOFs for the current finite element.
 
auto getNumberOfNodes () const
 Get the number of nodes in the current finite element.
 
EntityHandle getFEEntityHandle () const
 Get the entity handle of the current finite element.
 
MoFEMErrorCode getNodeData (const std::string field_name, VectorDouble &data, const bool reset_dofs=true)
 Get nodal data for a specific field.
 
- Public Member Functions inherited from MoFEM::BasicMethod
 BasicMethod ()
 Default constructor.
 
virtual ~BasicMethod ()=default
 Virtual destructor.
 
int getNinTheLoop () const
 Get current loop iteration index.
 
int getLoopSize () const
 Get total loop size.
 
auto getLoHiFERank () const
 Get processor rank range for finite element iteration.
 
auto getLoFERank () const
 Get lower processor rank for finite element iteration.
 
auto getHiFERank () const
 Get upper processor rank for finite element iteration.
 
unsigned int getFieldBitNumber (std::string field_name) const
 Get bit number for a specific field by name.
 
MoFEMErrorCode copyBasicMethod (const BasicMethod &basic)
 Copy data from another BasicMethod instance.
 
virtual MoFEMErrorCode preProcess ()
 Pre-processing function executed at loop initialization.
 
virtual MoFEMErrorCode operator() ()
 Main operator function executed for each loop iteration.
 
virtual MoFEMErrorCode postProcess ()
 Post-processing function executed at loop completion.
 
boost::weak_ptr< CacheTuplegetCacheWeakPtr () const
 Get the cache weak pointer object.
 
- Public Member Functions inherited from MoFEM::KspMethod
MoFEMErrorCode query_interface (boost::typeindex::type_index type_index, UnknownInterface **iface) const
 Query interface for type casting.
 
 KspMethod ()
 Default constructor.
 
virtual ~KspMethod ()=default
 Virtual destructor.
 
MoFEMErrorCode copyKsp (const KspMethod &ksp)
 Copy data from another KSP method.
 
- Public Member Functions inherited from MoFEM::PetscData
 PetscData ()
 Default constructor.
 
virtual ~PetscData ()=default
 Virtual destructor.
 
MoFEMErrorCode copyPetscData (const PetscData &petsc_data)
 Copy PETSc data from another instance.
 
- Public Member Functions inherited from MoFEM::UnknownInterface
template<class IFACE >
MoFEMErrorCode registerInterface (bool error_if_registration_failed=true)
 Register interface.
 
template<class IFACE >
MoFEMErrorCode getInterface (IFACE *&iface) const
 Get interface reference to pointer of interface.
 
template<class IFACE >
MoFEMErrorCode getInterface (IFACE **const iface) const
 Get interface pointer to pointer of interface.
 
template<class IFACE , typename boost::enable_if< boost::is_pointer< IFACE >, int >::type = 0>
IFACE getInterface () const
 Get interface pointer to pointer of interface.
 
template<class IFACE , typename boost::enable_if< boost::is_reference< IFACE >, int >::type = 0>
IFACE getInterface () const
 Get reference to interface.
 
template<class IFACE >
IFACE * getInterface () const
 Function returning pointer to interface.
 
virtual ~UnknownInterface ()=default
 
- Public Member Functions inherited from MoFEM::SnesMethod
MoFEMErrorCode query_interface (boost::typeindex::type_index type_index, UnknownInterface **iface) const
 Query interface for type casting.
 
 SnesMethod ()
 Default constructor.
 
virtual ~SnesMethod ()=default
 Virtual destructor.
 
MoFEMErrorCode copySnes (const SnesMethod &snes)
 Copy SNES data from another instance.
 
- Public Member Functions inherited from MoFEM::TSMethod
MoFEMErrorCode query_interface (boost::typeindex::type_index type_index, UnknownInterface **iface) const
 Query interface for type casting.
 
 TSMethod ()
 Default constructor.
 
virtual ~TSMethod ()=default
 Virtual destructor.
 
MoFEMErrorCode copyTs (const TSMethod &ts)
 Copy TS solver data from another instance.
 
- Public Member Functions inherited from MoFEM::TaoMethod
MoFEMErrorCode query_interface (boost::typeindex::type_index type_index, UnknownInterface **iface) const
 Query interface for type casting.
 
 TaoMethod ()
 Default constructor.
 
virtual ~TaoMethod ()=default
 Virtual destructor.
 
MoFEMErrorCode copyTao (const TaoMethod &tao)
 Copy TAO data from another instance.
 

Public Attributes

std::string feName
 Name of the finite element being processed.
 
boost::shared_ptr< const NumeredEntFiniteElementnumeredEntFiniteElementPtr
 Shared pointer to finite element database structure.
 
boost::function< bool(FEMethod *fe_method_ptr)> exeTestHook
 Test function to determine if element should be skipped.
 
- Public Attributes inherited from MoFEM::BasicMethod
int nInTheLoop
 Current index of processed method in the loop.
 
int loopSize
 Total number of methods to process in the loop.
 
std::pair< int, int > loHiFERank
 Processor rank range for distributed finite element iteration.
 
int rAnk
 Current processor rank in MPI communicator.
 
int sIze
 Total number of processors in MPI communicator.
 
const RefEntity_multiIndexrefinedEntitiesPtr
 Pointer to container of refined MoFEM DOF entities.
 
const RefElement_multiIndexrefinedFiniteElementsPtr
 Pointer to container of refined finite element entities.
 
const ProblemproblemPtr
 Raw pointer to current MoFEM problem instance.
 
const Field_multiIndexfieldsPtr
 Raw pointer to fields multi-index container.
 
const FieldEntity_multiIndexentitiesPtr
 Raw pointer to container of field entities.
 
const DofEntity_multiIndexdofsPtr
 Raw pointer to container of degree of freedom entities.
 
const FiniteElement_multiIndexfiniteElementsPtr
 Raw pointer to container of finite elements.
 
const EntFiniteElement_multiIndexfiniteElementsEntitiesPtr
 Raw pointer to container of finite element entities.
 
const FieldEntityEntFiniteElementAdjacencyMap_multiIndexadjacenciesPtr
 Raw pointer to container of adjacencies between DOFs and finite elements.
 
boost::function< MoFEMErrorCode()> preProcessHook
 Hook function for pre-processing operations.
 
boost::function< MoFEMErrorCode()> postProcessHook
 Hook function for post-processing operations.
 
boost::function< MoFEMErrorCode()> operatorHook
 Hook function for main operator execution.
 
boost::movelib::unique_ptr< boolvecAssembleSwitch
 Switch for vector assembly operations.
 
boost::movelib::unique_ptr< boolmatAssembleSwitch
 Switch for matrix assembly operations.
 
boost::weak_ptr< CacheTuplecacheWeakPtr
 Weak pointer to cached entity data.
 
- Public Attributes inherited from MoFEM::KspMethod
KSPContext ksp_ctx
 Current KSP computation context.
 
KSP ksp
 PETSc KSP linear solver object.
 
Vec & ksp_f
 Reference to residual vector in KSP context.
 
Mat & ksp_A
 Reference to system matrix in KSP context.
 
Mat & ksp_B
 Reference to preconditioner matrix in KSP context.
 
- Public Attributes inherited from MoFEM::PetscData
Switches data_ctx
 Current data context switches.
 
Vec f
 PETSc residual vector.
 
Mat A
 PETSc Jacobian matrix.
 
Mat B
 PETSc preconditioner matrix.
 
Vec x
 PETSc solution vector.
 
Vec dx
 PETSc solution increment vector.
 
Vec x_t
 PETSc first time derivative vector.
 
Vec x_tt
 PETSc second time derivative vector.
 
- Public Attributes inherited from MoFEM::SnesMethod
SNESContext snes_ctx
 Current SNES computation context.
 
SNES snes
 PETSc SNES nonlinear solver object.
 
Vec & snes_x
 Reference to current solution state vector.
 
Vec & snes_dx
 Reference to solution update/increment vector.
 
Vec & snes_f
 Reference to residual vector.
 
Mat & snes_A
 Reference to Jacobian matrix.
 
Mat & snes_B
 Reference to preconditioner of Jacobian matrix.
 
- Public Attributes inherited from MoFEM::TSMethod
TS ts
 PETSc time stepping solver object.
 
TSContext ts_ctx
 Current TS computation context.
 
PetscInt ts_step
 Current time step number.
 
PetscReal ts_a
 Shift parameter for U_t (see PETSc Time Solver documentation)
 
PetscReal ts_aa
 Shift parameter for U_tt (second time derivative)
 
PetscReal ts_t
 Current time value.
 
PetscReal ts_dt
 Current time step size.
 
Vec & ts_u
 Reference to current state vector U(t)
 
Vec & ts_u_t
 Reference to first time derivative of state vector dU/dt.
 
Vec & ts_u_tt
 Reference to second time derivative of state vector d²U/dt²
 
Vec & ts_F
 Reference to residual vector F(t,U,U_t,U_tt)
 
Mat & ts_A
 Reference to Jacobian matrix: dF/dU + a*dF/dU_t + aa*dF/dU_tt.
 
Mat & ts_B
 Reference to preconditioner matrix for ts_A.
 
- Public Attributes inherited from MoFEM::TaoMethod
TAOContext tao_ctx
 Current TAO computation context.
 
Tao tao
 PETSc TAO optimization solver object.
 
Vec & tao_x
 Reference to optimization variables vector.
 
Vec & tao_f
 Reference to gradient vector.
 
Mat & tao_A
 Reference to Hessian matrix.
 
Mat & tao_B
 Reference to preconditioner matrix for Hessian.
 

Additional Inherited Members

- Public Types inherited from MoFEM::KspMethod
enum  KSPContext { CTX_SETFUNCTION , CTX_OPERATORS , CTX_KSPNONE }
 Context enumeration for KSP solver phases. More...
 
- Public Types inherited from MoFEM::PetscData
enum  DataContext {
  CTX_SET_NONE = 0 , CTX_SET_F = 1 << 0 , CTX_SET_A = 1 << 1 , CTX_SET_B = 1 << 2 ,
  CTX_SET_X = 1 << 3 , CTX_SET_DX = 1 << 4 , CTX_SET_X_T = 1 << 5 , CTX_SET_X_TT = 1 << 6 ,
  CTX_SET_TIME = 1 << 7
}
 Enumeration for data context flags. More...
 
using Switches = std::bitset< 8 >
 Bitset type for context switches.
 
- Public Types inherited from MoFEM::SnesMethod
enum  SNESContext { CTX_SNESSETFUNCTION , CTX_SNESSETJACOBIAN , CTX_SNESNONE }
 Context enumeration for SNES solver phases. More...
 
- Public Types inherited from MoFEM::TSMethod
enum  TSContext {
  CTX_TSSETRHSFUNCTION , CTX_TSSETRHSJACOBIAN , CTX_TSSETIFUNCTION , CTX_TSSETIJACOBIAN ,
  CTX_TSTSMONITORSET , CTX_TSNONE
}
 Context enumeration for TS solver phases. More...
 
- Public Types inherited from MoFEM::TaoMethod
enum  TAOContext { CTX_TAO_OBJECTIVE , CTX_TAO_GRADIENT , CTX_TAO_HESSIAN , CTX_TAO_NONE }
 Context enumeration for TAO solver phases. More...
 
- Static Public Member Functions inherited from MoFEM::UnknownInterface
static MoFEMErrorCode getLibVersion (Version &version)
 Get library version.
 
static MoFEMErrorCode getFileVersion (moab::Interface &moab, Version &version)
 Get database major version.
 
static MoFEMErrorCode setFileVersion (moab::Interface &moab, Version version=Version(MoFEM_VERSION_MAJOR, MoFEM_VERSION_MINOR, MoFEM_VERSION_BUILD))
 Get database major version.
 
static MoFEMErrorCode getInterfaceVersion (Version &version)
 Get database major version.
 
- Static Public Attributes inherited from MoFEM::PetscData
static constexpr Switches CtxSetNone = PetscData::Switches(CTX_SET_NONE)
 No data switch.
 
static constexpr Switches CtxSetF = PetscData::Switches(CTX_SET_F)
 Residual vector switch.
 
static constexpr Switches CtxSetA = PetscData::Switches(CTX_SET_A)
 Jacobian matrix switch.
 
static constexpr Switches CtxSetB = PetscData::Switches(CTX_SET_B)
 Preconditioner matrix switch.
 
static constexpr Switches CtxSetX = PetscData::Switches(CTX_SET_X)
 Solution vector switch.
 
static constexpr Switches CtxSetDX = PetscData::Switches(CTX_SET_DX)
 Solution increment switch.
 
static constexpr Switches CtxSetX_T = PetscData::Switches(CTX_SET_X_T)
 First time derivative switch.
 
static constexpr Switches CtxSetX_TT = PetscData::Switches(CTX_SET_X_TT)
 Second time derivative switch.
 
static constexpr Switches CtxSetTime = PetscData::Switches(CTX_SET_TIME)
 Time value switch.
 

Detailed Description

Structure for user loop methods on finite elements.

This class provides a high-performance interface for finite element computations. It can be used to calculate stiffness matrices, residual vectors, load vectors, and other element-level quantities. While it is a low-level class, users seeking maximum speed and efficiency can use it directly.

This class is designed to work with Interface::loop_finite_elements, where the user-overloaded operator FEMethod::operator() is executed for each element in the problem. The class provides additional virtual methods (preProcess and postProcess) that can be overloaded by users for custom preprocessing and postprocessing operations.

Key features:

Examples
mofem/tutorials/cor-0to1/src/UnsaturatedFlow.hpp, and mofem/tutorials/cor-6/simple_elasticity.cpp.

Definition at line 664 of file LoopMethods.hpp.

Constructor & Destructor Documentation

◆ FEMethod()

MoFEM::FEMethod::FEMethod ( )
default

Member Function Documentation

◆ getColDofsPtr()

auto MoFEM::FEMethod::getColDofsPtr ( ) const
inline

Get pointer to column DOFs for the current finite element.

Returns a pointer to the container of column degrees of freedom associated with the current finite element. Column DOFs are used in matrix assembly for the column indices.

Returns
Pointer to column DOFs container

Definition at line 804 of file LoopMethods.hpp.

804 {
805 return numeredEntFiniteElementPtr->getColDofsPtr();
806 };
boost::shared_ptr< const NumeredEntFiniteElement > numeredEntFiniteElementPtr
Shared pointer to finite element database structure.

◆ getColFieldEnts()

auto & MoFEM::FEMethod::getColFieldEnts ( ) const
inline

Get reference to column field entities for the current finite element.

Returns
Reference to column field entities container

Definition at line 770 of file LoopMethods.hpp.

770 {
771 return numeredEntFiniteElementPtr->getColFieldEnts();
772 };

◆ getColFieldEntsPtr()

auto & MoFEM::FEMethod::getColFieldEntsPtr ( ) const
inline

Get shared pointer to column field entities for the current finite element.

Returns
Shared pointer to column field entities container

Definition at line 778 of file LoopMethods.hpp.

778 {
779 return numeredEntFiniteElementPtr->getColFieldEntsPtr();
780 };

◆ getDataDofsPtr()

auto MoFEM::FEMethod::getDataDofsPtr ( ) const
inline

Get pointer to DOF data for the current finite element.

Returns
Pointer to DOF data container

Definition at line 719 of file LoopMethods.hpp.

719 {
720 return numeredEntFiniteElementPtr->getDataDofsPtr();
721 };

◆ getDataFieldEnts()

const FieldEntity_vector_view & MoFEM::FEMethod::getDataFieldEnts ( ) const
inline

Get reference to data field entities for the current finite element.

Returns
Constant reference to field entities view

Definition at line 735 of file LoopMethods.hpp.

735 {
736 return numeredEntFiniteElementPtr->getDataFieldEnts();
737 }

◆ getDataFieldEntsPtr()

boost::shared_ptr< FieldEntity_vector_view > & MoFEM::FEMethod::getDataFieldEntsPtr ( ) const
inline

Get shared pointer to data field entities for the current finite element.

Returns
Shared pointer to field entities view (mutable access)

Definition at line 744 of file LoopMethods.hpp.

744 {
745 return const_cast<NumeredEntFiniteElement *>(
748 }
boost::shared_ptr< FieldEntity_vector_view > & getDataFieldEntsPtr() const
Get shared pointer to data field entities for the current finite element.

◆ getDataVectorDofsPtr()

auto MoFEM::FEMethod::getDataVectorDofsPtr ( ) const
inline

Get pointer to vector DOF data for the current finite element.

Returns
Pointer to vector DOF data container

Definition at line 727 of file LoopMethods.hpp.

727 {
728 return numeredEntFiniteElementPtr->getDataVectorDofsPtr();
729 };

◆ getFEEntityHandle()

EntityHandle FEMethod::getFEEntityHandle ( ) const
inline

Get the entity handle of the current finite element.

Returns the MOAB entity handle that uniquely identifies the current finite element in the mesh database.

Returns
MOAB entity handle of the finite element

Definition at line 848 of file LoopMethods.hpp.

848 {
849 return numeredEntFiniteElementPtr->getEnt();
850}

◆ getFEName()

auto MoFEM::FEMethod::getFEName ( ) const
inline

Get the name of the current finite element.

Returns the name string identifier of the finite element currently being processed in the loop.

Returns
String containing the finite element name

Definition at line 697 of file LoopMethods.hpp.

697{ return feName; }
std::string feName
Name of the finite element being processed.

◆ getNodeData()

MoFEMErrorCode FEMethod::getNodeData ( const std::string  field_name,
VectorDouble data,
const bool  reset_dofs = true 
)

Get nodal data for a specific field.

Extracts nodal values for the specified field from the current finite element and stores them in the provided data vector.

Parameters
field_nameName of the field to extract data for
dataVector to store the extracted nodal data
reset_dofsFlag to reset DOF values (default: true)
Returns
Error code indicating success or failure

Definition at line 193 of file LoopMethods.cpp.

195 {
197
198 // TODO: [CORE-60] Fix implementation not to use DOFs
199
200 auto get_nodes_field_data =
201 [&](boost::shared_ptr<FEDofEntity_multiIndex> &&dofs,
202 VectorDouble &nodes_data) {
204
205 auto bit_number = getFieldBitNumber(field_name);
206 auto &dofs_by_uid = dofs->get<Unique_mi_tag>();
207 auto dit =
208 dofs_by_uid.lower_bound(FieldEntity::getLocalUniqueIdCalculate(
209 bit_number, get_id_for_min_type<MBVERTEX>()));
210 auto hi_dit =
211 dofs_by_uid.upper_bound(FieldEntity::getLocalUniqueIdCalculate(
212 bit_number, get_id_for_max_type<MBVERTEX>()));
213
214 if (dit != hi_dit) {
215 auto &first_dof = **dit;
216 const int num_nodes = getNumberOfNodes();
217 const int nb_dof_idx = first_dof.getNbOfCoeffs();
218 const int max_nb_dofs = nb_dof_idx * num_nodes;
219 nodes_data.resize(max_nb_dofs, false);
220 nodes_data.clear();
221
222 for (; dit != hi_dit;) {
223 const auto &dof_ptr = *dit;
224 const auto &dof = *dof_ptr;
225 const auto &sn = *dof.getSideNumberPtr();
226 const int side_number = sn.side_number;
227
228 int pos = side_number * nb_dof_idx;
229 auto ent_filed_data_vec = dof.getEntFieldData();
230 for (int ii = 0; ii != nb_dof_idx; ++ii) {
231 nodes_data[pos] = ent_filed_data_vec[ii];
232 ++pos;
233 ++dit;
234 }
235 }
236
237 } else if (reset_dofs) {
238 nodes_data.resize(0, false);
239 }
240
242 };
243
244 return get_nodes_field_data(getDataDofsPtr(), data);
245
247}
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
UBlasVector< double > VectorDouble
Definition Types.hpp:68
constexpr auto field_name
unsigned int getFieldBitNumber(std::string field_name) const
Get bit number for a specific field by name.
auto getDataDofsPtr() const
Get pointer to DOF data for the current finite element.
auto getNumberOfNodes() const
Get the number of nodes in the current finite element.
UId getLocalUniqueIdCalculate()
Get the Local Unique Id Calculate object.

◆ getNumberOfNodes()

auto FEMethod::getNumberOfNodes ( ) const
inline

Get the number of nodes in the current finite element.

Returns the number of vertex nodes for the current finite element based on its entity type (e.g., 4 for tetrahedron, 8 for hexahedron).

Returns
Number of nodes in the finite element

Definition at line 844 of file LoopMethods.hpp.

844 {
845 return moab::CN::VerticesPerEntity(numeredEntFiniteElementPtr->getEntType());
846};

◆ getRowDofsPtr()

auto MoFEM::FEMethod::getRowDofsPtr ( ) const
inline

Get pointer to row DOFs for the current finite element.

Returns a pointer to the container of row degrees of freedom associated with the current finite element. Row DOFs are used in matrix assembly for the row indices.

Returns
Pointer to row DOFs container

Definition at line 791 of file LoopMethods.hpp.

791 {
792 return numeredEntFiniteElementPtr->getRowDofsPtr();
793 };

◆ getRowFieldEnts()

auto & MoFEM::FEMethod::getRowFieldEnts ( ) const
inline

Get reference to row field entities for the current finite element.

Returns
Reference to row field entities container

Definition at line 754 of file LoopMethods.hpp.

754 {
755 return numeredEntFiniteElementPtr->getRowFieldEnts();
756 };

◆ getRowFieldEntsPtr()

auto & MoFEM::FEMethod::getRowFieldEntsPtr ( ) const
inline

Get shared pointer to row field entities for the current finite element.

Returns
Shared pointer to row field entities container

Definition at line 762 of file LoopMethods.hpp.

762 {
763 return numeredEntFiniteElementPtr->getRowFieldEntsPtr();
764 };

◆ query_interface()

MoFEMErrorCode MoFEM::FEMethod::query_interface ( boost::typeindex::type_index  type_index,
UnknownInterface **  iface 
) const
inlinevirtual

Query interface for type casting.

Parameters
type_indexType information for interface querying
ifacePointer to interface object
Returns
Error code

Reimplemented from MoFEM::BasicMethod.

Definition at line 672 of file LoopMethods.hpp.

673 {
675 *iface = const_cast<FEMethod *>(this);
677 }
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...

Member Data Documentation

◆ exeTestHook

boost::function<bool(FEMethod *fe_method_ptr)> MoFEM::FEMethod::exeTestHook

Test function to determine if element should be skipped.

This optional hook function allows users to implement custom logic to determine whether a particular finite element should be skipped during the loop execution. If this function is set and returns false, the element will be skipped in MoFEM::Core::loop_finite_elements.

Note
This functionality is particularly useful for running elements only on specific bit levels or under certain conditions.
Parameters
fe_method_ptrPointer to the current FEMethod instance
Returns
true if element should be processed, false if it should be skipped

Definition at line 713 of file LoopMethods.hpp.

◆ feName

std::string MoFEM::FEMethod::feName

Name of the finite element being processed.

Definition at line 684 of file LoopMethods.hpp.

◆ numeredEntFiniteElementPtr

boost::shared_ptr<const NumeredEntFiniteElement> MoFEM::FEMethod::numeredEntFiniteElementPtr

Shared pointer to finite element database structure.

Examples
mofem/atom_tests/prism_elements_from_surface.cpp.

Definition at line 687 of file LoopMethods.hpp.


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