![]() |
v0.14.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 |
BasicMethod () | |
virtual | ~BasicMethod ()=default |
int | getNinTheLoop () const |
get number of evaluated element in the loop More... | |
int | getLoopSize () const |
get loop size More... | |
auto | getLoHiFERank () const |
Get lo and hi processor rank of iterated entities. More... | |
auto | getLoFERank () const |
Get upper rank in loop for iterating elements. More... | |
auto | getHiFERank () const |
Get upper rank in loop for iterating elements. More... | |
unsigned int | getFieldBitNumber (std::string field_name) const |
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... | |
boost::weak_ptr< CacheTuple > | getCacheWeakPtr () const |
Get the cache weak ptr object. More... | |
![]() | |
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 More... | |
![]() | |
MoFEMErrorCode | query_interface (boost::typeindex::type_index type_index, UnknownInterface **iface) const |
PetscData () | |
virtual | ~PetscData ()=default |
MoFEMErrorCode | copyPetscData (const PetscData &petsc_data) |
![]() | |
virtual MoFEMErrorCode | query_interface (boost::typeindex::type_index type_index, UnknownInterface **iface) const =0 |
template<class IFACE > | |
MoFEMErrorCode | registerInterface (bool error_if_registration_failed=true) |
Register 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 |
![]() | |
MoFEMErrorCode | query_interface (boost::typeindex::type_index type_index, UnknownInterface **iface) const |
SnesMethod () | |
virtual | ~SnesMethod ()=default |
MoFEMErrorCode | copySnes (const SnesMethod &snes) |
Copy snes data. More... | |
![]() | |
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. More... | |
Public Attributes | |
int | nInTheLoop |
number currently of processed method More... | |
int | loopSize |
local number oe methods to process More... | |
std::pair< int, int > | loHiFERank |
Llo and hi processor rank of iterated entities. More... | |
int | rAnk |
processor rank More... | |
int | sIze |
number of processors in communicator More... | |
const RefEntity_multiIndex * | refinedEntitiesPtr |
container of mofem dof entities More... | |
const RefElement_multiIndex * | refinedFiniteElementsPtr |
container of mofem finite element entities More... | |
const Problem * | problemPtr |
raw pointer to problem More... | |
const Field_multiIndex * | fieldsPtr |
raw pointer to fields container More... | |
const FieldEntity_multiIndex * | entitiesPtr |
raw pointer to container of field entities More... | |
const DofEntity_multiIndex * | dofsPtr |
raw pointer container of dofs More... | |
const FiniteElement_multiIndex * | finiteElementsPtr |
raw pointer to container finite elements More... | |
const EntFiniteElement_multiIndex * | finiteElementsEntitiesPtr |
const FieldEntityEntFiniteElementAdjacencyMap_multiIndex * | adjacenciesPtr |
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< bool > | vecAssembleSwitch |
boost::movelib::unique_ptr< bool > | matAssembleSwitch |
boost::weak_ptr< CacheTuple > | cacheWeakPtr |
![]() | |
KSPContext | ksp_ctx |
Context. More... | |
KSP | ksp |
KSP solver. More... | |
Vec & | ksp_f |
Mat & | ksp_A |
Mat & | ksp_B |
![]() | |
Switches | data_ctx |
Vec | f |
Mat | A |
Mat | B |
Vec | x |
Vec | x_t |
Vec | x_tt |
![]() | |
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... | |
![]() | |
TS | ts |
time solver More... | |
TSContext | ts_ctx |
PetscInt | ts_step |
time step number More... | |
PetscReal | ts_a |
shift for U_t (see PETSc Time Solver) More... | |
PetscReal | ts_aa |
shift for U_tt shift for U_tt More... | |
PetscReal | ts_t |
time More... | |
PetscReal | ts_dt |
time step size 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 | |
![]() | |
enum | KSPContext { CTX_SETFUNCTION , CTX_OPERATORS , CTX_KSPNONE } |
pass information about context of KSP/DM for with finite element is computed More... | |
![]() | |
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 > |
![]() | |
enum | SNESContext { CTX_SNESSETFUNCTION , CTX_SNESSETJACOBIAN , CTX_SNESNONE } |
![]() | |
enum | TSContext { CTX_TSSETRHSFUNCTION , CTX_TSSETRHSJACOBIAN , CTX_TSSETIFUNCTION , CTX_TSSETIJACOBIAN , CTX_TSTSMONITORSET , CTX_TSNONE } |
![]() | |
static MoFEMErrorCode | getLibVersion (Version &version) |
Get library version. More... | |
static MoFEMErrorCode | getFileVersion (moab::Interface &moab, Version &version) |
Get database major version. More... | |
static MoFEMErrorCode | setFileVersion (moab::Interface &moab, Version version=Version(MoFEM_VERSION_MAJOR, MoFEM_VERSION_MINOR, MoFEM_VERSION_BUILD)) |
Get database major version. More... | |
static MoFEMErrorCode | getInterfaceVersion (Version &version) |
Get database major version. More... | |
![]() | |
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) |
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 183 of file LoopMethods.hpp.
MoFEM::BasicMethod::BasicMethod | ( | ) |
Definition at line 109 of file LoopMethods.cpp.
|
virtualdefault |
MoFEMErrorCode MoFEM::BasicMethod::copyBasicMethod | ( | const BasicMethod & | basic | ) |
Copy data from other base method to this base method.
basic |
Definition at line 117 of file LoopMethods.cpp.
|
inline |
Get the cache weak ptr object.
Definition at line 344 of file LoopMethods.hpp.
|
inline |
Definition at line 270 of file LoopMethods.hpp.
|
inline |
Get upper rank in loop for iterating elements.
Definition at line 238 of file LoopMethods.hpp.
|
inline |
Get upper rank in loop for iterating elements.
Definition at line 231 of file LoopMethods.hpp.
|
inline |
Get lo and hi processor rank of iterated entities.
Definition at line 224 of file LoopMethods.hpp.
|
inline |
|
inline |
|
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 CountUp, CountDown, TestEntityMethod, MoFEM::Projection10NodeCoordsOnField, MoFEM::ProjectionFieldOn10NodeTet, MoFEM::ContactPrismElementForcesAndSourcesCore, MoFEM::EdgeElementForcesAndSourcesCore, MoFEM::FaceElementForcesAndSourcesCore, MoFEM::FatPrismElementForcesAndSourcesCore, MoFEM::FlatPrismElementForcesAndSourcesCore, MoFEM::ForcesAndSourcesCore, MoFEM::VertexElementForcesAndSourcesCore, MoFEM::VolumeElementForcesAndSourcesCore, MonitorPostProc, MonitorRestart, DirichletDisplacementBc, BcEntMethodDisp, BcEntMethodSpatial, DirichletDisplacementRemoveDofsBc, BasicFiniteElements::SaveVertexDofOnTag, NeumannForcesSurfaceComplexForLazy::MyTriangleSpatialFE, MonitorPostProc, FractureMechanics::CrackPropagation::PostProcVertexMethod, GelModule::Gel::MonitorPostProc, ContactPlasticity::MMonitor, ElectroPhysiology::Monitor, ElectroPhysiology::Monitor, ElecPhys::Monitor, StdRDOperators::Monitor, ReactionDiffusion::Monitor, UFOperators::Monitor, UFOperators2D::Monitor, MixTransport::UnsaturatedFlowElement::MonitorPostProc, ReactionDiffusionEquation::Monitor, MonitorIncompressible, PhotonDiffusion::Monitor, Monitor, and Monitor.
Definition at line 160 of file LoopMethods.cpp.
|
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 CountUp, CountDown, TestEntityMethod, MoFEM::Projection10NodeCoordsOnField, MoFEM::ForcesAndSourcesCore, MoFEM::PostProcBrokenMeshInMoabBase< EdgeElementForcesAndSourcesCore >, MoFEM::PostProcBrokenMeshInMoabBase< FaceElementForcesAndSourcesCore >, MoFEM::PostProcBrokenMeshInMoabBase< VolumeElementForcesAndSourcesCore >, CohesiveElement::AssembleRhsVectors, MonitorPostProc, MonitorRestart, DirichletDisplacementBc, BcEntMethodDisp, BcEntMethodSpatial, DirichletFixFieldAtEntitiesBc, DirichletDisplacementRemoveDofsBc, KelvinVoigtDamper::DamperFE, NonlinearElasticElement::MyVolumeFE, PostProcFatPrismOnRefinedMesh, PostProcFaceOnRefinedMesh, PostProcEdgeOnRefinedMesh, BasicFiniteElements::SaveVertexDofOnTag, Smoother::MyVolumeFE, MonitorPostProc, FractureMechanics::CrackPropagation::PostProcVertexMethod, GelModule::Gel::GelFE, GelModule::Gel::MonitorPostProc, ContactPlasticity::MMonitor, ElectroPhysiology::Monitor, ElectroPhysiology::Monitor, ElecPhys::Monitor, StdRDOperators::Monitor, ReactionDiffusion::Monitor, UFOperators::Monitor, UFOperators2D::Monitor, AssembleRhsVectors, AddLambdaVectorToFinternal, SolidShellModule::SolidShellPrismElement::SolidShell, SolidShellModule::SolidShellPrismElement::SolidShellError, SolidShellModule::SolidShellPrismElement::PostProcFatPrismOnTriangleOnRefinedMesh, DirichletDensityBc, Monitor, MixTransport::UnsaturatedFlowElement::MonitorPostProc, ApplyDirichletBc, ReactionDiffusionEquation::Monitor, MonitorIncompressible, PhotonDiffusion::Monitor, Monitor, Monitor, Monitor, Monitor, Monitor, and Monitor.
Definition at line 149 of file LoopMethods.cpp.
|
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 CountUp, CountDown, TestEntityMethod, MoFEM::Projection10NodeCoordsOnField, MoFEM::ProjectionFieldOn10NodeTet, MoFEM::ForcesAndSourcesCore, MoFEM::PostProcBrokenMeshInMoabBase< EdgeElementForcesAndSourcesCore >, MoFEM::PostProcBrokenMeshInMoabBase< FaceElementForcesAndSourcesCore >, MoFEM::PostProcBrokenMeshInMoabBase< VolumeElementForcesAndSourcesCore >, CohesiveElement::AssembleRhsVectors, MonitorPostProc, MonitorRestart, 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, FractureMechanics::CrackPropagation::PostProcVertexMethod, GelModule::Gel::GelFE, GelModule::Gel::MonitorPostProc, GroundSurfaceTemperature::PreProcess, GroundSurfaceTemperature::SolarRadiationPreProcessor, ContactPlasticity::MMonitor, ElectroPhysiology::Monitor, ElectroPhysiology::Monitor, ElecPhys::Monitor, StdRDOperators::Monitor, ReactionDiffusion::Monitor, UFOperators::Monitor, UFOperators2D::Monitor, AssembleRhsVectors, SolidShellModule::SolidShellPrismElement::SolidShell, SolidShellModule::SolidShellPrismElement::SolidShellError, MixTransport::UnsaturatedFlowElement::MonitorPostProc, ReactionDiffusionEquation::Monitor, MonitorIncompressible, PhotonDiffusion::Monitor, Monitor, and Monitor.
Definition at line 138 of file LoopMethods.cpp.
|
inlinevirtual |
Reimplemented from MoFEM::KspMethod.
Reimplemented in MoFEM::FEMethod, MoFEM::EntityMethod, and MoFEM::DofMethod.
Definition at line 185 of file LoopMethods.hpp.
const FieldEntityEntFiniteElementAdjacencyMap_multiIndex* MoFEM::BasicMethod::adjacenciesPtr |
raw pointer to container to adjacencies between dofs and finite elements
Definition at line 267 of file LoopMethods.hpp.
boost::weak_ptr<CacheTuple> MoFEM::BasicMethod::cacheWeakPtr |
Definition at line 351 of file LoopMethods.hpp.
const DofEntity_multiIndex* MoFEM::BasicMethod::dofsPtr |
raw pointer container of dofs
Definition at line 257 of file LoopMethods.hpp.
const FieldEntity_multiIndex* MoFEM::BasicMethod::entitiesPtr |
raw pointer to container of field entities
Definition at line 255 of file LoopMethods.hpp.
const Field_multiIndex* MoFEM::BasicMethod::fieldsPtr |
raw pointer to fields container
Definition at line 252 of file LoopMethods.hpp.
const EntFiniteElement_multiIndex* MoFEM::BasicMethod::finiteElementsEntitiesPtr |
raw pointer to container finite elements entities
Definition at line 263 of file LoopMethods.hpp.
const FiniteElement_multiIndex* MoFEM::BasicMethod::finiteElementsPtr |
raw pointer to container finite elements
Definition at line 260 of file LoopMethods.hpp.
std::pair<int, int> MoFEM::BasicMethod::loHiFERank |
Llo and hi processor rank of iterated entities.
Definition at line 217 of file LoopMethods.hpp.
int MoFEM::BasicMethod::loopSize |
local number oe methods to process
Definition at line 203 of file LoopMethods.hpp.
boost::movelib::unique_ptr<bool> MoFEM::BasicMethod::matAssembleSwitch |
Definition at line 349 of file LoopMethods.hpp.
int MoFEM::BasicMethod::nInTheLoop |
number currently of processed method
Definition at line 198 of file LoopMethods.hpp.
boost::function<MoFEMErrorCode()> MoFEM::BasicMethod::operatorHook |
Hook function for operator.
Definition at line 304 of file LoopMethods.hpp.
boost::function<MoFEMErrorCode()> MoFEM::BasicMethod::postProcessHook |
Hook function for post-processing.
Definition at line 299 of file LoopMethods.hpp.
boost::function<MoFEMErrorCode()> MoFEM::BasicMethod::preProcessHook |
Hook function for pre-processing.
Definition at line 294 of file LoopMethods.hpp.
const Problem* MoFEM::BasicMethod::problemPtr |
raw pointer to problem
Definition at line 250 of file LoopMethods.hpp.
int MoFEM::BasicMethod::rAnk |
processor rank
Definition at line 240 of file LoopMethods.hpp.
const RefEntity_multiIndex* MoFEM::BasicMethod::refinedEntitiesPtr |
container of mofem dof entities
Definition at line 245 of file LoopMethods.hpp.
const RefElement_multiIndex* MoFEM::BasicMethod::refinedFiniteElementsPtr |
container of mofem finite element entities
Definition at line 248 of file LoopMethods.hpp.
int MoFEM::BasicMethod::sIze |
number of processors in communicator
Definition at line 242 of file LoopMethods.hpp.
boost::movelib::unique_ptr<bool> MoFEM::BasicMethod::vecAssembleSwitch |
Definition at line 348 of file LoopMethods.hpp.