v0.14.0
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
 
 FEMethod ()=default
 
auto getFEName () const
 get finite element name More...
 
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 More...
 
int getLoopSize () const
 get loop size More...
 
auto getLoHiFERank () const
 Get lo and hi processor rank of iterated entities. More...
 
auto getLoFERank () const
 Get upper rank in loop for iterating elements. More...
 
auto getHiFERank () const
 Get upper rank in loop for iterating elements. More...
 
unsigned int getFieldBitNumber (std::string field_name) const
 
MoFEMErrorCode copyBasicMethod (const BasicMethod &basic)
 Copy data from other base method to this base method. More...
 
virtual MoFEMErrorCode preProcess ()
 function is run at the beginning of loop More...
 
virtual MoFEMErrorCode operator() ()
 function is run for every finite element More...
 
virtual MoFEMErrorCode postProcess ()
 function is run at the end of loop More...
 
boost::weak_ptr< CacheTuplegetCacheWeakPtr () const
 Get the cache weak ptr object. More...
 
- Public Member Functions inherited from MoFEM::KspMethod
 KspMethod ()
 
virtual ~KspMethod ()=default
 
MoFEMErrorCode copyKsp (const KspMethod &ksp)
 copy data form another method More...
 
- 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. More...
 
template<class IFACE >
MoFEMErrorCode getInterface (IFACE *&iface) const
 Get interface reference to pointer of interface. More...
 
template<class IFACE >
MoFEMErrorCode getInterface (IFACE **const iface) const
 Get interface pointer to pointer of interface. More...
 
template<class IFACE , typename boost::enable_if< boost::is_pointer< IFACE >, int >::type = 0>
IFACE getInterface () const
 Get interface pointer to pointer of interface. More...
 
template<class IFACE , typename boost::enable_if< boost::is_reference< IFACE >, int >::type = 0>
IFACE getInterface () const
 Get reference to interface. More...
 
template<class IFACE >
IFACE * getInterface () const
 Function returning pointer to interface. More...
 
virtual ~UnknownInterface ()=default
 
- Public Member Functions inherited from MoFEM::SnesMethod
 SnesMethod ()
 
virtual ~SnesMethod ()=default
 
MoFEMErrorCode copySnes (const SnesMethod &snes)
 Copy snes data. More...
 
- Public Member Functions inherited from MoFEM::TSMethod
 TSMethod ()
 
virtual ~TSMethod ()=default
 
MoFEMErrorCode copyTs (const TSMethod &ts)
 Copy TS solver data. More...
 

Public Attributes

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

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_X_T = 1 << 4, 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
}
 
- Static Public Member Functions inherited from MoFEM::UnknownInterface
static MoFEMErrorCode getLibVersion (Version &version)
 Get library version. More...
 
static MoFEMErrorCode getFileVersion (moab::Interface &moab, Version &version)
 Get database major version. More...
 
static MoFEMErrorCode setFileVersion (moab::Interface &moab, Version version=Version(MoFEM_VERSION_MAJOR, MoFEM_VERSION_MINOR, MoFEM_VERSION_BUILD))
 Get database major version. More...
 
static MoFEMErrorCode getInterfaceVersion (Version &version)
 Get database major version. More...
 
- 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 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
adolc_plasticity.cpp, child_and_parent.cpp, dm_build_partitioned_mesh.cpp, dynamic_first_order_con_law.cpp, free_surface.cpp, hanging_node_approx.cpp, heat_equation.cpp, level_set.cpp, nonlinear_dynamics.cpp, operators_tests.cpp, photon_diffusion.cpp, reaction_diffusion.cpp, seepage.cpp, shallow_wave.cpp, simple_elasticity.cpp, thermo_elastic.cpp, UnsaturatedFlow.hpp, and wave_equation.cpp.

Definition at line 369 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 443 of file LoopMethods.hpp.

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

◆ getColFieldEnts()

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

Definition at line 431 of file LoopMethods.hpp.

431  {
432  return numeredEntFiniteElementPtr->getColFieldEnts();
433  };

◆ getColFieldEntsPtr()

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

Definition at line 435 of file LoopMethods.hpp.

435  {
436  return numeredEntFiniteElementPtr->getColFieldEntsPtr();
437  };

◆ getDataDofsPtr()

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

Definition at line 404 of file LoopMethods.hpp.

404  {
405  return numeredEntFiniteElementPtr->getDataDofsPtr();
406  };

◆ getDataFieldEnts()

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

Definition at line 412 of file LoopMethods.hpp.

412  {
413  return numeredEntFiniteElementPtr->getDataFieldEnts();
414  }

◆ getDataFieldEntsPtr()

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

Definition at line 417 of file LoopMethods.hpp.

417  {
418  return const_cast<NumeredEntFiniteElement *>(
421  }

◆ getDataVectorDofsPtr()

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

Definition at line 408 of file LoopMethods.hpp.

408  {
409  return numeredEntFiniteElementPtr->getDataVectorDofsPtr();
410  };

◆ getFEEntityHandle()

EntityHandle FEMethod::getFEEntityHandle ( ) const
inline

Definition at line 460 of file LoopMethods.hpp.

460  {
461  return numeredEntFiniteElementPtr->getEnt();
462 }

◆ getFEName()

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

get finite element name

Returns
std::string

Definition at line 391 of file LoopMethods.hpp.

391 { return feName; }

◆ getNodeData()

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

Definition at line 172 of file LoopMethods.cpp.

174  {
176 
177  // TODO: [CORE-60] Fix implementation not to use DOFs
178 
179  auto get_nodes_field_data =
180  [&](boost::shared_ptr<FEDofEntity_multiIndex> &&dofs,
181  VectorDouble &nodes_data) {
183 
184  auto bit_number = getFieldBitNumber(field_name);
185  auto &dofs_by_uid = dofs->get<Unique_mi_tag>();
186  auto dit =
187  dofs_by_uid.lower_bound(FieldEntity::getLocalUniqueIdCalculate(
188  bit_number, get_id_for_min_type<MBVERTEX>()));
189  auto hi_dit =
190  dofs_by_uid.upper_bound(FieldEntity::getLocalUniqueIdCalculate(
191  bit_number, get_id_for_max_type<MBVERTEX>()));
192 
193  if (dit != hi_dit) {
194  auto &first_dof = **dit;
195  const int num_nodes = getNumberOfNodes();
196  const int nb_dof_idx = first_dof.getNbOfCoeffs();
197  const int max_nb_dofs = nb_dof_idx * num_nodes;
198  nodes_data.resize(max_nb_dofs, false);
199  nodes_data.clear();
200 
201  for (; dit != hi_dit;) {
202  const auto &dof_ptr = *dit;
203  const auto &dof = *dof_ptr;
204  const auto &sn = *dof.getSideNumberPtr();
205  const int side_number = sn.side_number;
206 
207  int pos = side_number * nb_dof_idx;
208  auto ent_filed_data_vec = dof.getEntFieldData();
209  for (int ii = 0; ii != nb_dof_idx; ++ii) {
210  nodes_data[pos] = ent_filed_data_vec[ii];
211  ++pos;
212  ++dit;
213  }
214  }
215 
216  } else if (reset_dofs) {
217  nodes_data.resize(0, false);
218  }
219 
221  };
222 
223  return get_nodes_field_data(getDataDofsPtr(), data);
224 
226 }

◆ getNumberOfNodes()

auto FEMethod::getNumberOfNodes ( ) const
inline

Definition at line 456 of file LoopMethods.hpp.

456  {
457  return moab::CN::VerticesPerEntity(numeredEntFiniteElementPtr->getEntType());
458 };

◆ getRowDofsPtr()

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

Definition at line 439 of file LoopMethods.hpp.

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

◆ getRowFieldEnts()

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

Definition at line 423 of file LoopMethods.hpp.

423  {
424  return numeredEntFiniteElementPtr->getRowFieldEnts();
425  };

◆ getRowFieldEntsPtr()

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

Definition at line 427 of file LoopMethods.hpp.

427  {
428  return numeredEntFiniteElementPtr->getRowFieldEntsPtr();
429  };

◆ query_interface()

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

Reimplemented from MoFEM::BasicMethod.

Definition at line 371 of file LoopMethods.hpp.

372  {
374  *iface = const_cast<FEMethod *>(this);
376  }

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 402 of file LoopMethods.hpp.

◆ feName

std::string MoFEM::FEMethod::feName

Name of finite element.

Definition at line 380 of file LoopMethods.hpp.

◆ numeredEntFiniteElementPtr

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

Pointer to finite element database structure

Examples
free_surface.cpp, and prism_elements_from_surface.cpp.

Definition at line 383 of file LoopMethods.hpp.


The documentation for this struct was generated from the following files:
MoFEMFunctionReturnHot
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:460
MoFEM::FEMethod::getDataDofsPtr
auto getDataDofsPtr() const
Definition: LoopMethods.hpp:404
MoFEM::FEMethod::feName
std::string feName
Name of finite element.
Definition: LoopMethods.hpp:380
FEMethod
MoFEM::FEMethod::getNumberOfNodes
auto getNumberOfNodes() const
Definition: LoopMethods.hpp:456
MoFEM::FieldEntity::getLocalUniqueIdCalculate
UId getLocalUniqueIdCalculate()
Get the Local Unique Id Calculate object.
Definition: FieldEntsMultiIndices.hpp:136
field_name
constexpr auto field_name
Definition: poisson_2d_homogeneous.cpp:13
MoFEM::BasicMethod::getFieldBitNumber
unsigned int getFieldBitNumber(std::string field_name) const
Definition: LoopMethods.hpp:270
MoFEM::Types::VectorDouble
UBlasVector< double > VectorDouble
Definition: Types.hpp:68
MoFEM::FEMethod::numeredEntFiniteElementPtr
boost::shared_ptr< const NumeredEntFiniteElement > numeredEntFiniteElementPtr
Definition: LoopMethods.hpp:383
MoFEMFunctionBeginHot
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:453
MoFEM::FEMethod::getDataFieldEntsPtr
boost::shared_ptr< FieldEntity_vector_view > & getDataFieldEntsPtr() const
Definition: LoopMethods.hpp:417
MoFEMFunctionReturn
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:429
MoFEMFunctionBegin
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:359