v0.9.1
Public Member Functions | Public Attributes | 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 ()=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

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

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

Constructor & Destructor Documentation

◆ BasicMethod()

MoFEM::BasicMethod::BasicMethod ( )

Definition at line 127 of file LoopMethods.cpp.

129  loopSize(0), rAnk(-1), sIze(-1), refinedEntitiesPtr(nullptr),
130  refinedFiniteElementsPtr(nullptr), problemPtr(nullptr),
131  fieldsPtr(nullptr), entitiesPtr(nullptr), dofsPtr(nullptr),
132  finiteElementsPtr(nullptr), finiteElementsEntitiesPtr(nullptr),
133  adjacenciesPtr(nullptr) {}
const DofEntity_multiIndex * dofsPtr
raw pointer container of dofs
const FieldEntityEntFiniteElementAdjacencyMap_multiIndex * adjacenciesPtr
const Field_multiIndex * fieldsPtr
raw pointer to fields container
const FiniteElement_multiIndex * finiteElementsPtr
raw pointer to container finite elements
int loopSize
local number oe methods to process
const FieldEntity_multiIndex * entitiesPtr
raw pointer to container of field entities
const RefElement_multiIndex * refinedFiniteElementsPtr
container of mofem finite element entities
int rAnk
processor rank
const Problem * problemPtr
raw pointer to problem
const RefEntity_multiIndex * refinedEntitiesPtr
container of mofem dof entities
int sIze
number of processors in communicator
int nInTheLoop
number currently of processed method
const EntFiniteElement_multiIndex * finiteElementsEntitiesPtr

◆ ~BasicMethod()

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

Member Function Documentation

◆ copyBasicMethod()

MoFEMErrorCode MoFEM::BasicMethod::copyBasicMethod ( const BasicMethod basic)

Copy data from other base method to this base method.

Parameters
basic
Returns
MoFEMErrorCode

Definition at line 135 of file LoopMethods.cpp.

135  {
137 
138  this->nInTheLoop = basic.nInTheLoop;
139  this->loopSize = basic.loopSize;
140  this->rAnk = basic.rAnk;
141  this->sIze = basic.sIze;
142  this->refinedEntitiesPtr = basic.refinedEntitiesPtr;
143  this->refinedFiniteElementsPtr = basic.refinedFiniteElementsPtr;
144  this->problemPtr = basic.problemPtr;
145  this->fieldsPtr = basic.fieldsPtr;
146  this->entitiesPtr = basic.entitiesPtr;
147  this->dofsPtr = basic.dofsPtr;
148  this->finiteElementsPtr = basic.finiteElementsPtr;
149  this->finiteElementsEntitiesPtr = basic.finiteElementsEntitiesPtr;
150  this->adjacenciesPtr = basic.adjacenciesPtr;
151 
153 }
const DofEntity_multiIndex * dofsPtr
raw pointer container of dofs
const FieldEntityEntFiniteElementAdjacencyMap_multiIndex * adjacenciesPtr
const Field_multiIndex * fieldsPtr
raw pointer to fields container
const FiniteElement_multiIndex * finiteElementsPtr
raw pointer to container finite elements
int loopSize
local number oe methods to process
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:507
const FieldEntity_multiIndex * entitiesPtr
raw pointer to container of field entities
const RefElement_multiIndex * refinedFiniteElementsPtr
container of mofem finite element entities
int rAnk
processor rank
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:514
const Problem * problemPtr
raw pointer to problem
const RefEntity_multiIndex * refinedEntitiesPtr
container of mofem dof entities
int sIze
number of processors in communicator
int nInTheLoop
number currently of processed method
const EntFiniteElement_multiIndex * finiteElementsEntitiesPtr

◆ getLoopSize()

int MoFEM::BasicMethod::getLoopSize ( ) const

get loop size

Definition at line 255 of file LoopMethods.hpp.

255 { return loopSize; }
int loopSize
local number oe methods to process

◆ getNinTheLoop()

int MoFEM::BasicMethod::getNinTheLoop ( ) const

get number of evaluated element in the loop

Definition at line 251 of file LoopMethods.hpp.

251 { return nInTheLoop; }
int nInTheLoop
number currently of processed method

◆ operator()()

MoFEMErrorCode MoFEM::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, ReactionDiffusion::Monitor, MixTransport::UnsaturatedFlowElement::MonitorPostProc, ElectroPhysiology::Monitor, ElecPhys::Monitor, FractureMechanics::CrackPropagation::PostProcVertexMethod, ElectroPhysiology::Monitor, StdRDOperators::Monitor, MoFEM::ForcesAndSourcesCore, UFOperators2D::Monitor, UFOperators::Monitor, MoFEM::VolumeElementForcesAndSourcesCoreSwitch< SWITCH >, MoFEM::FaceElementForcesAndSourcesCoreSwitch< SWITCH >, ReactionDiffusionEquation::Monitor, MonitorRestart, MoFEM::VolumeElementForcesAndSourcesCoreOnSideSwitch< SWITCH >, MoFEM::FatPrismElementForcesAndSourcesCore, MoFEM::ContactPrismElementForcesAndSourcesCore, MoFEM::ProjectionFieldOn10NodeTet, MoFEM::EdgeElementForcesAndSourcesCoreSwitch< SWITCH >, MoFEM::FaceElementForcesAndSourcesCoreOnSideSwitch< SWITCH >, MoFEM::FaceElementForcesAndSourcesCoreOnSideSwitch< FaceElementForcesAndSourcesCore::NO_HO_GEOMETRY >, MoFEM::FlatPrismElementForcesAndSourcesCore, NeumannForcesSurfaceComplexForLazy::MyTriangleSpatialFE, MonitorPostProc, MonitorPostProc, BasicFiniteElements::SaveVertexDofOnTag, MoFEM::VertexElementForcesAndSourcesCore, MoFEM::Projection10NodeCoordsOnField, CountDown, TestEntityMethod, and CountUp.

Definition at line 177 of file LoopMethods.cpp.

177  {
179  if (operatorHook) {
180  ierr = operatorHook();
181  CHKERRG(ierr);
182  } else {
183  SETERRQ(PETSC_COMM_SELF, MOFEM_OPERATION_UNSUCCESSFUL,
184  "should be implemented by user in derived class (operator)");
185  }
187 }
#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
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:514
boost::function< MoFEMErrorCode()> operatorHook
Hook function for operator.

◆ postProcess()

MoFEMErrorCode MoFEM::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, ReactionDiffusion::Monitor, ElectroPhysiology::Monitor, ElecPhys::Monitor, ElectroPhysiology::Monitor, FractureMechanics::CrackPropagation::PostProcVertexMethod, PostProcEdgeOnRefinedMesh, PostProcFatPrismOnRefinedMesh, StdRDOperators::Monitor, SolidShellModule::SolidShellPrismElement::PostProcFatPrismOnTriangleOnRefinedMesh, PostProcTemplateVolumeOnRefinedMesh< MoFEM::VolumeElementForcesAndSourcesCore >, MoFEM::ForcesAndSourcesCore, UFOperators2D::Monitor, UFOperators::Monitor, GelModule::Gel::GelFE, ApplyDirichletBc, SimpleContactProblem::PrintContactState, ReactionDiffusionEquation::Monitor, KelvinVoigtDamper::DamperFE, MonitorRestart, SolidShellModule::SolidShellPrismElement::SolidShellError, SolidShellModule::SolidShellPrismElement::SolidShell, AddLambdaVectorToFinternal, MonitorPostProc, DirichletFixFieldAtEntitiesBc, MoFEM::Projection10NodeCoordsOnField, CohesiveElement::AssembleRhsVectors, AssembleRhsVectors, Smoother::MyVolumeFE, MonitorPostProc, NonlinearElasticElement::MyVolumeFE, DirichletDisplacementBc, BasicFiniteElements::SaveVertexDofOnTag, CountDown, TestEntityMethod, and CountUp.

Definition at line 166 of file LoopMethods.cpp.

166  {
168  if (postProcessHook) {
169  ierr = postProcessHook();
170  CHKERRG(ierr);
171  } else {
172  SETERRQ(PETSC_COMM_SELF, MOFEM_OPERATION_UNSUCCESSFUL,
173  "should be implemented by user in derived class (postProcess)");
174  }
176 }
#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
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:514
boost::function< MoFEMErrorCode()> postProcessHook
Hook function for post-processing.

◆ preProcess()

MoFEMErrorCode MoFEM::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, ReactionDiffusion::Monitor, MixTransport::UnsaturatedFlowElement::MonitorPostProc, ElectroPhysiology::Monitor, ElecPhys::Monitor, EdgeSlidingConstrains::MyEdgeFE, ElectroPhysiology::Monitor, FractureMechanics::CrackPropagation::PostProcVertexMethod, PostProcEdgeOnRefinedMesh, PostProcFatPrismOnRefinedMesh, StdRDOperators::Monitor, PostProcTemplateVolumeOnRefinedMesh< MoFEM::VolumeElementForcesAndSourcesCore >, MoFEM::ForcesAndSourcesCore, UFOperators2D::Monitor, UFOperators::Monitor, GelModule::Gel::GelFE, SurfaceSlidingConstrains::MyTriangleFE, GroundSurfaceTemperature::SolarRadiationPreProcessor, GroundSurfaceTemperature::PreProcess, SimpleContactProblem::PrintContactState, ReactionDiffusionEquation::Monitor, KelvinVoigtDamper::DamperFE, MonitorRestart, SolidShellModule::SolidShellPrismElement::SolidShellError, SolidShellModule::SolidShellPrismElement::SolidShell, MoFEM::ProjectionFieldOn10NodeTet, DirichletFixFieldAtEntitiesBc, NeumannForcesSurfaceComplexForLazy::MyTriangleSpatialFE, CohesiveElement::AssembleRhsVectors, AssembleRhsVectors, MonitorPostProc, MonitorPostProc, Smoother::MyVolumeFE, NonlinearElasticElement::MyVolumeFE, DirichletDisplacementBc, MoFEM::Projection10NodeCoordsOnField, CountDown, BasicFiniteElements::SaveVertexDofOnTag, TestEntityMethod, and CountUp.

Definition at line 155 of file LoopMethods.cpp.

155  {
157  if (preProcessHook) {
158  ierr = preProcessHook();
159  CHKERRG(ierr);
160  } else {
161  SETERRQ(PETSC_COMM_SELF, MOFEM_OPERATION_UNSUCCESSFUL,
162  "should be implemented by user in derived class (preProcess)");
163  }
165 }
#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
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:514
boost::function< MoFEMErrorCode()> preProcessHook
Hook function for pre-processing.

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

228  {
229  if (uuid == IDD_MOFEMBasicMethod) {
230  *iface = const_cast<BasicMethod *>(this);
232  }
233  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "unknown interface");
234  }
static const MOFEMuuid IDD_MOFEMBasicMethod
Definition: LoopMethods.hpp:31
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:514

Member Data Documentation

◆ adjacenciesPtr

const FieldEntityEntFiniteElementAdjacencyMap_multiIndex* MoFEM::BasicMethod::adjacenciesPtr

raw pointer to container to adjacencies between dofs and finite elements

Definition at line 284 of file LoopMethods.hpp.

◆ dofsPtr

const DofEntity_multiIndex* MoFEM::BasicMethod::dofsPtr

raw pointer container of dofs

Definition at line 274 of file LoopMethods.hpp.

◆ entitiesPtr

const FieldEntity_multiIndex* MoFEM::BasicMethod::entitiesPtr

raw pointer to container of field entities

Definition at line 272 of file LoopMethods.hpp.

◆ fieldsPtr

const Field_multiIndex* MoFEM::BasicMethod::fieldsPtr

raw pointer to fields container

Definition at line 269 of file LoopMethods.hpp.

◆ finiteElementsEntitiesPtr

const EntFiniteElement_multiIndex* MoFEM::BasicMethod::finiteElementsEntitiesPtr

raw pointer to container finite elements entities

Definition at line 280 of file LoopMethods.hpp.

◆ finiteElementsPtr

const FiniteElement_multiIndex* MoFEM::BasicMethod::finiteElementsPtr

raw pointer to container finite elements

Definition at line 277 of file LoopMethods.hpp.

◆ loopSize

int MoFEM::BasicMethod::loopSize

local number oe methods to process

Definition at line 247 of file LoopMethods.hpp.

◆ matAssembleSwitch

boost::movelib::unique_ptr<bool> MoFEM::BasicMethod::matAssembleSwitch

Definition at line 338 of file LoopMethods.hpp.

◆ nInTheLoop

int MoFEM::BasicMethod::nInTheLoop

number currently of processed method

Definition at line 242 of file LoopMethods.hpp.

◆ operatorHook

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

Hook function for operator.

Definition at line 308 of file LoopMethods.hpp.

◆ postProcessHook

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

Hook function for post-processing.

Definition at line 303 of file LoopMethods.hpp.

◆ preProcessHook

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

Hook function for pre-processing.

Definition at line 298 of file LoopMethods.hpp.

◆ problemPtr

const Problem* MoFEM::BasicMethod::problemPtr

raw pointer to problem

Definition at line 267 of file LoopMethods.hpp.

◆ rAnk

int MoFEM::BasicMethod::rAnk

processor rank

Definition at line 257 of file LoopMethods.hpp.

◆ refinedEntitiesPtr

const RefEntity_multiIndex* MoFEM::BasicMethod::refinedEntitiesPtr

container of mofem dof entities

Definition at line 262 of file LoopMethods.hpp.

◆ refinedFiniteElementsPtr

const RefElement_multiIndex* MoFEM::BasicMethod::refinedFiniteElementsPtr

container of mofem finite element entities

Definition at line 265 of file LoopMethods.hpp.

◆ sIze

int MoFEM::BasicMethod::sIze

number of processors in communicator

Definition at line 259 of file LoopMethods.hpp.

◆ vecAssembleSwitch

boost::movelib::unique_ptr<bool> MoFEM::BasicMethod::vecAssembleSwitch

Definition at line 337 of file LoopMethods.hpp.


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