v0.15.0
Loading...
Searching...
No Matches
MoFEM::BasicMethod Struct Reference

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

#include "src/interfaces/LoopMethods.hpp"

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

Public Member Functions

MoFEMErrorCode query_interface (boost::typeindex::type_index type_index, UnknownInterface **iface) const
 
 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

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

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

Constructor & Destructor Documentation

◆ BasicMethod()

MoFEM::BasicMethod::BasicMethod ( )

Definition at line 130 of file LoopMethods.cpp.

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

◆ ~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 138 of file LoopMethods.cpp.

138 {
140
141 this->nInTheLoop = basic.nInTheLoop;
142 this->loopSize = basic.loopSize;
143 this->rAnk = basic.rAnk;
144 this->sIze = basic.sIze;
145 this->refinedEntitiesPtr = basic.refinedEntitiesPtr;
146 this->refinedFiniteElementsPtr = basic.refinedFiniteElementsPtr;
147 this->problemPtr = basic.problemPtr;
148 this->fieldsPtr = basic.fieldsPtr;
149 this->entitiesPtr = basic.entitiesPtr;
150 this->dofsPtr = basic.dofsPtr;
151 this->finiteElementsPtr = basic.finiteElementsPtr;
152 this->finiteElementsEntitiesPtr = basic.finiteElementsEntitiesPtr;
153 this->adjacenciesPtr = basic.adjacenciesPtr;
154 this->cacheWeakPtr = basic.cacheWeakPtr;
155
157}
#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 ...
boost::weak_ptr< CacheTuple > cacheWeakPtr

◆ getCacheWeakPtr()

boost::weak_ptr< CacheTuple > MoFEM::BasicMethod::getCacheWeakPtr ( ) const
inline

Get the cache weak ptr object.

Note
This store problem information on entities about DOFs. Each problem store different information. If you iterate over finite elements in preprocessor of TS solve element, us TS cache in the loop. Otherwise you will create undetermined behaviour or segmentation error. This is necessary compromise over bug resilience for memory saving and performance.
Returns
boost::weak_ptr<CacheTuple>
Examples
adolc_plasticity.cpp, dynamic_first_order_con_law.cpp, and free_surface.cpp.

Definition at line 379 of file LoopMethods.hpp.

379 {
380 return cacheWeakPtr;
381 }

◆ getFieldBitNumber()

unsigned int MoFEM::BasicMethod::getFieldBitNumber ( std::string field_name) const
inline

Definition at line 305 of file LoopMethods.hpp.

305 {
306 if (fieldsPtr) {
307 auto field_it = fieldsPtr->get<FieldName_mi_tag>().find(field_name);
308 if (field_it != fieldsPtr->get<FieldName_mi_tag>().end())
309 return (*field_it)->getBitNumber();
310 else
311 return BITFEID_SIZE;
312 } else {
313 THROW_MESSAGE("Pointer to fields multi-index is not set");
314 return BITFEID_SIZE;
315 }
316 }
#define BITFEID_SIZE
max number of finite elements
#define THROW_MESSAGE(msg)
Throw MoFEM exception.
Field_multiIndex::index< FieldName_mi_tag >::type::iterator field_it
constexpr auto field_name

◆ getHiFERank()

auto MoFEM::BasicMethod::getHiFERank ( ) const
inline

Get upper rank in loop for iterating elements.

Returns
loHiFERank.first

Definition at line 273 of file LoopMethods.hpp.

273{ return loHiFERank.second; }
std::pair< int, int > loHiFERank
Llo and hi processor rank of iterated entities.

◆ getLoFERank()

auto MoFEM::BasicMethod::getLoFERank ( ) const
inline

Get upper rank in loop for iterating elements.

Returns
loHiFERank.first

Definition at line 266 of file LoopMethods.hpp.

266{ return loHiFERank.first; }

◆ getLoHiFERank()

auto MoFEM::BasicMethod::getLoHiFERank ( ) const
inline

Get lo and hi processor rank of iterated entities.

Returns
raturn std::pair<int, int> loHiFERank

Definition at line 259 of file LoopMethods.hpp.

259{ return loHiFERank; }

◆ getLoopSize()

int MoFEM::BasicMethod::getLoopSize ( ) const
inline

get loop size

Definition at line 246 of file LoopMethods.hpp.

246{ return loopSize; }

◆ getNinTheLoop()

int MoFEM::BasicMethod::getNinTheLoop ( ) const
inline

get number of evaluated element in the loop

Definition at line 242 of file LoopMethods.hpp.

242{ return nInTheLoop; }

◆ 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 BasicFiniteElements::SaveVertexDofOnTag, BasicFiniteElements::SaveVertexDofOnTag, BcEntMethodDisp, BcEntMethodSpatial, CountDown, CountUp, DirichletDisplacementBc, DirichletDisplacementRemoveDofsBc, EshelbianPlasticity::EshelbianMonitor, MixTransport::UnsaturatedFlowElement::MonitorPostProc, MoFEM::ContactPrismElementForcesAndSourcesCore, MoFEM::EdgeElementForcesAndSourcesCore, MoFEM::FaceElementForcesAndSourcesCore, MoFEM::FatPrismElementForcesAndSourcesCore, MoFEM::FlatPrismElementForcesAndSourcesCore, MoFEM::ForcesAndSourcesCore, MoFEM::PipelineManager::MeshsetFE, MoFEM::Projection10NodeCoordsOnField, MoFEM::ProjectionFieldOn10NodeTet, MoFEM::VertexElementForcesAndSourcesCore, MoFEM::VolumeElementForcesAndSourcesCore, MoFEM::VolumeElementForcesAndSourcesCoreOnSide, Monitor, Monitor, MonitorIncompressible, MonitorPostProc, MonitorPostProc, MonitorRestart, NeumannForcesSurfaceComplexForLazy::MyTriangleSpatialFE, PhotonDiffusion::Monitor, ReactionDiffusionEquation::Monitor, and TestEntityMethod.

Definition at line 181 of file LoopMethods.cpp.

181 {
183 if (operatorHook) {
184 ierr = operatorHook();
185 CHKERRG(ierr);
186 } else {
187 SETERRQ(PETSC_COMM_SELF, MOFEM_OPERATION_UNSUCCESSFUL,
188 "should be implemented by user in derived class (operator)");
189 }
191}
#define CHKERRG(n)
Check error code of MoFEM/MOAB/PETSc function.
@ MOFEM_OPERATION_UNSUCCESSFUL
Definition definitions.h:34
static MoFEMErrorCodeGeneric< PetscErrorCode > ierr
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 ApplyDirichletBc, BasicFiniteElements::SaveVertexDofOnTag, BasicFiniteElements::SaveVertexDofOnTag, BcEntMethodDisp, BcEntMethodSpatial, CohesiveElement::AssembleRhsVectors, ConvectiveMassElement::MyVolumeFE, ConvectiveMassElement::UpdateAndControl, CountDown, CountUp, DirichletDisplacementBc, DirichletDisplacementRemoveDofsBc, DirichletFixFieldAtEntitiesBc, EshelbianPlasticity::EshelbianMonitor, KelvinVoigtDamper::DamperFE, MixTransport::UnsaturatedFlowElement::MonitorPostProc, MoFEM::ForcesAndSourcesCore, MoFEM::PipelineManager::MeshsetFE, MoFEM::PostProcBrokenMeshInMoabBase< EdgeElementForcesAndSourcesCore >, MoFEM::PostProcBrokenMeshInMoabBase< FaceElementForcesAndSourcesCore >, MoFEM::PostProcBrokenMeshInMoabBase< VolumeElementForcesAndSourcesCore >, MoFEM::Projection10NodeCoordsOnField, Monitor, Monitor, Monitor, Monitor, Monitor, Monitor, Monitor, Monitor, MonitorIncompressible, MonitorPostProc, MonitorPostProc, MonitorRestart, NonlinearElasticElement::MyVolumeFE, PhotonDiffusion::Monitor, PostProcEdgeOnRefinedMesh, PostProcFaceOnRefinedMesh, PostProcFatPrismOnRefinedMesh, PostProcTemplateVolumeOnRefinedMesh< MoFEM::VolumeElementForcesAndSourcesCore >, ReactionDiffusionEquation::Monitor, Smoother::MyVolumeFE, Smoother::MyVolumeFE, and TestEntityMethod.

Definition at line 170 of file LoopMethods.cpp.

170 {
172 if (postProcessHook) {
174 CHKERRG(ierr);
175 } else {
176 SETERRQ(PETSC_COMM_SELF, MOFEM_OPERATION_UNSUCCESSFUL,
177 "should be implemented by user in derived class (postProcess)");
178 }
180}
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 BasicFiniteElements::SaveVertexDofOnTag, BasicFiniteElements::SaveVertexDofOnTag, BcEntMethodDisp, BcEntMethodSpatial, CohesiveElement::AssembleRhsVectors, ConvectiveMassElement::MyVolumeFE, ConvectiveMassElement::ShellResidualElement, ConvectiveMassElement::UpdateAndControl, CountDown, CountUp, DirichletDisplacementBc, DirichletDisplacementRemoveDofsBc, DirichletFixFieldAtEntitiesBc, EdgeSlidingConstrains::MyEdgeFE, EdgeSlidingConstrains::MyEdgeFE, EshelbianPlasticity::EshelbianMonitor, FluidPressure::MyTriangleFE, KelvinVoigtDamper::DamperFE, MixTransport::UnsaturatedFlowElement::MonitorPostProc, MoFEM::ForcesAndSourcesCore, MoFEM::PipelineManager::MeshsetFE, MoFEM::PostProcBrokenMeshInMoabBase< EdgeElementForcesAndSourcesCore >, MoFEM::PostProcBrokenMeshInMoabBase< FaceElementForcesAndSourcesCore >, MoFEM::PostProcBrokenMeshInMoabBase< VolumeElementForcesAndSourcesCore >, MoFEM::Projection10NodeCoordsOnField, MoFEM::ProjectionFieldOn10NodeTet, Monitor, Monitor, MonitorIncompressible, MonitorPostProc, MonitorPostProc, MonitorRestart, NeumannForcesSurfaceComplexForLazy::MyTriangleSpatialFE, NonlinearElasticElement::MyVolumeFE, PhotonDiffusion::Monitor, PostProcEdgeOnRefinedMesh, PostProcFaceOnRefinedMesh, PostProcFatPrismOnRefinedMesh, PostProcTemplateVolumeOnRefinedMesh< MoFEM::VolumeElementForcesAndSourcesCore >, ReactionDiffusionEquation::Monitor, Smoother::MyVolumeFE, Smoother::MyVolumeFE, SurfaceSlidingConstrains::MyTriangleFE, SurfaceSlidingConstrains::MyTriangleFE, and TestEntityMethod.

Definition at line 159 of file LoopMethods.cpp.

159 {
161 if (preProcessHook) {
163 CHKERRG(ierr);
164 } else {
165 SETERRQ(PETSC_COMM_SELF, MOFEM_OPERATION_UNSUCCESSFUL,
166 "should be implemented by user in derived class (preProcess)");
167 }
169}
boost::function< MoFEMErrorCode()> preProcessHook
Hook function for pre-processing.

◆ query_interface()

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

Implements MoFEM::UnknownInterface.

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

Definition at line 220 of file LoopMethods.hpp.

221 {
223 *iface = const_cast<BasicMethod *>(this);
225 }

Member Data Documentation

◆ adjacenciesPtr

const FieldEntityEntFiniteElementAdjacencyMap_multiIndex* MoFEM::BasicMethod::adjacenciesPtr

raw pointer to container to adjacencies between dofs and finite elements

Definition at line 302 of file LoopMethods.hpp.

◆ cacheWeakPtr

boost::weak_ptr<CacheTuple> MoFEM::BasicMethod::cacheWeakPtr

Definition at line 386 of file LoopMethods.hpp.

◆ dofsPtr

const DofEntity_multiIndex* MoFEM::BasicMethod::dofsPtr

raw pointer container of dofs

Definition at line 292 of file LoopMethods.hpp.

◆ entitiesPtr

const FieldEntity_multiIndex* MoFEM::BasicMethod::entitiesPtr

raw pointer to container of field entities

Definition at line 290 of file LoopMethods.hpp.

◆ fieldsPtr

const Field_multiIndex* MoFEM::BasicMethod::fieldsPtr

raw pointer to fields container

Definition at line 287 of file LoopMethods.hpp.

◆ finiteElementsEntitiesPtr

const EntFiniteElement_multiIndex* MoFEM::BasicMethod::finiteElementsEntitiesPtr

raw pointer to container finite elements entities

Definition at line 298 of file LoopMethods.hpp.

◆ finiteElementsPtr

const FiniteElement_multiIndex* MoFEM::BasicMethod::finiteElementsPtr

raw pointer to container finite elements

Definition at line 295 of file LoopMethods.hpp.

◆ loHiFERank

std::pair<int, int> MoFEM::BasicMethod::loHiFERank

Llo and hi processor rank of iterated entities.

Definition at line 252 of file LoopMethods.hpp.

◆ loopSize

int MoFEM::BasicMethod::loopSize

local number oe methods to process

Definition at line 238 of file LoopMethods.hpp.

◆ matAssembleSwitch

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

Definition at line 384 of file LoopMethods.hpp.

◆ nInTheLoop

int MoFEM::BasicMethod::nInTheLoop

number currently of processed method

Definition at line 233 of file LoopMethods.hpp.

◆ operatorHook

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

Hook function for operator.

Definition at line 339 of file LoopMethods.hpp.

◆ postProcessHook

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

Hook function for post-processing.

Definition at line 334 of file LoopMethods.hpp.

◆ preProcessHook

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

Hook function for pre-processing.

Definition at line 329 of file LoopMethods.hpp.

◆ problemPtr

const Problem* MoFEM::BasicMethod::problemPtr

raw pointer to problem

Examples
simple_elasticity.cpp.

Definition at line 285 of file LoopMethods.hpp.

◆ rAnk

int MoFEM::BasicMethod::rAnk

processor rank

Definition at line 275 of file LoopMethods.hpp.

◆ refinedEntitiesPtr

const RefEntity_multiIndex* MoFEM::BasicMethod::refinedEntitiesPtr

container of mofem dof entities

Definition at line 280 of file LoopMethods.hpp.

◆ refinedFiniteElementsPtr

const RefElement_multiIndex* MoFEM::BasicMethod::refinedFiniteElementsPtr

container of mofem finite element entities

Definition at line 283 of file LoopMethods.hpp.

◆ sIze

int MoFEM::BasicMethod::sIze

number of processors in communicator

Definition at line 277 of file LoopMethods.hpp.

◆ vecAssembleSwitch

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

Definition at line 383 of file LoopMethods.hpp.


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