v0.14.0
photon_diffusion.cpp
Go to the documentation of this file.
1 /**
2  * \file photon_diffusion.cpp
3  * \example photon_diffusion.cpp
4  *
5  **/
6 
7 #include <stdlib.h>
8 #include <cmath>
10 
11 using namespace MoFEM;
12 
13 static char help[] = "...\n\n";
14 
15 template <int DIM> struct ElementsAndOps {};
16 
17 //! [Define dimension]
18 constexpr int SPACE_DIM = 3; //< Space dimension of problem, mesh
19 //! [Define dimension]
20 
26 using PostProcFaceEle =
28 
30 
32 
38  PETSC>::LinearForm<GAUSS>::OpBaseTimesScalar<1>;
42  PETSC>::LinearForm<GAUSS>::OpSource<1, 1>;
43 
47  PETSC>::LinearForm<GAUSS>::OpBaseTimesScalar<1>;
49  PETSC>::LinearForm<GAUSS>::OpSource<1, 1>;
50 
51 const double n = 1.44; ///< refractive index of diffusive medium
52 const double c = 30.; ///< speed of light (cm/ns)
53 const double v = c / n; ///< phase velocity of light in medium (cm/ns)
54 const double inv_v = 1. / v;
55 
56 double mu_a; ///< absorption coefficient (cm^-1)
57 double mu_sp; ///< scattering coefficient (cm^-1)
58 double D;
59 double A;
60 double h;
61 
62 PetscBool from_initial = PETSC_TRUE;
63 PetscBool output_volume = PETSC_FALSE;
64 PetscBool output_camera = PETSC_FALSE;
65 
66 int order = 2;
68 
69 char init_data_file_name[255] = "init_file.dat";
70 int numHoLevels = 1;
71 
72 
73 struct PhotonDiffusion {
74 public:
76 
77  // Declaration of the main function to run analysis
78  MoFEMErrorCode runProgram();
79 
80 private:
81  // Declaration of other main functions called in runProgram()
82  MoFEMErrorCode readMesh();
83  MoFEMErrorCode createCommonData();
84  MoFEMErrorCode setupProblem();
85  MoFEMErrorCode setIntegrationRules();
86  MoFEMErrorCode initialCondition();
87  MoFEMErrorCode boundaryCondition();
88  MoFEMErrorCode assembleSystem();
89  MoFEMErrorCode solveSystem();
90  MoFEMErrorCode outputResults();
91 
92  // Main interfaces
93  MoFEM::Interface &mField;
94 
95  // Object to mark boundary entities for the assembling of domain elements
96  boost::shared_ptr<std::vector<unsigned char>> boundaryMarker;
97 
98  boost::shared_ptr<FEMethod> domainLhsFEPtr;
99  boost::shared_ptr<FEMethod> boundaryLhsFEPtr;
100  boost::shared_ptr<FEMethod> boundaryRhsFEPtr;
101 
102  struct CommonData {
103  boost::shared_ptr<VectorDouble> approxVals;
105 
106  enum VecElements {
107  VALUES_INTEG = 0,
108  LAST_ELEMENT
109  };
110  };
111 
112  boost::shared_ptr<CommonData> commonDataPtr;
113 
114  struct OpCameraInteg : public BoundaryEleOp {
115  boost::shared_ptr<CommonData> commonDataPtr;
116  OpCameraInteg(boost::shared_ptr<CommonData> common_data_ptr)
117  : BoundaryEleOp("PHOTON_FLUENCE_RATE", OPROW),
118  commonDataPtr(common_data_ptr) {
119  std::fill(&doEntities[MBVERTEX], &doEntities[MBMAXTYPE], false);
120  doEntities[MBTRI] = doEntities[MBQUAD] = true;
121  }
122  MoFEMErrorCode doWork(int side, EntityType type, EntData &data);
123  };
124 
126 
127  boost::shared_ptr<VolSideFe> sideOpFe;
128 
129  OpGetScalarFieldGradientValuesOnSkin(boost::shared_ptr<VolSideFe> side_fe)
130  : BoundaryEleOp("PHOTON_FLUENCE_RATE", OPROW), sideOpFe(side_fe) {}
131 
132  MoFEMErrorCode doWork(int side, EntityType type,
135  if (type != MBVERTEX)
137  CHKERR loopSideVolumes("dFE", *sideOpFe);
139  }
140  };
141 
142  struct Monitor : public FEMethod {
143 
144  Monitor(SmartPetscObj<DM> dm, boost::shared_ptr<PostProcEle> post_proc,
145  boost::shared_ptr<PostProcFaceEle> skin_post_proc,
146  boost::shared_ptr<BoundaryEle> skin_post_proc_integ,
147  boost::shared_ptr<CommonData> common_data_ptr)
148  : dM(dm), postProc(post_proc), skinPostProc(skin_post_proc),
149  skinPostProcInteg(skin_post_proc_integ),
150  commonDataPtr(common_data_ptr){};
151 
155  }
159  }
160 
163 
164  CHKERR VecZeroEntries(commonDataPtr->petscVec);
165  CHKERR VecGhostUpdateBegin(commonDataPtr->petscVec, INSERT_VALUES,
166  SCATTER_FORWARD);
167  CHKERR VecGhostUpdateEnd(commonDataPtr->petscVec, INSERT_VALUES,
168  SCATTER_FORWARD);
169  CHKERR DMoFEMLoopFiniteElements(dM, "CAMERA_FE", skinPostProcInteg);
170  CHKERR VecAssemblyBegin(commonDataPtr->petscVec);
171  CHKERR VecAssemblyEnd(commonDataPtr->petscVec);
172  CHKERR VecGhostUpdateBegin(commonDataPtr->petscVec, ADD_VALUES,
173  SCATTER_REVERSE);
174  CHKERR VecGhostUpdateEnd(commonDataPtr->petscVec, ADD_VALUES,
175  SCATTER_REVERSE);
176  CHKERR VecGhostUpdateBegin(commonDataPtr->petscVec, INSERT_VALUES,
177  SCATTER_FORWARD);
178  CHKERR VecGhostUpdateEnd(commonDataPtr->petscVec, INSERT_VALUES,
179  SCATTER_FORWARD);
180  const double *array;
181  CHKERR VecGetArrayRead(commonDataPtr->petscVec, &array);
182  MOFEM_LOG("PHOTON", Sev::inform) << "Fluence rate integral: " << array[0];
183 
184  if (ts_step % save_every_nth_step == 0) {
185  if (output_volume) {
186  CHKERR DMoFEMLoopFiniteElements(dM, "dFE", postProc);
187  CHKERR postProc->writeFile("out_volume_" +
188  boost::lexical_cast<std::string>(ts_step) +
189  ".h5m");
190  }
191  if (output_camera && skinPostProc) {
192  CHKERR DMoFEMLoopFiniteElements(dM, "CAMERA_FE", skinPostProc);
193  CHKERR skinPostProc->writeFile(
194  "out_camera_" + boost::lexical_cast<std::string>(ts_step) +
195  ".h5m");
196  }
197  }
199  }
200 
201  private:
203  boost::shared_ptr<PostProcEle> postProc;
204  boost::shared_ptr<PostProcFaceEle> skinPostProc;
205  boost::shared_ptr<BoundaryEle> skinPostProcInteg;
206  boost::shared_ptr<CommonData> commonDataPtr;
207  };
208 };
209 
210 PhotonDiffusion::PhotonDiffusion(MoFEM::Interface &m_field) : mField(m_field) {}
211 
214 
215  auto *simple = mField.getInterface<Simple>();
217  CHKERR simple->getOptions();
218  CHKERR simple->loadFile();
219 
221 }
222 
225  commonDataPtr = boost::make_shared<CommonData>();
226  PetscInt ghosts[1] = {0};
227  if (!mField.get_comm_rank())
228  commonDataPtr->petscVec =
229  createGhostVector(mField.get_comm(), 1, 1, 0, ghosts);
230  else
231  commonDataPtr->petscVec =
232  createGhostVector(mField.get_comm(), 0, 1, 1, ghosts);
233  commonDataPtr->approxVals = boost::make_shared<VectorDouble>();
235 }
236 
239 
240  auto *simple = mField.getInterface<Simple>();
241  CHKERR simple->addDomainField("PHOTON_FLUENCE_RATE", H1,
243  CHKERR simple->addBoundaryField("PHOTON_FLUENCE_RATE", H1,
245 
246  CHKERR PetscOptionsGetString(PETSC_NULL, "", "-initial_file",
247  init_data_file_name, 255, PETSC_NULL);
248 
249  CHKERR PetscOptionsGetBool(PETSC_NULL, "", "-from_initial", &from_initial,
250  PETSC_NULL);
251  CHKERR PetscOptionsGetBool(PETSC_NULL, "", "-output_volume", &output_volume,
252  PETSC_NULL);
253  CHKERR PetscOptionsGetBool(PETSC_NULL, "", "-output_camera", &output_camera,
254  PETSC_NULL);
255  CHKERR PetscOptionsGetScalar(PETSC_NULL, "", "-mu_a", &mu_a, PETSC_NULL);
256  CHKERR PetscOptionsGetScalar(PETSC_NULL, "", "-mu_sp", &mu_sp, PETSC_NULL);
257  CHKERR PetscOptionsGetScalar(PETSC_NULL, "", "-coef_A", &A, PETSC_NULL);
258 
259  CHKERR PetscOptionsGetInt(PETSC_NULL, "", "-order", &order, PETSC_NULL);
260  CHKERR PetscOptionsGetInt(PETSC_NULL, "", "-save_step", &save_every_nth_step,
261  PETSC_NULL);
262 
263  h = 0.5 / A;
264  D = 1. / (3. * (mu_a + mu_sp));
265 
266  MOFEM_LOG("PHOTON", Sev::inform) << "Refractive index: " << n;
267  MOFEM_LOG("PHOTON", Sev::inform) << "Speed of light (cm/ns): " << c;
268  MOFEM_LOG("PHOTON", Sev::inform) << "Phase velocity in medium (cm/ns): " << v;
269  MOFEM_LOG("PHOTON", Sev::inform) << "Inverse velocity : " << inv_v;
270  MOFEM_LOG("PHOTON", Sev::inform)
271  << "Absorption coefficient (cm^-1): " << mu_a;
272  MOFEM_LOG("PHOTON", Sev::inform)
273  << "Scattering coefficient (cm^-1): " << mu_sp;
274  MOFEM_LOG("PHOTON", Sev::inform) << "Diffusion coefficient D : " << D;
275  MOFEM_LOG("PHOTON", Sev::inform) << "Coefficient A : " << A;
276  MOFEM_LOG("PHOTON", Sev::inform) << "Coefficient h : " << h;
277 
278  MOFEM_LOG("PHOTON", Sev::inform) << "Approximation order: " << order;
279  MOFEM_LOG("PHOTON", Sev::inform) << "Save step: " << save_every_nth_step;
280 
281  CHKERR simple->setFieldOrder("PHOTON_FLUENCE_RATE", order);
282 
283  auto set_camera_skin_fe = [&]() {
285  auto meshset_mng = mField.getInterface<MeshsetsManager>();
286 
287  Range camera_surface;
288  const std::string block_name = "CAM";
289  bool add_fe = false;
290 
292  if (bit->getName().compare(0, block_name.size(), block_name) == 0) {
293  MOFEM_LOG("PHOTON", Sev::inform) << "Found CAM block";
294  CHKERR mField.get_moab().get_entities_by_dimension(
295  bit->getMeshset(), 2, camera_surface, true);
296  add_fe = true;
297  }
298  }
299 
300  MOFEM_LOG("PHOTON", Sev::noisy) << "CAM block entities:\n"
301  << camera_surface;
302 
303  if (add_fe) {
304  CHKERR mField.add_finite_element("CAMERA_FE");
306  "PHOTON_FLUENCE_RATE");
308  "CAMERA_FE");
309  }
311  };
312 
313  auto my_simple_set_up = [&]() {
315  CHKERR simple->defineFiniteElements();
316  CHKERR simple->defineProblem(PETSC_TRUE);
317  CHKERR simple->buildFields();
318  CHKERR simple->buildFiniteElements();
319 
320  if (mField.check_finite_element("CAMERA_FE")) {
321  CHKERR mField.build_finite_elements("CAMERA_FE");
322  CHKERR DMMoFEMAddElement(simple->getDM(), "CAMERA_FE");
323  }
324 
325  CHKERR simple->buildProblem();
327  };
328 
329  CHKERR set_camera_skin_fe();
330  CHKERR my_simple_set_up();
331 
333 }
334 
337 
338  auto integration_rule = [](int o_row, int o_col, int approx_order) {
339  return 2 * approx_order;
340  };
341 
342  auto *pipeline_mng = mField.getInterface<PipelineManager>();
343  CHKERR pipeline_mng->setDomainRhsIntegrationRule(integration_rule);
344  CHKERR pipeline_mng->setDomainLhsIntegrationRule(integration_rule);
345  CHKERR pipeline_mng->setBoundaryLhsIntegrationRule(integration_rule);
346  CHKERR pipeline_mng->setBoundaryRhsIntegrationRule(integration_rule);
347 
349 }
350 
353 
355 }
356 
359  auto bc_mng = mField.getInterface<BcManager>();
360  auto *simple = mField.getInterface<Simple>();
361  CHKERR bc_mng->pushMarkDOFsOnEntities(simple->getProblemName(), "EXT",
362  "PHOTON_FLUENCE_RATE", 0, 0, false);
363 
364  // Get boundary edges marked in block named "BOUNDARY_CONDITION"
365  Range boundary_ents;
367  std::string entity_name = it->getName();
368  if (entity_name.compare(0, 3, "INT") == 0) {
369  CHKERR it->getMeshsetIdEntitiesByDimension(mField.get_moab(), 1,
370  boundary_ents, true);
371  }
372  }
373  // Add vertices to boundary entities
374  Range boundary_verts;
375  CHKERR mField.get_moab().get_connectivity(boundary_ents, boundary_verts,
376  true);
377  boundary_ents.merge(boundary_verts);
378 
379  // Remove DOFs as homogeneous boundary condition is used
380  CHKERR mField.getInterface<ProblemsManager>()->removeDofsOnEntities(
381  simple->getProblemName(), "PHOTON_FLUENCE_RATE", boundary_ents);
382 
384 }
385 
388 
389  auto bc_mng = mField.getInterface<BcManager>();
390  auto &bc_map = bc_mng->getBcMapByBlockName();
391 
392  auto add_domain_base_ops = [&](auto &pipeline) {
393  auto jac_ptr = boost::make_shared<MatrixDouble>();
394  auto inv_jac_ptr = boost::make_shared<MatrixDouble>();
395  auto det_ptr = boost::make_shared<VectorDouble>();
396  pipeline.push_back(new OpCalculateHOJac<3>(jac_ptr));
397  pipeline.push_back(new OpInvertMatrix<3>(jac_ptr, det_ptr, inv_jac_ptr));
398  pipeline.push_back(new OpSetHOInvJacToScalarBases<3>(H1, inv_jac_ptr));
399  pipeline.push_back(new OpSetHOWeights(det_ptr));
400  };
401 
402  auto add_domain_lhs_ops = [&](auto &pipeline) {
403  pipeline.push_back(new OpDomainGradGrad(
404  "PHOTON_FLUENCE_RATE", "PHOTON_FLUENCE_RATE",
405  [](double, double, double) -> double { return D; }));
406  auto get_mass_coefficient = [&](const double, const double, const double) {
407  return inv_v * domainLhsFEPtr->ts_a + mu_a;
408  };
409  pipeline.push_back(new OpDomainMass(
410  "PHOTON_FLUENCE_RATE", "PHOTON_FLUENCE_RATE", get_mass_coefficient));
411  };
412 
413  auto add_domain_rhs_ops = [&](auto &pipeline) {
414  auto grad_u_at_gauss_pts = boost::make_shared<MatrixDouble>();
415  auto u_at_gauss_pts = boost::make_shared<VectorDouble>();
416  auto dot_u_at_gauss_pts = boost::make_shared<VectorDouble>();
417  pipeline.push_back(new OpCalculateScalarFieldGradient<SPACE_DIM>(
418  "PHOTON_FLUENCE_RATE", grad_u_at_gauss_pts));
419  pipeline.push_back(new OpCalculateScalarFieldValues("PHOTON_FLUENCE_RATE",
420  u_at_gauss_pts));
421  pipeline.push_back(new OpCalculateScalarFieldValuesDot(
422  "PHOTON_FLUENCE_RATE", dot_u_at_gauss_pts));
423  pipeline.push_back(new OpDomainGradTimesVec(
424  "PHOTON_FLUENCE_RATE", grad_u_at_gauss_pts,
425  [](double, double, double) -> double { return D; }));
426  pipeline.push_back(new OpDomainTimesScalarField(
427  "PHOTON_FLUENCE_RATE", dot_u_at_gauss_pts,
428  [](const double, const double, const double) { return inv_v; }));
429  pipeline.push_back(new OpDomainTimesScalarField(
430  "PHOTON_FLUENCE_RATE", u_at_gauss_pts,
431  [](const double, const double, const double) { return mu_a; }));
432  };
433 
434  auto add_boundary_base_ops = [&](auto &pipeline) {
435  pipeline.push_back(new OpSetHOWeightsOnFace());
436  };
437 
438  auto add_boundary_lhs_ops = [&](auto &pipeline) {
439  for (auto b : bc_map) {
440  if (std::regex_match(b.first, std::regex("(.*)EXT(.*)"))) {
441  pipeline.push_back(new OpBoundaryMass(
442  "PHOTON_FLUENCE_RATE", "PHOTON_FLUENCE_RATE",
443 
444  [](const double, const double, const double) { return h; },
445 
446  b.second->getBcEntsPtr()));
447  }
448  }
449  };
450 
451  auto add_boundary_rhs_ops = [&](auto &pipeline) {
452  auto u_at_gauss_pts = boost::make_shared<VectorDouble>();
453  pipeline.push_back(new OpCalculateScalarFieldValues("PHOTON_FLUENCE_RATE",
454  u_at_gauss_pts));
455  for (auto b : bc_map) {
456  if (std::regex_match(b.first, std::regex("(.*)EXT(.*)"))) {
457  pipeline.push_back(new OpBoundaryTimeScalarField(
458  "PHOTON_FLUENCE_RATE", u_at_gauss_pts,
459 
460  [](const double, const double, const double) { return h; },
461 
462  b.second->getBcEntsPtr()));
463  }
464  }
465  };
466 
467  auto pipeline_mng = mField.getInterface<PipelineManager>();
468  add_domain_base_ops(pipeline_mng->getOpDomainLhsPipeline());
469  add_domain_base_ops(pipeline_mng->getOpDomainRhsPipeline());
470  add_domain_lhs_ops(pipeline_mng->getOpDomainLhsPipeline());
471  add_domain_rhs_ops(pipeline_mng->getOpDomainRhsPipeline());
472 
473  add_boundary_base_ops(pipeline_mng->getOpBoundaryLhsPipeline());
474  add_boundary_base_ops(pipeline_mng->getOpBoundaryRhsPipeline());
475  add_boundary_lhs_ops(pipeline_mng->getOpBoundaryLhsPipeline());
476  add_boundary_rhs_ops(pipeline_mng->getOpBoundaryRhsPipeline());
477 
478  domainLhsFEPtr = pipeline_mng->getDomainLhsFE();
479  boundaryLhsFEPtr = pipeline_mng->getBoundaryLhsFE();
480  boundaryRhsFEPtr = pipeline_mng->getBoundaryRhsFE();
481 
483 }
484 
487 
488  auto *simple = mField.getInterface<Simple>();
489  auto *pipeline_mng = mField.getInterface<PipelineManager>();
490 
491  auto create_post_process_element = [&]() {
492  auto post_froc_fe = boost::make_shared<PostProcEle>(mField);
493  auto u_ptr = boost::make_shared<VectorDouble>();
494  auto grad_ptr = boost::make_shared<MatrixDouble>();
495  post_froc_fe->getOpPtrVector().push_back(
496  new OpCalculateScalarFieldValues("PHOTON_FLUENCE_RATE", u_ptr));
497  post_froc_fe->getOpPtrVector().push_back(
498  new OpCalculateScalarFieldGradient<SPACE_DIM>("PHOTON_FLUENCE_RATE",
499  grad_ptr));
500  post_froc_fe->getOpPtrVector().push_back(new OpPPMap(
501  post_froc_fe->getPostProcMesh(), post_froc_fe->getMapGaussPts(),
502  {{"PHOTON_FLUENCE_RATE", u_ptr}},
503  {{"GRAD_PHOTON_FLUENCE_RATE", grad_ptr}}, {}, {}));
504  return post_froc_fe;
505  };
506 
507  auto create_post_process_camera_element = [&]() {
508  if (mField.check_finite_element("CAMERA_FE")) {
509 
510  auto post_proc_skin = boost::make_shared<PostProcFaceEle>(mField);
511 
512  auto u_ptr = boost::make_shared<VectorDouble>();
513  auto grad_ptr = boost::make_shared<MatrixDouble>();
514 
516  mField, simple->getDomainFEName(), SPACE_DIM);
517 
518  auto jac_ptr = boost::make_shared<MatrixDouble>();
519  auto inv_jac_ptr = boost::make_shared<MatrixDouble>();
520  auto det_ptr = boost::make_shared<VectorDouble>();
521 
522  // push operators to side element
523  op_loop_side->getOpPtrVector().push_back(
524  new OpCalculateHOJac<SPACE_DIM>(jac_ptr));
525  op_loop_side->getOpPtrVector().push_back(
526  new OpInvertMatrix<SPACE_DIM>(jac_ptr, det_ptr, inv_jac_ptr));
527  op_loop_side->getOpPtrVector().push_back(
528  new OpSetHOInvJacToScalarBases<SPACE_DIM>(H1, inv_jac_ptr));
529  op_loop_side->getOpPtrVector().push_back(
530  new OpCalculateScalarFieldValues("PHOTON_FLUENCE_RATE", u_ptr));
531  op_loop_side->getOpPtrVector().push_back(
532  new OpCalculateScalarFieldGradient<SPACE_DIM>("PHOTON_FLUENCE_RATE", grad_ptr));
533  // push op to boundary element
534  post_proc_skin->getOpPtrVector().push_back(op_loop_side);
535 
536  post_proc_skin->getOpPtrVector().push_back(new OpPPMap(
537  post_proc_skin->getPostProcMesh(), post_proc_skin->getMapGaussPts(),
538  {{"PHOTON_FLUENCE_RATE", u_ptr}},
539  {{"GRAD_PHOTON_FLUENCE_RATE", grad_ptr}}, {}, {}));
540 
541  return post_proc_skin;
542  } else {
543  return boost::shared_ptr<PostProcFaceEle>();
544  }
545  };
546 
547  auto create_post_process_integ_camera_element = [&]() {
548  if (mField.check_finite_element("CAMERA_FE")) {
549 
550  auto post_proc_integ_skin = boost::make_shared<BoundaryEle>(mField);
551  post_proc_integ_skin->getOpPtrVector().push_back(
552  new OpSetHOWeightsOnFace());
553  post_proc_integ_skin->getOpPtrVector().push_back(
554  new OpCalculateScalarFieldValues("PHOTON_FLUENCE_RATE",
555  commonDataPtr->approxVals));
556  post_proc_integ_skin->getOpPtrVector().push_back(
557  new OpCameraInteg(commonDataPtr));
558 
559  return post_proc_integ_skin;
560  } else {
561  return boost::shared_ptr<BoundaryEle>();
562  }
563  };
564 
565  auto set_time_monitor = [&](auto dm, auto solver) {
567  boost::shared_ptr<Monitor> monitor_ptr(new Monitor(
568  dm, create_post_process_element(), create_post_process_camera_element(),
569  create_post_process_integ_camera_element(), commonDataPtr));
570  boost::shared_ptr<ForcesAndSourcesCore> null;
571  CHKERR DMMoFEMTSSetMonitor(dm, solver, simple->getDomainFEName(),
572  monitor_ptr, null, null);
574  };
575 
576  auto dm = simple->getDM();
577  auto X = createDMVector(dm);
578 
579  if (from_initial) {
580 
581  MOFEM_LOG("PHOTON", Sev::inform) << "reading vector in binary from file "
582  << init_data_file_name << " ...";
583  PetscViewer viewer;
584  PetscViewerBinaryOpen(PETSC_COMM_WORLD, init_data_file_name, FILE_MODE_READ,
585  &viewer);
586  VecLoad(X, viewer);
587 
588  CHKERR DMoFEMMeshToLocalVector(dm, X, INSERT_VALUES, SCATTER_REVERSE);
589  }
590 
591  auto solver = pipeline_mng->createTS();
592 
593  CHKERR TSSetSolution(solver, X);
594  CHKERR set_time_monitor(dm, solver);
595  CHKERR TSSetSolution(solver, X);
596  CHKERR TSSetFromOptions(solver);
597  CHKERR TSSetUp(solver);
598  CHKERR TSSolve(solver, NULL);
599 
600  CHKERR VecGhostUpdateBegin(X, INSERT_VALUES, SCATTER_FORWARD);
601  CHKERR VecGhostUpdateEnd(X, INSERT_VALUES, SCATTER_FORWARD);
602  CHKERR DMoFEMMeshToLocalVector(dm, X, INSERT_VALUES, SCATTER_REVERSE);
603 
605 }
606 
609 
610  // Processes to set output results are integrated in solveSystem()
611 
613 }
614 
617 
618  CHKERR readMesh();
627 
629 }
630 
632  EntData &data) {
634  const int nb_integration_pts = getGaussPts().size2();
635  const double area = getMeasure();
636  auto t_w = getFTensor0IntegrationWeight();
637  auto t_val = getFTensor0FromVec(*(commonDataPtr->approxVals));
638 
639  double values_integ = 0;
640 
641  for (int gg = 0; gg != nb_integration_pts; ++gg) {
642  const double alpha = t_w * area;
643 
644  values_integ += alpha * t_val;
645 
646  ++t_w;
647  ++t_val;
648  }
649 
650  constexpr std::array<int, 1> indices = {CommonData::VALUES_INTEG};
651  std::array<double, 1> values;
652  values[0] = values_integ;
653  CHKERR VecSetValues(commonDataPtr->petscVec, 1, indices.data(), values.data(),
654  ADD_VALUES);
656 }
657 
658 int main(int argc, char *argv[]) {
659 
660  // Initialisation of MoFEM/PETSc and MOAB data structures
661  const char param_file[] = "param_file.petsc";
662  MoFEM::Core::Initialize(&argc, &argv, param_file, help);
663 
664  // Add logging channel for example
665  auto core_log = logging::core::get();
666  core_log->add_sink(
667  LogManager::createSink(LogManager::getStrmWorld(), "PHOTON"));
668  LogManager::setLog("PHOTON");
669  MOFEM_LOG_TAG("PHOTON", "photon_diffusion")
670 
671  // Error handling
672  try {
673  // Register MoFEM discrete manager in PETSc
674  DMType dm_name = "DMMOFEM";
675  CHKERR DMRegister_MoFEM(dm_name);
676 
677  // Create MOAB instance
678  moab::Core mb_instance; // mesh database
679  moab::Interface &moab = mb_instance; // mesh database interface
680 
681  // Create MoFEM instance
682  MoFEM::Core core(moab); // finite element database
683  MoFEM::Interface &m_field = core; // finite element interface
684 
685  // Run the main analysis
686  PhotonDiffusion heat_problem(m_field);
687  CHKERR heat_problem.runProgram();
688  }
689  CATCH_ERRORS;
690 
691  // Finish work: cleaning memory, getting statistics, etc.
693 
694  return 0;
695 }
n
const double n
refractive index of diffusive medium
Definition: photon_diffusion.cpp:51
MoFEMFunctionReturnHot
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:460
MoFEM::UnknownInterface::getInterface
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.
Definition: UnknownInterface.hpp:93
MoFEM::EntitiesFieldData::EntData
Data on single entity (This is passed as argument to DataOperator::doWork)
Definition: EntitiesFieldData.hpp:128
PhotonDiffusion::Monitor::postProc
boost::shared_ptr< PostProcEle > postProc
Definition: photon_diffusion.cpp:203
BoundaryEle
ElementsAndOps< SPACE_DIM >::BoundaryEle BoundaryEle
Definition: child_and_parent.cpp:39
PhotonDiffusion::boundaryRhsFEPtr
boost::shared_ptr< FEMethod > boundaryRhsFEPtr
Definition: photon_diffusion.cpp:100
MoFEM::CoreTmp< 0 >
Core (interface) class.
Definition: Core.hpp:82
H1
@ H1
continuous field
Definition: definitions.h:85
main
int main(int argc, char *argv[])
Definition: photon_diffusion.cpp:658
help
static char help[]
Definition: photon_diffusion.cpp:13
PhotonDiffusion::OpGetScalarFieldGradientValuesOnSkin
Definition: photon_diffusion.cpp:125
OpDomainGradTimesVec
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::LinearForm< GAUSS >::OpGradTimesTensor< 1, 1, SPACE_DIM > OpDomainGradTimesVec
Definition: initial_diffusion.cpp:34
PhotonDiffusion::Monitor::skinPostProc
boost::shared_ptr< PostProcFaceEle > skinPostProc
Definition: photon_diffusion.cpp:204
order
int order
Definition: photon_diffusion.cpp:66
A
double A
Definition: photon_diffusion.cpp:59
MoFEM::ProblemsManager
Problem manager is used to build and partition problems.
Definition: ProblemsManager.hpp:21
MoFEM::FEMethod
structure for User Loop Methods on finite elements
Definition: LoopMethods.hpp:369
PhotonDiffusion::commonDataPtr
boost::shared_ptr< CommonData > commonDataPtr
Definition: photon_diffusion.cpp:112
PhotonDiffusion
Definition: initial_diffusion.cpp:74
MoFEM::CoreInterface::get_comm
virtual MPI_Comm & get_comm() const =0
PhotonDiffusion::boundaryLhsFEPtr
boost::shared_ptr< FEMethod > boundaryLhsFEPtr
Definition: photon_diffusion.cpp:99
OpDomainGradGrad
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::BiLinearForm< GAUSS >::OpGradGrad< 1, 1, SPACE_DIM > OpDomainGradGrad
Definition: photon_diffusion.cpp:36
PhotonDiffusion::CommonData::VecElements
VecElements
Definition: photon_diffusion.cpp:106
MoFEM::Exceptions::MoFEMErrorCode
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
Definition: Exceptions.hpp:56
OpBoundaryTimeScalarField
FormsIntegrators< BoundaryEleOp >::Assembly< PETSC >::LinearForm< GAUSS >::OpBaseTimesScalar< 1 > OpBoundaryTimeScalarField
Definition: photon_diffusion.cpp:47
MoFEM::CoreInterface::get_comm_rank
virtual int get_comm_rank() const =0
MoFEM::OpCalculateScalarFieldValuesDot
OpCalculateScalarFieldValuesFromPetscVecImpl< PetscData::CTX_SET_X_T > OpCalculateScalarFieldValuesDot
Definition: UserDataOperators.hpp:273
OpDomainGradTimesVec
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::LinearForm< GAUSS >::OpGradTimesTensor< 1, 1, SPACE_DIM > OpDomainGradTimesVec
Definition: photon_diffusion.cpp:40
MoFEM::PETSC
@ PETSC
Definition: FormsIntegrators.hpp:105
MoFEM::PipelineManager
PipelineManager interface.
Definition: PipelineManager.hpp:24
PhotonDiffusion::setIntegrationRules
MoFEMErrorCode setIntegrationRules()
Definition: initial_diffusion.cpp:236
MoFEM::DMoFEMMeshToLocalVector
PetscErrorCode DMoFEMMeshToLocalVector(DM dm, Vec l, InsertMode mode, ScatterMode scatter_mode)
set local (or ghosted) vector values on mesh for partition only
Definition: DMMoFEM.cpp:523
BasicFiniteElements.hpp
MoFEM::CoreTmp< 0 >::Finalize
static MoFEMErrorCode Finalize()
Checks for options to be called at the conclusion of the program.
Definition: Core.cpp:112
PhotonDiffusion::Monitor::commonDataPtr
boost::shared_ptr< CommonData > commonDataPtr
Definition: photon_diffusion.cpp:206
PhotonDiffusion::Monitor::preProcess
MoFEMErrorCode preProcess()
function is run at the beginning of loop
Definition: photon_diffusion.cpp:152
MoFEM::Simple
Simple interface for fast problem set-up.
Definition: Simple.hpp:27
PhotonDiffusion::PhotonDiffusion
PhotonDiffusion(MoFEM::Interface &m_field)
Definition: initial_diffusion.cpp:138
OpDomainMass
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::BiLinearForm< GAUSS >::OpMass< 1, FIELD_DIM > OpDomainMass
Definition: child_and_parent.cpp:53
MoFEM::VecSetValues
MoFEMErrorCode VecSetValues(Vec V, const EntitiesFieldData::EntData &data, const double *ptr, InsertMode iora)
Assemble PETSc vector.
Definition: EntitiesFieldData.hpp:1579
PhotonDiffusion::readMesh
MoFEMErrorCode readMesh()
Definition: initial_diffusion.cpp:140
MoFEM::DMMoFEMAddElement
PetscErrorCode DMMoFEMAddElement(DM dm, std::string fe_name)
add element to dm
Definition: DMMoFEM.cpp:497
numHoLevels
int numHoLevels
Definition: photon_diffusion.cpp:70
MoFEM::DeprecatedCoreInterface
Deprecated interface functions.
Definition: DeprecatedCoreInterface.hpp:16
PhotonDiffusion::domainLhsFEPtr
boost::shared_ptr< FEMethod > domainLhsFEPtr
Definition: photon_diffusion.cpp:98
OpBoundaryMass
FormsIntegrators< BoundaryEleOp >::Assembly< PETSC >::BiLinearForm< GAUSS >::OpMass< 1, SPACE_DIM > OpBoundaryMass
[Only used with Hencky/nonlinear material]
Definition: seepage.cpp:71
MoFEM::OpCalculateScalarFieldGradient
Get field gradients at integration pts for scalar filed rank 0, i.e. vector field.
Definition: UserDataOperators.hpp:1293
MoFEM::OpCalculateHOJac
Definition: HODataOperators.hpp:267
MoFEM::Interface
DeprecatedCoreInterface Interface
Definition: Interface.hpp:2002
PhotonDiffusion::OpGetScalarFieldGradientValuesOnSkin::doWork
MoFEMErrorCode doWork(int side, EntityType type, DataForcesAndSourcesCore::EntData &data)
Definition: photon_diffusion.cpp:132
PhotonDiffusion::setupProblem
MoFEMErrorCode setupProblem()
Definition: initial_diffusion.cpp:151
MoFEM::PostProcBrokenMeshInMoab
Definition: PostProcBrokenMeshInMoabBase.hpp:667
output_volume
PetscBool output_volume
Definition: photon_diffusion.cpp:63
save_every_nth_step
int save_every_nth_step
Definition: photon_diffusion.cpp:67
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:548
PhotonDiffusion::OpGetScalarFieldGradientValuesOnSkin::OpGetScalarFieldGradientValuesOnSkin
OpGetScalarFieldGradientValuesOnSkin(boost::shared_ptr< VolSideFe > side_fe)
Definition: photon_diffusion.cpp:129
MoFEM::createDMVector
auto createDMVector(DM dm)
Get smart vector from DM.
Definition: DMMoFEM.hpp:1099
MoFEM::CoreInterface::add_finite_element
virtual MoFEMErrorCode add_finite_element(const std::string &fe_name, enum MoFEMTypes bh=MF_EXCL, int verb=DEFAULT_VERBOSITY)=0
add finite element
MoFEM::createGhostVector
auto createGhostVector(MPI_Comm comm, PetscInt n, PetscInt N, PetscInt nghost, const PetscInt ghosts[])
Create smart ghost vector.
Definition: PetscSmartObj.hpp:179
from_initial
PetscBool from_initial
Definition: photon_diffusion.cpp:62
PhotonDiffusion::Monitor::skinPostProcInteg
boost::shared_ptr< BoundaryEle > skinPostProcInteg
Definition: photon_diffusion.cpp:205
MoFEM::CoreInterface::get_moab
virtual moab::Interface & get_moab()=0
PhotonDiffusion::Monitor::operator()
MoFEMErrorCode operator()()
function is run for every finite element
Definition: photon_diffusion.cpp:156
MoFEM
implementation of Data Operators for Forces and Sources
Definition: Common.hpp:10
mu_sp
double mu_sp
scattering coefficient (cm^-1)
Definition: photon_diffusion.cpp:57
OpDomainTimesScalarField
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::LinearForm< GAUSS >::OpBaseTimesScalar< 1 > OpDomainTimesScalarField
Definition: initial_diffusion.cpp:32
MoFEM::BcManager
Simple interface for fast problem set-up.
Definition: BcManager.hpp:25
BoundaryEleOp
mu_a
double mu_a
absorption coefficient (cm^-1)
Definition: photon_diffusion.cpp:56
PhotonDiffusion::Monitor::postProcess
MoFEMErrorCode postProcess()
function is run at the end of loop
Definition: photon_diffusion.cpp:161
PhotonDiffusion::mField
MoFEM::Interface & mField
Definition: initial_diffusion.cpp:135
MoFEM::CoreInterface::build_finite_elements
virtual MoFEMErrorCode build_finite_elements(int verb=DEFAULT_VERBOSITY)=0
Build finite elements.
PhotonDiffusion::initialCondition
MoFEMErrorCode initialCondition()
Definition: initial_diffusion.cpp:252
simple
void simple(double P1[], double P2[], double P3[], double c[], const int N)
Definition: acoustic.cpp:69
PhotonDiffusion::CommonData::petscVec
SmartPetscObj< Vec > petscVec
Definition: photon_diffusion.cpp:104
MoFEM::OpSetHOInvJacToScalarBases
Set inverse jacobian to base functions.
Definition: HODataOperators.hpp:73
double
bit
auto bit
set bit
Definition: hanging_node_approx.cpp:75
OpDomainTimesScalarField
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::LinearForm< GAUSS >::OpBaseTimesScalar< 1 > OpDomainTimesScalarField
Definition: photon_diffusion.cpp:38
convert.type
type
Definition: convert.py:64
MoFEM::BcManager::getBcMapByBlockName
BcMapByBlockName & getBcMapByBlockName()
Get the bc map.
Definition: BcManager.hpp:207
OpBoundaryMass
FormsIntegrators< BoundaryEleOp >::Assembly< PETSC >::BiLinearForm< GAUSS >::OpMass< 1, 1 > OpBoundaryMass
Definition: photon_diffusion.cpp:45
OpBoundarySource
FormsIntegrators< EdgeEleOp >::Assembly< PETSC >::LinearForm< GAUSS >::OpSource< 1, 1 > OpBoundarySource
Definition: helmholtz.cpp:31
PhotonDiffusion::Monitor::dM
SmartPetscObj< DM > dM
Definition: photon_diffusion.cpp:202
MoFEM::FormsIntegrators::Assembly
Assembly methods.
Definition: FormsIntegrators.hpp:312
OpPPMap
OpPostProcMapInMoab< SPACE_DIM, SPACE_DIM > OpPPMap
Definition: photon_diffusion.cpp:29
PhotonDiffusion::CommonData
Definition: photon_diffusion.cpp:102
MoFEM::getFTensor0FromVec
static auto getFTensor0FromVec(ublas::vector< T, A > &data)
Get tensor rank 0 (scalar) form data vector.
Definition: Templates.hpp:135
h
double h
Definition: photon_diffusion.cpp:60
MoFEM::OpCalculateScalarFieldValues
Get value at integration points for scalar field.
Definition: UserDataOperators.hpp:82
MoFEM::DMRegister_MoFEM
PetscErrorCode DMRegister_MoFEM(const char sname[])
Register MoFEM problem.
Definition: DMMoFEM.cpp:43
D
double D
Definition: photon_diffusion.cpp:58
MoFEM::PostProcBrokenMeshInMoab< VolumeElementForcesAndSourcesCore >
Definition: PostProcBrokenMeshInMoabBase.hpp:670
PhotonDiffusion::OpCameraInteg::commonDataPtr
boost::shared_ptr< CommonData > commonDataPtr
Definition: photon_diffusion.cpp:115
MoFEM::OpSetHOWeightsOnFace
Modify integration weights on face to take in account higher-order geometry.
Definition: HODataOperators.hpp:122
MoFEM::FaceElementForcesAndSourcesCore
Face finite element.
Definition: FaceElementForcesAndSourcesCore.hpp:23
PhotonDiffusion::Monitor::Monitor
Monitor(SmartPetscObj< DM > dm, boost::shared_ptr< PostProcEle > post_proc, boost::shared_ptr< PostProcFaceEle > skin_post_proc, boost::shared_ptr< BoundaryEle > skin_post_proc_integ, boost::shared_ptr< CommonData > common_data_ptr)
Definition: photon_diffusion.cpp:144
MoFEM::OpCalculateHOJac< 3 >
Definition: HODataOperators.hpp:269
PhotonDiffusion::runProgram
MoFEMErrorCode runProgram()
Definition: initial_diffusion.cpp:383
PhotonDiffusion::createCommonData
MoFEMErrorCode createCommonData()
Definition: photon_diffusion.cpp:223
PhotonDiffusion::OpCameraInteg
Definition: photon_diffusion.cpp:114
MOFEM_LOG_TAG
#define MOFEM_LOG_TAG(channel, tag)
Tag channel.
Definition: LogManager.hpp:339
PhotonDiffusion::CommonData::VALUES_INTEG
@ VALUES_INTEG
Definition: photon_diffusion.cpp:107
MoFEM::CoreInterface::add_ents_to_finite_element_by_dim
virtual MoFEMErrorCode add_ents_to_finite_element_by_dim(const EntityHandle entities, const int dim, const std::string &name, const bool recursive=true)=0
add entities to finite element
PhotonDiffusion::CommonData::approxVals
boost::shared_ptr< VectorDouble > approxVals
Definition: photon_diffusion.cpp:103
MoFEM::VolumeElementForcesAndSourcesCore
Volume finite element base.
Definition: VolumeElementForcesAndSourcesCore.hpp:26
PhotonDiffusion::solveSystem
MoFEMErrorCode solveSystem()
Definition: initial_diffusion.cpp:318
BiLinearForm
SPACE_DIM
constexpr int SPACE_DIM
[Define dimension]
Definition: photon_diffusion.cpp:18
OpGradTimesTensor
FormsIntegrators< DomainEleOp >::Assembly< A >::LinearForm< I >::OpGradTimesTensor< 1, FIELD_DIM, SPACE_DIM > OpGradTimesTensor
Definition: operators_tests.cpp:48
PhotonDiffusion::boundaryCondition
MoFEMErrorCode boundaryCondition()
Definition: initial_diffusion.cpp:257
c
const double c
speed of light (cm/ns)
Definition: photon_diffusion.cpp:52
EntData
EntitiesFieldData::EntData EntData
Definition: child_and_parent.cpp:37
v
const double v
phase velocity of light in medium (cm/ns)
Definition: photon_diffusion.cpp:53
inv_v
const double inv_v
Definition: photon_diffusion.cpp:54
integration_rule
auto integration_rule
Definition: free_surface.cpp:185
VolSideFe
VolumeElementForcesAndSourcesCoreOnSide VolSideFe
Definition: mortar_contact_thermal.cpp:29
ElementsAndOps
Definition: child_and_parent.cpp:18
Range
output_camera
PetscBool output_camera
Definition: photon_diffusion.cpp:64
DomainEleOp
PhotonDiffusion::boundaryMarker
boost::shared_ptr< std::vector< unsigned char > > boundaryMarker
Definition: photon_diffusion.cpp:96
PhotonDiffusion::assembleSystem
MoFEMErrorCode assembleSystem()
Definition: initial_diffusion.cpp:284
MoFEM::CoreTmp< 0 >::Initialize
static MoFEMErrorCode Initialize(int *argc, char ***args, const char file[], const char help[])
Initializes the MoFEM database PETSc, MOAB and MPI.
Definition: Core.cpp:72
PhotonDiffusion::Monitor
Definition: photon_diffusion.cpp:142
MOFEM_LOG
#define MOFEM_LOG(channel, severity)
Log.
Definition: LogManager.hpp:308
OpDomainGradGrad
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::BiLinearForm< GAUSS >::OpGradGrad< 1, 1, SPACE_DIM > OpDomainGradGrad
Definition: helmholtz.cpp:25
CATCH_ERRORS
#define CATCH_ERRORS
Catch errors.
Definition: definitions.h:385
_IT_CUBITMESHSETS_BY_SET_TYPE_FOR_LOOP_
#define _IT_CUBITMESHSETS_BY_SET_TYPE_FOR_LOOP_(MESHSET_MANAGER, CUBITBCTYPE, IT)
Iterator that loops over a specific Cubit MeshSet having a particular BC meshset in a moFEM field.
Definition: MeshsetsManager.hpp:71
MoFEM::Core
CoreTmp< 0 > Core
Definition: Core.hpp:1140
PhotonDiffusion::OpCameraInteg::doWork
MoFEMErrorCode doWork(int side, EntityType type, EntData &data)
Definition: photon_diffusion.cpp:631
UserDataOperator
ForcesAndSourcesCore::UserDataOperator UserDataOperator
Definition: HookeElement.hpp:75
PhotonDiffusion::OpCameraInteg::OpCameraInteg
OpCameraInteg(boost::shared_ptr< CommonData > common_data_ptr)
Definition: photon_diffusion.cpp:116
BLOCKSET
@ BLOCKSET
Definition: definitions.h:161
Monitor
[Push operators to pipeline]
Definition: adolc_plasticity.cpp:333
approx_order
int approx_order
Definition: test_broken_space.cpp:50
MoFEM::OpSetHOWeights
Set inverse jacobian to base functions.
Definition: HODataOperators.hpp:159
MoFEM::DMMoFEMTSSetMonitor
PetscErrorCode DMMoFEMTSSetMonitor(DM dm, TS ts, const std::string fe_name, boost::shared_ptr< MoFEM::FEMethod > method, boost::shared_ptr< MoFEM::BasicMethod > pre_only, boost::shared_ptr< MoFEM::BasicMethod > post_only)
Set Monitor To TS solver.
Definition: DMMoFEM.cpp:1056
AINSWORTH_LEGENDRE_BASE
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base .
Definition: definitions.h:60
MoFEM::PetscOptionsGetString
PetscErrorCode PetscOptionsGetString(PetscOptions *, const char pre[], const char name[], char str[], size_t size, PetscBool *set)
Definition: DeprecatedPetsc.hpp:172
MoFEM::MeshsetsManager
Interface for managing meshsets containing materials and boundary conditions.
Definition: MeshsetsManager.hpp:104
MoFEM::OpInvertMatrix
Definition: UserDataOperators.hpp:3249
OpDomainMass
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::BiLinearForm< GAUSS >::OpMass< 1, 1 > OpDomainMass
Definition: photon_diffusion.cpp:34
init_data_file_name
char init_data_file_name[255]
Definition: photon_diffusion.cpp:69
MoFEM::CoreInterface::modify_finite_element_add_field_data
virtual MoFEMErrorCode modify_finite_element_add_field_data(const std::string &fe_name, const std::string name_filed)=0
set finite element field data
MoFEM::PetscOptionsGetScalar
PetscErrorCode PetscOptionsGetScalar(PetscOptions *, const char pre[], const char name[], PetscScalar *dval, PetscBool *set)
Definition: DeprecatedPetsc.hpp:162
MoFEM::CoreInterface::check_finite_element
virtual bool check_finite_element(const std::string &name) const =0
Check if finite element is in database.
MoFEM::VolumeElementForcesAndSourcesCoreOnSide
Base volume element used to integrate on skeleton.
Definition: VolumeElementForcesAndSourcesCoreOnSide.hpp:22
DomainEle
ElementsAndOps< SPACE_DIM >::DomainEle DomainEle
Definition: child_and_parent.cpp:34
MoFEM::SmartPetscObj< Vec >
OpDomainSource
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::LinearForm< GAUSS >::OpSource< 1, FIELD_DIM > OpDomainSource
Definition: child_and_parent.cpp:55
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
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
PhotonDiffusion::OpGetScalarFieldGradientValuesOnSkin::sideOpFe
boost::shared_ptr< VolSideFe > sideOpFe
Definition: photon_diffusion.cpp:127
MoFEM::OpLoopSide
Element used to execute operators on side of the element.
Definition: ForcesAndSourcesCore.hpp:1290
MoFEMFunctionBegin
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:359
MoFEM::PostProcBrokenMeshInMoab< FaceElementForcesAndSourcesCore >
Definition: PostProcBrokenMeshInMoabBase.hpp:677
MoFEM::PetscOptionsGetBool
PetscErrorCode PetscOptionsGetBool(PetscOptions *, const char pre[], const char name[], PetscBool *bval, PetscBool *set)
Definition: DeprecatedPetsc.hpp:182
MoFEM::OpPostProcMapInMoab
Post post-proc data at points from hash maps.
Definition: PostProcBrokenMeshInMoabBase.hpp:698
PhotonDiffusion::outputResults
MoFEMErrorCode outputResults()
Definition: initial_diffusion.cpp:349