![]() |
v0.15.0 |
Data structure to exchange data between MoFEM and User Loop Methods. More...
#include "src/interfaces/LoopMethods.hpp"
Public Member Functions | |
| MoFEMErrorCode | query_interface (boost::typeindex::type_index type_index, UnknownInterface **iface) const |
| Query interface for type casting. | |
| BasicMethod () | |
| Default constructor. | |
| virtual | ~BasicMethod ()=default |
| Virtual destructor. | |
| int | getNinTheLoop () const |
| Get current loop iteration index. | |
| int | getLoopSize () const |
| Get total loop size. | |
| auto | getLoHiFERank () const |
| Get processor rank range for finite element iteration. | |
| auto | getLoFERank () const |
| Get lower processor rank for finite element iteration. | |
| auto | getHiFERank () const |
| Get upper processor rank for finite element iteration. | |
| unsigned int | getFieldBitNumber (std::string field_name) const |
| Get bit number for a specific field by name. | |
| MoFEMErrorCode | copyBasicMethod (const BasicMethod &basic) |
| Copy data from another BasicMethod instance. | |
| virtual MoFEMErrorCode | preProcess () |
| Pre-processing function executed at loop initialization. | |
| virtual MoFEMErrorCode | operator() () |
| Main operator function executed for each loop iteration. | |
| virtual MoFEMErrorCode | postProcess () |
| Post-processing function executed at loop completion. | |
| boost::weak_ptr< CacheTuple > | getCacheWeakPtr () const |
| Get the cache weak pointer object. | |
Public Member Functions inherited from MoFEM::KspMethod | |
| MoFEMErrorCode | query_interface (boost::typeindex::type_index type_index, UnknownInterface **iface) const |
| Query interface for type casting. | |
| KspMethod () | |
| Default constructor. | |
| virtual | ~KspMethod ()=default |
| Virtual destructor. | |
| MoFEMErrorCode | copyKsp (const KspMethod &ksp) |
| Copy data from another KSP method. | |
Public Member Functions inherited from MoFEM::PetscData | |
| PetscData () | |
| Default constructor. | |
| virtual | ~PetscData ()=default |
| Virtual destructor. | |
| MoFEMErrorCode | copyPetscData (const PetscData &petsc_data) |
| Copy PETSc data from another instance. | |
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 |
| Query interface for type casting. | |
| SnesMethod () | |
| Default constructor. | |
| virtual | ~SnesMethod ()=default |
| Virtual destructor. | |
| MoFEMErrorCode | copySnes (const SnesMethod &snes) |
| Copy SNES data from another instance. | |
Public Member Functions inherited from MoFEM::TSMethod | |
| MoFEMErrorCode | query_interface (boost::typeindex::type_index type_index, UnknownInterface **iface) const |
| Query interface for type casting. | |
| TSMethod () | |
| Default constructor. | |
| virtual | ~TSMethod ()=default |
| Virtual destructor. | |
| MoFEMErrorCode | copyTs (const TSMethod &ts) |
| Copy TS solver data from another instance. | |
Public Member Functions inherited from MoFEM::TaoMethod | |
| MoFEMErrorCode | query_interface (boost::typeindex::type_index type_index, UnknownInterface **iface) const |
| Query interface for type casting. | |
| TaoMethod () | |
| Default constructor. | |
| virtual | ~TaoMethod ()=default |
| Virtual destructor. | |
| MoFEMErrorCode | copyTao (const TaoMethod &tao) |
| Copy TAO data from another instance. | |
Public Attributes | |
| int | nInTheLoop |
| Current index of processed method in the loop. | |
| int | loopSize |
| Total number of methods to process in the loop. | |
| std::pair< int, int > | loHiFERank |
| Processor rank range for distributed finite element iteration. | |
| int | rAnk |
| Current processor rank in MPI communicator. | |
| int | sIze |
| Total number of processors in MPI communicator. | |
| const RefEntity_multiIndex * | refinedEntitiesPtr |
| Pointer to container of refined MoFEM DOF entities. | |
| const RefElement_multiIndex * | refinedFiniteElementsPtr |
| Pointer to container of refined finite element entities. | |
| const Problem * | problemPtr |
| Raw pointer to current MoFEM problem instance. | |
| const Field_multiIndex * | fieldsPtr |
| Raw pointer to fields multi-index container. | |
| const FieldEntity_multiIndex * | entitiesPtr |
| Raw pointer to container of field entities. | |
| const DofEntity_multiIndex * | dofsPtr |
| Raw pointer to container of degree of freedom entities. | |
| const FiniteElement_multiIndex * | finiteElementsPtr |
| Raw pointer to container of finite elements. | |
| const EntFiniteElement_multiIndex * | finiteElementsEntitiesPtr |
| Raw pointer to container of finite element entities. | |
| const FieldEntityEntFiniteElementAdjacencyMap_multiIndex * | adjacenciesPtr |
| Raw pointer to container of adjacencies between DOFs and finite elements. | |
| boost::function< MoFEMErrorCode()> | preProcessHook |
| Hook function for pre-processing operations. | |
| boost::function< MoFEMErrorCode()> | postProcessHook |
| Hook function for post-processing operations. | |
| boost::function< MoFEMErrorCode()> | operatorHook |
| Hook function for main operator execution. | |
| boost::movelib::unique_ptr< bool > | vecAssembleSwitch |
| Switch for vector assembly operations. | |
| boost::movelib::unique_ptr< bool > | matAssembleSwitch |
| Switch for matrix assembly operations. | |
| boost::weak_ptr< CacheTuple > | cacheWeakPtr |
| Weak pointer to cached entity data. | |
Public Attributes inherited from MoFEM::KspMethod | |
| KSPContext | ksp_ctx |
| Current KSP computation context. | |
| KSP | ksp |
| PETSc KSP linear solver object. | |
| Vec & | ksp_f |
| Reference to residual vector in KSP context. | |
| Mat & | ksp_A |
| Reference to system matrix in KSP context. | |
| Mat & | ksp_B |
| Reference to preconditioner matrix in KSP context. | |
Public Attributes inherited from MoFEM::PetscData | |
| Switches | data_ctx |
| Current data context switches. | |
| Vec | f |
| PETSc residual vector. | |
| Mat | A |
| PETSc Jacobian matrix. | |
| Mat | B |
| PETSc preconditioner matrix. | |
| Vec | x |
| PETSc solution vector. | |
| Vec | dx |
| PETSc solution increment vector. | |
| Vec | x_t |
| PETSc first time derivative vector. | |
| Vec | x_tt |
| PETSc second time derivative vector. | |
Public Attributes inherited from MoFEM::SnesMethod | |
| SNESContext | snes_ctx |
| Current SNES computation context. | |
| SNES | snes |
| PETSc SNES nonlinear solver object. | |
| Vec & | snes_x |
| Reference to current solution state vector. | |
| Vec & | snes_dx |
| Reference to solution update/increment vector. | |
| Vec & | snes_f |
| Reference to residual vector. | |
| Mat & | snes_A |
| Reference to Jacobian matrix. | |
| Mat & | snes_B |
| Reference to preconditioner of Jacobian matrix. | |
Public Attributes inherited from MoFEM::TSMethod | |
| TS | ts |
| PETSc time stepping solver object. | |
| TSContext | ts_ctx |
| Current TS computation context. | |
| PetscInt | ts_step |
| Current time step number. | |
| PetscReal | ts_a |
| Shift parameter for U_t (see PETSc Time Solver documentation) | |
| PetscReal | ts_aa |
| Shift parameter for U_tt (second time derivative) | |
| PetscReal | ts_t |
| Current time value. | |
| PetscReal | ts_dt |
| Current time step size. | |
| Vec & | ts_u |
| Reference to current state vector U(t) | |
| Vec & | ts_u_t |
| Reference to first time derivative of state vector dU/dt. | |
| Vec & | ts_u_tt |
| Reference to second time derivative of state vector d²U/dt² | |
| Vec & | ts_F |
| Reference to residual vector F(t,U,U_t,U_tt) | |
| Mat & | ts_A |
| Reference to Jacobian matrix: dF/dU + a*dF/dU_t + aa*dF/dU_tt. | |
| Mat & | ts_B |
| Reference to preconditioner matrix for ts_A. | |
Public Attributes inherited from MoFEM::TaoMethod | |
| TAOContext | tao_ctx |
| Current TAO computation context. | |
| Tao | tao |
| PETSc TAO optimization solver object. | |
| Vec & | tao_x |
| Reference to optimization variables vector. | |
| Vec & | tao_f |
| Reference to gradient vector. | |
| Mat & | tao_A |
| Reference to Hessian matrix. | |
| Mat & | tao_B |
| Reference to preconditioner matrix for Hessian. | |
Additional Inherited Members | |
Public Types inherited from MoFEM::KspMethod | |
| enum | KSPContext { CTX_SETFUNCTION , CTX_OPERATORS , CTX_KSPNONE } |
| Context enumeration for KSP solver phases. 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 } |
| Enumeration for data context flags. More... | |
| using | Switches = std::bitset< 8 > |
| Bitset type for context switches. | |
Public Types inherited from MoFEM::SnesMethod | |
| enum | SNESContext { CTX_SNESSETFUNCTION , CTX_SNESSETJACOBIAN , CTX_SNESNONE } |
| Context enumeration for SNES solver phases. More... | |
Public Types inherited from MoFEM::TSMethod | |
| enum | TSContext { CTX_TSSETRHSFUNCTION , CTX_TSSETRHSJACOBIAN , CTX_TSSETIFUNCTION , CTX_TSSETIJACOBIAN , CTX_TSTSMONITORSET , CTX_TSNONE } |
| Context enumeration for TS solver phases. More... | |
Public Types inherited from MoFEM::TaoMethod | |
| enum | TAOContext { CTX_TAO_OBJECTIVE , CTX_TAO_GRADIENT , CTX_TAO_HESSIAN , CTX_TAO_NONE } |
| Context enumeration for TAO solver phases. More... | |
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) |
| No data switch. | |
| static constexpr Switches | CtxSetF = PetscData::Switches(CTX_SET_F) |
| Residual vector switch. | |
| static constexpr Switches | CtxSetA = PetscData::Switches(CTX_SET_A) |
| Jacobian matrix switch. | |
| static constexpr Switches | CtxSetB = PetscData::Switches(CTX_SET_B) |
| Preconditioner matrix switch. | |
| static constexpr Switches | CtxSetX = PetscData::Switches(CTX_SET_X) |
| Solution vector switch. | |
| static constexpr Switches | CtxSetDX = PetscData::Switches(CTX_SET_DX) |
| Solution increment switch. | |
| static constexpr Switches | CtxSetX_T = PetscData::Switches(CTX_SET_X_T) |
| First time derivative switch. | |
| static constexpr Switches | CtxSetX_TT = PetscData::Switches(CTX_SET_X_TT) |
| Second time derivative switch. | |
| static constexpr Switches | CtxSetTime = PetscData::Switches(CTX_SET_TIME) |
| Time value switch. | |
Data structure to exchange data between MoFEM and User Loop Methods.
This is a comprehensive base class that combines all PETSc solver contexts (KSP, SNES, TS, TAO) and provides fundamental data exchange capabilities between MoFEM and user-defined functions. It stores information about multi-indices, problem databases, and provides access to finite element data structures for loop-based computations.
Definition at line 393 of file LoopMethods.hpp.
| MoFEM::BasicMethod::BasicMethod | ( | ) |
Default constructor.
Definition at line 130 of file LoopMethods.cpp.
|
virtualdefault |
Virtual destructor.
| MoFEMErrorCode MoFEM::BasicMethod::copyBasicMethod | ( | const BasicMethod & | basic | ) |
Copy data from another BasicMethod instance.
This function copies all relevant data from another BasicMethod instance including solver contexts, database pointers, and configuration settings.
| basic | Source BasicMethod to copy from |
Definition at line 138 of file LoopMethods.cpp.
|
inline |
Get the cache weak pointer object.
Returns a weak pointer to the cache tuple that stores problem information on entities about DOFs. Each problem stores different information.
Definition at line 633 of file LoopMethods.hpp.
|
inline |
Get bit number for a specific field by name.
This function looks up a field by name in the fields container and returns its bit number, which is used for field identification and bit operations in MoFEM's field management system.
| field_name | Name of the field to look up |
| MoFEM | exception if fieldsPtr is not set |
Definition at line 524 of file LoopMethods.hpp.
|
inline |
Get upper processor rank for finite element iteration.
Returns the highest processor rank in the range of processors that will process finite elements in the current loop.
Definition at line 483 of file LoopMethods.hpp.
|
inline |
Get lower processor rank for finite element iteration.
Returns the lowest processor rank in the range of processors that will process finite elements in the current loop.
Definition at line 473 of file LoopMethods.hpp.
|
inline |
Get processor rank range for finite element iteration.
Returns the pair containing the lowest and highest processor ranks involved in the current finite element loop iteration.
Definition at line 463 of file LoopMethods.hpp.
|
inline |
Get total loop size.
Definition at line 444 of file LoopMethods.hpp.
|
inline |
Get current loop iteration index.
Definition at line 438 of file LoopMethods.hpp.
|
virtual |
Main operator function executed for each loop iteration.
This virtual function is called for every item (finite element, entity, etc.) in the loop. It is the core computation function typically used for:
Reimplemented in CountUp, CountDown, TestEntityMethod, MoFEM::Projection10NodeCoordsOnField, MoFEM::ProjectionFieldOn10NodeTet, MoFEM::ContactPrismElementForcesAndSourcesCore, MoFEM::EdgeElementForcesAndSourcesCore, MoFEM::FaceElementForcesAndSourcesCore, MoFEM::FatPrismElementForcesAndSourcesCore, MoFEM::FlatPrismElementForcesAndSourcesCore, MoFEM::ForcesAndSourcesCore, MoFEM::VertexElementForcesAndSourcesCore, MoFEM::VolumeElementForcesAndSourcesCore, MoFEM::VolumeElementForcesAndSourcesCoreOnSide, MoFEM::PipelineManager::MeshsetFE, BasicFiniteElements::SaveVertexDofOnTag, MixTransport::UnsaturatedFlowElement::MonitorPostProc, ReactionDiffusionEquation::Monitor, MonitorIncompressible, PhotonDiffusion::Monitor, Monitor, Monitor, MonitorPostProc, MonitorRestart, DirichletDisplacementBc, BcEntMethodDisp, BcEntMethodSpatial, DirichletDisplacementRemoveDofsBc, BasicFiniteElements::SaveVertexDofOnTag, NeumannForcesSurfaceComplexForLazy::MyTriangleSpatialFE, and MonitorPostProc.
Definition at line 181 of file LoopMethods.cpp.
|
virtual |
Post-processing function executed at loop completion.
This virtual function is called once at the end of the loop. It is typically used for:
Example of iterating over DOFs:
Reimplemented in CountUp, CountDown, TestEntityMethod, MoFEM::Projection10NodeCoordsOnField, MoFEM::ForcesAndSourcesCore, MoFEM::PipelineManager::MeshsetFE, MoFEM::PostProcBrokenMeshInMoabBase< EdgeElementForcesAndSourcesCore >, MoFEM::PostProcBrokenMeshInMoabBase< FaceElementForcesAndSourcesCore >, MoFEM::PostProcBrokenMeshInMoabBase< VolumeElementForcesAndSourcesCore >, BasicFiniteElements::SaveVertexDofOnTag, Smoother::MyVolumeFE, Monitor, MixTransport::UnsaturatedFlowElement::MonitorPostProc, ApplyDirichletBc, ReactionDiffusionEquation::Monitor, MonitorIncompressible, PhotonDiffusion::Monitor, Monitor, Monitor, Monitor, Monitor, Monitor, Monitor, Monitor, CohesiveElement::AssembleRhsVectors, MonitorPostProc, MonitorRestart, ConvectiveMassElement::MyVolumeFE, ConvectiveMassElement::UpdateAndControl, DirichletDisplacementBc, BcEntMethodDisp, BcEntMethodSpatial, DirichletFixFieldAtEntitiesBc, DirichletDisplacementRemoveDofsBc, KelvinVoigtDamper::DamperFE, NonlinearElasticElement::MyVolumeFE, PostProcFatPrismOnRefinedMesh, PostProcFaceOnRefinedMesh, PostProcEdgeOnRefinedMesh, BasicFiniteElements::SaveVertexDofOnTag, Smoother::MyVolumeFE, MonitorPostProc, EshelbianPlasticity::ContactTree, and Monitor.
Definition at line 170 of file LoopMethods.cpp.
|
virtual |
Pre-processing function executed at loop initialization.
This virtual function is called once at the beginning of the loop. It is typically used for:
Reimplemented in CountUp, CountDown, TestEntityMethod, MoFEM::Projection10NodeCoordsOnField, MoFEM::ProjectionFieldOn10NodeTet, MoFEM::ForcesAndSourcesCore, MoFEM::PipelineManager::MeshsetFE, MoFEM::PostProcBrokenMeshInMoabBase< EdgeElementForcesAndSourcesCore >, MoFEM::PostProcBrokenMeshInMoabBase< FaceElementForcesAndSourcesCore >, MoFEM::PostProcBrokenMeshInMoabBase< VolumeElementForcesAndSourcesCore >, BasicFiniteElements::SaveVertexDofOnTag, Smoother::MyVolumeFE, SurfaceSlidingConstrains::MyTriangleFE, EdgeSlidingConstrains::MyEdgeFE, MixTransport::UnsaturatedFlowElement::MonitorPostProc, ReactionDiffusionEquation::Monitor, MonitorIncompressible, PhotonDiffusion::Monitor, Monitor, Monitor, CohesiveElement::AssembleRhsVectors, MonitorPostProc, MonitorRestart, ConvectiveMassElement::MyVolumeFE, ConvectiveMassElement::UpdateAndControl, ConvectiveMassElement::ShellResidualElement, DirichletDisplacementBc, BcEntMethodDisp, BcEntMethodSpatial, DirichletFixFieldAtEntitiesBc, DirichletDisplacementRemoveDofsBc, FluidPressure::MyTriangleFE, KelvinVoigtDamper::DamperFE, NonlinearElasticElement::MyVolumeFE, PostProcFatPrismOnRefinedMesh, PostProcFaceOnRefinedMesh, PostProcEdgeOnRefinedMesh, BasicFiniteElements::SaveVertexDofOnTag, Smoother::MyVolumeFE, NeumannForcesSurfaceComplexForLazy::MyTriangleSpatialFE, SurfaceSlidingConstrains::MyTriangleFE, EdgeSlidingConstrains::MyEdgeFE, MonitorPostProc, and EshelbianPlasticity::ContactTree.
Definition at line 159 of file LoopMethods.cpp.
|
inlinevirtual |
Query interface for type casting.
| type_index | Type information for interface querying |
| iface | Pointer to interface object |
Implements MoFEM::UnknownInterface.
Reimplemented in MoFEM::FEMethod, MoFEM::EntityMethod, and MoFEM::DofMethod.
Definition at line 401 of file LoopMethods.hpp.
| const FieldEntityEntFiniteElementAdjacencyMap_multiIndex* MoFEM::BasicMethod::adjacenciesPtr |
Raw pointer to container of adjacencies between DOFs and finite elements.
Definition at line 511 of file LoopMethods.hpp.
| boost::weak_ptr<CacheTuple> MoFEM::BasicMethod::cacheWeakPtr |
Weak pointer to cached entity data.
Definition at line 640 of file LoopMethods.hpp.
| const DofEntity_multiIndex* MoFEM::BasicMethod::dofsPtr |
Raw pointer to container of degree of freedom entities.
Definition at line 502 of file LoopMethods.hpp.
| const FieldEntity_multiIndex* MoFEM::BasicMethod::entitiesPtr |
Raw pointer to container of field entities.
Definition at line 500 of file LoopMethods.hpp.
| const Field_multiIndex* MoFEM::BasicMethod::fieldsPtr |
Raw pointer to fields multi-index container.
Definition at line 497 of file LoopMethods.hpp.
| const EntFiniteElement_multiIndex* MoFEM::BasicMethod::finiteElementsEntitiesPtr |
Raw pointer to container of finite element entities.
Definition at line 508 of file LoopMethods.hpp.
| const FiniteElement_multiIndex* MoFEM::BasicMethod::finiteElementsPtr |
Raw pointer to container of finite elements.
Definition at line 505 of file LoopMethods.hpp.
| std::pair<int, int> MoFEM::BasicMethod::loHiFERank |
Processor rank range for distributed finite element iteration.
This pair stores the low and high processor ranks that define the range of processors involved in the current finite element iteration. Used for parallel processing and load balancing across MPI processes.
Definition at line 453 of file LoopMethods.hpp.
| int MoFEM::BasicMethod::loopSize |
Total number of methods to process in the loop.
This value represents the total count of items (elements, entities, etc.) that will be processed in the current loop execution.
Definition at line 432 of file LoopMethods.hpp.
| boost::movelib::unique_ptr<bool> MoFEM::BasicMethod::matAssembleSwitch |
Switch for matrix assembly operations.
Definition at line 638 of file LoopMethods.hpp.
| int MoFEM::BasicMethod::nInTheLoop |
Current index of processed method in the loop.
This counter tracks which method/element is currently being processed in the loop iteration, useful for debugging and progress monitoring.
Definition at line 424 of file LoopMethods.hpp.
| boost::function<MoFEMErrorCode()> MoFEM::BasicMethod::operatorHook |
Hook function for main operator execution.
User-defined function that can be executed during the main operator call, providing additional customization points in the computation.
Definition at line 570 of file LoopMethods.hpp.
| boost::function<MoFEMErrorCode()> MoFEM::BasicMethod::postProcessHook |
Hook function for post-processing operations.
User-defined function that is executed after the main loop completes. Can be used for cleanup, finalization, or custom post-processing.
Definition at line 562 of file LoopMethods.hpp.
| boost::function<MoFEMErrorCode()> MoFEM::BasicMethod::preProcessHook |
Hook function for pre-processing operations.
User-defined function that is executed before the main loop begins. Can be used for initialization, setup operations, or custom preprocessing.
Definition at line 554 of file LoopMethods.hpp.
| const Problem* MoFEM::BasicMethod::problemPtr |
Raw pointer to current MoFEM problem instance.
Definition at line 495 of file LoopMethods.hpp.
| int MoFEM::BasicMethod::rAnk |
Current processor rank in MPI communicator.
Definition at line 485 of file LoopMethods.hpp.
| const RefEntity_multiIndex* MoFEM::BasicMethod::refinedEntitiesPtr |
Pointer to container of refined MoFEM DOF entities.
Definition at line 490 of file LoopMethods.hpp.
| const RefElement_multiIndex* MoFEM::BasicMethod::refinedFiniteElementsPtr |
Pointer to container of refined finite element entities.
Definition at line 493 of file LoopMethods.hpp.
| int MoFEM::BasicMethod::sIze |
Total number of processors in MPI communicator.
Definition at line 487 of file LoopMethods.hpp.
| boost::movelib::unique_ptr<bool> MoFEM::BasicMethod::vecAssembleSwitch |
Switch for vector assembly operations.
Definition at line 637 of file LoopMethods.hpp.