v0.15.0
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
9template <int DIM> struct Monitor : public FEMethod {
10
12 SmartPetscObj<DM> dm,
13 std::pair<boost::shared_ptr<PostProcEle>,
14 boost::shared_ptr<SkinPostProcEle>>
15 pair_post_proc_fe,
16 boost::shared_ptr<DomainEle> reaction_fe,
17 std::tuple<SmartPetscObj<Vec>, SmartPetscObj<VecScatter>> ux_scatter,
18 std::tuple<SmartPetscObj<Vec>, SmartPetscObj<VecScatter>> uy_scatter,
19 std::tuple<SmartPetscObj<Vec>, SmartPetscObj<VecScatter>> uz_scatter,
20 std::array<double, 3> pass_field_eval_coords,
21 boost::shared_ptr<SetPtsData> pass_field_eval_data,
22 boost::shared_ptr<std::map<std::string, boost::shared_ptr<VectorDouble>>>
23 scalar_field_ptrs,
24 boost::shared_ptr<std::map<std::string, boost::shared_ptr<MatrixDouble>>>
25 vec_field_ptrs,
26 boost::shared_ptr<std::map<std::string, boost::shared_ptr<MatrixDouble>>>
27 sym_tensor_field_ptrs,
28 boost::shared_ptr<std::map<std::string, boost::shared_ptr<MatrixDouble>>>
29 tensor_field_ptrs)
30 : dM(dm), reactionFe(reaction_fe), uXScatter(ux_scatter),
31 uYScatter(uy_scatter), uZScatter(uz_scatter),
32 fieldEvalCoords(pass_field_eval_coords),
33 fieldEvalData(pass_field_eval_data), scalarFieldPtrs(scalar_field_ptrs),
34 vecFieldPtrs(vec_field_ptrs), symTensorFieldPtrs(sym_tensor_field_ptrs),
35 tensorFieldPtrs(tensor_field_ptrs) {
36 postProcFe = pair_post_proc_fe.first;
37 skinPostProcFe = pair_post_proc_fe.second;
38 };
39
40 MoFEMErrorCode preProcess() { return 0; }
41 MoFEMErrorCode operator()() { return 0; }
42
43private:
44 SmartPetscObj<DM> dM;
45 boost::shared_ptr<PostProcEle> postProcFe;
46 boost::shared_ptr<SkinPostProcEle> skinPostProcFe;
47 boost::shared_ptr<DomainEle> reactionFe;
48 std::tuple<SmartPetscObj<Vec>, SmartPetscObj<VecScatter>> uXScatter;
49 std::tuple<SmartPetscObj<Vec>, SmartPetscObj<VecScatter>> uYScatter;
50 std::tuple<SmartPetscObj<Vec>, SmartPetscObj<VecScatter>> uZScatter;
51 std::array<double, 3> fieldEvalCoords;
52 boost::shared_ptr<SetPtsData> fieldEvalData;
53 boost::shared_ptr<std::map<std::string, boost::shared_ptr<VectorDouble>>>
55 boost::shared_ptr<std::map<std::string, boost::shared_ptr<MatrixDouble>>>
57 boost::shared_ptr<std::map<std::string, boost::shared_ptr<MatrixDouble>>>
59 boost::shared_ptr<std::map<std::string, boost::shared_ptr<MatrixDouble>>>
61
62protected:
63 MoFEMErrorCode postProcess();
64};
65
66template <int DIM> MoFEMErrorCode Monitor<DIM>::postProcess() {
68
69 MoFEM::Interface *m_field_ptr;
70 CHKERR DMoFEMGetInterfacePtr(dM, &m_field_ptr);
71 auto *simple = m_field_ptr->getInterface<Simple>();
72
73 if (do_eval_field) {
74
75 CHKERR m_field_ptr->getInterface<FieldEvaluatorInterface>()
76 ->evalFEAtThePoint<DIM>(
77 fieldEvalCoords.data(), 1e-12, simple->getProblemName(),
78 simple->getDomainFEName(), fieldEvalData,
79 m_field_ptr->get_comm_rank(), m_field_ptr->get_comm_rank(),
80 getCacheWeakPtr().lock(), MF_EXIST, QUIET);
81
82 auto process_scalar_field =
83 [](const std::string label,
84 const boost::shared_ptr<VectorDouble> scalarFieldPtr) {
85 if (scalarFieldPtr->size()) {
86 auto t_scalar_holder = getFTensor0FromVec(*scalarFieldPtr);
87
88 MOFEM_LOG("SYNC", Sev::inform)
89 << "For " << label << " field, " << t_scalar_holder;
90 }
91 };
92
93 auto process_vector_field =
94 [](const std::string label,
95 const boost::shared_ptr<MatrixDouble> vecFieldPtr) {
96 if (vecFieldPtr->size1()) {
97 auto t_vec_holder = getFTensor1FromMat<DIM>(*vecFieldPtr);
98
99 std::string vec_holder_string = "";
100 for (int i = 0; i < DIM; i++) {
101 vec_holder_string +=
102 " " + boost::lexical_cast<std::string>(t_vec_holder(i));
103 }
104
105 MOFEM_LOG("SYNC", Sev::inform)
106 << "For " << label << " field," << vec_holder_string;
107 }
108 };
109
110 auto process_sym_tensor_field =
111 [](const std::string label,
112 const boost::shared_ptr<MatrixDouble> symTensorFieldPtr) {
113 if (symTensorFieldPtr->size1()) {
114 auto t_sym_tensor_holder =
115 getFTensor2SymmetricFromMat<DIM>(*symTensorFieldPtr);
116
117 std::string sym_tensor_holder_string = "";
118 for (int i = 0; i < DIM; i++) {
119 for (int j = i; j < DIM; j++) {
120 sym_tensor_holder_string +=
121 ", entry " + boost::lexical_cast<std::string>(i) +
122 boost::lexical_cast<std::string>(j) + " = ";
123 sym_tensor_holder_string +=
124 boost::lexical_cast<std::string>(t_sym_tensor_holder(i, j));
125 }
126 }
127
128 MOFEM_LOG("SYNC", Sev::inform)
129 << "For " << label << " field" << sym_tensor_holder_string;
130 }
131 };
132
133 auto process_tensor_field =
134 [](const std::string label,
135 const boost::shared_ptr<MatrixDouble> tensorFieldPtr) {
136 if (tensorFieldPtr->size1()) {
137 auto t_tensor_holder =
138 getFTensor2FromMat<DIM, DIM>(*tensorFieldPtr);
139
140 std::string tensor_holder_string = "";
141 for (int i = 0; i < DIM; i++) {
142 for (int j = 0; j < DIM; j++) {
143 tensor_holder_string +=
144 ", entry " + boost::lexical_cast<std::string>(i) +
145 boost::lexical_cast<std::string>(j) + " = ";
146 tensor_holder_string +=
147 boost::lexical_cast<std::string>(t_tensor_holder(i, j));
148 }
149 }
150
151 MOFEM_LOG("SYNC", Sev::inform)
152 << "For " << label << " field" << tensor_holder_string;
153 }
154 };
155
156 auto process_fields = [](auto &fieldPtr, auto process_field) {
157 for (const auto &pair : fieldPtr) {
158 const std::string &label = pair.first;
159 const boost::shared_ptr<MatrixDouble> &ptr = pair.second;
160
161 if (ptr) {
162 process_field(label, ptr);
163 }
164 }
165 };
166
167 if (scalarFieldPtrs) {
168 for (const auto &pair : *scalarFieldPtrs) {
169 const std::string &label = pair.first;
170 const boost::shared_ptr<VectorDouble> &ptr = pair.second;
171
172 if (ptr) {
173 process_scalar_field(label, ptr);
174 }
175 }
176 }
177
178 if (vecFieldPtrs) {
179 for (const auto &pair : *vecFieldPtrs) {
180 const std::string &label = pair.first;
181 const boost::shared_ptr<MatrixDouble> &ptr = pair.second;
182
183 if (ptr) {
184 process_vector_field(label, ptr);
185 }
186 }
187 }
188
189 if (symTensorFieldPtrs) {
190 for (const auto &pair : *symTensorFieldPtrs) {
191 const std::string &label = pair.first;
192 const boost::shared_ptr<MatrixDouble> &ptr = pair.second;
193
194 if (ptr) {
195 process_sym_tensor_field(label, ptr);
196 }
197 }
198 }
199
200 if (tensorFieldPtrs) {
201 for (const auto &pair : *tensorFieldPtrs) {
202 const std::string &label = pair.first;
203 const boost::shared_ptr<MatrixDouble> &ptr = pair.second;
204
205 if (ptr) {
206 process_tensor_field(label, ptr);
207 }
208 }
209 }
210 }
211
212 MOFEM_LOG_SEVERITY_SYNC(m_field_ptr->get_comm(), Sev::inform);
213
214 auto make_vtk = [&]() {
216 if (postProcFe) {
217 CHKERR DMoFEMLoopFiniteElements(dM, "dFE", postProcFe, getCacheWeakPtr());
218 CHKERR postProcFe->writeFile(
219 "out_plastic_" + boost::lexical_cast<std::string>(ts_step) + ".h5m");
220 }
221 if (skinPostProcFe) {
222 CHKERR DMoFEMLoopFiniteElements(dM, "bFE", skinPostProcFe,
223 getCacheWeakPtr());
224 CHKERR skinPostProcFe->writeFile(
225 "out_skin_plastic_" + boost::lexical_cast<std::string>(ts_step) +
226 ".h5m");
227 }
229 };
230
231 auto calculate_reaction = [&]() {
233 auto r = createDMVector(dM);
234 reactionFe->f = r;
235 CHKERR VecZeroEntries(r);
236 CHKERR DMoFEMLoopFiniteElements(dM, "dFE", reactionFe, getCacheWeakPtr());
237
238#ifndef NDEBUG
239 auto post_proc_residual = [&](auto dm, auto f_res, auto out_name) {
241 auto post_proc_fe =
242 boost::make_shared<PostProcBrokenMeshInMoab<DomainEle>>(*m_field_ptr);
243 using OpPPMap = OpPostProcMapInMoab<DIM, DIM>;
244 auto u_vec = boost::make_shared<MatrixDouble>();
245 post_proc_fe->getOpPtrVector().push_back(
246 new OpCalculateVectorFieldValues<DIM>("U", u_vec, f_res));
247 post_proc_fe->getOpPtrVector().push_back(
248
249 new OpPPMap(
250
251 post_proc_fe->getPostProcMesh(), post_proc_fe->getMapGaussPts(),
252
253 {},
254
255 {{"RES", u_vec}},
256
257 {}, {})
258
259 );
260
261 CHKERR DMoFEMLoopFiniteElements(dM, "dFE", post_proc_fe);
262 post_proc_fe->writeFile("res.h5m");
264 };
265
266 CHKERR post_proc_residual(dM, r, "reaction");
267#endif // NDEBUG
268
270 };
271
272 auto print_max_min = [&](auto &tuple, const std::string msg) {
274 CHKERR VecScatterBegin(std::get<1>(tuple), ts_u, std::get<0>(tuple),
275 INSERT_VALUES, SCATTER_FORWARD);
276 CHKERR VecScatterEnd(std::get<1>(tuple), ts_u, std::get<0>(tuple),
277 INSERT_VALUES, SCATTER_FORWARD);
278 double max, min;
279 CHKERR VecMax(std::get<0>(tuple), PETSC_NULLPTR, &max);
280 CHKERR VecMin(std::get<0>(tuple), PETSC_NULLPTR, &min);
281 MOFEM_LOG_C("PLASTICITY", Sev::inform, "%s time %3.4e min %3.4e max %3.4e",
282 msg.c_str(), ts_t, min, max);
284 };
285
286 int se = 1;
287 CHKERR PetscOptionsGetInt(PETSC_NULLPTR, "", "-save_every", &se, PETSC_NULLPTR);
288 if (!(ts_step % se)) {
289 CHKERR make_vtk();
290 }
291 if (reactionFe)
292 CHKERR calculate_reaction();
293 CHKERR print_max_min(uXScatter, "Ux");
294 CHKERR print_max_min(uYScatter, "Uy");
295 if (DIM == 3)
296 CHKERR print_max_min(uZScatter, "Uz");
297
299}
300
301}; // namespace PlasticOps
#define MOFEM_LOG_SEVERITY_SYNC(comm, severity)
Synchronise "SYNC" on curtain severity level.
#define MOFEM_LOG_C(channel, severity, format,...)
void simple(double P1[], double P2[], double P3[], double c[], const int N)
Definition acoustic.cpp:69
@ QUIET
@ MF_EXIST
#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()
#define CHKERR
Inline error check.
#define MOFEM_LOG(channel, severity)
Log.
FTensor::Index< 'i', SPACE_DIM > i
FTensor::Index< 'j', 3 > j
OpPostProcMapInMoab< SPACE_DIM, SPACE_DIM > OpPPMap
PetscBool do_eval_field
Evaluate field.
Definition plastic.cpp:119
virtual MPI_Comm & get_comm() const =0
virtual int get_comm_rank() const =0
Deprecated interface functions.
Post post-proc data at points from hash maps.
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.
boost::shared_ptr< std::map< std::string, boost::shared_ptr< MatrixDouble > > > symTensorFieldPtrs
std::tuple< SmartPetscObj< Vec >, SmartPetscObj< VecScatter > > uXScatter
boost::shared_ptr< SetPtsData > fieldEvalData
boost::shared_ptr< std::map< std::string, boost::shared_ptr< MatrixDouble > > > vecFieldPtrs
MoFEMErrorCode operator()()
boost::shared_ptr< DomainEle > reactionFe
MoFEMErrorCode preProcess()
boost::shared_ptr< std::map< std::string, boost::shared_ptr< VectorDouble > > > scalarFieldPtrs
boost::shared_ptr< SkinPostProcEle > skinPostProcFe
SmartPetscObj< DM > dM
std::tuple< SmartPetscObj< Vec >, SmartPetscObj< VecScatter > > uYScatter
Monitor(SmartPetscObj< DM > dm, std::pair< boost::shared_ptr< PostProcEle >, boost::shared_ptr< SkinPostProcEle > > pair_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, std::array< double, 3 > pass_field_eval_coords, boost::shared_ptr< SetPtsData > pass_field_eval_data, boost::shared_ptr< std::map< std::string, boost::shared_ptr< VectorDouble > > > scalar_field_ptrs, boost::shared_ptr< std::map< std::string, boost::shared_ptr< MatrixDouble > > > vec_field_ptrs, boost::shared_ptr< std::map< std::string, boost::shared_ptr< MatrixDouble > > > sym_tensor_field_ptrs, boost::shared_ptr< std::map< std::string, boost::shared_ptr< MatrixDouble > > > tensor_field_ptrs)
MoFEMErrorCode postProcess()
boost::shared_ptr< std::map< std::string, boost::shared_ptr< MatrixDouble > > > tensorFieldPtrs
std::tuple< SmartPetscObj< Vec >, SmartPetscObj< VecScatter > > uZScatter
std::array< double, 3 > fieldEvalCoords
boost::shared_ptr< PostProcEle > postProcFe