v0.14.0
PlasticOpsMonitor.hpp
Go to the documentation of this file.
1 
2 
3 /** \file PlasticOpsMonitor.hpp
4  * \example PlasticOpsMonitor.hpp
5  */
6 
7 namespace PlasticOps {
8 
9 template <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, DIM> 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 
43 private:
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, DIM> 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 
62 protected:
64 };
65 
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_NULL, &max);
280  CHKERR VecMin(std::get<0>(tuple), PETSC_NULL, &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_NULL, "", "-save_every", &se, PETSC_NULL);
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
MoFEM::UnknownInterface::getInterface
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.
Definition: UnknownInterface.hpp:93
PlasticOps::Monitor::fieldEvalData
boost::shared_ptr< SetPtsData > fieldEvalData
Definition: PlasticOpsMonitor.hpp:52
MOFEM_LOG_SEVERITY_SYNC
#define MOFEM_LOG_SEVERITY_SYNC(comm, severity)
Synchronise "SYNC" on curtain severity level.
Definition: LogManager.hpp:352
PlasticOps::Monitor::operator()
MoFEMErrorCode operator()()
Definition: PlasticOpsMonitor.hpp:41
PlasticOps::Monitor::preProcess
MoFEMErrorCode preProcess()
Definition: PlasticOpsMonitor.hpp:40
PlasticOps::Monitor::uYScatter
std::tuple< SmartPetscObj< Vec >, SmartPetscObj< VecScatter > > uYScatter
Definition: PlasticOpsMonitor.hpp:49
PlasticOps::Monitor::scalarFieldPtrs
boost::shared_ptr< std::map< std::string, boost::shared_ptr< VectorDouble > > > scalarFieldPtrs
Definition: PlasticOpsMonitor.hpp:54
MoFEM::CoreInterface::get_comm
virtual MPI_Comm & get_comm() const =0
do_eval_field
PetscBool do_eval_field
Evaluate field.
Definition: plastic.cpp:119
MoFEM::Exceptions::MoFEMErrorCode
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
Definition: Exceptions.hpp:56
PlasticOps::Monitor::symTensorFieldPtrs
boost::shared_ptr< std::map< std::string, boost::shared_ptr< MatrixDouble > > > symTensorFieldPtrs
Definition: PlasticOpsMonitor.hpp:58
MoFEM::CoreInterface::get_comm_rank
virtual int get_comm_rank() const =0
PlasticOps::Monitor::uXScatter
std::tuple< SmartPetscObj< Vec >, SmartPetscObj< VecScatter > > uXScatter
Definition: PlasticOpsMonitor.hpp:48
PlasticOps::Monitor::postProcFe
boost::shared_ptr< PostProcEle > postProcFe
Definition: PlasticOpsMonitor.hpp:45
sdf.r
int r
Definition: sdf.py:8
MoFEM::DeprecatedCoreInterface
Deprecated interface functions.
Definition: DeprecatedCoreInterface.hpp:16
PlasticOps::Monitor::skinPostProcFe
boost::shared_ptr< SkinPostProcEle > skinPostProcFe
Definition: PlasticOpsMonitor.hpp:46
FEMethod
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:548
MoFEM::createDMVector
auto createDMVector(DM dm)
Get smart vector from DM.
Definition: DMMoFEM.hpp:1099
PlasticOps::Monitor::tensorFieldPtrs
boost::shared_ptr< std::map< std::string, boost::shared_ptr< MatrixDouble > > > tensorFieldPtrs
Definition: PlasticOpsMonitor.hpp:60
MOFEM_LOG_C
#define MOFEM_LOG_C(channel, severity, format,...)
Definition: LogManager.hpp:311
simple
void simple(double P1[], double P2[], double P3[], double c[], const int N)
Definition: acoustic.cpp:69
PlasticOps::Monitor::fieldEvalCoords
std::array< double, DIM > fieldEvalCoords
Definition: PlasticOpsMonitor.hpp:51
OpPPMap
OpPostProcMapInMoab< SPACE_DIM, SPACE_DIM > OpPPMap
Definition: photon_diffusion.cpp:29
MoFEM::getFTensor0FromVec
static auto getFTensor0FromVec(ublas::vector< T, A > &data)
Get tensor rank 0 (scalar) form data vector.
Definition: Templates.hpp:135
PlasticOps::Monitor::vecFieldPtrs
boost::shared_ptr< std::map< std::string, boost::shared_ptr< MatrixDouble > > > vecFieldPtrs
Definition: PlasticOpsMonitor.hpp:56
PlasticOps::Monitor
Definition: PlasticOpsMonitor.hpp:9
i
FTensor::Index< 'i', SPACE_DIM > i
Definition: hcurl_divergence_operator_2d.cpp:27
PlasticOps::Monitor::postProcess
MoFEMErrorCode postProcess()
Definition: PlasticOpsMonitor.hpp:66
PlasticOps
Definition: plastic.cpp:182
PlasticOps::Monitor::uZScatter
std::tuple< SmartPetscObj< Vec >, SmartPetscObj< VecScatter > > uZScatter
Definition: PlasticOpsMonitor.hpp:50
MOFEM_LOG
#define MOFEM_LOG(channel, severity)
Log.
Definition: LogManager.hpp:308
MoFEM::DMoFEMGetInterfacePtr
PetscErrorCode DMoFEMGetInterfacePtr(DM dm, MoFEM::Interface **m_field_ptr)
Get pointer to MoFEM::Interface.
Definition: DMMoFEM.cpp:414
j
FTensor::Index< 'j', 3 > j
Definition: matrix_function.cpp:19
PlasticOps::Monitor::Monitor
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, DIM > 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)
Definition: PlasticOpsMonitor.hpp:11
QUIET
@ QUIET
Definition: definitions.h:221
PlasticOps::Monitor::dM
SmartPetscObj< DM > dM
Definition: PlasticOpsMonitor.hpp:44
MoFEM::DMoFEMLoopFiniteElements
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:586
MF_EXIST
@ MF_EXIST
Definition: definitions.h:113
MoFEM::PetscOptionsGetInt
PetscErrorCode PetscOptionsGetInt(PetscOptions *, const char pre[], const char name[], PetscInt *ivalue, PetscBool *set)
Definition: DeprecatedPetsc.hpp:142
MoFEMFunctionReturn
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:429
PlasticOps::Monitor::reactionFe
boost::shared_ptr< DomainEle > reactionFe
Definition: PlasticOpsMonitor.hpp:47
MoFEMFunctionBegin
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:359
MoFEM::OpPostProcMapInMoab
Post post-proc data at points from hash maps.
Definition: PostProcBrokenMeshInMoabBase.hpp:709