v0.13.2
Loading...
Searching...
No Matches
Public Member Functions | Private Attributes | List of all members
PlasticOps::Monitor Struct Reference

#include <users_modules/tutorials/adv-0/src/PlasticOpsMonitor.hpp>

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

Public Member Functions

 Monitor (SmartPetscObj< DM > &dm, boost::shared_ptr< PostProcEle > post_proc_fe, boost::shared_ptr< DomainEle > reaction_fe, std::tuple< SmartPetscObj< Vec >, SmartPetscObj< VecScatter > > ux_scatter, std::tuple< SmartPetscObj< Vec >, SmartPetscObj< VecScatter > > uy_scatter, std::tuple< SmartPetscObj< Vec >, SmartPetscObj< VecScatter > > uz_scatter)
 
MoFEMErrorCode preProcess ()
 
MoFEMErrorCode operator() ()
 
MoFEMErrorCode postProcess ()
 

Private Attributes

SmartPetscObj< DM > dM
 
boost::shared_ptr< PostProcElepostProcFe
 
boost::shared_ptr< DomainElereactionFe
 
std::tuple< SmartPetscObj< Vec >, SmartPetscObj< VecScatter > > uXScatter
 
std::tuple< SmartPetscObj< Vec >, SmartPetscObj< VecScatter > > uYScatter
 
std::tuple< SmartPetscObj< Vec >, SmartPetscObj< VecScatter > > uZScatter
 

Detailed Description

Definition at line 9 of file PlasticOpsMonitor.hpp.

Constructor & Destructor Documentation

◆ Monitor()

PlasticOps::Monitor::Monitor ( SmartPetscObj< DM > &  dm,
boost::shared_ptr< PostProcEle post_proc_fe,
boost::shared_ptr< DomainEle reaction_fe,
std::tuple< SmartPetscObj< Vec >, SmartPetscObj< VecScatter > >  ux_scatter,
std::tuple< SmartPetscObj< Vec >, SmartPetscObj< VecScatter > >  uy_scatter,
std::tuple< SmartPetscObj< Vec >, SmartPetscObj< VecScatter > >  uz_scatter 
)
inline

Definition at line 11 of file PlasticOpsMonitor.hpp.

16 : dM(dm), postProcFe(post_proc_fe), reactionFe(reaction_fe),
17 uXScatter(ux_scatter), uYScatter(uy_scatter), uZScatter(uz_scatter){};
SmartPetscObj< DM > dM
std::tuple< SmartPetscObj< Vec >, SmartPetscObj< VecScatter > > uXScatter
std::tuple< SmartPetscObj< Vec >, SmartPetscObj< VecScatter > > uZScatter
boost::shared_ptr< PostProcEle > postProcFe
std::tuple< SmartPetscObj< Vec >, SmartPetscObj< VecScatter > > uYScatter
boost::shared_ptr< DomainEle > reactionFe

Member Function Documentation

◆ operator()()

MoFEMErrorCode PlasticOps::Monitor::operator() ( )
inline
Examples
PlasticOpsMonitor.hpp.

Definition at line 20 of file PlasticOpsMonitor.hpp.

20{ return 0; }

◆ postProcess()

MoFEMErrorCode PlasticOps::Monitor::postProcess ( )
inline
Examples
PlasticOpsMonitor.hpp.

Definition at line 22 of file PlasticOpsMonitor.hpp.

22 {
24
25 auto make_vtk = [&]() {
27 CHKERR DMoFEMLoopFiniteElements(dM, "dFE", postProcFe, getCacheWeakPtr());
28 CHKERR postProcFe->writeFile(
29 "out_plastic_" + boost::lexical_cast<std::string>(ts_step) + ".h5m");
31 };
32
33 auto calculate_reaction = [&]() {
35 auto r = createDMVector(dM);
36 reactionFe->f = r;
37 CHKERR VecZeroEntries(r);
39 CHKERR VecGhostUpdateBegin(r, ADD_VALUES, SCATTER_REVERSE);
40 CHKERR VecGhostUpdateEnd(r, ADD_VALUES, SCATTER_REVERSE);
41 CHKERR VecAssemblyBegin(r);
42 CHKERR VecAssemblyEnd(r);
43
44 double sum;
45 CHKERR VecSum(r, &sum);
46 MOFEM_LOG_C("EXAMPLE", Sev::inform, "reaction time %3.4e %3.4e", ts_t,
47 sum);
48
50 };
51
52 auto print_max_min = [&](auto &tuple, const std::string msg) {
54 CHKERR VecScatterBegin(std::get<1>(tuple), ts_u, std::get<0>(tuple),
55 INSERT_VALUES, SCATTER_FORWARD);
56 CHKERR VecScatterEnd(std::get<1>(tuple), ts_u, std::get<0>(tuple),
57 INSERT_VALUES, SCATTER_FORWARD);
58 double max, min;
59 CHKERR VecMax(std::get<0>(tuple), PETSC_NULL, &max);
60 CHKERR VecMin(std::get<0>(tuple), PETSC_NULL, &min);
61 MOFEM_LOG_C("EXAMPLE", Sev::inform, "%s time %3.4e min %3.4e max %3.4e",
62 msg.c_str(), ts_t, min, max);
64 };
65
66 CHKERR make_vtk();
67 if (reactionFe)
68 CHKERR calculate_reaction();
69 CHKERR print_max_min(uXScatter, "Ux");
70 CHKERR print_max_min(uYScatter, "Uy");
71 if constexpr (SPACE_DIM == 3)
72 CHKERR print_max_min(uZScatter, "Uz");
73
75 }
#define MOFEM_LOG_C(channel, severity, format,...)
Definition: LogManager.hpp:311
constexpr int SPACE_DIM
#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
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
auto createDMVector(DM dm)
Get smart vector from DM.
Definition: DMMoFEM.hpp:1003
int r
Definition: sdf.py:5

◆ preProcess()

MoFEMErrorCode PlasticOps::Monitor::preProcess ( )
inline
Examples
PlasticOpsMonitor.hpp.

Definition at line 19 of file PlasticOpsMonitor.hpp.

19{ return 0; }

Member Data Documentation

◆ dM

SmartPetscObj<DM> PlasticOps::Monitor::dM
private
Examples
PlasticOpsMonitor.hpp.

Definition at line 78 of file PlasticOpsMonitor.hpp.

◆ postProcFe

boost::shared_ptr<PostProcEle> PlasticOps::Monitor::postProcFe
private
Examples
PlasticOpsMonitor.hpp.

Definition at line 79 of file PlasticOpsMonitor.hpp.

◆ reactionFe

boost::shared_ptr<DomainEle> PlasticOps::Monitor::reactionFe
private
Examples
PlasticOpsMonitor.hpp.

Definition at line 80 of file PlasticOpsMonitor.hpp.

◆ uXScatter

std::tuple<SmartPetscObj<Vec>, SmartPetscObj<VecScatter> > PlasticOps::Monitor::uXScatter
private
Examples
PlasticOpsMonitor.hpp.

Definition at line 81 of file PlasticOpsMonitor.hpp.

◆ uYScatter

std::tuple<SmartPetscObj<Vec>, SmartPetscObj<VecScatter> > PlasticOps::Monitor::uYScatter
private
Examples
PlasticOpsMonitor.hpp.

Definition at line 82 of file PlasticOpsMonitor.hpp.

◆ uZScatter

std::tuple<SmartPetscObj<Vec>, SmartPetscObj<VecScatter> > PlasticOps::Monitor::uZScatter
private
Examples
PlasticOpsMonitor.hpp.

Definition at line 83 of file PlasticOpsMonitor.hpp.


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