v0.15.0
Loading...
Searching...
No Matches
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
 
 FEMethod ()=default
 
auto getFEName () const
 get finite element name
 
auto getDataDofsPtr () const
 
auto getDataVectorDofsPtr () const
 
const FieldEntity_vector_viewgetDataFieldEnts () const
 
boost::shared_ptr< FieldEntity_vector_view > & getDataFieldEntsPtr () const
 
auto & getRowFieldEnts () const
 
auto & getRowFieldEntsPtr () const
 
auto & getColFieldEnts () const
 
auto & getColFieldEntsPtr () const
 
auto getRowDofsPtr () const
 
auto getColDofsPtr () const
 
auto getNumberOfNodes () const
 
EntityHandle getFEEntityHandle () const
 
MoFEMErrorCode getNodeData (const std::string field_name, VectorDouble &data, const bool reset_dofs=true)
 
- Public Member Functions inherited from MoFEM::BasicMethod
 BasicMethod ()
 
virtual ~BasicMethod ()=default
 
int getNinTheLoop () const
 get number of evaluated element in the loop
 
int getLoopSize () const
 get loop size
 
auto getLoHiFERank () const
 Get lo and hi processor rank of iterated entities.
 
auto getLoFERank () const
 Get upper rank in loop for iterating elements.
 
auto getHiFERank () const
 Get upper rank in loop for iterating elements.
 
unsigned int getFieldBitNumber (std::string field_name) const
 
MoFEMErrorCode copyBasicMethod (const BasicMethod &basic)
 Copy data from other base method to this base method.
 
virtual MoFEMErrorCode preProcess ()
 function is run at the beginning of loop
 
virtual MoFEMErrorCode operator() ()
 function is run for every finite element
 
virtual MoFEMErrorCode postProcess ()
 function is run at the end of loop
 
boost::weak_ptr< CacheTuplegetCacheWeakPtr () const
 Get the cache weak ptr object.
 
- Public Member Functions inherited from MoFEM::KspMethod
MoFEMErrorCode query_interface (boost::typeindex::type_index type_index, UnknownInterface **iface) const
 
 KspMethod ()
 
virtual ~KspMethod ()=default
 
MoFEMErrorCode copyKsp (const KspMethod &ksp)
 copy data form another method
 
- Public Member Functions inherited from MoFEM::PetscData
 PetscData ()
 
virtual ~PetscData ()=default
 
MoFEMErrorCode copyPetscData (const PetscData &petsc_data)
 
- 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
 
 SnesMethod ()
 
virtual ~SnesMethod ()=default
 
MoFEMErrorCode copySnes (const SnesMethod &snes)
 Copy snes data.
 
- Public Member Functions inherited from MoFEM::TSMethod
MoFEMErrorCode query_interface (boost::typeindex::type_index type_index, UnknownInterface **iface) const
 
 TSMethod ()
 
virtual ~TSMethod ()=default
 
MoFEMErrorCode copyTs (const TSMethod &ts)
 Copy TS solver data.
 
- Public Member Functions inherited from MoFEM::TaoMethod
MoFEMErrorCode query_interface (boost::typeindex::type_index type_index, UnknownInterface **iface) const
 
 TaoMethod ()
 
virtual ~TaoMethod ()=default
 
MoFEMErrorCode copyTao (const TaoMethod &tao)
 Copy TAO data.
 

Public Attributes

std::string feName
 Name of finite element.
 
boost::shared_ptr< const NumeredEntFiniteElementnumeredEntFiniteElementPtr
 
boost::function< bool(FEMethod *fe_method_ptr)> exeTestHook
 Tet if element to skip element.
 
- Public Attributes inherited from MoFEM::BasicMethod
int nInTheLoop
 number currently of processed method
 
int loopSize
 local number oe methods to process
 
std::pair< int, int > loHiFERank
 Llo and hi processor rank of iterated entities.
 
int rAnk
 processor rank
 
int sIze
 number of processors in communicator
 
const RefEntity_multiIndexrefinedEntitiesPtr
 container of mofem dof entities
 
const RefElement_multiIndexrefinedFiniteElementsPtr
 container of mofem finite element entities
 
const ProblemproblemPtr
 raw pointer to problem
 
const Field_multiIndexfieldsPtr
 raw pointer to fields container
 
const FieldEntity_multiIndexentitiesPtr
 raw pointer to container of field entities
 
const DofEntity_multiIndexdofsPtr
 raw pointer container of dofs
 
const FiniteElement_multiIndexfiniteElementsPtr
 raw pointer to container finite elements
 
const EntFiniteElement_multiIndexfiniteElementsEntitiesPtr
 
const FieldEntityEntFiniteElementAdjacencyMap_multiIndexadjacenciesPtr
 
boost::function< MoFEMErrorCode()> preProcessHook
 Hook function for pre-processing.
 
boost::function< MoFEMErrorCode()> postProcessHook
 Hook function for post-processing.
 
boost::function< MoFEMErrorCode()> operatorHook
 Hook function for operator.
 
boost::movelib::unique_ptr< boolvecAssembleSwitch
 
boost::movelib::unique_ptr< boolmatAssembleSwitch
 
boost::weak_ptr< CacheTuplecacheWeakPtr
 
- Public Attributes inherited from MoFEM::KspMethod
KSPContext ksp_ctx
 Context.
 
KSP ksp
 KSP solver.
 
Vec & ksp_f
 
Mat & ksp_A
 
Mat & ksp_B
 
- Public Attributes inherited from MoFEM::PetscData
Switches data_ctx
 
Vec f
 
Mat A
 
Mat B
 
Vec x
 
Vec dx
 
Vec x_t
 
Vec x_tt
 
- Public Attributes inherited from MoFEM::SnesMethod
SNESContext snes_ctx
 
SNES snes
 snes solver
 
Vec & snes_x
 state vector
 
Vec & snes_dx
 solution update
 
Vec & snes_f
 residual
 
Mat & snes_A
 jacobian matrix
 
Mat & snes_B
 preconditioner of jacobian matrix
 
- Public Attributes inherited from MoFEM::TSMethod
TS ts
 time solver
 
TSContext ts_ctx
 
PetscInt ts_step
 time step number
 
PetscReal ts_a
 shift for U_t (see PETSc Time Solver)
 
PetscReal ts_aa
 shift for U_tt shift for U_tt
 
PetscReal ts_t
 time
 
PetscReal ts_dt
 time step size
 
Vec & ts_u
 state vector
 
Vec & ts_u_t
 time derivative of state vector
 
Vec & ts_u_tt
 second time derivative of state vector
 
Vec & ts_F
 residual vector
 
Mat & ts_A
 
Mat & ts_B
 Preconditioner for ts_A.
 
- Public Attributes inherited from MoFEM::TaoMethod
TAOContext tao_ctx
 
Tao tao
 tao solver
 
Vec & tao_x
 
Vec & tao_f
 state vector
 
Mat & tao_A
 gradient vector
 
Mat & tao_B
 hessian matrix
 

Additional Inherited Members

- Public Types inherited from MoFEM::KspMethod
enum  KSPContext { CTX_SETFUNCTION , CTX_OPERATORS , CTX_KSPNONE }
 pass information about context of KSP/DM for with finite element is computed 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
}
 
using Switches = std::bitset<8>
 
- Public Types inherited from MoFEM::SnesMethod
enum  SNESContext { CTX_SNESSETFUNCTION , CTX_SNESSETJACOBIAN , CTX_SNESNONE }
 
- Public Types inherited from MoFEM::TSMethod
enum  TSContext {
  CTX_TSSETRHSFUNCTION , CTX_TSSETRHSJACOBIAN , CTX_TSSETIFUNCTION , CTX_TSSETIJACOBIAN ,
  CTX_TSTSMONITORSET , CTX_TSNONE
}
 
- Public Types inherited from MoFEM::TaoMethod
enum  TAOContext { CTX_TAO_OBJECTIVE , CTX_TAO_GRADIENT , CTX_TAO_HESSIAN , CTX_TAO_NONE }
 
- 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)
 
static constexpr Switches CtxSetF = PetscData::Switches(CTX_SET_F)
 
static constexpr Switches CtxSetA = PetscData::Switches(CTX_SET_A)
 
static constexpr Switches CtxSetB = PetscData::Switches(CTX_SET_B)
 
static constexpr Switches CtxSetX = PetscData::Switches(CTX_SET_X)
 
static constexpr Switches CtxSetDX = PetscData::Switches(CTX_SET_DX)
 
static constexpr Switches CtxSetX_T = PetscData::Switches(CTX_SET_X_T)
 
static constexpr Switches CtxSetX_TT = PetscData::Switches(CTX_SET_X_TT)
 
static constexpr Switches CtxSetTime = PetscData::Switches(CTX_SET_TIME)
 

Detailed Description

structure for User Loop Methods on finite elements

It can be used to calculate stiffness matrices, residuals, load vectors etc. It is low level class however in some class users looking for speed and efficiency, can use it directly.

This class is used with Interface::loop_finite_elements, where user overloaded operator FEMethod::operator() is executed for each element in the problem. Class have to additional methods which are overloaded by user, FEMethod::preProcess() and FEMethod::postProcess() executed at beginning and end of the loop over problem elements, respectively.

Examples
UnsaturatedFlow.hpp, and simple_elasticity.cpp.

Definition at line 404 of file LoopMethods.hpp.

Constructor & Destructor Documentation

◆ FEMethod()

MoFEM::FEMethod::FEMethod ( )
default

Member Function Documentation

◆ getColDofsPtr()

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

Definition at line 478 of file LoopMethods.hpp.

478 {
479 return numeredEntFiniteElementPtr->getColDofsPtr();
480 };
boost::shared_ptr< const NumeredEntFiniteElement > numeredEntFiniteElementPtr

◆ getColFieldEnts()

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

Definition at line 466 of file LoopMethods.hpp.

466 {
467 return numeredEntFiniteElementPtr->getColFieldEnts();
468 };

◆ getColFieldEntsPtr()

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

Definition at line 470 of file LoopMethods.hpp.

470 {
471 return numeredEntFiniteElementPtr->getColFieldEntsPtr();
472 };

◆ getDataDofsPtr()

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

Definition at line 439 of file LoopMethods.hpp.

439 {
440 return numeredEntFiniteElementPtr->getDataDofsPtr();
441 };

◆ getDataFieldEnts()

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

Definition at line 447 of file LoopMethods.hpp.

447 {
448 return numeredEntFiniteElementPtr->getDataFieldEnts();
449 }

◆ getDataFieldEntsPtr()

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

Definition at line 452 of file LoopMethods.hpp.

452 {
453 return const_cast<NumeredEntFiniteElement *>(
455 ->getDataFieldEntsPtr();
456 }

◆ getDataVectorDofsPtr()

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

Definition at line 443 of file LoopMethods.hpp.

443 {
444 return numeredEntFiniteElementPtr->getDataVectorDofsPtr();
445 };

◆ getFEEntityHandle()

EntityHandle FEMethod::getFEEntityHandle ( ) const
inline

Definition at line 495 of file LoopMethods.hpp.

495 {
496 return numeredEntFiniteElementPtr->getEnt();
497}

◆ getFEName()

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

get finite element name

Returns
std::string

Definition at line 426 of file LoopMethods.hpp.

426{ return feName; }
std::string feName
Name of finite element.

◆ getNodeData()

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

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()
EntityHandle get_id_for_max_type()
EntityHandle get_id_for_min_type()
constexpr auto field_name
unsigned int getFieldBitNumber(std::string field_name) const
auto getDataDofsPtr() const
auto getNumberOfNodes() const
UId getLocalUniqueIdCalculate()
Get the Local Unique Id Calculate object.

◆ getNumberOfNodes()

auto FEMethod::getNumberOfNodes ( ) const
inline

Definition at line 491 of file LoopMethods.hpp.

491 {
492 return moab::CN::VerticesPerEntity(numeredEntFiniteElementPtr->getEntType());
493};

◆ getRowDofsPtr()

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

Definition at line 474 of file LoopMethods.hpp.

474 {
475 return numeredEntFiniteElementPtr->getRowDofsPtr();
476 };

◆ getRowFieldEnts()

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

Definition at line 458 of file LoopMethods.hpp.

458 {
459 return numeredEntFiniteElementPtr->getRowFieldEnts();
460 };

◆ getRowFieldEntsPtr()

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

Definition at line 462 of file LoopMethods.hpp.

462 {
463 return numeredEntFiniteElementPtr->getRowFieldEntsPtr();
464 };

◆ query_interface()

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

Reimplemented from MoFEM::BasicMethod.

Definition at line 406 of file LoopMethods.hpp.

407 {
409 *iface = const_cast<FEMethod *>(this);
411 }
#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

Tet if element to skip element.

If is set and return false elemnent us skiped in MoFEM::Core::loop_finite_elements

Note
That functionality is used to run elements on particular bit levels

Definition at line 437 of file LoopMethods.hpp.

◆ feName

std::string MoFEM::FEMethod::feName

Name of finite element.

Definition at line 415 of file LoopMethods.hpp.

◆ numeredEntFiniteElementPtr

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

Pointer to finite element database structure

Examples
prism_elements_from_surface.cpp.

Definition at line 418 of file LoopMethods.hpp.


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