v0.14.0
hanging_node_approx.cpp
Go to the documentation of this file.
1 /**
2  * \example hanging_node_approx.cpp
3  *
4  * Testing approximation with hanging nodes.
5  *
6  */
7 
8 #include <MoFEM.hpp>
9 
10 using namespace MoFEM;
11 
12 static char help[] = "...\n\n";
13 
14 constexpr char FIELD_NAME[] = "U";
15 constexpr int FIELD_DIM = 1;
16 constexpr int SPACE_DIM = 2;
17 constexpr int nb_ref_levels = 3; ///< Three levels of refinement
18 
20 using DomainParentEle =
22 using BoundaryEle =
24 using BoundaryParentEle =
28 
30 
31 template <int FIELD_DIM> struct ApproxFieldFunction;
32 template <int FIELD_DIM> struct ApproxFieldFunctionDerivative;
33 
34 /**
35  * @brief third order polynomial used for testing
36  *
37  */
38 template <> struct ApproxFieldFunction<1> {
39  auto operator()(const double x, const double y, const double z) {
40  return x * x + y * y + x * y * y + x * x * y;
41  }
42 };
43 
44 /**
45  * @brief third order polynomial used for testing
46  *
47  */
48 template <> struct ApproxFieldFunctionDerivative<1> {
49  auto operator()(const double x, const double y, const double z) {
50  // x * x + y * y + x * y * y + x * x * y
51 
52  return FTensor::Tensor1<double, SPACE_DIM>{2 * x + y * y + 2 * x * y,
53  2 * y + 2 * x * y + x * x};
54  }
55 };
56 
57 /**
58  * @brief evaluate mass matrix
59  *
60  */
63 
64 /**
65  * @brief evaluate source, i.e. rhs vector
66  *
67  */
69  PETSC>::LinearForm<GAUSS>::OpSource<1, FIELD_DIM>;
70 
71 /**
72  * @brief set bit
73  *
74  */
75 auto bit = [](auto l) { return BitRefLevel().set(l); };
76 
77 /**
78  * @brief set bit to marker
79  *
80  * Marker is used to mark field entities on skin on which we have hanging nodes
81  */
82 auto marker = [](auto l) {
83  return BitRefLevel().set(BITREFLEVEL_SIZE - 1 - l);
84 };
85 
86 /**
87  * @brief set levels of projection operators, which project field data from
88  * parent entities, to child, up to to level, i.e. last mesh refinement.
89  *
90  */
91 template <typename PARENT_FE>
93  boost::shared_ptr<FEMethod> &fe_top,
95  int verbosity, LogManager::SeverityLevel sev) {
96 
97  BitRefLevel bit_marker;
98  for (auto l = 1; l <= nb_ref_levels; ++l)
99  bit_marker |= marker(l);
100 
101  /**
102  * @brief Collect data from parent elements to child
103  */
104  boost::function<void(boost::shared_ptr<ForcesAndSourcesCore>, int)>
105  add_parent_level =
106  [&](boost::shared_ptr<ForcesAndSourcesCore> parent_fe_pt, int level) {
107 
108  // Evaluate if not last parent element
109  if (level > 0) {
110 
111  // Create domain parent FE
112  auto fe_ptr_current = boost::shared_ptr<ForcesAndSourcesCore>(
113  new PARENT_FE(m_field));
114  if (op == DomainEleOp::OPSPACE) {
115  // Push base function
116  if (typeid(PARENT_FE) == typeid(DomainParentEle))
118  fe_ptr_current->getOpPtrVector(), {H1});
119  }
120 
121  // Call next level
122  add_parent_level(
123  boost::dynamic_pointer_cast<ForcesAndSourcesCore>(
124  fe_ptr_current),
125  level - 1);
126 
127  // Add data to curent fe level
128  if (op == DomainEleOp::OPSPACE) {
129 
130  // Only base
131  parent_fe_pt->getOpPtrVector().push_back(
132 
133  new OpAddParentEntData(
134 
135  H1, op, fe_ptr_current,
136 
137  BitRefLevel().set(), bit(0).flip(),
138 
139  bit_marker, BitRefLevel().set(),
140 
141  verbosity, sev));
142 
143  } else {
144 
145  // Filed data
146  parent_fe_pt->getOpPtrVector().push_back(
147 
148  new OpAddParentEntData(
149 
150  FIELD_NAME, op, fe_ptr_current,
151 
152  BitRefLevel().set(), bit(0).flip(),
153 
154  bit_marker, BitRefLevel().set(),
155 
156  verbosity, sev));
157  }
158  }
159  };
160 
161  add_parent_level(boost::dynamic_pointer_cast<ForcesAndSourcesCore>(fe_top),
162  nb_ref_levels);
163 };
164 
165 /**
166  * @brief lambda function used to select elements on which finite element
167  * pipelines are executed.
168  *
169  * @note childs elements on pipeline, retrieve data from parents using operators
170  * pushed by \ref set_parent_dofs
171  *
172  */
173 auto test_bit_child = [](FEMethod *fe_ptr) {
174  return fe_ptr->numeredEntFiniteElementPtr->getBitRefLevel().test(
175  nb_ref_levels);
176 };
177 
178 struct AtomTest {
179 
180  AtomTest(MoFEM::Interface &m_field) : mField(m_field) {}
181 
182  MoFEMErrorCode runProblem();
183 
184 private:
185  MoFEM::Interface &mField;
186  Simple *simpleInterface;
187 
188  static ApproxFieldFunction<FIELD_DIM> approxFunction;
190 
191  /**
192  * @brief red mesh and randomly refine three times
193  *
194  * @return MoFEMErrorCode
195  */
196  MoFEMErrorCode readMesh();
197 
198  /**
199  * @brief add field, and set up problem
200  *
201  * @return MoFEMErrorCode
202  */
203  MoFEMErrorCode setupProblem();
204  MoFEMErrorCode assembleSystem();
205  MoFEMErrorCode solveSystem();
206  MoFEMErrorCode checkResults();
207  MoFEMErrorCode printResults();
208  struct CommonData {
209  boost::shared_ptr<VectorDouble> approxVals;
210  boost::shared_ptr<MatrixDouble> divApproxVals;
211  SmartPetscObj<Vec> L2Vec;
212  SmartPetscObj<Vec> resVec;
213  };
214 
215  template <int FIELD_DIM> struct OpError;
216 
217  template <int FIELD_DIM> struct OpErrorSkel;
218 };
219 
224 template <> struct AtomTest::OpError<1> : public DomainEleOp {
225  boost::shared_ptr<CommonData> commonDataPtr;
226  OpError(boost::shared_ptr<CommonData> &common_data_ptr)
227  : DomainEleOp(FIELD_NAME, OPROW), commonDataPtr(common_data_ptr) {}
228  MoFEMErrorCode doWork(int side, EntityType type, EntData &data) {
230 
231  if (const size_t nb_dofs = data.getIndices().size()) {
232 
234 
235  const int nb_integration_pts = getGaussPts().size2();
236  auto t_w = getFTensor0IntegrationWeight();
237  auto t_val = getFTensor0FromVec(*(commonDataPtr->approxVals));
238  auto t_grad_val =
239  getFTensor1FromMat<SPACE_DIM>(*(commonDataPtr->divApproxVals));
240  auto t_coords = getFTensor1CoordsAtGaussPts();
241 
242  VectorDouble nf(nb_dofs, false);
243  nf.clear();
244 
245  const double volume = getMeasure();
246 
247  auto t_row_base = data.getFTensor0N();
248  double error = 0;
249  for (int gg = 0; gg != nb_integration_pts; ++gg) {
250 
251  const double alpha = t_w * volume;
252  double diff = t_val - AtomTest::approxFunction(t_coords(0), t_coords(1),
253  t_coords(2));
254 
255  auto t_grad_diff =
256  AtomTest::divApproxFunction(t_coords(0), t_coords(1), t_coords(2));
257  t_grad_diff(i) -= t_grad_val(i);
258 
259  error += alpha * (pow(diff, 2) + t_grad_diff(i) * t_grad_diff(i));
260 
261  for (size_t r = 0; r != nb_dofs; ++r) {
262  nf[r] += alpha * t_row_base * diff;
263  ++t_row_base;
264  }
265 
266  ++t_w;
267  ++t_val;
268  ++t_grad_val;
269  ++t_coords;
270  }
271 
272  const int index = 0;
273  CHKERR VecSetValue(commonDataPtr->L2Vec, index, error, ADD_VALUES);
274  CHKERR VecSetValues(commonDataPtr->resVec, data, &nf[0], ADD_VALUES);
275  }
276 
278  }
279 };
280 
281 template <> struct AtomTest::OpErrorSkel<1> : public BoundaryEleOp {
282  boost::shared_ptr<CommonData> commonDataPtr;
283  OpErrorSkel(boost::shared_ptr<CommonData> &common_data_ptr)
284  : BoundaryEleOp(H1, OPSPACE), commonDataPtr(common_data_ptr) {}
285  MoFEMErrorCode doWork(int side, EntityType type, EntData &data) {
287 
289 
290  const int nb_integration_pts = getGaussPts().size2();
291  auto t_w = getFTensor0IntegrationWeight();
292  auto t_val = getFTensor0FromVec(*(commonDataPtr->approxVals));
293  auto t_coords = getFTensor1CoordsAtGaussPts();
294 
295  const double volume = getMeasure();
296 
297  double error2 = 0;
298  for (int gg = 0; gg != nb_integration_pts; ++gg) {
299 
300  const double alpha = t_w * volume;
301  double diff = t_val - AtomTest::approxFunction(t_coords(0), t_coords(1),
302  t_coords(2));
303  error2 += alpha * (pow(diff, 2));
304 
305  ++t_w;
306  ++t_val;
307  ++t_coords;
308  }
309 
310  MOFEM_LOG("SELF", Sev::verbose) << "Boundary error " << sqrt(error2);
311 
312  constexpr double eps = 1e-8;
313  if (sqrt(error2) > eps)
314  SETERRQ1(PETSC_COMM_SELF, MOFEM_ATOM_TEST_INVALID,
315  "Error on boundary = %6.4e", sqrt(error2));
316 
318  }
319 };
320 
321 //! [Run programme]
324  CHKERR readMesh();
325  CHKERR setupProblem();
326  CHKERR assembleSystem();
327  CHKERR solveSystem();
328  CHKERR checkResults();
329  CHKERR printResults();
331 }
332 //! [Run programme]
333 
334 //! [Read mesh]
336  BitRefManager *bit_mng = mField.getInterface<BitRefManager>();
337  ParallelComm *pcomm =
338  ParallelComm::get_pcomm(&mField.get_moab(), MYPCOMM_INDEX);
339  Skinner skin(&mField.get_moab());
341 
342  CHKERR mField.getInterface(simpleInterface);
343  CHKERR simpleInterface->getOptions();
344  CHKERR simpleInterface->loadFile();
345 
346  MOFEM_LOG("WORLD", Sev::verbose) << "Dim " << simpleInterface->getDim();
347 
348  auto &moab = mField.get_moab();
349 
350  Range level0_ents;
351  CHKERR mField.getInterface<BitRefManager>()->getEntitiesByDimAndRefLevel(
352  bit(0), BitRefLevel().set(), SPACE_DIM, level0_ents);
353  Range level0_skin;
354  CHKERR skin.find_skin(0, level0_ents, false, level0_skin);
355  CHKERR pcomm->filter_pstatus(level0_skin,
356  PSTATUS_SHARED | PSTATUS_MULTISHARED,
357  PSTATUS_NOT, -1, nullptr);
358 
359  auto refine_mesh = [&](auto l) {
361 
362  auto refine = mField.getInterface<MeshRefinement>();
363 
364  auto meshset_level0_ptr = get_temp_meshset_ptr(moab);
365  CHKERR bit_mng->getEntitiesByDimAndRefLevel(bit(l - 1), BitRefLevel().set(),
366  SPACE_DIM, *meshset_level0_ptr);
367 
368  // random mesh refinement
369  auto meshset_ref_edges_ptr = get_temp_meshset_ptr(moab);
370 
371  Range els;
372  CHKERR moab.get_entities_by_dimension(*meshset_level0_ptr, SPACE_DIM, els);
373  CHKERR bit_mng->filterEntitiesByRefLevel(bit(l - 1), bit(l - 1), els);
374 
375  Range ele_to_refine;
376 
377  if (l == 1) {
378  int ii = 0;
379  for (auto t : els) {
380  if ((ii % 2)) {
381  ele_to_refine.insert(t);
382  std::vector<EntityHandle> adj_edges;
383  CHKERR mField.get_moab().get_adjacencies(&t, 1, SPACE_DIM - 1, false,
384  adj_edges);
385  CHKERR moab.add_entities(*meshset_ref_edges_ptr, &*adj_edges.begin(),
386  adj_edges.size());
387  }
388  ++ii;
389  }
390  } else {
391  Range level_skin;
392  CHKERR skin.find_skin(0, els, false, level_skin);
393  CHKERR pcomm->filter_pstatus(level_skin,
394  PSTATUS_SHARED | PSTATUS_MULTISHARED,
395  PSTATUS_NOT, -1, nullptr);
396  level_skin = subtract(level_skin, level0_skin);
397  Range adj;
398  CHKERR mField.get_moab().get_adjacencies(level_skin, SPACE_DIM, false,
399  adj, moab::Interface::UNION);
400  els = subtract(els, adj);
401  ele_to_refine.merge(els);
402  Range adj_edges;
403  CHKERR mField.get_moab().get_adjacencies(
404  els, SPACE_DIM - 1, false, adj_edges, moab::Interface::UNION);
405  CHKERR moab.add_entities(*meshset_ref_edges_ptr, adj_edges);
406  }
407 
408  CHKERR refine->addVerticesInTheMiddleOfEdges(*meshset_ref_edges_ptr, bit(l),
409  false, VERBOSE);
410  CHKERR refine->refineTrisHangingNodes(*meshset_level0_ptr, bit(l), VERBOSE);
411  CHKERR bit_mng->updateRangeByChildren(level0_skin, level0_skin);
413  simpleInterface->getBoundaryMeshSet(), bit(l - 1), BitRefLevel(),
414  bit(l), BitRefLevel(), simpleInterface->getBoundaryMeshSet(), MBMAXTYPE,
415  true);
416 
417  CHKERR bit_mng->writeBitLevelByDim(
418  bit(l), BitRefLevel().set(), SPACE_DIM,
419  (boost::lexical_cast<std::string>(l) + "_ref_mesh.vtk").c_str(), "VTK",
420  "");
421  CHKERR bit_mng->writeBitLevelByDim(
422  bit(l), bit(l), MBTRI,
423  (boost::lexical_cast<std::string>(l) + "_only_ref_mesh.vtk").c_str(),
424  "VTK", "");
425 
427  };
428 
429  auto mark_skins = [&](auto l, auto m) {
431  Range ents;
433  ents);
434  Range level_skin;
435  CHKERR skin.find_skin(0, ents, false, level_skin);
436  CHKERR pcomm->filter_pstatus(level_skin,
437  PSTATUS_SHARED | PSTATUS_MULTISHARED,
438  PSTATUS_NOT, -1, nullptr);
439  level_skin = subtract(level_skin, level0_skin);
440  CHKERR mField.get_moab().get_adjacencies(level_skin, 0, false, level_skin,
441  moab::Interface::UNION);
442  CHKERR bit_mng->addBitRefLevel(level_skin, marker(m));
444  };
445 
446  BitRefLevel bit_sum;
447  for (auto l = 0; l != nb_ref_levels; ++l) {
448  CHKERR refine_mesh(l + 1);
449  CHKERR mark_skins(l, l + 1);
450  CHKERR mark_skins(l + 1, l + 1);
451  bit_sum |= bit(l);
452  }
453  bit_sum |= bit(nb_ref_levels);
454 
455  simpleInterface->getBitRefLevel() = bit_sum;
456  simpleInterface->getBitRefLevelMask() = BitRefLevel().set();
457 
459 }
460 //! [Read mesh]
461 
462 //! [Set up problem]
465  // Add field
466  CHKERR simpleInterface->addDomainField(FIELD_NAME, H1,
468  CHKERR simpleInterface->addBoundaryField(FIELD_NAME, H1,
470 
471  constexpr int order = 3;
472  CHKERR simpleInterface->setFieldOrder(FIELD_NAME, order);
473 
474  // Simple interface will resolve adjacency to DOFs of parent of the element.
475  // Using that information MAtrixManager allocate appropriately size of
476  // matrix.
477  simpleInterface->getParentAdjacencies() = true;
478 
479  CHKERR simpleInterface->setUp();
480 
481  ProblemsManager *prb_mng = mField.getInterface<ProblemsManager>();
482 
483  // remove obsolete DOFs from problem
484 
485  for (int l = 0; l != nb_ref_levels; ++l) {
486  CHKERR prb_mng->removeDofsOnEntities(simpleInterface->getProblemName(),
487  FIELD_NAME, bit(l), bit(l));
488  CHKERR prb_mng->removeDofsOnEntities(simpleInterface->getProblemName(),
489  FIELD_NAME, marker(l + 1),
490  bit(l).flip());
491  }
493 }
494 //! [Set up problem]
495 
496 //! [Push operators to pipeline]
499  PipelineManager *pipeline_mng = mField.getInterface<PipelineManager>();
500 
501  auto rule = [](int, int, int p) -> int { return 2 * p + 1; };
502 
503  CHKERR pipeline_mng->setDomainLhsIntegrationRule(rule);
504  CHKERR pipeline_mng->setDomainRhsIntegrationRule(rule);
505 
506  pipeline_mng->getDomainLhsFE()->exeTestHook = test_bit_child;
507  pipeline_mng->getDomainRhsFE()->exeTestHook = test_bit_child;
508 
509  auto beta = [](const double, const double, const double) { return 1; };
510  set_parent_dofs<DomainParentEle>(mField, pipeline_mng->getDomainLhsFE(),
511  DomainEleOp::OPSPACE, QUIET, Sev::noisy);
512  set_parent_dofs<DomainParentEle>(mField, pipeline_mng->getDomainLhsFE(),
513  DomainEleOp::OPROW, QUIET, Sev::noisy);
514  set_parent_dofs<DomainParentEle>(mField, pipeline_mng->getDomainLhsFE(),
515  DomainEleOp::OPCOL, QUIET, Sev::noisy);
516  pipeline_mng->getOpDomainLhsPipeline().push_back(
517  new OpDomainMass(FIELD_NAME, FIELD_NAME, beta));
518 
519  auto field_op_row = new ForcesAndSourcesCore::UserDataOperator(
520  FIELD_NAME, DomainEleOp::OPROW);
521  field_op_row->doWorkRhsHook = [](DataOperator *op_ptr, int side,
522  EntityType type,
525  if (type == MBENTITYSET) {
526 
527  MOFEM_LOG("SELF", Sev::verbose)
528  << "ROW: side/type: " << side << "/" << CN::EntityTypeName(type)
529  << " op space/base: " << FieldSpaceNames[data.getSpace()] << "/"
530  << ApproximationBaseNames[data.getBase()] << " DOFs "
531  << data.getIndices() << " nb base functions " << data.getN().size2()
532  << " nb base functions integration points " << data.getN().size1();
533 
534  for (auto &field_ent : data.getFieldEntities()) {
535  MOFEM_LOG("SELF", Sev::verbose)
536  << "\t" << CN::EntityTypeName(field_ent->getEntType());
537  }
538  }
540  };
541 
542  set_parent_dofs<DomainParentEle>(mField, pipeline_mng->getDomainRhsFE(),
543  DomainEleOp::OPSPACE, NOISY, Sev::verbose);
544  set_parent_dofs<DomainParentEle>(mField, pipeline_mng->getDomainRhsFE(),
545  DomainEleOp::OPROW, NOISY, Sev::noisy);
546  pipeline_mng->getOpDomainRhsPipeline().push_back(field_op_row);
547 
548  pipeline_mng->getOpDomainRhsPipeline().push_back(
549  new OpDomainSource(FIELD_NAME, approxFunction));
550 
552 }
553 //! [Push operators to pipeline]
554 
555 //! [Solve]
558  PipelineManager *pipeline_mng = mField.getInterface<PipelineManager>();
559 
560  MOFEM_LOG("WORLD", Sev::inform) << "Solve problem";
561 
562  auto solver = pipeline_mng->createKSP();
563  CHKERR KSPSetFromOptions(solver);
564  CHKERR KSPSetUp(solver);
565 
566  auto dm = simpleInterface->getDM();
567  auto D = createDMVector(dm);
568  auto F = vectorDuplicate(D);
569 
570  CHKERR KSPSolve(solver, F, D);
571  CHKERR VecGhostUpdateBegin(D, INSERT_VALUES, SCATTER_FORWARD);
572  CHKERR VecGhostUpdateEnd(D, INSERT_VALUES, SCATTER_FORWARD);
573  CHKERR DMoFEMMeshToLocalVector(dm, D, INSERT_VALUES, SCATTER_REVERSE);
575 }
576 
577 //! [Check results]
580  PipelineManager *pipeline_mng = mField.getInterface<PipelineManager>();
581  pipeline_mng->getDomainLhsFE().reset();
582  pipeline_mng->getDomainRhsFE().reset();
583  pipeline_mng->getBoundaryRhsFE().reset();
584 
585  auto rule = [](int, int, int p) -> int { return 2 * p + 1; };
586  CHKERR pipeline_mng->setDomainRhsIntegrationRule(rule);
587  CHKERR pipeline_mng->setBoundaryRhsIntegrationRule(rule);
588  pipeline_mng->getDomainRhsFE()->exeTestHook = test_bit_child;
589  pipeline_mng->getBoundaryRhsFE()->exeTestHook = test_bit_child;
590 
591  auto common_data_ptr = boost::make_shared<CommonData>();
592  common_data_ptr->resVec = createDMVector(simpleInterface->getDM());
593  common_data_ptr->L2Vec = createVectorMPI(
594  mField.get_comm(), (!mField.get_comm_rank()) ? 1 : 0, 1);
595  common_data_ptr->approxVals = boost::make_shared<VectorDouble>();
596  common_data_ptr->divApproxVals = boost::make_shared<MatrixDouble>();
597 
599  pipeline_mng->getOpDomainRhsPipeline(), {H1});
600  set_parent_dofs<DomainParentEle>(mField, pipeline_mng->getDomainRhsFE(),
601  DomainEleOp::OPSPACE, QUIET, Sev::noisy);
602  set_parent_dofs<DomainParentEle>(mField, pipeline_mng->getDomainRhsFE(),
603  DomainEleOp::OPROW, VERBOSE, Sev::noisy);
604 
605  pipeline_mng->getOpDomainRhsPipeline().push_back(
607  FIELD_NAME, common_data_ptr->divApproxVals));
608  pipeline_mng->getOpDomainRhsPipeline().push_back(
610  common_data_ptr->approxVals));
611 
612  pipeline_mng->getOpDomainRhsPipeline().push_back(
613  new OpError<FIELD_DIM>(common_data_ptr));
614 
615  set_parent_dofs<BoundaryParentEle>(mField, pipeline_mng->getBoundaryRhsFE(),
616  BoundaryEleOp::OPSPACE, QUIET, Sev::noisy);
617  set_parent_dofs<BoundaryParentEle>(mField, pipeline_mng->getBoundaryRhsFE(),
618  BoundaryEleOp::OPROW, VERBOSE, Sev::noisy);
619  pipeline_mng->getOpBoundaryRhsPipeline().push_back(
621  common_data_ptr->approxVals));
622  pipeline_mng->getOpBoundaryRhsPipeline().push_back(
623  new OpErrorSkel<FIELD_DIM>(common_data_ptr));
624 
625  CHKERR VecZeroEntries(common_data_ptr->L2Vec);
626  CHKERR VecZeroEntries(common_data_ptr->resVec);
627 
628  CHKERR pipeline_mng->loopFiniteElements();
629 
630  CHKERR VecAssemblyBegin(common_data_ptr->L2Vec);
631  CHKERR VecAssemblyEnd(common_data_ptr->L2Vec);
632  CHKERR VecAssemblyBegin(common_data_ptr->resVec);
633  CHKERR VecAssemblyEnd(common_data_ptr->resVec);
634  double nrm2;
635  CHKERR VecNorm(common_data_ptr->resVec, NORM_2, &nrm2);
636  const double *array;
637  CHKERR VecGetArrayRead(common_data_ptr->L2Vec, &array);
638  MOFEM_LOG_C("WORLD", Sev::inform, "Error %6.4e Vec norm %6.4e\n",
639  std::sqrt(array[0]), nrm2);
640 
641  constexpr double eps = 1e-8;
642  if (nrm2 > eps)
643  SETERRQ1(PETSC_COMM_SELF, MOFEM_ATOM_TEST_INVALID,
644  "Not converged solution err = %6.4e", nrm2);
645  if (std::sqrt(array[0]) > eps)
646  SETERRQ1(PETSC_COMM_SELF, MOFEM_ATOM_TEST_INVALID,
647  "Error in approximation err = %6.4e", std::sqrt(array[0]));
648 
649  CHKERR VecRestoreArrayRead(common_data_ptr->L2Vec, &array);
650 
652 }
653 //! [Check results]
654 
657  PipelineManager *pipeline_mng = mField.getInterface<PipelineManager>();
658  pipeline_mng->getDomainLhsFE().reset();
659  pipeline_mng->getDomainRhsFE().reset();
660 
661  auto rule = [](int, int, int p) -> int { return -1; };
662  CHKERR pipeline_mng->setDomainRhsIntegrationRule(rule);
663 
664  static_cast<ForcesAndSourcesCore *>(pipeline_mng->getDomainRhsFE().get())
665  ->setRuleHook = [](ForcesAndSourcesCore *fe_raw_ptr, int order_row,
666  int order_col, int order_data) -> MoFEMErrorCode {
668  fe_raw_ptr->gaussPts.resize(3, 3);
669  fe_raw_ptr->gaussPts(0, 0) = 0;
670  fe_raw_ptr->gaussPts(1, 0) = 0;
671  fe_raw_ptr->gaussPts(2, 0) = 0;
672  fe_raw_ptr->gaussPts(0, 1) = 1;
673  fe_raw_ptr->gaussPts(1, 1) = 0;
674  fe_raw_ptr->gaussPts(2, 1) = 0;
675  fe_raw_ptr->gaussPts(0, 2) = 0;
676  fe_raw_ptr->gaussPts(1, 2) = 1;
677  fe_raw_ptr->gaussPts(2, 2) = 0;
679  };
680 
681  auto field_op_row = new ForcesAndSourcesCore::UserDataOperator(
682  FIELD_NAME, DomainEleOp::OPROW);
683 
684  auto approx_vals = boost::make_shared<VectorDouble>();
685 
686  auto &moab = mField.get_moab();
687  Tag th;
688  double def_val[] = {0};
689  CHKERR moab.tag_get_handle("FIELD", 1, MB_TYPE_DOUBLE, th,
690  MB_TAG_CREAT | MB_TAG_SPARSE, &def_val);
691 
692  field_op_row->doWorkRhsHook = [&](DataOperator *base_op_ptr, int side,
693  EntityType type,
696  if (type == MBVERTEX) {
697  auto op_ptr =
699  base_op_ptr);
700  auto t_field = getFTensor0FromVec(*approx_vals);
701  auto nb_gauss_pts = op_ptr->getGaussPts().size2();
702  if (nb_gauss_pts != 3)
703  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
704  "Should be three guass pts.");
705  auto conn = op_ptr->getConn();
706  for (auto gg = 0; gg != nb_gauss_pts; ++gg) {
707  const double v = t_field;
708  CHKERR moab.tag_set_data(th, &conn[gg], 1, &v);
709  ++t_field;
710  }
711  }
713  };
714 
715  set_parent_dofs<DomainParentEle>(mField, pipeline_mng->getDomainRhsFE(),
716  DomainEleOp::OPSPACE, VERBOSE, Sev::noisy);
717  set_parent_dofs<DomainParentEle>(mField, pipeline_mng->getDomainRhsFE(),
718  DomainEleOp::OPROW, VERBOSE, Sev::noisy);
719  pipeline_mng->getOpDomainRhsPipeline().push_back(
720  new OpCalculateScalarFieldValues(FIELD_NAME, approx_vals));
721  pipeline_mng->getOpDomainRhsPipeline().push_back(field_op_row);
722  CHKERR pipeline_mng->loopFiniteElements();
723 
724  CHKERR mField.getInterface<BitRefManager>()->writeBitLevelByType(
725  bit(nb_ref_levels), BitRefLevel().set(), MBTRI, "out.vtk", "VTK", "");
726 
728 }
729 
730 int main(int argc, char *argv[]) {
731 
732  // Initialisation of MoFEM/PETSc and MOAB data structures
733  MoFEM::Core::Initialize(&argc, &argv, NULL, help);
734 
735  try {
736 
737  //! [Register MoFEM discrete manager in PETSc]
738  DMType dm_name = "DMMOFEM";
739  CHKERR DMRegister_MoFEM(dm_name);
740  //! [Register MoFEM discrete manager in PETSc
741 
742  //! [Create MoAB]
743  moab::Core mb_instance; ///< mesh database
744  moab::Interface &moab = mb_instance; ///< mesh database interface
745  //! [Create MoAB]
746 
747  //! [Create MoFEM]
748  MoFEM::Core core(moab); ///< finite element database
749  MoFEM::Interface &m_field = core; ///< finite element database insterface
750  //! [Create MoFEM]
751 
752  //! [AtomTest]
753  AtomTest ex(m_field);
754  CHKERR ex.runProblem();
755  //! [AtomTest]
756  }
757  CATCH_ERRORS;
758 
760 }
AtomTest
Definition: child_and_parent.cpp:57
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
OpDomainSource
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::LinearForm< GAUSS >::OpSource< 1, FIELD_DIM > OpDomainSource
evaluate source, i.e. rhs vector
Definition: hanging_node_approx.cpp:69
MoFEM::EntitiesFieldData::EntData
Data on single entity (This is passed as argument to DataOperator::doWork)
Definition: EntitiesFieldData.hpp:128
BoundaryEle
ElementsAndOps< SPACE_DIM >::BoundaryEle BoundaryEle
Definition: child_and_parent.cpp:39
MoFEM::PipelineManager::getDomainRhsFE
boost::shared_ptr< FEMethod > & getDomainRhsFE()
Definition: PipelineManager.hpp:405
MYPCOMM_INDEX
#define MYPCOMM_INDEX
default communicator number PCOMM
Definition: definitions.h:228
MoFEM::CoreTmp< 0 >
Core (interface) class.
Definition: Core.hpp:82
H1
@ H1
continuous field
Definition: definitions.h:85
OpError
Definition: initial_diffusion.cpp:61
MoFEM::DataOperator
base operator to do operations at Gauss Pt. level
Definition: DataOperators.hpp:24
FTensor::Tensor1
Definition: Tensor1_value.hpp:8
MoFEM::ProblemsManager
Problem manager is used to build and partition problems.
Definition: ProblemsManager.hpp:21
NOISY
@ NOISY
Definition: definitions.h:224
MoFEM::FEMethod
structure for User Loop Methods on finite elements
Definition: LoopMethods.hpp:369
MoFEM::PipelineManager::ElementsAndOpsByDim
Definition: PipelineManager.hpp:38
ApproximationBaseNames
const static char *const ApproximationBaseNames[]
Definition: definitions.h:72
MoFEM::PipelineManager::loopFiniteElements
MoFEMErrorCode loopFiniteElements(SmartPetscObj< DM > dm=nullptr)
Iterate finite elements.
Definition: PipelineManager.cpp:19
AtomTest::CommonData::divApproxVals
boost::shared_ptr< MatrixDouble > divApproxVals
Definition: hanging_node_approx.cpp:210
MoFEM::Exceptions::MoFEMErrorCode
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
Definition: Exceptions.hpp:56
EntData
EntitiesFieldData::EntData EntData
Definition: hanging_node_approx.cpp:29
MoFEM::PipelineManager::getOpDomainRhsPipeline
boost::ptr_deque< UserDataOperator > & getOpDomainRhsPipeline()
Get the Op Domain Rhs Pipeline object.
Definition: PipelineManager.hpp:773
MoFEM::th
Tag th
Definition: Projection10NodeCoordsOnField.cpp:122
MoFEM::ForcesAndSourcesCore::UserDataOperator::OpType
OpType
Controls loop over entities on element.
Definition: ForcesAndSourcesCore.hpp:566
MoFEM::PETSC
@ PETSC
Definition: FormsIntegrators.hpp:105
MoFEM::PipelineManager
PipelineManager interface.
Definition: PipelineManager.hpp:24
MoFEM.hpp
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
MoFEM::PipelineManager::getBoundaryRhsFE
boost::shared_ptr< FEMethod > & getBoundaryRhsFE()
Definition: PipelineManager.hpp:413
MoFEM::CoreTmp< 0 >::Finalize
static MoFEMErrorCode Finalize()
Checks for options to be called at the conclusion of the program.
Definition: Core.cpp:112
AtomTest::assembleSystem
MoFEMErrorCode assembleSystem()
[Push operators to pipeline]
Definition: child_and_parent.cpp:290
nb_ref_levels
constexpr int nb_ref_levels
Three levels of refinement.
Definition: hanging_node_approx.cpp:17
test_bit_child
auto test_bit_child
lambda function used to select elements on which finite element pipelines are executed.
Definition: hanging_node_approx.cpp:173
MoFEM::Simple
Simple interface for fast problem set-up.
Definition: Simple.hpp:27
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:1589
sdf.r
int r
Definition: sdf.py:8
SPACE_DIM
constexpr int SPACE_DIM
Definition: hanging_node_approx.cpp:16
MoFEM::BitRefManager::filterEntitiesByRefLevel
MoFEMErrorCode filterEntitiesByRefLevel(const BitRefLevel bit, const BitRefLevel mask, Range &ents, int verb=QUIET) const
filter entities by bit ref level
Definition: BitRefManager.cpp:746
order
constexpr int order
Definition: dg_projection.cpp:18
MoFEM::DeprecatedCoreInterface
Deprecated interface functions.
Definition: DeprecatedCoreInterface.hpp:16
MoFEM::OpCalculateScalarFieldGradient
Get field gradients at integration pts for scalar filed rank 0, i.e. vector field.
Definition: UserDataOperators.hpp:1293
help
static char help[]
Definition: hanging_node_approx.cpp:12
MoFEM::Interface
DeprecatedCoreInterface Interface
Definition: Interface.hpp:2010
MoFEM::EntitiesFieldData::EntData::getFTensor0N
FTensor::Tensor0< FTensor::PackPtr< double *, 1 > > getFTensor0N(const FieldApproximationBase base)
Get base function as Tensor0.
Definition: EntitiesFieldData.hpp:1502
AtomTest::OpError< 1 >::OpError
OpError(boost::shared_ptr< CommonData > &common_data_ptr)
Definition: hanging_node_approx.cpp:226
AtomTest::OpError
Operator to evaluate errors.
Definition: child_and_parent.cpp:82
VERBOSE
@ VERBOSE
Definition: definitions.h:222
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:548
AtomTest::OpError< 1 >::doWork
MoFEMErrorCode doWork(int side, EntityType type, EntData &data)
Definition: hanging_node_approx.cpp:228
MoFEM::createDMVector
auto createDMVector(DM dm)
Get smart vector from DM.
Definition: DMMoFEM.hpp:1099
MoFEM
implementation of Data Operators for Forces and Sources
Definition: Common.hpp:10
AtomTest::OpErrorSkel< 1 >::doWork
MoFEMErrorCode doWork(int side, EntityType type, EntData &data)
Definition: hanging_node_approx.cpp:285
BoundaryEleOp
MoFEM::BitRefManager::getEntitiesByDimAndRefLevel
MoFEMErrorCode getEntitiesByDimAndRefLevel(const BitRefLevel bit, const BitRefLevel mask, const int dim, const EntityHandle meshset, int verb=0) const
add all ents from ref level given by bit to meshset
Definition: BitRefManager.cpp:822
MoFEM::FaceElementForcesAndSourcesCore::UserDataOperator
default operator for TRI element
Definition: FaceElementForcesAndSourcesCore.hpp:94
MoFEM::MeshRefinement
Mesh refinement interface.
Definition: MeshRefinement.hpp:26
MOFEM_LOG_C
#define MOFEM_LOG_C(channel, severity, format,...)
Definition: LogManager.hpp:311
MoFEM::get_temp_meshset_ptr
auto get_temp_meshset_ptr(moab::Interface &moab)
Create smart pointer to temporary meshset.
Definition: Templates.hpp:1886
FIELD_DIM
constexpr int FIELD_DIM
Definition: hanging_node_approx.cpp:15
double
bit
auto bit
set bit
Definition: hanging_node_approx.cpp:75
MoFEM::ProblemsManager::removeDofsOnEntities
MoFEMErrorCode removeDofsOnEntities(const std::string problem_name, const std::string field_name, const Range ents, const int lo_coeff=0, const int hi_coeff=MAX_DOFS_ON_ENTITY, const int lo_order=0, const int hi_order=100, int verb=VERBOSE, const bool debug=false)
Remove DOFs from problem.
Definition: ProblemsManager.cpp:2954
convert.type
type
Definition: convert.py:64
ApproxFieldFunction< 1 >::operator()
auto operator()(const double x, const double y, const double z)
Definition: hanging_node_approx.cpp:39
MoFEM::FormsIntegrators::Assembly
Assembly methods.
Definition: FormsIntegrators.hpp:317
MoFEM::getFTensor0FromVec
static auto getFTensor0FromVec(ublas::vector< T, A > &data)
Get tensor rank 0 (scalar) form data vector.
Definition: Templates.hpp:135
DomainParentEle
ElementsAndOps< SPACE_DIM >::DomianParentEle DomainParentEle
Definition: child_and_parent.cpp:35
AtomTest::divApproxFunction
static ApproxFieldFunctionDerivative< FIELD_DIM > divApproxFunction
Definition: hanging_node_approx.cpp:189
MoFEM::OpCalculateScalarFieldValues
Get value at integration points for scalar field.
Definition: UserDataOperators.hpp:82
MoFEM::PipelineManager::createKSP
SmartPetscObj< KSP > createKSP(SmartPetscObj< DM > dm=nullptr)
Create KSP (linear) solver.
Definition: PipelineManager.cpp:49
MoFEM::EntitiesFieldData::EntData::getIndices
const VectorInt & getIndices() const
Get global indices of dofs on entity.
Definition: EntitiesFieldData.hpp:1214
MoFEM::BitRefManager::updateMeshsetByEntitiesChildren
MoFEMErrorCode updateMeshsetByEntitiesChildren(const EntityHandle parent, const BitRefLevel &parent_bit, const BitRefLevel &parent_mask, const BitRefLevel &child_bit, const BitRefLevel &child_mask, const EntityHandle child, EntityType child_type, const bool recursive=false, int verb=0)
Get child entities form meshset containing parent entities.
Definition: BitRefManager.cpp:1029
MoFEM::DMRegister_MoFEM
PetscErrorCode DMRegister_MoFEM(const char sname[])
Register MoFEM problem.
Definition: DMMoFEM.cpp:43
MoFEM::BitRefManager::addBitRefLevel
MoFEMErrorCode addBitRefLevel(const Range &ents, const BitRefLevel &bit, int verb=QUIET) const
add bit ref level to ref entity
Definition: BitRefManager.cpp:554
AtomTest::OpErrorSkel< 1 >::commonDataPtr
boost::shared_ptr< CommonData > commonDataPtr
Definition: hanging_node_approx.cpp:282
AtomTest::OpErrorSkel
Definition: hanging_node_approx.cpp:217
ApproxFieldFunctionDerivative< 1 >::operator()
auto operator()(const double x, const double y, const double z)
Definition: hanging_node_approx.cpp:49
MoFEM::BitRefManager::writeBitLevelByDim
MoFEMErrorCode writeBitLevelByDim(const BitRefLevel bit, const BitRefLevel mask, const int dim, const char *file_name, const char *file_type, const char *options, const bool check_for_empty=true) const
Write bit ref level to file.
Definition: BitRefManager.cpp:660
MoFEM::PipelineManager::setDomainRhsIntegrationRule
MoFEMErrorCode setDomainRhsIntegrationRule(RuleHookFun rule)
Definition: PipelineManager.hpp:530
AtomTest::CommonData
Collected data use d by operator to evaluate errors for the test.
Definition: child_and_parent.cpp:76
AtomTest::AtomTest
AtomTest(MoFEM::Interface &m_field)
Definition: hanging_node_approx.cpp:180
MoFEM::LogManager::SeverityLevel
SeverityLevel
Severity levels.
Definition: LogManager.hpp:33
AtomTest::approxFunction
static ApproxFieldFunction< FIELD_DIM > approxFunction
Definition: child_and_parent.cpp:67
MoFEM::PipelineManager::getOpDomainLhsPipeline
boost::ptr_deque< UserDataOperator > & getOpDomainLhsPipeline()
Get the Op Domain Lhs Pipeline object.
Definition: PipelineManager.hpp:749
t
constexpr double t
plate stiffness
Definition: plate.cpp:58
MoFEM::ForcesAndSourcesCore::UserDataOperator
friend class UserDataOperator
Definition: ForcesAndSourcesCore.hpp:482
i
FTensor::Index< 'i', SPACE_DIM > i
Definition: hcurl_divergence_operator_2d.cpp:27
BiLinearForm
FTensor::Index< 'i', SPACE_DIM >
AtomTest::solveSystem
MoFEMErrorCode solveSystem()
[Push operators to pipeline]
Definition: child_and_parent.cpp:377
MoFEM::AddHOOps
Add operators pushing bases from local to physical configuration.
Definition: HODataOperators.hpp:413
MoFEM::PipelineManager::setDomainLhsIntegrationRule
MoFEMErrorCode setDomainLhsIntegrationRule(RuleHookFun rule)
Definition: PipelineManager.hpp:503
MoFEM::BitRefManager::updateRangeByChildren
MoFEMErrorCode updateRangeByChildren(const Range &parent, Range &child, MoFEMTypes bh=MF_ZERO)
Update range by childrens.
Definition: BitRefManager.cpp:1144
main
int main(int argc, char *argv[])
Definition: hanging_node_approx.cpp:730
MoFEM::ForcesAndSourcesCore
structure to get information form mofem into EntitiesFieldData
Definition: ForcesAndSourcesCore.hpp:22
v
const double v
phase velocity of light in medium (cm/ns)
Definition: initial_diffusion.cpp:40
Range
MoFEM::PipelineManager::getDomainLhsFE
boost::shared_ptr< FEMethod > & getDomainLhsFE()
Definition: PipelineManager.hpp:401
DomainEleOp
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
FieldSpaceNames
const static char *const FieldSpaceNames[]
Definition: definitions.h:92
MOFEM_LOG
#define MOFEM_LOG(channel, severity)
Log.
Definition: LogManager.hpp:308
MoFEM::vectorDuplicate
SmartPetscObj< Vec > vectorDuplicate(Vec vec)
Create duplicate vector of smart vector.
Definition: PetscSmartObj.hpp:221
CATCH_ERRORS
#define CATCH_ERRORS
Catch errors.
Definition: definitions.h:385
MoFEM::Core
CoreTmp< 0 > Core
Definition: Core.hpp:1148
ApproxFieldFunctionDerivative
Definition: hanging_node_approx.cpp:32
UserDataOperator
ForcesAndSourcesCore::UserDataOperator UserDataOperator
Definition: HookeElement.hpp:75
eps
static const double eps
Definition: check_base_functions_derivatives_on_tet.cpp:11
BITREFLEVEL_SIZE
#define BITREFLEVEL_SIZE
max number of refinements
Definition: definitions.h:232
AtomTest::OpErrorSkel< 1 >::OpErrorSkel
OpErrorSkel(boost::shared_ptr< CommonData > &common_data_ptr)
Definition: hanging_node_approx.cpp:283
AINSWORTH_LEGENDRE_BASE
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base .
Definition: definitions.h:60
MOFEM_DATA_INCONSISTENCY
@ MOFEM_DATA_INCONSISTENCY
Definition: definitions.h:31
MoFEM::ForcesAndSourcesCore::gaussPts
MatrixDouble gaussPts
Matrix of integration points.
Definition: ForcesAndSourcesCore.hpp:109
AtomTest::checkResults
MoFEMErrorCode checkResults()
[Check results]
Definition: dg_projection.cpp:213
MoFEM::createVectorMPI
auto createVectorMPI(MPI_Comm comm, PetscInt n, PetscInt N)
Create MPI Vector.
Definition: PetscSmartObj.hpp:202
OpDomainMass
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::BiLinearForm< GAUSS >::OpMass< 1, FIELD_DIM > OpDomainMass
evaluate mass matrix
Definition: hanging_node_approx.cpp:62
MoFEM::PipelineManager::getOpBoundaryRhsPipeline
boost::ptr_deque< UserDataOperator > & getOpBoundaryRhsPipeline()
Get the Op Boundary Rhs Pipeline object.
Definition: PipelineManager.hpp:821
MoFEM::Types::VectorDouble
UBlasVector< double > VectorDouble
Definition: Types.hpp:68
AtomTest::readMesh
MoFEMErrorCode readMesh()
[Run programme]
Definition: child_and_parent.cpp:208
ReactionDiffusionEquation::D
const double D
diffusivity
Definition: reaction_diffusion.cpp:20
MoFEM::PipelineManager::setBoundaryRhsIntegrationRule
MoFEMErrorCode setBoundaryRhsIntegrationRule(RuleHookFun rule)
Definition: PipelineManager.hpp:584
m
FTensor::Index< 'm', 3 > m
Definition: shallow_wave.cpp:80
MoFEM::BitRefManager
Managing BitRefLevels.
Definition: BitRefManager.hpp:21
MoFEM::Types::BitRefLevel
std::bitset< BITREFLEVEL_SIZE > BitRefLevel
Bit structure attached to each entity identifying to what mesh entity is attached.
Definition: Types.hpp:40
MOFEM_ATOM_TEST_INVALID
@ MOFEM_ATOM_TEST_INVALID
Definition: definitions.h:40
marker
auto marker
set bit to marker
Definition: hanging_node_approx.cpp:82
BoundaryParentEle
ElementsAndOps< SPACE_DIM >::BoundaryParentEle BoundaryParentEle
Definition: child_and_parent.cpp:41
MoFEMFunctionBeginHot
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:453
QUIET
@ QUIET
Definition: definitions.h:221
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
set_parent_dofs
auto set_parent_dofs(MoFEM::Interface &m_field, boost::shared_ptr< FEMethod > &fe_top, ForcesAndSourcesCore::UserDataOperator::OpType op, int verbosity, LogManager::SeverityLevel sev)
set levels of projection operators, which project field data from parent entities,...
Definition: hanging_node_approx.cpp:92
convert.int
int
Definition: convert.py:64
MoFEMFunctionReturn
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:429
MoFEM::OpAddParentEntData
Operator to project base functions from parent entity to child.
Definition: MeshProjectionDataOperators.hpp:66
AtomTest::printResults
MoFEMErrorCode printResults()
[Check results]
Definition: hanging_node_approx.cpp:655
AtomTest::setupProblem
MoFEMErrorCode setupProblem()
[Read mesh]
Definition: child_and_parent.cpp:268
AtomTest::runProblem
MoFEMErrorCode runProblem()
[Run programme]
Definition: child_and_parent.cpp:188
MoFEMFunctionBegin
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:359
l
FTensor::Index< 'l', 3 > l
Definition: matrix_function.cpp:21
FIELD_NAME
constexpr char FIELD_NAME[]
Definition: hanging_node_approx.cpp:14
ApproxFieldFunction
Definition: child_and_parent.cpp:43
F
@ F
Definition: free_surface.cpp:394