v0.13.2
Loading...
Searching...
No Matches
PlasticOpsMonitor.hpp
Go to the documentation of this file.
1
2
3/** \file PlasticOpsMonitor.hpp
4 * \example PlasticOpsMonitor.hpp
5 */
6
7namespace PlasticOps {
8
9struct Monitor : public FEMethod {
10
11 Monitor(SmartPetscObj<DM> &dm, boost::shared_ptr<PostProcEle> post_proc_fe,
12 boost::shared_ptr<DomainEle> reaction_fe,
13 std::tuple<SmartPetscObj<Vec>, SmartPetscObj<VecScatter>> ux_scatter,
14 std::tuple<SmartPetscObj<Vec>, SmartPetscObj<VecScatter>> uy_scatter,
15 std::tuple<SmartPetscObj<Vec>, SmartPetscObj<VecScatter>> uz_scatter)
16 : dM(dm), postProcFe(post_proc_fe), reactionFe(reaction_fe),
17 uXScatter(ux_scatter), uYScatter(uy_scatter), uZScatter(uz_scatter){};
18
19 MoFEMErrorCode preProcess() { return 0; }
20 MoFEMErrorCode operator()() { return 0; }
21
22 MoFEMErrorCode postProcess() {
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);
38 CHKERR DMoFEMLoopFiniteElements(dM, "dFE", reactionFe);
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 }
76
77private:
78 SmartPetscObj<DM> dM;
79 boost::shared_ptr<PostProcEle> postProcFe;
80 boost::shared_ptr<DomainEle> reactionFe;
81 std::tuple<SmartPetscObj<Vec>, SmartPetscObj<VecScatter>> uXScatter;
82 std::tuple<SmartPetscObj<Vec>, SmartPetscObj<VecScatter>> uYScatter;
83 std::tuple<SmartPetscObj<Vec>, SmartPetscObj<VecScatter>> uZScatter;
84};
85
86}; // namespace PlasticOps
#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
SmartPetscObj< DM > dM
std::tuple< SmartPetscObj< Vec >, SmartPetscObj< VecScatter > > uXScatter
std::tuple< SmartPetscObj< Vec >, SmartPetscObj< VecScatter > > uZScatter
MoFEMErrorCode preProcess()
MoFEMErrorCode postProcess()
boost::shared_ptr< PostProcEle > postProcFe
std::tuple< SmartPetscObj< Vec >, SmartPetscObj< VecScatter > > uYScatter
MoFEMErrorCode operator()()
boost::shared_ptr< DomainEle > reactionFe
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)