v0.9.1
Public Member Functions | Public Attributes | List of all members
MoFEM::FEMethod Struct Reference

structure for User Loop Methods on finite elementsIt 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. More...

#include <src/interfaces/LoopMethods.hpp>

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

Public Member Functions

MoFEMErrorCode query_interface (const MOFEMuuid &uuid, UnknownInterface **iface) const
 
 FEMethod ()
 
MoFEMErrorCode getNumberOfNodes (int &num_nodes) const
 Get number of DOFs on element. More...
 
EntityHandle getFEEntityHandle () const
 
MoFEMErrorCode getNodeData (const std::string field_name, VectorDouble &data, const bool reset_dofs=true)
 
template<class MULTIINDEX >
MULTIINDEX::iterator get_begin (const MULTIINDEX &index, const std::string &field_name, const EntityType type) const
 
template<class MULTIINDEX >
MULTIINDEX::iterator get_end (const MULTIINDEX &index, const std::string &field_name, const EntityType type) const
 
template<class MULTIINDEX >
MULTIINDEX::iterator get_begin (const MULTIINDEX &index, const std::string &field_name) const
 
template<class MULTIINDEX >
MULTIINDEX::iterator get_end (const MULTIINDEX &index, const std::string &field_name) const
 
- 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...
 
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...
 
- 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
 
- Public Member Functions inherited from MoFEM::UnknownInterface
template<class IFACE >
MoFEMErrorCode registerInterface (const MOFEMuuid &uuid, bool error_if_registration_failed=true)
 Register interface. More...
 
template<class IFACE , bool VERIFY = false>
MoFEMErrorCode getInterface (const MOFEMuuid &uuid, IFACE *&iface) const
 Get interface by uuid and return reference to pointer of interface. More...
 
template<class IFACE >
MoFEMErrorCode getInterface (IFACE *&iface) const
 Get interface refernce 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
 
virtual MoFEMErrorCode getLibVersion (Version &version) const
 Get library version. More...
 
virtual const MoFEMErrorCode getFileVersion (moab::Interface &moab, Version &version) const
 Get database major version. More...
 
virtual MoFEMErrorCode getInterfaceVersion (Version &version) const
 Get database major version. More...
 
template<>
MoFEMErrorCode getInterface (const MOFEMuuid &uuid, UnknownInterface *&iface) const
 
- Public Member Functions inherited from MoFEM::SnesMethod
 SnesMethod ()
 
virtual ~SnesMethod ()=default
 
MoFEMErrorCode copySnes (const SnesMethod &snes)
 Copy snes data. More...
 
DEPRECATED MoFEMErrorCode setSnesCtx (SNESContext ctx)
 
- Public Member Functions inherited from MoFEM::TSMethod
 TSMethod ()
 
virtual ~TSMethod ()=default
 
MoFEMErrorCode copyTs (const TSMethod &ts)
 Copy TS solver data. More...
 
DEPRECATED MoFEMErrorCode setTsCtx (TSContext ctx)
 

Public Attributes

std::string feName
 Name of finite element. More...
 
boost::shared_ptr< const NumeredEntFiniteElementnumeredEntFiniteElementPtr
 
boost::shared_ptr< const FENumeredDofEntity_multiIndexrowPtr
 Pointer to finite element rows dofs view. More...
 
boost::shared_ptr< const FENumeredDofEntity_multiIndexcolPtr
 Pointer to finite element columns dofs view. More...
 
boost::shared_ptr< const FEDofEntity_multiIndexdataPtr
 Pointer to finite element data dofs. More...
 
boost::shared_ptr< const FieldEntity_vector_viewrowFieldEntsPtr
 Pointer to finite element field entities row view. More...
 
boost::shared_ptr< const FieldEntity_vector_viewcolFieldEntsPtr
 Pointer to finite element field entities column view. More...
 
boost::shared_ptr< const FieldEntity_multiIndex_spaceType_viewdataFieldEntsPtr
 Pointer to finite element field entities data view. More...
 
- Public Attributes inherited from MoFEM::BasicMethod
int nInTheLoop
 number currently of processed method More...
 
int loopSize
 local number oe methods to process 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
 
- 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 More...
 
PetscReal ts_a
 shift for U_tt (see PETSc Time Solver) More...
 
PetscReal ts_v
 shift for U_t shift for U_t More...
 
PetscReal ts_t
 time 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 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)
 
- Protected Member Functions inherited from MoFEM::UnknownInterface
boost::typeindex::type_index getClassIdx (const MOFEMuuid &uid) const
 Get type name for interface Id. More...
 
MOFEMuuid getUId (const boost::typeindex::type_index &class_idx) const
 Get interface Id for class name. More...
 

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
dm_build_partitioned_mesh.cpp, EshelbianPlasticity.cpp, reaction_diffusion_equation.cpp, simple_elasticity.cpp, and UnsaturatedFlow.hpp.

Definition at line 356 of file LoopMethods.hpp.

Constructor & Destructor Documentation

◆ FEMethod()

FEMethod::FEMethod ( )

Definition at line 190 of file LoopMethods.cpp.

Member Function Documentation

◆ get_begin() [1/2]

template<class MULTIINDEX >
MULTIINDEX::iterator MoFEM::FEMethod::get_begin ( const MULTIINDEX &  index,
const std::string &  field_name,
const EntityType  type 
) const

Definition at line 401 of file LoopMethods.hpp.

403  {
404  return index.lower_bound(boost::make_tuple(field_name, type));
405  }

◆ get_begin() [2/2]

template<class MULTIINDEX >
MULTIINDEX::iterator MoFEM::FEMethod::get_begin ( const MULTIINDEX &  index,
const std::string &  field_name 
) const

Definition at line 444 of file LoopMethods.hpp.

445  {
446  return index.lower_bound(field_name);
447  }

◆ get_end() [1/2]

template<class MULTIINDEX >
MULTIINDEX::iterator MoFEM::FEMethod::get_end ( const MULTIINDEX &  index,
const std::string &  field_name,
const EntityType  type 
) const

Definition at line 407 of file LoopMethods.hpp.

409  {
410  return index.upper_bound(boost::make_tuple(field_name, type));
411  }

◆ get_end() [2/2]

template<class MULTIINDEX >
MULTIINDEX::iterator MoFEM::FEMethod::get_end ( const MULTIINDEX &  index,
const std::string &  field_name 
) const

Definition at line 449 of file LoopMethods.hpp.

450  {
451  return index.upper_bound(field_name);
452  }

◆ getFEEntityHandle()

EntityHandle FEMethod::getFEEntityHandle ( ) const

Definition at line 485 of file LoopMethods.hpp.

485  {
486  return numeredEntFiniteElementPtr->getEnt();
487 }
boost::shared_ptr< const NumeredEntFiniteElement > numeredEntFiniteElementPtr

◆ getNodeData()

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

Definition at line 232 of file LoopMethods.cpp.

234  {
236 
237  auto get_nodes_field_data = [&](FEDofEntity_multiIndex &dofs,
238  VectorDouble &nodes_data) {
240 
241  auto &dofs_by_name_and_type = dofs.get<Composite_Name_And_Type_mi_tag>();
242 
243  auto tuple = boost::make_tuple(field_name, MBVERTEX);
244  auto dit = dofs_by_name_and_type.lower_bound(tuple);
245  auto hi_dit = dofs.get<Composite_Name_And_Type_mi_tag>().upper_bound(tuple);
246 
247  if (dit != hi_dit) {
248  auto &first_dof = **dit;
249  int num_nodes;
250  CHKERR getNumberOfNodes(num_nodes);
251  const int nb_dof_idx = first_dof.getNbOfCoeffs();
252  const int max_nb_dofs = nb_dof_idx * num_nodes;
253  nodes_data.resize(max_nb_dofs, false);
254  nodes_data.clear();
255 
256  for (; dit != hi_dit;) {
257  const auto &dof_ptr = *dit;
258  const auto &dof = *dof_ptr;
259  const auto &sn = *dof.sideNumberPtr;
260  const int side_number = sn.side_number;
261 
262  int pos = side_number * nb_dof_idx;
263  auto ent_filed_data_vec = dof.getEntFieldData();
264  for (int ii = 0; ii != nb_dof_idx; ++ii) {
265  nodes_data[pos] = ent_filed_data_vec[ii];
266  ++pos;
267  ++dit;
268  }
269  }
270 
271  } else if(reset_dofs) {
272  nodes_data.resize(0, false);
273  }
274 
276  };
277 
278  return get_nodes_field_data(const_cast<FEDofEntity_multiIndex &>(
279  numeredEntFiniteElementPtr->getDataDofs()),
280  data);
281 
283 }
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:483
boost::shared_ptr< const NumeredEntFiniteElement > numeredEntFiniteElementPtr
#define CHKERR
Inline error check.
Definition: definitions.h:602
multi_index_container< boost::shared_ptr< FEDofEntity >, indexed_by< ordered_unique< tag< Unique_mi_tag >, const_mem_fun< FEDofEntity::interface_type_DofEntity, const UId &, &FEDofEntity::getGlobalUniqueId > >, ordered_non_unique< tag< Ent_mi_tag >, const_mem_fun< FEDofEntity::interface_type_DofEntity, EntityHandle, &FEDofEntity::getEnt > >, ordered_non_unique< tag< FieldName_mi_tag >, const_mem_fun< FEDofEntity::interface_type_Field, boost::string_ref, &FEDofEntity::getNameRef > >, ordered_non_unique< tag< EntType_mi_tag >, const_mem_fun< FEDofEntity::interface_type_RefEntity, EntityType, &FEDofEntity::getEntType > >, ordered_non_unique< tag< Composite_Name_And_Type_mi_tag >, composite_key< FEDofEntity, const_mem_fun< FEDofEntity::interface_type_Field, boost::string_ref, &FEDofEntity::getNameRef >, const_mem_fun< FEDofEntity::interface_type_RefEntity, EntityType, &FEDofEntity::getEntType > > >, ordered_non_unique< tag< Composite_Name_And_Ent_mi_tag >, composite_key< FEDofEntity, const_mem_fun< FEDofEntity::interface_type_Field, boost::string_ref, &FEDofEntity::getNameRef >, const_mem_fun< FEDofEntity::interface_type_DofEntity, EntityHandle, &FEDofEntity::getEnt > > >, ordered_non_unique< tag< Composite_Name_Type_And_Side_Number_mi_tag >, composite_key< FEDofEntity, const_mem_fun< FEDofEntity::interface_type_Field, boost::string_ref, &FEDofEntity::getNameRef >, const_mem_fun< FEDofEntity::interface_type_RefEntity, EntityType, &FEDofEntity::getEntType >, KeyFromKey< member< SideNumber, char, &SideNumber::side_number >, member< FEDofEntity::BaseFEEntity, boost::shared_ptr< SideNumber >, &FEDofEntity::sideNumberPtr > > > > > > FEDofEntity_multiIndex
MultiIndex container keeps FEDofEntity.
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:413
ublas::vector< double, DoubleAllocator > VectorDouble
Definition: Types.hpp:73
MoFEMErrorCode getNumberOfNodes(int &num_nodes) const
Get number of DOFs on element.

◆ getNumberOfNodes()

MoFEMErrorCode FEMethod::getNumberOfNodes ( int &  num_nodes) const

Get number of DOFs on element.

Definition at line 198 of file LoopMethods.cpp.

198  {
200 
201  EntityHandle handle = numeredEntFiniteElementPtr->getEnt();
202  if (handle) {
203  switch (static_cast<EntityType>(handle >> MB_ID_WIDTH)) {
204  case MBVERTEX:
205  num_nodes = 1;
206  break;
207  case MBEDGE:
208  num_nodes = 2;
209  break;
210  case MBTRI:
211  num_nodes = 3;
212  break;
213  case MBQUAD:
214  num_nodes = 4;
215  break;
216  case MBTET:
217  num_nodes = 4;
218  break;
219  case MBPRISM:
220  num_nodes = 6;
221  break;
222  default:
223  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "not implemented");
224  }
225  } else {
226  num_nodes = 0;
227  }
228 
230 }
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:507
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:514
boost::shared_ptr< const NumeredEntFiniteElement > numeredEntFiniteElementPtr
#define MB_ID_WIDTH
Definition: definitions.h:295

◆ query_interface()

MoFEMErrorCode MoFEM::FEMethod::query_interface ( const MOFEMuuid uuid,
UnknownInterface **  iface 
) const
virtual

Reimplemented from MoFEM::BasicMethod.

Definition at line 358 of file LoopMethods.hpp.

359  {
361  if (uuid == IDD_MOFEMFEMethod) {
362  *iface = const_cast<FEMethod *>(this);
364  }
365 
366  ierr = query_interface(uuid, iface);
367  CHKERRG(ierr);
369  }
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:507
#define CHKERRG(n)
Check error code of MoFEM/MOAB/PETSc function.
Definition: definitions.h:550
MoFEMErrorCode query_interface(const MOFEMuuid &uuid, UnknownInterface **iface) const
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:514
static const MOFEMuuid IDD_MOFEMFEMethod
Definition: LoopMethods.hpp:33

Member Data Documentation

◆ colFieldEntsPtr

boost::shared_ptr<const FieldEntity_vector_view> MoFEM::FEMethod::colFieldEntsPtr

Pointer to finite element field entities column view.

Definition at line 388 of file LoopMethods.hpp.

◆ colPtr

boost::shared_ptr<const FENumeredDofEntity_multiIndex> MoFEM::FEMethod::colPtr

Pointer to finite element columns dofs view.

Definition at line 381 of file LoopMethods.hpp.

◆ dataFieldEntsPtr

boost::shared_ptr<const FieldEntity_multiIndex_spaceType_view> MoFEM::FEMethod::dataFieldEntsPtr

Pointer to finite element field entities data view.

Definition at line 390 of file LoopMethods.hpp.

◆ dataPtr

boost::shared_ptr<const FEDofEntity_multiIndex> MoFEM::FEMethod::dataPtr

Pointer to finite element data dofs.

Definition at line 383 of file LoopMethods.hpp.

◆ feName

std::string MoFEM::FEMethod::feName

Name of finite element.

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

◆ rowFieldEntsPtr

boost::shared_ptr<const FieldEntity_vector_view> MoFEM::FEMethod::rowFieldEntsPtr

Pointer to finite element field entities row view.

Definition at line 386 of file LoopMethods.hpp.

◆ rowPtr

boost::shared_ptr<const FENumeredDofEntity_multiIndex> MoFEM::FEMethod::rowPtr

Pointer to finite element rows dofs view.

Definition at line 379 of file LoopMethods.hpp.


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