v0.14.0
Loading...
Searching...
No Matches
Public Member Functions | Public Attributes | Static Public Attributes | Private Attributes | List of all members
Monitor Struct Reference

[Push operators to pipeline] More...

Inheritance diagram for Monitor:
[legend]
Collaboration diagram for Monitor:
[legend]

Public Member Functions

 Monitor (SmartPetscObj< DM > dm, MoFEM::Interface &m_field, boost::shared_ptr< PostProcEle > post_proc, boost::shared_ptr< PostProcFaceEle > post_proc_bdry, boost::shared_ptr< MatrixDouble > velocity_field_ptr, boost::shared_ptr< MatrixDouble > x2_field_ptr, boost::shared_ptr< MatrixDouble > geometry_field_ptr, std::array< double, 3 > pass_field_eval_coords, boost::shared_ptr< SetPtsData > pass_field_eval_data)
 
MoFEMErrorCode postProcess ()
 function is run at the end of loop More...
 
 Monitor (SmartPetscObj< DM > dm, boost::shared_ptr< PostProcEle > post_proc)
 
MoFEMErrorCode preProcess ()
 function is run at the beginning of loop More...
 
MoFEMErrorCode operator() ()
 function is run for every finite element More...
 
MoFEMErrorCode postProcess ()
 function is run at the end of loop More...
 
 Monitor (SmartPetscObj< DM > dm, boost::shared_ptr< PostProcEle > post_proc)
 
MoFEMErrorCode preProcess ()
 function is run at the beginning of loop More...
 
MoFEMErrorCode operator() ()
 function is run for every finite element More...
 
MoFEMErrorCode postProcess ()
 function is run at the end of loop More...
 
 Monitor (SmartPetscObj< DM > dm, boost::shared_ptr< PostProcEleBdy > post_proc)
 
MoFEMErrorCode postProcess ()
 function is run at the end of loop More...
 
 Monitor (SmartPetscObj< DM > dm, boost::shared_ptr< PostProcEle > post_proc)
 
MoFEMErrorCode postProcess ()
 function is run at the end of loop More...
 
 Monitor (SmartPetscObj< DM > dm, boost::shared_ptr< PostProcEle > post_proc)
 
MoFEMErrorCode postProcess ()
 function is run at the end of loop More...
 
 Monitor (SmartPetscObj< DM > dm, boost::shared_ptr< moab::Core > post_proc_mesh, boost::shared_ptr< PostProcEleDomainCont > post_proc, boost::shared_ptr< PostProcEleBdyCont > post_proc_edge, std::pair< boost::shared_ptr< BoundaryEle >, boost::shared_ptr< VectorDouble > > p)
 
MoFEMErrorCode postProcess ()
 function is run at the end of loop More...
 
- Public Member Functions inherited from MoFEM::FEMethod
MoFEMErrorCode query_interface (boost::typeindex::type_index type_index, UnknownInterface **iface) const
 
 FEMethod ()=default
 
auto getFEName () const
 get finite element name More...
 
auto getDataDofsPtr () const
 
auto getDataVectorDofsPtr () const
 
const FieldEntity_vector_viewgetDataFieldEnts () const
 
boost::shared_ptr< FieldEntity_vector_view > & getDataFieldEntsPtr () const
 
auto & getRowFieldEnts () const
 
auto & getRowFieldEntsPtr () const
 
auto & getColFieldEnts () const
 
auto & getColFieldEntsPtr () const
 
auto getRowDofsPtr () const
 
auto getColDofsPtr () const
 
auto getNumberOfNodes () const
 
EntityHandle getFEEntityHandle () const
 
MoFEMErrorCode getNodeData (const std::string field_name, VectorDouble &data, const bool reset_dofs=true)
 
- Public Member Functions inherited from MoFEM::BasicMethod
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< CacheTuplegetCacheWeakPtr () const
 Get the cache weak ptr object. More...
 
- 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 More...
 
- Public Member Functions inherited from MoFEM::PetscData
MoFEMErrorCode query_interface (boost::typeindex::type_index type_index, UnknownInterface **iface) const
 
 PetscData ()
 
virtual ~PetscData ()=default
 
MoFEMErrorCode copyPetscData (const PetscData &petsc_data)
 
- Public Member Functions inherited from MoFEM::UnknownInterface
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
 
- 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. More...
 
- 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. More...
 

Public Attributes

MoFEM::InterfacemField
 
- Public Attributes inherited from MoFEM::FEMethod
std::string feName
 Name of finite element. More...
 
boost::shared_ptr< const NumeredEntFiniteElementnumeredEntFiniteElementPtr
 
boost::function< bool(FEMethod *fe_method_ptr)> exeTestHook
 Tet if element to skip element. More...
 
- Public Attributes inherited from MoFEM::BasicMethod
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_multiIndexrefinedEntitiesPtr
 container of mofem dof entities More...
 
const RefElement_multiIndexrefinedFiniteElementsPtr
 container of mofem finite element entities More...
 
const ProblemproblemPtr
 raw pointer to problem More...
 
const Field_multiIndexfieldsPtr
 raw pointer to fields container More...
 
const FieldEntity_multiIndexentitiesPtr
 raw pointer to container of field entities More...
 
const DofEntity_multiIndexdofsPtr
 raw pointer container of dofs More...
 
const FiniteElement_multiIndexfiniteElementsPtr
 raw pointer to container finite elements More...
 
const EntFiniteElement_multiIndexfiniteElementsEntitiesPtr
 
const FieldEntityEntFiniteElementAdjacencyMap_multiIndexadjacenciesPtr
 
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< boolvecAssembleSwitch
 
boost::movelib::unique_ptr< boolmatAssembleSwitch
 
boost::weak_ptr< CacheTuplecacheWeakPtr
 
- Public Attributes inherited from MoFEM::KspMethod
KSPContext ksp_ctx
 Context. More...
 
KSP ksp
 KSP solver. More...
 
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 x_t
 
Vec x_tt
 
- Public Attributes inherited from MoFEM::SnesMethod
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...
 
- Public Attributes inherited from MoFEM::TSMethod
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...
 

Static Public Attributes

static constexpr int saveEveryNthStep = 1
 
- 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 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)
 

Private Attributes

SmartPetscObj< DM > dM
 
boost::shared_ptr< PostProcElepostProc
 
boost::shared_ptr< PostProcFaceElepostProcBdy
 
boost::shared_ptr< MatrixDouble > velocityFieldPtr
 
boost::shared_ptr< MatrixDouble > x2FieldPtr
 
boost::shared_ptr< MatrixDouble > geometryFieldPtr
 
std::array< double, 3 > fieldEvalCoords
 
boost::shared_ptr< SetPtsDatafieldEvalData
 
boost::shared_ptr< PostProcEleBdypostProc
 
boost::shared_ptr< moab::Core > postProcMesh
 
boost::shared_ptr< PostProcEleDomainContpostProc
 
boost::shared_ptr< PostProcEleBdyContpostProcEdge
 
boost::shared_ptr< BoundaryEleliftFE
 
boost::shared_ptr< VectorDouble > liftVec
 

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_X_T = 1 << 4 , 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
}
 
- Static Public Member Functions inherited from MoFEM::UnknownInterface
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...
 

Detailed Description

[Push operators to pipeline]

[Push operators to pip]

Monitor solution.

Monitor solution

This functions is called by TS solver at the end of each step. It is used to output results to the hard drive.

This functions is called by TS solver at the end of each step. It is used to output results to the hard drive.

This functions is called by TS solver at the end of each step. It is used to output results to the hard drive.

Monitor solution

This functions is called by TS kso at the end of each step. It is used

Examples
PlasticOpsMonitor.hpp, dynamic_first_order_con_law.cpp, free_surface.cpp, heat_equation.cpp, nonlinear_elastic.cpp, photon_diffusion.cpp, plastic.cpp, reaction_diffusion.cpp, shallow_wave.cpp, and wave_equation.cpp.

Definition at line 783 of file dynamic_first_order_con_law.cpp.

Constructor & Destructor Documentation

◆ Monitor() [1/7]

Monitor::Monitor ( SmartPetscObj< DM >  dm,
MoFEM::Interface m_field,
boost::shared_ptr< PostProcEle post_proc,
boost::shared_ptr< PostProcFaceEle post_proc_bdry,
boost::shared_ptr< MatrixDouble >  velocity_field_ptr,
boost::shared_ptr< MatrixDouble >  x2_field_ptr,
boost::shared_ptr< MatrixDouble >  geometry_field_ptr,
std::array< double, 3 >  pass_field_eval_coords,
boost::shared_ptr< SetPtsData pass_field_eval_data 
)
inline

Definition at line 785 of file dynamic_first_order_con_law.cpp.

793 : dM(dm), mField(m_field), postProc(post_proc),
794 postProcBdy(post_proc_bdry), velocityFieldPtr(velocity_field_ptr),
795 x2FieldPtr(x2_field_ptr), geometryFieldPtr(geometry_field_ptr),
796 fieldEvalCoords(pass_field_eval_coords),
797 fieldEvalData(pass_field_eval_data){};
boost::shared_ptr< PostProcFaceEle > postProcBdy
std::array< double, 3 > fieldEvalCoords
MoFEM::Interface & mField
boost::shared_ptr< MatrixDouble > geometryFieldPtr
SmartPetscObj< DM > dM
boost::shared_ptr< MatrixDouble > velocityFieldPtr
boost::shared_ptr< SetPtsData > fieldEvalData
boost::shared_ptr< MatrixDouble > x2FieldPtr
boost::shared_ptr< PostProcEle > postProc

◆ Monitor() [2/7]

Monitor::Monitor ( SmartPetscObj< DM >  dm,
boost::shared_ptr< PostProcEle post_proc 
)
inline

Definition at line 68 of file heat_equation.cpp.

69 : dM(dm), postProc(post_proc){};

◆ Monitor() [3/7]

Monitor::Monitor ( SmartPetscObj< DM >  dm,
boost::shared_ptr< PostProcEle post_proc 
)
inline

Definition at line 65 of file wave_equation.cpp.

66 : dM(dm), postProc(post_proc){};

◆ Monitor() [4/7]

Monitor::Monitor ( SmartPetscObj< DM >  dm,
boost::shared_ptr< PostProcEleBdy post_proc 
)
inline

Definition at line 215 of file nonlinear_elastic.cpp.

216 : dM(dm), postProc(post_proc){};

◆ Monitor() [5/7]

Monitor::Monitor ( SmartPetscObj< DM >  dm,
boost::shared_ptr< PostProcEle post_proc 
)
inline

Definition at line 313 of file nonlinear_dynamic_elastic.cpp.

314 : dM(dm), postProc(post_proc){};

◆ Monitor() [6/7]

Monitor::Monitor ( SmartPetscObj< DM >  dm,
boost::shared_ptr< PostProcEle post_proc 
)
inline

Definition at line 731 of file shallow_wave.cpp.

732 : dM(dm), postProc(post_proc){};

◆ Monitor() [7/7]

Monitor::Monitor ( SmartPetscObj< DM >  dm,
boost::shared_ptr< moab::Core >  post_proc_mesh,
boost::shared_ptr< PostProcEleDomainCont post_proc,
boost::shared_ptr< PostProcEleBdyCont post_proc_edge,
std::pair< boost::shared_ptr< BoundaryEle >, boost::shared_ptr< VectorDouble > >  p 
)
inline

Definition at line 1990 of file free_surface.cpp.

1996 : dM(dm), postProcMesh(post_proc_mesh), postProc(post_proc),
1997 postProcEdge(post_proc_edge), liftFE(p.first), liftVec(p.second) {}
static Index< 'p', 3 > p
boost::shared_ptr< moab::Core > postProcMesh
boost::shared_ptr< VectorDouble > liftVec
boost::shared_ptr< BoundaryEle > liftFE
boost::shared_ptr< PostProcEleBdyCont > postProcEdge

Member Function Documentation

◆ operator()() [1/2]

MoFEMErrorCode Monitor::operator() ( )
inlinevirtual

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 from MoFEM::BasicMethod.

Examples
heat_equation.cpp, and wave_equation.cpp.

Definition at line 72 of file heat_equation.cpp.

72{ return 0; }

◆ operator()() [2/2]

MoFEMErrorCode Monitor::operator() ( )
inlinevirtual

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 from MoFEM::BasicMethod.

Definition at line 69 of file wave_equation.cpp.

69{ return 0; }

◆ postProcess() [1/7]

MoFEMErrorCode Monitor::postProcess ( )
inlinevirtual

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 from MoFEM::BasicMethod.

Examples
dynamic_first_order_con_law.cpp, free_surface.cpp, heat_equation.cpp, nonlinear_elastic.cpp, shallow_wave.cpp, and wave_equation.cpp.

Definition at line 798 of file dynamic_first_order_con_law.cpp.

798 {
800
801 // cerr << "wagawaga\n";
802 auto *simple = mField.getInterface<Simple>();
803
804 if (SPACE_DIM == 3) {
805 CHKERR mField.getInterface<FieldEvaluatorInterface>()->evalFEAtThePoint3D(
806 fieldEvalCoords.data(), 1e-12, simple->getProblemName(),
807 simple->getDomainFEName(), fieldEvalData, mField.get_comm_rank(),
809 } else {
810 CHKERR mField.getInterface<FieldEvaluatorInterface>()->evalFEAtThePoint2D(
811 fieldEvalCoords.data(), 1e-12, simple->getProblemName(),
812 simple->getDomainFEName(), fieldEvalData, mField.get_comm_rank(),
814 }
815
816 if (velocityFieldPtr->size1()) {
817 auto t_vel = getFTensor1FromMat<SPACE_DIM>(*velocityFieldPtr);
818 auto t_x2_field = getFTensor1FromMat<SPACE_DIM>(*x2FieldPtr);
819 auto t_geom = getFTensor1FromMat<SPACE_DIM>(*geometryFieldPtr);
820
821 double u_x = t_x2_field(0) - t_geom(0);
822 double u_y = t_x2_field(1) - t_geom(1);
823 double u_z = t_x2_field(2) - t_geom(2);
824
825 MOFEM_LOG("SYNC", Sev::inform)
826 << "Velocities x: " << t_vel(0) << " y: " << t_vel(1)
827 << " z: " << t_vel(2) << "\n";
828 MOFEM_LOG("SYNC", Sev::inform) << "Displacement x: " << u_x
829 << " y: " << u_y << " z: " << u_z << "\n";
830 }
831
833 std::regex((boost::format("%s(.*)") % "Data_Vertex").str()))) {
834 Range ents;
835 mField.get_moab().get_entities_by_dimension(m->getMeshset(), 0, ents,
836 true);
837 auto print_vets = [](boost::shared_ptr<FieldEntity> ent_ptr) {
839 if (!(ent_ptr->getPStatus() & PSTATUS_NOT_OWNED)) {
840 MOFEM_LOG("SYNC", Sev::inform)
841 << "Velocities: " << ent_ptr->getEntFieldData()[0] << " "
842 << ent_ptr->getEntFieldData()[1] << " "
843 << ent_ptr->getEntFieldData()[2] << "\n";
844 }
846 };
847 CHKERR mField.getInterface<FieldBlas>()->fieldLambdaOnEntities(
848 print_vets, "V", &ents);
849 }
851
852 PetscBool print_volume = PETSC_FALSE;
853 CHKERR PetscOptionsGetBool(PETSC_NULL, "", "-print_volume", &print_volume,
854 PETSC_NULL);
855
856 PetscBool print_skin = PETSC_FALSE;
857 CHKERR PetscOptionsGetBool(PETSC_NULL, "", "-print_skin", &print_skin,
858 PETSC_NULL);
859
860 int save_every_nth_step = 1;
861 CHKERR PetscOptionsGetInt(PETSC_NULL, "", "-save_step",
862 &save_every_nth_step, PETSC_NULL);
863 if (ts_step % save_every_nth_step == 0) {
864
865 if (print_volume) {
867 CHKERR postProc->writeFile(
868 "out_step_" + boost::lexical_cast<std::string>(ts_step) + ".h5m");
869 }
870
871 if (print_skin) {
873 CHKERR postProcBdy->writeFile(
874 "out_boundary_" + boost::lexical_cast<std::string>(ts_step) +
875 ".h5m");
876 }
877 }
879 }
#define MOFEM_LOG_SEVERITY_SYNC(comm, severity)
Synchronise "SYNC" on curtain severity level.
Definition: LogManager.hpp:352
void simple(double P1[], double P2[], double P3[], double c[], const int N)
Definition: acoustic.cpp:69
@ QUIET
Definition: definitions.h:208
@ MF_EXIST
Definition: definitions.h:100
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:346
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:416
#define CHKERR
Inline error check.
Definition: definitions.h:535
FTensor::Index< 'm', SPACE_DIM > m
constexpr int SPACE_DIM
PetscErrorCode DMoFEMLoopFiniteElements(DM dm, const char fe_name[], MoFEM::FEMethod *method, CacheTupleWeakPtr cache_ptr=CacheTupleSharedPtr())
Executes FEMethod for finite elements in DM.
Definition: DMMoFEM.cpp:572
#define MOFEM_LOG(channel, severity)
Log.
Definition: LogManager.hpp:308
MoFEMErrorCode getCubitMeshsetPtr(const int ms_id, const CubitBCType cubit_bc_type, const CubitMeshSets **cubit_meshset_ptr) const
get cubit meshset
PetscErrorCode PetscOptionsGetInt(PetscOptions *, const char pre[], const char name[], PetscInt *ivalue, PetscBool *set)
PetscErrorCode PetscOptionsGetBool(PetscOptions *, const char pre[], const char name[], PetscBool *bval, PetscBool *set)
int save_every_nth_step
boost::weak_ptr< CacheTuple > getCacheWeakPtr() const
Get the cache weak ptr object.
virtual moab::Interface & get_moab()=0
virtual MPI_Comm & get_comm() const =0
virtual int get_comm_rank() const =0
Basic algebra on fields.
Definition: FieldBlas.hpp:21
Field evaluator interface.
Interface for managing meshsets containing materials and boundary conditions.
Simple interface for fast problem set-up.
Definition: Simple.hpp:27
PetscInt ts_step
time step number
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface refernce to pointer of interface.

◆ postProcess() [2/7]

MoFEMErrorCode Monitor::postProcess ( )
inlinevirtual

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 from MoFEM::BasicMethod.

Definition at line 76 of file heat_equation.cpp.

76 {
78 if (ts_step % saveEveryNthStep == 0) {
80 CHKERR postProc->writeFile(
81 "out_level_" + boost::lexical_cast<std::string>(ts_step) + ".h5m");
82 }
84 }
static constexpr int saveEveryNthStep

◆ postProcess() [3/7]

MoFEMErrorCode Monitor::postProcess ( )
inlinevirtual

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 from MoFEM::BasicMethod.

Definition at line 73 of file wave_equation.cpp.

73 {
75 if (ts_step % saveEveryNthStep == 0) {
77 CHKERR postProc->writeFile(
78 "out_level_" + boost::lexical_cast<std::string>(ts_step) + ".h5m");
79 }
81 }

◆ postProcess() [4/7]

MoFEMErrorCode Monitor::postProcess ( )
inlinevirtual

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 from MoFEM::BasicMethod.

Definition at line 217 of file nonlinear_elastic.cpp.

217 {
219 constexpr int save_every_nth_step = 1;
220 if (ts_step % save_every_nth_step == 0) {
222 CHKERR postProc->writeFile(
223 "out_step_" + boost::lexical_cast<std::string>(ts_step) + ".h5m");
224 }
226 }

◆ postProcess() [5/7]

MoFEMErrorCode Monitor::postProcess ( )
inlinevirtual

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 from MoFEM::BasicMethod.

Definition at line 315 of file nonlinear_dynamic_elastic.cpp.

315 {
317 constexpr int save_every_nth_step = 1;
318 if (ts_step % save_every_nth_step == 0) {
320 CHKERR postProc->writeFile(
321 "out_step_" + boost::lexical_cast<std::string>(ts_step) + ".h5m");
322 }
324 }

◆ postProcess() [6/7]

MoFEMErrorCode Monitor::postProcess ( )
inlinevirtual

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 from MoFEM::BasicMethod.

Definition at line 733 of file shallow_wave.cpp.

733 {
735 constexpr int save_every_nth_step = 50;
736 if (ts_step % save_every_nth_step == 0) {
738 CHKERR postProc->writeFile(
739 "out_step_" + boost::lexical_cast<std::string>(ts_step) + ".h5m");
740 MOFEM_LOG("SW", Sev::verbose)
741 << "writing vector in binary to vector.dat ...";
742 PetscViewer viewer;
743 PetscViewerBinaryOpen(PETSC_COMM_WORLD, "vector.dat", FILE_MODE_WRITE,
744 &viewer);
745 VecView(ts_u, viewer);
746 PetscViewerDestroy(&viewer);
747 }
749 }
Vec & ts_u
state vector

◆ postProcess() [7/7]

MoFEMErrorCode Monitor::postProcess ( )
inlinevirtual

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 from MoFEM::BasicMethod.

Definition at line 1998 of file free_surface.cpp.

1998 {
2000
2001 MOFEM_LOG("FS", Sev::verbose) << "Monitor";
2002 constexpr int save_every_nth_step = 1;
2003 if (ts_step % save_every_nth_step == 0) {
2004 MOFEM_LOG("FS", Sev::verbose) << "Mesh pre proc";
2005 MoFEM::Interface *m_field_ptr;
2006 CHKERR DMoFEMGetInterfacePtr(dM, &m_field_ptr);
2007 auto post_proc_begin =
2008 boost::make_shared<PostProcBrokenMeshInMoabBaseBegin>(*m_field_ptr,
2009 postProcMesh);
2010 auto post_proc_end = boost::make_shared<PostProcBrokenMeshInMoabBaseEnd>(
2011 *m_field_ptr, postProcMesh);
2012 CHKERR DMoFEMPreProcessFiniteElements(dM, post_proc_begin->getFEMethod());
2014 this->getCacheWeakPtr());
2016 this->getCacheWeakPtr());
2017 CHKERR DMoFEMPostProcessFiniteElements(dM, post_proc_end->getFEMethod());
2018 CHKERR post_proc_end->writeFile(
2019 "out_step_" + boost::lexical_cast<std::string>(ts_step) + ".h5m");
2020 MOFEM_LOG("FS", Sev::verbose) << "Mesh pre proc done";
2021 }
2022
2023 liftVec->resize(SPACE_DIM, false);
2024 liftVec->clear();
2026 MPI_Allreduce(MPI_IN_PLACE, &(*liftVec)[0], SPACE_DIM, MPI_DOUBLE, MPI_SUM,
2027 MPI_COMM_WORLD);
2028 MOFEM_LOG("FS", Sev::inform)
2029 << "Step " << ts_step << " time " << ts_t
2030 << " lift vec x: " << (*liftVec)[0] << " y: " << (*liftVec)[1];
2031
2033 }
PetscErrorCode DMoFEMPostProcessFiniteElements(DM dm, MoFEM::FEMethod *method)
execute finite element method for each element in dm (problem)
Definition: DMMoFEM.cpp:542
PetscErrorCode DMoFEMGetInterfacePtr(DM dm, MoFEM::Interface **m_field_ptr)
Get pointer to MoFEM::Interface.
Definition: DMMoFEM.cpp:400
PetscErrorCode DMoFEMPreProcessFiniteElements(DM dm, MoFEM::FEMethod *method)
execute finite element method for each element in dm (problem)
Definition: DMMoFEM.cpp:532
Deprecated interface functions.
PetscReal ts_t
time

◆ preProcess() [1/2]

MoFEMErrorCode Monitor::preProcess ( )
inlinevirtual

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 from MoFEM::BasicMethod.

Examples
heat_equation.cpp, and wave_equation.cpp.

Definition at line 71 of file heat_equation.cpp.

71{ return 0; }

◆ preProcess() [2/2]

MoFEMErrorCode Monitor::preProcess ( )
inlinevirtual

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 from MoFEM::BasicMethod.

Definition at line 68 of file wave_equation.cpp.

68{ return 0; }

Member Data Documentation

◆ dM

SmartPetscObj< DM > Monitor::dM
private

◆ fieldEvalCoords

std::array<double, 3> Monitor::fieldEvalCoords
private

◆ fieldEvalData

boost::shared_ptr<SetPtsData> Monitor::fieldEvalData
private

◆ geometryFieldPtr

boost::shared_ptr<MatrixDouble> Monitor::geometryFieldPtr
private

◆ liftFE

boost::shared_ptr<BoundaryEle> Monitor::liftFE
private
Examples
free_surface.cpp.

Definition at line 2040 of file free_surface.cpp.

◆ liftVec

boost::shared_ptr<VectorDouble> Monitor::liftVec
private
Examples
free_surface.cpp.

Definition at line 2041 of file free_surface.cpp.

◆ mField

MoFEM::Interface& Monitor::mField

◆ postProc [1/3]

boost::shared_ptr< PostProcEle > Monitor::postProc
private

◆ postProc [2/3]

boost::shared_ptr<PostProcEleBdy> Monitor::postProc
private

Definition at line 230 of file nonlinear_elastic.cpp.

◆ postProc [3/3]

boost::shared_ptr<PostProcEleDomainCont> Monitor::postProc
private

Definition at line 2038 of file free_surface.cpp.

◆ postProcBdy

boost::shared_ptr<PostProcFaceEle> Monitor::postProcBdy
private

◆ postProcEdge

boost::shared_ptr<PostProcEleBdyCont> Monitor::postProcEdge
private
Examples
free_surface.cpp.

Definition at line 2039 of file free_surface.cpp.

◆ postProcMesh

boost::shared_ptr<moab::Core> Monitor::postProcMesh
private
Examples
free_surface.cpp.

Definition at line 2037 of file free_surface.cpp.

◆ saveEveryNthStep

static constexpr int Monitor::saveEveryNthStep = 1
staticconstexpr
Examples
heat_equation.cpp, and wave_equation.cpp.

Definition at line 74 of file heat_equation.cpp.

◆ velocityFieldPtr

boost::shared_ptr<MatrixDouble> Monitor::velocityFieldPtr
private

◆ x2FieldPtr

boost::shared_ptr<MatrixDouble> Monitor::x2FieldPtr
private

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