![]() |
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 |
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< CacheTuple > | getCacheWeakPtr () const |
Get the cache weak ptr object. | |
![]() | |
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 | |
![]() | |
PetscData () | |
virtual | ~PetscData ()=default |
MoFEMErrorCode | copyPetscData (const PetscData &petsc_data) |
![]() | |
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 |
![]() | |
MoFEMErrorCode | query_interface (boost::typeindex::type_index type_index, UnknownInterface **iface) const |
SnesMethod () | |
virtual | ~SnesMethod ()=default |
MoFEMErrorCode | copySnes (const SnesMethod &snes) |
Copy snes data. | |
![]() | |
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. | |
![]() | |
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_multiIndex * | refinedEntitiesPtr |
container of mofem dof entities | |
const RefElement_multiIndex * | refinedFiniteElementsPtr |
container of mofem finite element entities | |
const Problem * | problemPtr |
raw pointer to problem | |
const Field_multiIndex * | fieldsPtr |
raw pointer to fields container | |
const FieldEntity_multiIndex * | entitiesPtr |
raw pointer to container of field entities | |
const DofEntity_multiIndex * | dofsPtr |
raw pointer container of dofs | |
const FiniteElement_multiIndex * | finiteElementsPtr |
raw pointer to container finite elements | |
const EntFiniteElement_multiIndex * | finiteElementsEntitiesPtr |
const FieldEntityEntFiniteElementAdjacencyMap_multiIndex * | adjacenciesPtr |
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< bool > | vecAssembleSwitch |
boost::movelib::unique_ptr< bool > | matAssembleSwitch |
boost::weak_ptr< CacheTuple > | cacheWeakPtr |
![]() | |
KSPContext | ksp_ctx |
Context. | |
KSP | ksp |
KSP solver. | |
Vec & | ksp_f |
Mat & | ksp_A |
Mat & | ksp_B |
![]() | |
Switches | data_ctx |
Vec | f |
Mat | A |
Mat | B |
Vec | x |
Vec | dx |
Vec | x_t |
Vec | x_tt |
![]() | |
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 | |
![]() | |
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. | |
![]() | |
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 | |
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.
MoFEM::BasicMethod::BasicMethod | ( | ) |
Definition at line 130 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 138 of file LoopMethods.cpp.
|
inline |
Get the cache weak ptr object.
Definition at line 379 of file LoopMethods.hpp.
|
inline |
Definition at line 305 of file LoopMethods.hpp.
|
inline |
Get upper rank in loop for iterating elements.
Definition at line 273 of file LoopMethods.hpp.
|
inline |
Get upper rank in loop for iterating elements.
Definition at line 266 of file LoopMethods.hpp.
|
inline |
Get lo and hi processor rank of iterated entities.
Definition at line 259 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, 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 |
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::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, and EshelbianPlasticity::ContactTree.
Definition at line 170 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::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 |
Implements MoFEM::UnknownInterface.
Reimplemented in MoFEM::FEMethod, MoFEM::EntityMethod, and MoFEM::DofMethod.
Definition at line 220 of file LoopMethods.hpp.
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.
boost::weak_ptr<CacheTuple> MoFEM::BasicMethod::cacheWeakPtr |
Definition at line 386 of file LoopMethods.hpp.
const DofEntity_multiIndex* MoFEM::BasicMethod::dofsPtr |
raw pointer container of dofs
Definition at line 292 of file LoopMethods.hpp.
const FieldEntity_multiIndex* MoFEM::BasicMethod::entitiesPtr |
raw pointer to container of field entities
Definition at line 290 of file LoopMethods.hpp.
const Field_multiIndex* MoFEM::BasicMethod::fieldsPtr |
raw pointer to fields container
Definition at line 287 of file LoopMethods.hpp.
const EntFiniteElement_multiIndex* MoFEM::BasicMethod::finiteElementsEntitiesPtr |
raw pointer to container finite elements entities
Definition at line 298 of file LoopMethods.hpp.
const FiniteElement_multiIndex* MoFEM::BasicMethod::finiteElementsPtr |
raw pointer to container finite elements
Definition at line 295 of file LoopMethods.hpp.
std::pair<int, int> MoFEM::BasicMethod::loHiFERank |
Llo and hi processor rank of iterated entities.
Definition at line 252 of file LoopMethods.hpp.
int MoFEM::BasicMethod::loopSize |
local number oe methods to process
Definition at line 238 of file LoopMethods.hpp.
boost::movelib::unique_ptr<bool> MoFEM::BasicMethod::matAssembleSwitch |
Definition at line 384 of file LoopMethods.hpp.
int MoFEM::BasicMethod::nInTheLoop |
number currently of processed method
Definition at line 233 of file LoopMethods.hpp.
boost::function<MoFEMErrorCode()> MoFEM::BasicMethod::operatorHook |
Hook function for operator.
Definition at line 339 of file LoopMethods.hpp.
boost::function<MoFEMErrorCode()> MoFEM::BasicMethod::postProcessHook |
Hook function for post-processing.
Definition at line 334 of file LoopMethods.hpp.
boost::function<MoFEMErrorCode()> MoFEM::BasicMethod::preProcessHook |
Hook function for pre-processing.
Definition at line 329 of file LoopMethods.hpp.
const Problem* MoFEM::BasicMethod::problemPtr |
raw pointer to problem
Definition at line 285 of file LoopMethods.hpp.
int MoFEM::BasicMethod::rAnk |
processor rank
Definition at line 275 of file LoopMethods.hpp.
const RefEntity_multiIndex* MoFEM::BasicMethod::refinedEntitiesPtr |
container of mofem dof entities
Definition at line 280 of file LoopMethods.hpp.
const RefElement_multiIndex* MoFEM::BasicMethod::refinedFiniteElementsPtr |
container of mofem finite element entities
Definition at line 283 of file LoopMethods.hpp.
int MoFEM::BasicMethod::sIze |
number of processors in communicator
Definition at line 277 of file LoopMethods.hpp.
boost::movelib::unique_ptr<bool> MoFEM::BasicMethod::vecAssembleSwitch |
Definition at line 383 of file LoopMethods.hpp.