v0.15.0
Loading...
Searching...
No Matches
Public Member Functions | Public Attributes | List of all members
MonitorPostProc Struct Reference
Inheritance diagram for MonitorPostProc:
[legend]
Collaboration diagram for MonitorPostProc:
[legend]

Public Member Functions

 MonitorPostProc (MoFEM::Interface &m_field, std::map< int, NonlinearElasticElement::BlockData > &set_of_blocks, NonlinearElasticElement::MyVolumeFE &fe_elastic_energy, ConvectiveMassElement::MyVolumeFE &fe_kinetic_energy)
 
MoFEMErrorCode preProcess ()
 Pre-processing function executed at loop initialization.
 
MoFEMErrorCode operator() ()
 Main operator function executed for each loop iteration.
 
MoFEMErrorCode postProcess ()
 Post-processing function executed at loop completion.
 
 MonitorPostProc (MoFEM::Interface &m_field)
 
MoFEMErrorCode preProcess ()
 Pre-processing function executed at loop initialization.
 
MoFEMErrorCode operator() ()
 Main operator function executed for each loop iteration.
 
MoFEMErrorCode postProcess ()
 Post-processing function executed at loop completion.
 
- Public Member Functions inherited from MoFEM::FEMethod
MoFEMErrorCode query_interface (boost::typeindex::type_index type_index, UnknownInterface **iface) const
 Query interface for type casting.
 
 FEMethod ()=default
 Default constructor.
 
auto getFEName () const
 Get the name of the current finite element.
 
auto getDataDofsPtr () const
 Get pointer to DOF data for the current finite element.
 
auto getDataVectorDofsPtr () const
 Get pointer to vector DOF data for the current finite element.
 
const FieldEntity_vector_viewgetDataFieldEnts () const
 Get reference to data field entities for the current finite element.
 
boost::shared_ptr< FieldEntity_vector_view > & getDataFieldEntsPtr () const
 Get shared pointer to data field entities for the current finite element.
 
auto & getRowFieldEnts () const
 Get reference to row field entities for the current finite element.
 
auto & getRowFieldEntsPtr () const
 Get shared pointer to row field entities for the current finite element.
 
auto & getColFieldEnts () const
 Get reference to column field entities for the current finite element.
 
auto & getColFieldEntsPtr () const
 Get shared pointer to column field entities for the current finite element.
 
auto getRowDofsPtr () const
 Get pointer to row DOFs for the current finite element.
 
auto getColDofsPtr () const
 Get pointer to column DOFs for the current finite element.
 
auto getNumberOfNodes () const
 Get the number of nodes in the current finite element.
 
EntityHandle getFEEntityHandle () const
 Get the entity handle of the current finite element.
 
MoFEMErrorCode getNodeData (const std::string field_name, VectorDouble &data, const bool reset_dofs=true)
 Get nodal data for a specific field.
 
- Public Member Functions inherited from MoFEM::BasicMethod
 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.
 
boost::weak_ptr< CacheTuplegetCacheWeakPtr () 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

MoFEM::InterfacemField
 
PostProcVolumeOnRefinedMesh postProc
 
std::map< int, NonlinearElasticElement::BlockData > & setOfBlocks
 
NonlinearElasticElement::MyVolumeFEfeElasticEnergy
 calculate elastic energy
 
ConvectiveMassElement::MyVolumeFEfeKineticEnergy
 calculate elastic energy
 
bool iNit
 
int pRT
 
int * step
 
PostProcFaceOnRefinedMesh skinPostProc
 
PetscBool saveSkin
 
- Public Attributes inherited from MoFEM::FEMethod
std::string feName
 Name of the finite element being processed.
 
boost::shared_ptr< const NumeredEntFiniteElementnumeredEntFiniteElementPtr
 Shared pointer to finite element database structure.
 
boost::function< bool(FEMethod *fe_method_ptr)> exeTestHook
 Test function to determine if element should be skipped.
 
- Public Attributes inherited from MoFEM::BasicMethod
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_multiIndexrefinedEntitiesPtr
 Pointer to container of refined MoFEM DOF entities.
 
const RefElement_multiIndexrefinedFiniteElementsPtr
 Pointer to container of refined finite element entities.
 
const ProblemproblemPtr
 Raw pointer to current MoFEM problem instance.
 
const Field_multiIndexfieldsPtr
 Raw pointer to fields multi-index container.
 
const FieldEntity_multiIndexentitiesPtr
 Raw pointer to container of field entities.
 
const DofEntity_multiIndexdofsPtr
 Raw pointer to container of degree of freedom entities.
 
const FiniteElement_multiIndexfiniteElementsPtr
 Raw pointer to container of finite elements.
 
const EntFiniteElement_multiIndexfiniteElementsEntitiesPtr
 Raw pointer to container of finite element entities.
 
const FieldEntityEntFiniteElementAdjacencyMap_multiIndexadjacenciesPtr
 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< boolvecAssembleSwitch
 Switch for vector assembly operations.
 
boost::movelib::unique_ptr< boolmatAssembleSwitch
 Switch for matrix assembly operations.
 
boost::weak_ptr< CacheTuplecacheWeakPtr
 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.
 

Detailed Description

Examples
mofem/tutorials/cor-0to1/src/UnsaturatedFlow.hpp, mofem/users_modules/basic_finite_elements/nonlinear_elasticity/nonlinear_dynamics.cpp, mofem/users_modules/bone_remodelling/src/Remodeling.hpp, and mofem/users_modules/bone_remodelling/src/impl/Remodeling.cpp.

Definition at line 30 of file nonlinear_dynamics.cpp.

Constructor & Destructor Documentation

◆ MonitorPostProc() [1/2]

MonitorPostProc::MonitorPostProc ( MoFEM::Interface m_field,
std::map< int, NonlinearElasticElement::BlockData > &  set_of_blocks,
NonlinearElasticElement::MyVolumeFE fe_elastic_energy,
ConvectiveMassElement::MyVolumeFE fe_kinetic_energy 
)
inline

Definition at line 45 of file nonlinear_dynamics.cpp.

50 : FEMethod(), mField(m_field), postProc(m_field),
51 setOfBlocks(set_of_blocks), feElasticEnergy(fe_elastic_energy),
52 feKineticEnergy(fe_kinetic_energy), iNit(false) {
53
54 double def_t_val = 0;
55 const EntityHandle root_meshset = mField.get_moab().get_root_set();
56
57 Tag th_step;
58 rval = m_field.get_moab().tag_get_handle(
59 "_TsStep_", 1, MB_TYPE_INTEGER, th_step,
60 MB_TAG_CREAT | MB_TAG_EXCL | MB_TAG_MESH, &def_t_val);
61 if (rval == MB_ALREADY_ALLOCATED) {
62 MOAB_THROW(m_field.get_moab().tag_get_by_ptr(th_step, &root_meshset, 1,
63 (const void **)&step));
64 } else {
65 MOAB_THROW(m_field.get_moab().tag_set_data(th_step, &root_meshset, 1,
66 &def_t_val));
67 MOAB_THROW(m_field.get_moab().tag_get_by_ptr(th_step, &root_meshset, 1,
68 (const void **)&step));
69 }
70
71 PetscBool flg = PETSC_TRUE;
72 CHK_THROW_MESSAGE(PetscOptionsGetInt(PETSC_NULLPTR, PETSC_NULLPTR,
73 "-my_output_prt", &pRT, &flg),
74 "Can not get option");
75 if (flg != PETSC_TRUE) {
76 pRT = 10;
77 }
78 }
#define MOAB_THROW(err)
Check error code of MoAB function and throw MoFEM exception.
#define CHK_THROW_MESSAGE(err, msg)
Check and throw MoFEM exception.
static MoFEMErrorCodeGeneric< moab::ErrorCode > rval
PetscErrorCode PetscOptionsGetInt(PetscOptions *, const char pre[], const char name[], PetscInt *ivalue, PetscBool *set)
virtual moab::Interface & get_moab()=0
FEMethod()=default
Default constructor.
NonlinearElasticElement::MyVolumeFE & feElasticEnergy
calculate elastic energy
MoFEM::Interface & mField
std::map< int, NonlinearElasticElement::BlockData > & setOfBlocks
ConvectiveMassElement::MyVolumeFE & feKineticEnergy
calculate elastic energy
PostProcVolumeOnRefinedMesh postProc

◆ MonitorPostProc() [2/2]

MonitorPostProc::MonitorPostProc ( MoFEM::Interface m_field)
inline

Definition at line 60 of file thermal_unsteady.cpp.

61 : FEMethod(), mField(m_field), postProc(m_field), skinPostProc(m_field),
62 iNit(false) {
63
64 PetscBool flg = PETSC_TRUE;
65 CHKERR PetscOptionsGetInt(PETSC_NULLPTR, PETSC_NULLPTR, "-my_output_prt", &pRT,
66 &flg);
67 CHKERRABORT(PETSC_COMM_WORLD, ierr);
68 if (flg != PETSC_TRUE) {
69 pRT = 1;
70 }
71 saveSkin = PETSC_TRUE;
72 CHKERR PetscOptionsGetBool(PETSC_NULLPTR, PETSC_NULLPTR, "-my_save_skin",
73 &saveSkin, PETSC_NULLPTR);
74 }
#define CHKERR
Inline error check.
static MoFEMErrorCodeGeneric< PetscErrorCode > ierr
PetscErrorCode PetscOptionsGetBool(PetscOptions *, const char pre[], const char name[], PetscBool *bval, PetscBool *set)
PostProcFaceOnRefinedMesh skinPostProc

Member Function Documentation

◆ operator()() [1/2]

MoFEMErrorCode MonitorPostProc::operator() ( )
inlinevirtual

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:

  • Calculating element local matrices and vectors
  • Assembling contributions to global system
  • Element-level post-processing operations
Returns
Error code indicating success or failure

Reimplemented from MoFEM::BasicMethod.

Examples
mofem/users_modules/basic_finite_elements/nonlinear_elasticity/nonlinear_dynamics.cpp.

Definition at line 128 of file nonlinear_dynamics.cpp.

128 {
131 }
#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 ...

◆ operator()() [2/2]

MoFEMErrorCode MonitorPostProc::operator() ( )
inlinevirtual

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:

  • Calculating element local matrices and vectors
  • Assembling contributions to global system
  • Element-level post-processing operations
Returns
Error code indicating success or failure

Reimplemented from MoFEM::BasicMethod.

Definition at line 81 of file thermal_unsteady.cpp.

◆ postProcess() [1/2]

MoFEMErrorCode MonitorPostProc::postProcess ( )
inlinevirtual

Post-processing function executed at loop completion.

This virtual function is called once at the end of the loop. It is typically used for:

  • Assembling matrices and vectors to global system
  • Computing global variables (e.g., total internal energy)
  • Finalizing computation results
  • Cleanup operations

Example of iterating over DOFs:

for(_IT_GET_FEROW_BY_NAME_DOFS_FOR_LOOP_(this,"DISPLACEMENT",it)) {
// Process each DOF
}
Returns
Error code indicating success or failure

Reimplemented from MoFEM::BasicMethod.

Examples
mofem/users_modules/basic_finite_elements/nonlinear_elasticity/nonlinear_dynamics.cpp.

Definition at line 133 of file nonlinear_dynamics.cpp.

◆ postProcess() [2/2]

MoFEMErrorCode MonitorPostProc::postProcess ( )
inlinevirtual

Post-processing function executed at loop completion.

This virtual function is called once at the end of the loop. It is typically used for:

  • Assembling matrices and vectors to global system
  • Computing global variables (e.g., total internal energy)
  • Finalizing computation results
  • Cleanup operations

Example of iterating over DOFs:

for(_IT_GET_FEROW_BY_NAME_DOFS_FOR_LOOP_(this,"DISPLACEMENT",it)) {
// Process each DOF
}
Returns
Error code indicating success or failure

Reimplemented from MoFEM::BasicMethod.

Definition at line 86 of file thermal_unsteady.cpp.

86 {
88
89 if (!iNit) {
90 CHKERR addHOOpsVol("MESH_NODE_POSITIONS", postProc, true, false, false,
91 false);
96 CHKERR postProc.addFieldValuesPostProc("MESH_NODE_POSITIONS");
97
98 CHKERR addHOOpsFace3D("MESH_NODE_POSITIONS", skinPostProc, false, false);
101
102 iNit = true;
103 }
104 int step;
105 CHKERR TSGetTimeStepNumber(ts, &step);
106
107 if (pRT && (step) % pRT == 0) {
108 // CHKERR mField.loop_finite_elements("DMTHERMAL","THERMAL_FE",postProc);
109 // std::ostringstream sss;
110 // sss << "out_thermal_" << step << ".h5m";
111 // CHKERR postProc.writeFile(sss.str().c_str());
112 if (saveSkin) {
113 CHKERR mField.loop_finite_elements("DMTHERMAL", "POST_PROC_SKIN",
115 std::ostringstream sss;
116 sss << "out_skin_" << step << ".h5m";
117 CHKERR skinPostProc.writeFile(sss.str().c_str());
118 }
119 }
121 }
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
MoFEMErrorCode addFieldValuesPostProc(const std::string field_name, Vec v=PETSC_NULLPTR)
Add operator to post-process L2, H1, Hdiv, Hcurl field value.
MoFEMErrorCode addFieldValuesGradientPostProc(const std::string field_name, Vec v=PETSC_NULLPTR)
Add operator to post-process L2 or H1 field gradient.
MoFEMErrorCode writeFile(const std::string file_name, const char *file_type="MOAB", const char *file_options="PARALLEL=WRITE_PART")
wrote results in (MOAB) format, use "file_name.h5m"
virtual MoFEMErrorCode loop_finite_elements(const std::string problem_name, const std::string &fe_name, FEMethod &method, boost::shared_ptr< NumeredEntFiniteElement_multiIndex > fe_ptr=nullptr, MoFEMTypes bh=MF_EXIST, CacheTupleWeakPtr cache_ptr=CacheTupleSharedPtr(), int verb=DEFAULT_VERBOSITY)=0
Make a loop over finite elements.
MoFEMErrorCode addHOOpsFace3D(const std::string field, E &e, bool hcurl, bool hdiv)
MoFEMErrorCode addHOOpsVol(const std::string field, E &e, bool h1, bool hcurl, bool hdiv, bool l2)
TS ts
PETSc time stepping solver object.
MoFEMErrorCode generateReferenceElementMesh()
MoFEMErrorCode generateReferenceElementMesh()
Generate reference mesh on single element.

◆ preProcess() [1/2]

MoFEMErrorCode MonitorPostProc::preProcess ( )
inlinevirtual

Pre-processing function executed at loop initialization.

This virtual function is called once at the beginning of the loop. It is typically used for:

  • Zeroing matrices and vectors
  • Computing shape functions on reference elements
  • Preprocessing boundary conditions
  • Setting up initial computation state
Returns
Error code indicating success or failure

Reimplemented from MoFEM::BasicMethod.

Examples
mofem/users_modules/basic_finite_elements/nonlinear_elasticity/nonlinear_dynamics.cpp.

Definition at line 80 of file nonlinear_dynamics.cpp.

80 {
82
83 if (!iNit) {
85 if(mField.check_field("MESH_NODE_POSITIONS"))
86 CHKERR addHOOpsVol("MESH_NODE_POSITIONS", postProc, true, false, false,
87 false);
90 CHKERR postProc.addFieldValuesPostProc("MESH_NODE_POSITIONS");
92
93 std::map<int, NonlinearElasticElement::BlockData>::iterator sit =
94 setOfBlocks.begin();
95 for (; sit != setOfBlocks.end(); sit++) {
96 postProc.getOpPtrVector().push_back(new PostProcStress(
98 sit->second, postProc.commonData, true));
99 }
100
101 iNit = true;
102 }
103
104 if ((*step) % pRT == 0) {
105 CHKERR mField.loop_finite_elements("DYNAMICS", "MASS_ELEMENT", postProc);
106 std::ostringstream sss;
107 sss << "out_values_" << (*step) << ".h5m";
108 CHKERR postProc.writeFile(sss.str().c_str());
109 }
110
113 CHKERR mField.loop_finite_elements("DYNAMICS", "ELASTIC", feElasticEnergy);
114
116 CHKERR mField.loop_finite_elements("DYNAMICS", "MASS_ELEMENT",
118 double E = feElasticEnergy.eNergy;
119 double T = feKineticEnergy.eNergy;
121 "DYNAMIC", Sev::inform,
122 "%d Time %3.2e Elastic energy %3.2e Kinetic Energy %3.2e Total %3.2e\n",
123 ts_step, ts_t, E, T, E + T);
124
126 }
#define MOFEM_LOG_C(channel, severity, format,...)
virtual bool check_field(const std::string &name) const =0
check if field is in database
@ CTX_SNESNONE
No specific SNES context.
SNESContext snes_ctx
Current SNES computation context.
@ CTX_TSNONE
No specific TS context.
PetscReal ts_t
Current time value.
TSContext ts_ctx
Current TS computation context.
PetscInt ts_step
Current time step number.
std::vector< EntityHandle > mapGaussPts

◆ preProcess() [2/2]

MoFEMErrorCode MonitorPostProc::preProcess ( )
inlinevirtual

Pre-processing function executed at loop initialization.

This virtual function is called once at the beginning of the loop. It is typically used for:

  • Zeroing matrices and vectors
  • Computing shape functions on reference elements
  • Preprocessing boundary conditions
  • Setting up initial computation state
Returns
Error code indicating success or failure

Reimplemented from MoFEM::BasicMethod.

Definition at line 76 of file thermal_unsteady.cpp.

Member Data Documentation

◆ feElasticEnergy

NonlinearElasticElement::MyVolumeFE& MonitorPostProc::feElasticEnergy

◆ feKineticEnergy

ConvectiveMassElement::MyVolumeFE& MonitorPostProc::feKineticEnergy

◆ iNit

bool MonitorPostProc::iNit

◆ mField

MoFEM::Interface & MonitorPostProc::mField

◆ postProc

PostProcVolumeOnRefinedMesh MonitorPostProc::postProc

◆ pRT

int MonitorPostProc::pRT

◆ saveSkin

PetscBool MonitorPostProc::saveSkin

Definition at line 58 of file thermal_unsteady.cpp.

◆ setOfBlocks

std::map<int, NonlinearElasticElement::BlockData>& MonitorPostProc::setOfBlocks

◆ skinPostProc

PostProcFaceOnRefinedMesh MonitorPostProc::skinPostProc

Definition at line 54 of file thermal_unsteady.cpp.

◆ step

int* MonitorPostProc::step

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