v0.15.5
Loading...
Searching...
No Matches
Public Member Functions | Protected Member Functions | Private Attributes | List of all members
PlasticOps::Monitor< DIM > Struct Template Reference

#include "tutorials/adv-0_plasticity/src/PlasticOpsMonitor.hpp"

Inheritance diagram for PlasticOps::Monitor< DIM >:
[legend]
Collaboration diagram for PlasticOps::Monitor< DIM >:
[legend]

Public Member Functions

 Monitor (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)
 
virtual ~Monitor ()=default
 
MoFEMErrorCode preProcess ()
 
MoFEMErrorCode operator() ()
 
 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 preProcess ()
 
MoFEMErrorCode operator() ()
 

Protected Member Functions

MoFEMErrorCode postProcess ()
 
MoFEMErrorCode postProcess ()
 

Private Attributes

DM dM
 
boost::shared_ptr< PostProcElepostProcFe
 
boost::shared_ptr< SkinPostProcEleskinPostProcFe
 
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
 
std::array< double, 3 > fieldEvalCoords
 
boost::shared_ptr< SetPtsDatafieldEvalData
 
boost::shared_ptr< std::map< std::string, boost::shared_ptr< VectorDouble > > > scalarFieldPtrs
 
boost::shared_ptr< std::map< std::string, boost::shared_ptr< MatrixDouble > > > vecFieldPtrs
 
boost::shared_ptr< std::map< std::string, boost::shared_ptr< MatrixDouble > > > symTensorFieldPtrs
 
boost::shared_ptr< std::map< std::string, boost::shared_ptr< MatrixDouble > > > tensorFieldPtrs
 
SmartPetscObj< DM > dM
 

Detailed Description

template<int DIM>
struct PlasticOps::Monitor< DIM >

Definition at line 9 of file PlasticOpsMonitor.hpp.

Constructor & Destructor Documentation

◆ Monitor() [1/2]

template<int DIM>
PlasticOps::Monitor< DIM >::Monitor ( 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 
)
inline

Definition at line 11 of file PlasticOpsMonitor.hpp.

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 };
std::tuple< SmartPetscObj< Vec >, SmartPetscObj< VecScatter > > uYScatter
std::array< double, 3 > fieldEvalCoords
std::tuple< SmartPetscObj< Vec >, SmartPetscObj< VecScatter > > uZScatter
boost::shared_ptr< SetPtsData > fieldEvalData
boost::shared_ptr< PostProcEle > postProcFe
boost::shared_ptr< std::map< std::string, boost::shared_ptr< MatrixDouble > > > symTensorFieldPtrs
boost::shared_ptr< DomainEle > reactionFe
boost::shared_ptr< std::map< std::string, boost::shared_ptr< MatrixDouble > > > tensorFieldPtrs
boost::shared_ptr< SkinPostProcEle > skinPostProcFe
boost::shared_ptr< std::map< std::string, boost::shared_ptr< MatrixDouble > > > vecFieldPtrs
boost::shared_ptr< std::map< std::string, boost::shared_ptr< VectorDouble > > > scalarFieldPtrs
std::tuple< SmartPetscObj< Vec >, SmartPetscObj< VecScatter > > uXScatter

◆ ~Monitor()

template<int DIM>
virtual PlasticOps::Monitor< DIM >::~Monitor ( )
virtualdefault

◆ Monitor() [2/2]

template<int DIM>
PlasticOps::Monitor< DIM >::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 
)
inline

Definition at line 11 of file PlasticOpsMonitor.hpp.

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 };

Member Function Documentation

◆ operator()() [1/2]

template<int DIM>
MoFEMErrorCode PlasticOps::Monitor< DIM >::operator() ( )
inline

◆ operator()() [2/2]

template<int DIM>
MoFEMErrorCode PlasticOps::Monitor< DIM >::operator() ( )
inline

Definition at line 41 of file PlasticOpsMonitor.hpp.

41{ return 0; }

◆ postProcess() [1/2]

template<int DIM>
MoFEMErrorCode Monitor::postProcess ( )
protected
Examples
PlasticOpsMonitor.hpp, and mofem/tutorials/adv-0_plasticity/src/PlasticOpsMonitor.hpp.

Definition at line 68 of file PlasticOpsMonitor.hpp.

68 {
70
71 MoFEM::Interface *m_field_ptr;
72 CHKERR DMoFEMGetInterfacePtr(dM, &m_field_ptr);
73 auto *simple = m_field_ptr->getInterface<Simple>();
74
75 if (do_eval_field) {
76
77 CHKERR m_field_ptr->getInterface<FieldEvaluatorInterface>()
78 ->evalFEAtThePoint<DIM>(
79 fieldEvalCoords.data(), 1e-12, simple->getProblemName(),
80 simple->getDomainFEName(), fieldEvalData,
81 m_field_ptr->get_comm_rank(), m_field_ptr->get_comm_rank(),
82 getCacheWeakPtr().lock(), MF_EXIST, QUIET);
83
84 auto process_scalar_field =
85 [](const std::string label,
86 const boost::shared_ptr<VectorDouble> scalarFieldPtr) {
87 if (scalarFieldPtr->size()) {
88 auto t_scalar_holder = getFTensor0FromVec(*scalarFieldPtr);
89
90 MOFEM_LOG("SYNC", Sev::inform)
91 << "For " << label << " field, " << t_scalar_holder;
92 }
93 };
94
95 auto process_vector_field =
96 [](const std::string label,
97 const boost::shared_ptr<MatrixDouble> vecFieldPtr) {
98 if (vecFieldPtr->size1()) {
99 auto t_vec_holder = getFTensor1FromMat<DIM>(*vecFieldPtr);
100
101 std::string vec_holder_string = "";
102 for (int i = 0; i < DIM; i++) {
103 vec_holder_string +=
104 " " + boost::lexical_cast<std::string>(t_vec_holder(i));
105 }
106
107 MOFEM_LOG("SYNC", Sev::inform)
108 << "For " << label << " field," << vec_holder_string;
109 }
110 };
111
112 auto process_sym_tensor_field =
113 [](const std::string label,
114 const boost::shared_ptr<MatrixDouble> symTensorFieldPtr) {
115 if (symTensorFieldPtr->size1()) {
116 auto t_sym_tensor_holder =
117 getFTensor2SymmetricFromMat<DIM>(*symTensorFieldPtr);
118
119 std::string sym_tensor_holder_string = "";
120 for (int i = 0; i < DIM; i++) {
121 for (int j = i; j < DIM; j++) {
122 sym_tensor_holder_string +=
123 ", entry " + boost::lexical_cast<std::string>(i) +
124 boost::lexical_cast<std::string>(j) + " = ";
125 sym_tensor_holder_string +=
126 boost::lexical_cast<std::string>(t_sym_tensor_holder(i, j));
127 }
128 }
129
130 MOFEM_LOG("SYNC", Sev::inform)
131 << "For " << label << " field" << sym_tensor_holder_string;
132 }
133 };
134
135 auto process_tensor_field =
136 [](const std::string label,
137 const boost::shared_ptr<MatrixDouble> tensorFieldPtr) {
138 if (tensorFieldPtr->size1()) {
139 auto t_tensor_holder =
140 getFTensor2FromMat<DIM, DIM>(*tensorFieldPtr);
141
142 std::string tensor_holder_string = "";
143 for (int i = 0; i < DIM; i++) {
144 for (int j = 0; j < DIM; j++) {
145 tensor_holder_string +=
146 ", entry " + boost::lexical_cast<std::string>(i) +
147 boost::lexical_cast<std::string>(j) + " = ";
148 tensor_holder_string +=
149 boost::lexical_cast<std::string>(t_tensor_holder(i, j));
150 }
151 }
152
153 MOFEM_LOG("SYNC", Sev::inform)
154 << "For " << label << " field" << tensor_holder_string;
155 }
156 };
157
158 auto process_fields = [](auto &fieldPtr, auto process_field) {
159 for (const auto &pair : fieldPtr) {
160 const std::string &label = pair.first;
161 const boost::shared_ptr<MatrixDouble> &ptr = pair.second;
162
163 if (ptr) {
164 process_field(label, ptr);
165 }
166 }
167 };
168
169 if (scalarFieldPtrs) {
170 for (const auto &pair : *scalarFieldPtrs) {
171 const std::string &label = pair.first;
172 const boost::shared_ptr<VectorDouble> &ptr = pair.second;
173
174 if (ptr) {
175 process_scalar_field(label, ptr);
176 }
177 }
178 }
179
180 if (vecFieldPtrs) {
181 for (const auto &pair : *vecFieldPtrs) {
182 const std::string &label = pair.first;
183 const boost::shared_ptr<MatrixDouble> &ptr = pair.second;
184
185 if (ptr) {
186 process_vector_field(label, ptr);
187 }
188 }
189 }
190
191 if (symTensorFieldPtrs) {
192 for (const auto &pair : *symTensorFieldPtrs) {
193 const std::string &label = pair.first;
194 const boost::shared_ptr<MatrixDouble> &ptr = pair.second;
195
196 if (ptr) {
197 process_sym_tensor_field(label, ptr);
198 }
199 }
200 }
201
202 if (tensorFieldPtrs) {
203 for (const auto &pair : *tensorFieldPtrs) {
204 const std::string &label = pair.first;
205 const boost::shared_ptr<MatrixDouble> &ptr = pair.second;
206
207 if (ptr) {
208 process_tensor_field(label, ptr);
209 }
210 }
211 }
212 }
213
214 MOFEM_LOG_SEVERITY_SYNC(m_field_ptr->get_comm(), Sev::inform);
215
216 auto make_vtk = [&]() {
218 if (postProcFe) {
219 CHKERR DMoFEMLoopFiniteElements(dM, "dFE", postProcFe, getCacheWeakPtr());
220 CHKERR postProcFe->writeFile(
221 "out_plastic_" + boost::lexical_cast<std::string>(ts_step) + ".h5m");
222 }
223 if (skinPostProcFe) {
225 getCacheWeakPtr());
226 CHKERR skinPostProcFe->writeFile(
227 "out_skin_plastic_" + boost::lexical_cast<std::string>(ts_step) +
228 ".h5m");
229 }
231 };
232
233 auto calculate_reaction = [&]() {
235 auto r = createDMVector(dM);
236 reactionFe->f = r;
237 CHKERR VecZeroEntries(r);
238 CHKERR DMoFEMLoopFiniteElements(dM, "dFE", reactionFe, getCacheWeakPtr());
239
240#ifndef NDEBUG
241 auto post_proc_residual = [&](auto dm, auto f_res, auto out_name) {
243 auto post_proc_fe =
244 boost::make_shared<PostProcBrokenMeshInMoab<DomainEle>>(*m_field_ptr);
245 using OpPPMap = OpPostProcMapInMoab<DIM, DIM>;
246 auto u_vec = boost::make_shared<MatrixDouble>();
247 post_proc_fe->getOpPtrVector().push_back(
248 new OpCalculateVectorFieldValues<DIM>("U", u_vec, f_res));
249 post_proc_fe->getOpPtrVector().push_back(
250
251 new OpPPMap(
252
253 post_proc_fe->getPostProcMesh(), post_proc_fe->getMapGaussPts(),
254
255 {},
256
257 {{"RES", u_vec}},
258
259 {}, {})
260
261 );
262
263 CHKERR DMoFEMLoopFiniteElements(dM, "dFE", post_proc_fe);
264 post_proc_fe->writeFile("res.h5m");
266 };
267
268 CHKERR post_proc_residual(dM, r, "reaction");
269#endif // NDEBUG
270
272 };
273
274 auto print_max_min = [&](auto &tuple, const std::string msg) {
276 CHKERR VecScatterBegin(std::get<1>(tuple), ts_u, std::get<0>(tuple),
277 INSERT_VALUES, SCATTER_FORWARD);
278 CHKERR VecScatterEnd(std::get<1>(tuple), ts_u, std::get<0>(tuple),
279 INSERT_VALUES, SCATTER_FORWARD);
280 double max, min;
281 CHKERR VecMax(std::get<0>(tuple), PETSC_NULLPTR, &max);
282 CHKERR VecMin(std::get<0>(tuple), PETSC_NULLPTR, &min);
283 MOFEM_LOG_C("PLASTICITY", Sev::inform, "%s time %3.4e min %3.4e max %3.4e",
284 msg.c_str(), ts_t, min, max);
286 };
287
288 int se = 1;
289 CHKERR PetscOptionsGetInt(PETSC_NULLPTR, "", "-save_every", &se, PETSC_NULLPTR);
290 if (!(ts_step % se)) {
291 CHKERR make_vtk();
292 }
293 if (reactionFe)
294 CHKERR calculate_reaction();
295 CHKERR print_max_min(uXScatter, "Ux");
296 CHKERR print_max_min(uYScatter, "Uy");
297 if (DIM == 3)
298 CHKERR print_max_min(uZScatter, "Uz");
299
301}
#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.
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:576
auto createDMVector(DM dm, RowColData rc=RowColData::COL)
Get smart vector from DM.
Definition DMMoFEM.hpp:1237
PetscErrorCode DMoFEMGetInterfacePtr(DM dm, MoFEM::Interface **m_field_ptr)
Get pointer to MoFEM::Interface.
Definition DMMoFEM.cpp:410
#define MOFEM_LOG(channel, severity)
Log.
FTensor::Index< 'i', SPACE_DIM > i
FTensor::Index< 'j', 3 > j
PetscErrorCode PetscOptionsGetInt(PetscOptions *, const char pre[], const char name[], PetscInt *ivalue, PetscBool *set)
int r
Definition sdf.py:205
OpPostProcMapInMoab< SPACE_DIM, SPACE_DIM > OpPPMap
virtual MPI_Comm & get_comm() const =0
virtual int get_comm_rank() const =0
Deprecated interface functions.
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.
PetscBool do_eval_field
Evaluate field.
Definition plastic.cpp:119

◆ postProcess() [2/2]

template<int DIM>
MoFEMErrorCode PlasticOps::Monitor< DIM >::postProcess ( )
protected

◆ preProcess() [1/2]

template<int DIM>
MoFEMErrorCode PlasticOps::Monitor< DIM >::preProcess ( )
inline

◆ preProcess() [2/2]

template<int DIM>
MoFEMErrorCode PlasticOps::Monitor< DIM >::preProcess ( )
inline

Definition at line 40 of file PlasticOpsMonitor.hpp.

40{ return 0; }

Member Data Documentation

◆ dM [1/2]

template<int DIM>
DM PlasticOps::Monitor< DIM >::dM
private

◆ dM [2/2]

template<int DIM>
SmartPetscObj<DM> PlasticOps::Monitor< DIM >::dM
private

Definition at line 44 of file PlasticOpsMonitor.hpp.

◆ fieldEvalCoords

template<int DIM>
std::array< double, 3 > Monitor::fieldEvalCoords
private

◆ fieldEvalData

template<int DIM>
boost::shared_ptr< SetPtsData > Monitor::fieldEvalData
private

◆ postProcFe

template<int DIM>
boost::shared_ptr< PostProcEle > Monitor::postProcFe
private

◆ reactionFe

template<int DIM>
boost::shared_ptr< DomainEle > Monitor::reactionFe
private

◆ scalarFieldPtrs

template<int DIM>
boost::shared_ptr< std::map< std::string, boost::shared_ptr< VectorDouble > > > Monitor::scalarFieldPtrs
private

◆ skinPostProcFe

template<int DIM>
boost::shared_ptr< SkinPostProcEle > Monitor::skinPostProcFe
private

◆ symTensorFieldPtrs

template<int DIM>
boost::shared_ptr< std::map< std::string, boost::shared_ptr< MatrixDouble > > > Monitor::symTensorFieldPtrs
private

◆ tensorFieldPtrs

template<int DIM>
boost::shared_ptr< std::map< std::string, boost::shared_ptr< MatrixDouble > > > Monitor::tensorFieldPtrs
private

◆ uXScatter

template<int DIM>
std::tuple< SmartPetscObj< Vec >, SmartPetscObj< VecScatter > > Monitor::uXScatter
private

◆ uYScatter

template<int DIM>
std::tuple< SmartPetscObj< Vec >, SmartPetscObj< VecScatter > > Monitor::uYScatter
private

◆ uZScatter

template<int DIM>
std::tuple< SmartPetscObj< Vec >, SmartPetscObj< VecScatter > > Monitor::uZScatter
private

◆ vecFieldPtrs

template<int DIM>
boost::shared_ptr< std::map< std::string, boost::shared_ptr< MatrixDouble > > > Monitor::vecFieldPtrs
private

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