v0.8.16
Public Member Functions | Public Attributes | Private Member Functions | List of all members
MoFEM::BasicMethod Struct Reference

Data structure to exchange data between mofem and User Loop Methods.It allows to exchange data between MoFEM and user functions. It stores information about multi-indices. More...

#include <src/interfaces/LoopMethods.hpp>

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

Public Member Functions

MoFEMErrorCode query_interface (const MOFEMuuid &uuid, UnknownInterface **iface) const
 
 BasicMethod ()
 
virtual ~BasicMethod ()
 
int getNinTheLoop () const
 get number of evaluated element in the loop More...
 
int getLoopSize () const
 get loop size More...
 
MoFEMErrorCode copyBasicMethod (const BasicMethod &basic)
 
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 ()
 
MoFEMErrorCode setKspCtx (const KSPContext &ctx)
 set operator type More...
 
MoFEMErrorCode setKsp (KSP ksp)
 set solver More...
 
MoFEMErrorCode copyKsp (const KspMethod &ksp)
 copy data form another method More...
 
- 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 ()
 
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 ()
 
MoFEMErrorCode setSnesCtx (const SNESContext &ctx)
 Set SNES context. More...
 
MoFEMErrorCode setSnes (SNES snes)
 Set SNES instance. More...
 
MoFEMErrorCode copySnes (const SnesMethod &snes)
 Copy snes data. More...
 
- Public Member Functions inherited from MoFEM::TSMethod
 TSMethod ()
 
virtual ~TSMethod ()
 
MoFEMErrorCode setTsCtx (const TSContext &ctx)
 Set Ts context. More...
 
MoFEMErrorCode copyTs (const TSMethod &ts)
 Copy TS solver data. More...
 
MoFEMErrorCode setTs (TS _ts)
 Set TS solver. More...
 

Public Attributes

int nInTheLoop
 
int loopSize
 
int rAnk
 
int sIze
 
const RefEntity_multiIndexrefinedEntitiesPtr
 
const RefElement_multiIndexrefinedFiniteElementsPtr
 
const ProblemproblemPtr
 
const Field_multiIndexfieldsPtr
 
const FieldEntity_multiIndexentitiesPtr
 
const DofEntity_multiIndexdofsPtr
 
const FiniteElement_multiIndexfiniteElementsPtr
 
const EntFiniteElement_multiIndexfiniteElementsEntitiesPtr
 
const FieldEntityEntFiniteElementAdjacencyMap_multiIndexadjacenciesPtr
 
boost::function< MoFEMErrorCode()> preProcessHook
 
boost::function< MoFEMErrorCode()> postProcessHook
 
boost::function< MoFEMErrorCode()> operatorHook
 
- Public Attributes inherited from MoFEM::KspMethod
KSPContext ksp_ctx
 Context. More...
 
KSP ksp
 KSP solver. More...
 
Vec ksp_f
 the right hand side vector More...
 
Mat ksp_A
 matrix More...
 
Mat ksp_B
 preconditioner matrix More...
 
- Public Attributes inherited from MoFEM::SnesMethod
SNESContext snes_ctx
 
SNES snes
 
Vec snes_x
 
Vec snes_f
 
Mat snes_A
 
Mat snes_B
 
- Public Attributes inherited from MoFEM::TSMethod
TSContext ts_ctx
 
TS ts
 
Vec ts_u
 
Vec ts_u_t
 
Vec ts_F
 
Mat ts_A
 
Mat ts_B
 
PetscInt ts_step
 
PetscReal ts_a
 
PetscReal ts_t
 

Private Member Functions

void iNit ()
 

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::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
}
 
- 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

Data structure to exchange data between mofem and User Loop Methods.

It allows to exchange data between MoFEM and user functions. It stores information about multi-indices.

Definition at line 207 of file LoopMethods.hpp.

Constructor & Destructor Documentation

◆ BasicMethod()

BasicMethod::BasicMethod ( )

Definition at line 91 of file LoopMethods.cpp.

94  problemPtr(NULL), fieldsPtr(NULL), entitiesPtr(NULL), dofsPtr(NULL),
96  adjacenciesPtr(NULL), preProcessHook(NULL), postProcessHook(NULL),
97  operatorHook(NULL) {}
const DofEntity_multiIndex * dofsPtr
const FieldEntityEntFiniteElementAdjacencyMap_multiIndex * adjacenciesPtr
const Field_multiIndex * fieldsPtr
const FiniteElement_multiIndex * finiteElementsPtr
const FieldEntity_multiIndex * entitiesPtr
const RefElement_multiIndex * refinedFiniteElementsPtr
boost::function< MoFEMErrorCode()> postProcessHook
const Problem * problemPtr
const RefEntity_multiIndex * refinedEntitiesPtr
boost::function< MoFEMErrorCode()> operatorHook
const EntFiniteElement_multiIndex * finiteElementsEntitiesPtr
boost::function< MoFEMErrorCode()> preProcessHook

◆ ~BasicMethod()

virtual MoFEM::BasicMethod::~BasicMethod ( )
virtual

Definition at line 219 of file LoopMethods.hpp.

219 {};

Member Function Documentation

◆ copyBasicMethod()

MoFEMErrorCode BasicMethod::copyBasicMethod ( const BasicMethod basic)

Definition at line 99 of file LoopMethods.cpp.

99  {
101 
102  this->nInTheLoop = basic.nInTheLoop;
103  this->loopSize = basic.loopSize;
104  this->rAnk = basic.rAnk;
105  this->sIze = basic.sIze;
106  this->refinedEntitiesPtr = basic.refinedEntitiesPtr;
107  this->refinedFiniteElementsPtr = basic.refinedFiniteElementsPtr;
108  this->problemPtr = basic.problemPtr;
109  this->fieldsPtr = basic.fieldsPtr;
110  this->entitiesPtr = basic.entitiesPtr;
111  this->dofsPtr = basic.dofsPtr;
112  this->finiteElementsPtr = basic.finiteElementsPtr;
113  this->finiteElementsEntitiesPtr = basic.finiteElementsEntitiesPtr;
114  this->adjacenciesPtr = basic.adjacenciesPtr;
115 
117 }
const DofEntity_multiIndex * dofsPtr
const FieldEntityEntFiniteElementAdjacencyMap_multiIndex * adjacenciesPtr
const Field_multiIndex * fieldsPtr
const FiniteElement_multiIndex * finiteElementsPtr
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:483
const FieldEntity_multiIndex * entitiesPtr
const RefElement_multiIndex * refinedFiniteElementsPtr
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return() ...
Definition: definitions.h:490
const Problem * problemPtr
const RefEntity_multiIndex * refinedEntitiesPtr
const EntFiniteElement_multiIndex * finiteElementsEntitiesPtr

◆ getLoopSize()

int MoFEM::BasicMethod::getLoopSize ( ) const

get loop size

Definition at line 230 of file LoopMethods.hpp.

230 { return loopSize; }

◆ getNinTheLoop()

int MoFEM::BasicMethod::getNinTheLoop ( ) const

get number of evaluated element in the loop

Definition at line 226 of file LoopMethods.hpp.

226 { return nInTheLoop; }

◆ iNit()

void MoFEM::BasicMethod::iNit ( )
private

◆ operator()()

MoFEMErrorCode BasicMethod::operator() ( )
virtual

function is run for every finite element

It is used to calculate element local matrices and assembly. It can be used for post-processing.

Reimplemented in GelModule::Gel::MonitorPostProc, MixTransport::UnsaturatedFlowElement::MonitorPostProc, FractureMechanics::CrackPropagation::PostProcVertexMethod, MoFEM::ForcesAndSourcesCore, MoFEM::FaceElementForcesAndSourcesCore, MoFEM::VolumeElementForcesAndSourcesCore, MoFEM::FatPrismElementForcesAndSourcesCore, MoFEM::FlatPrismElementForcesAndSourcesCore, MonitorRestart, MoFEM::ProjectionFieldOn10NodeTet, MoFEM::EdgeElementForcesAndSourcesCore, MonitorPostProc, NeummanForcesSurfaceComplexForLazy::MyTriangleSpatialFE, MonitorPostProc, MoFEM::VertexElementForcesAndSourcesCore, BasicFiniteElements::SaveVertexDofOnTag, and MoFEM::Projection10NodeCoordsOnField.

Definition at line 141 of file LoopMethods.cpp.

141  {
143  if (operatorHook) {
144  ierr = operatorHook();
145  CHKERRG(ierr);
146  } else {
147  SETERRQ(PETSC_COMM_SELF, MOFEM_OPERATION_UNSUCCESSFUL,
148  "should be implemented by user in derived class (operator)");
149  }
151 }
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:483
#define CHKERRG(n)
Check error code of MoFEM/MOAB/PETSc function.
Definition: definitions.h:526
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return() ...
Definition: definitions.h:490
static MoFEMErrorCodeGeneric< PetscErrorCode > ierr
Definition: Common.hpp:80
boost::function< MoFEMErrorCode()> operatorHook

◆ postProcess()

MoFEMErrorCode BasicMethod::postProcess ( )
virtual

function is run at the end of loop

It is used to assembly matrices and vectors, calculating global variables, f.e. total internal energy, ect.

Iterating over dofs: Example1 iterating over dofs in row by name of the field for(IT_GET_FEROW_BY_NAME_DOFS_FOR_LOOP(this,"DISPLACEMENT",it)) { ... }

Reimplemented in GelModule::Gel::MonitorPostProc, MixTransport::UnsaturatedFlowElement::MonitorPostProc, PostProcFaceOnRefinedMesh, PostProcFatPrismOnRefinedMesh, FractureMechanics::CrackPropagation::PostProcVertexMethod, MoFEM::ForcesAndSourcesCore, SolidShellModule::SolidShellPrismElement::PostProcFatPrismOnTriangleOnRefinedMesh, PostProcTemplateVolumeOnRefinedMesh< MoFEM::VolumeElementForcesAndSourcesCore >, GelModule::Gel::GelFE, MoFEM::FatPrismElementForcesAndSourcesCore, ApplyDirichletBc, KelvinVoigtDamper::DamperFE, MonitorRestart, SolidShellModule::SolidShellPrismElement::SolidShellError, SolidShellModule::SolidShellPrismElement::SolidShell, DirichletFixFieldAtEntitiesBc, AddLambdaVectorToFinternal, MonitorPostProc, MoFEM::Projection10NodeCoordsOnField, AssembleRhsVectors, CohesiveElement::AssembleRhsVectors, Smoother::MyVolumeFE, MonitorPostProc, NonlinearElasticElement::MyVolumeFE, DirichletDisplacementBc, and BasicFiniteElements::SaveVertexDofOnTag.

Definition at line 130 of file LoopMethods.cpp.

130  {
132  if (postProcessHook) {
133  ierr = postProcessHook();
134  CHKERRG(ierr);
135  } else {
136  SETERRQ(PETSC_COMM_SELF, MOFEM_OPERATION_UNSUCCESSFUL,
137  "should be implemented by user in derived class (postProcess)");
138  }
140 }
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:483
#define CHKERRG(n)
Check error code of MoFEM/MOAB/PETSc function.
Definition: definitions.h:526
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return() ...
Definition: definitions.h:490
boost::function< MoFEMErrorCode()> postProcessHook
static MoFEMErrorCodeGeneric< PetscErrorCode > ierr
Definition: Common.hpp:80

◆ preProcess()

MoFEMErrorCode BasicMethod::preProcess ( )
virtual

function is run at the beginning of loop

It is used to zeroing matrices and vectors, calculation of shape functions on reference element, preprocessing boundary conditions, etc.

Reimplemented in GelModule::Gel::MonitorPostProc, MixTransport::UnsaturatedFlowElement::MonitorPostProc, PostProcFaceOnRefinedMesh, PostProcFatPrismOnRefinedMesh, FractureMechanics::CrackPropagation::PostProcVertexMethod, MoFEM::ForcesAndSourcesCore, PostProcTemplateVolumeOnRefinedMesh< MoFEM::VolumeElementForcesAndSourcesCore >, GelModule::Gel::GelFE, MoFEM::FatPrismElementForcesAndSourcesCore, GroundSurfaceTemerature::SolarRadiationPreProcessor, GroundSurfaceTemerature::PreProcess, KelvinVoigtDamper::DamperFE, MonitorRestart, SolidShellModule::SolidShellPrismElement::SolidShellError, SolidShellModule::SolidShellPrismElement::SolidShell, MoFEM::ProjectionFieldOn10NodeTet, DirichletFixFieldAtEntitiesBc, NeummanForcesSurfaceComplexForLazy::MyTriangleSpatialFE, AssembleRhsVectors, CohesiveElement::AssembleRhsVectors, MonitorPostProc, MonitorPostProc, Smoother::MyVolumeFE, NonlinearElasticElement::MyVolumeFE, DirichletDisplacementBc, MoFEM::Projection10NodeCoordsOnField, FluidPressure::MyTriangleFE, and BasicFiniteElements::SaveVertexDofOnTag.

Definition at line 119 of file LoopMethods.cpp.

119  {
121  if (preProcessHook) {
122  ierr = preProcessHook();
123  CHKERRG(ierr);
124  } else {
125  SETERRQ(PETSC_COMM_SELF, MOFEM_OPERATION_UNSUCCESSFUL,
126  "should be implemented by user in derived class (preProcess)");
127  }
129 }
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:483
#define CHKERRG(n)
Check error code of MoFEM/MOAB/PETSc function.
Definition: definitions.h:526
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return() ...
Definition: definitions.h:490
static MoFEMErrorCodeGeneric< PetscErrorCode > ierr
Definition: Common.hpp:80
boost::function< MoFEMErrorCode()> preProcessHook

◆ query_interface()

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

Reimplemented from MoFEM::KspMethod.

Reimplemented in MoFEM::DofMethod, MoFEM::EntityMethod, and MoFEM::FEMethod.

Definition at line 209 of file LoopMethods.hpp.

210  {
211  if (uuid == IDD_MOFEMBasicMethod) {
212  *iface = const_cast<BasicMethod *>(this);
214  }
215  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "unknown interface");
216  }
static const MOFEMuuid IDD_MOFEMBasicMethod
Definition: LoopMethods.hpp:29
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return() ...
Definition: definitions.h:490

Member Data Documentation

◆ adjacenciesPtr

const FieldEntityEntFiniteElementAdjacencyMap_multiIndex* MoFEM::BasicMethod::adjacenciesPtr

Definition at line 241 of file LoopMethods.hpp.

◆ dofsPtr

const DofEntity_multiIndex* MoFEM::BasicMethod::dofsPtr

Definition at line 238 of file LoopMethods.hpp.

◆ entitiesPtr

const FieldEntity_multiIndex* MoFEM::BasicMethod::entitiesPtr

Definition at line 237 of file LoopMethods.hpp.

◆ fieldsPtr

const Field_multiIndex* MoFEM::BasicMethod::fieldsPtr

Definition at line 236 of file LoopMethods.hpp.

◆ finiteElementsEntitiesPtr

const EntFiniteElement_multiIndex* MoFEM::BasicMethod::finiteElementsEntitiesPtr

Definition at line 240 of file LoopMethods.hpp.

◆ finiteElementsPtr

const FiniteElement_multiIndex* MoFEM::BasicMethod::finiteElementsPtr

Definition at line 239 of file LoopMethods.hpp.

◆ loopSize

int MoFEM::BasicMethod::loopSize

Definition at line 222 of file LoopMethods.hpp.

◆ nInTheLoop

int MoFEM::BasicMethod::nInTheLoop

Definition at line 219 of file LoopMethods.hpp.

◆ operatorHook

boost::function<MoFEMErrorCode()> MoFEM::BasicMethod::operatorHook

Hook function for operator

Definition at line 258 of file LoopMethods.hpp.

◆ postProcessHook

boost::function<MoFEMErrorCode()> MoFEM::BasicMethod::postProcessHook

Hook function for post-processing

Definition at line 253 of file LoopMethods.hpp.

◆ preProcessHook

boost::function<MoFEMErrorCode()> MoFEM::BasicMethod::preProcessHook

Hook function for pre-processing

Definition at line 248 of file LoopMethods.hpp.

◆ problemPtr

const Problem* MoFEM::BasicMethod::problemPtr

Definition at line 235 of file LoopMethods.hpp.

◆ rAnk

int MoFEM::BasicMethod::rAnk

Definition at line 232 of file LoopMethods.hpp.

◆ refinedEntitiesPtr

const RefEntity_multiIndex* MoFEM::BasicMethod::refinedEntitiesPtr

Definition at line 233 of file LoopMethods.hpp.

◆ refinedFiniteElementsPtr

const RefElement_multiIndex* MoFEM::BasicMethod::refinedFiniteElementsPtr

Definition at line 234 of file LoopMethods.hpp.

◆ sIze

int MoFEM::BasicMethod::sIze

Definition at line 232 of file LoopMethods.hpp.


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