v0.14.0
plot_base.cpp
Go to the documentation of this file.
1 /**
2  * \file plot_base.cpp
3  * \example plot_base.cpp
4  *
5  * Utility for plotting base functions for different spaces, polynomial bases
6  */
7 
8 #include <MoFEM.hpp>
9 
10 using namespace MoFEM;
11 
12 static char help[] = "...\n\n";
13 
14 
15 
16 template <int DIM> struct ElementsAndOps {};
17 
18 template <> struct ElementsAndOps<2> {
20 };
21 
22 template <> struct ElementsAndOps<3> {
24 };
25 
26 constexpr int SPACE_DIM =
27  EXECUTABLE_DIMENSION; //< Space dimension of problem, mesh
28 
33 
34 struct MyPostProc : public PostProcEle {
36 
37  MoFEMErrorCode generateReferenceElementMesh();
38  MoFEMErrorCode setGaussPts(int order);
39 
40  MoFEMErrorCode preProcess();
41  MoFEMErrorCode postProcess();
42 
43 protected:
44  ublas::matrix<int> refEleMap;
46 };
47 
48 struct Example {
49 
50  Example(MoFEM::Interface &m_field) : mField(m_field) {}
51 
52  MoFEMErrorCode runProblem();
53 
54 private:
55  MoFEM::Interface &mField;
56  Simple *simpleInterface;
57 
58  MoFEMErrorCode readMesh();
59  MoFEMErrorCode setupProblem();
60  MoFEMErrorCode setIntegrationRules();
61  MoFEMErrorCode createCommonData();
62  MoFEMErrorCode boundaryCondition();
63  MoFEMErrorCode assembleSystem();
64  MoFEMErrorCode solveSystem();
65  MoFEMErrorCode outputResults();
66  MoFEMErrorCode checkResults();
67 
70 };
71 
72 //! [Run programme]
75  CHKERR readMesh();
76  CHKERR setupProblem();
77  CHKERR setIntegrationRules();
78  CHKERR createCommonData();
79  CHKERR boundaryCondition();
80  CHKERR assembleSystem();
81  CHKERR solveSystem();
82  CHKERR outputResults();
83  CHKERR checkResults();
85 }
86 //! [Run programme]
87 
88 //! [Read mesh]
91 
92  PetscBool load_file = PETSC_FALSE;
93  CHKERR PetscOptionsGetBool(PETSC_NULL, "", "-load_file", &load_file,
94  PETSC_NULL);
95 
96  if (load_file == PETSC_FALSE) {
97 
98  auto &moab = mField.get_moab();
99 
100  if (SPACE_DIM == 3) {
101 
102  // create one tet
103  double tet_coords[] = {0, 0, 0, 1, 0, 0, 0, 1, 0, 0, 0, 1};
104  EntityHandle nodes[4];
105  for (int nn = 0; nn < 4; nn++) {
106  CHKERR moab.create_vertex(&tet_coords[3 * nn], nodes[nn]);
107  }
108  EntityHandle tet;
109  CHKERR moab.create_element(MBTET, nodes, 4, tet);
110  Range adj;
111  for (auto d : {1, 2})
112  CHKERR moab.get_adjacencies(&tet, 1, d, true, adj);
113  }
114 
115  if (SPACE_DIM == 2) {
116 
117  // create one triangle
118  double tri_coords[] = {0, 0, 0, 1, 0, 0, 0, 1, 0};
119  EntityHandle nodes[3];
120  for (int nn = 0; nn < 3; nn++) {
121  CHKERR moab.create_vertex(&tri_coords[3 * nn], nodes[nn]);
122  }
123  EntityHandle tri;
124  CHKERR moab.create_element(MBTRI, nodes, 3, tri);
125  Range adj;
126  CHKERR moab.get_adjacencies(&tri, 1, 1, true, adj);
127  }
128 
129  CHKERR mField.rebuild_database();
130  CHKERR mField.getInterface(simpleInterface);
131  simpleInterface->setDim(SPACE_DIM);
132 
133  // Add all elements to database
134  CHKERR mField.getInterface<BitRefManager>()->setBitRefLevelByDim(
135  0, SPACE_DIM, simpleInterface->getBitRefLevel());
136 
137  } else {
138 
139  CHKERR mField.getInterface(simpleInterface);
140  CHKERR simpleInterface->getOptions();
141  CHKERR simpleInterface->loadFile();
142  }
143 
145 }
146 //! [Read mesh]
147 
148 //! [Set up problem]
151  // Add field
152 
153  // Declare elements
154  enum bases { AINSWORTH, AINSWORTH_LOBATTO, DEMKOWICZ, BERNSTEIN, LASBASETOP };
155  const char *list_bases[] = {"ainsworth", "ainsworth_lobatto", "demkowicz",
156  "bernstein"};
157 
158  PetscBool flg;
159  PetscInt choice_base_value = AINSWORTH;
160  CHKERR PetscOptionsGetEList(PETSC_NULL, NULL, "-base", list_bases, LASBASETOP,
161  &choice_base_value, &flg);
162  if (flg != PETSC_TRUE)
163  SETERRQ(PETSC_COMM_SELF, MOFEM_IMPOSSIBLE_CASE, "base not set");
165  if (choice_base_value == AINSWORTH)
167  if (choice_base_value == AINSWORTH_LOBATTO)
168  base = AINSWORTH_LOBATTO_BASE;
169  else if (choice_base_value == DEMKOWICZ)
170  base = DEMKOWICZ_JACOBI_BASE;
171  else if (choice_base_value == BERNSTEIN)
173 
174  const char *list_continuity[] = {"continuous", "discontinuous"};
175  PetscInt choice_continuity_value = CONTINUOUS;
176  CHKERR PetscOptionsGetEList(PETSC_NULL, NULL, "-continuity", list_continuity,
177  LASTCONTINUITY, &choice_continuity_value, &flg);
178 
179  FieldContinuity continuity;
180  if (choice_continuity_value == CONTINUOUS)
181  continuity = CONTINUOUS;
182  else if (choice_continuity_value == DISCONTINUOUS)
183  continuity = DISCONTINUOUS;
184  else
185  SETERRQ(PETSC_COMM_SELF, MOFEM_IMPOSSIBLE_CASE, "continuity not set");
186 
187  enum spaces { H1SPACE, L2SPACE, HCURLSPACE, HDIVSPACE, LASBASETSPACE };
188  const char *list_spaces[] = {"h1", "l2", "hcurl", "hdiv"};
189  PetscInt choice_space_value = H1SPACE;
190  CHKERR PetscOptionsGetEList(PETSC_NULL, NULL, "-space", list_spaces,
191  LASBASETSPACE, &choice_space_value, &flg);
192  if (flg != PETSC_TRUE)
193  SETERRQ(PETSC_COMM_SELF, MOFEM_IMPOSSIBLE_CASE, "space not set");
194  space = H1;
195  if (choice_space_value == H1SPACE)
196  space = H1;
197  else if (choice_space_value == L2SPACE)
198  space = L2;
199  else if (choice_space_value == HCURLSPACE)
200  space = HCURL;
201  else if (choice_space_value == HDIVSPACE)
202  space = HDIV;
203 
204  AinsworthOrderHooks::broken_nbfacetri_edge_hdiv = [](int p) { return p; };
205  AinsworthOrderHooks::broken_nbfacetri_face_hdiv = [](int p) { return p; };
206  AinsworthOrderHooks::broken_nbvolumetet_edge_hdiv = [](int p) { return p; };
207  AinsworthOrderHooks::broken_nbvolumetet_face_hdiv = [](int p) { return p; };
208  AinsworthOrderHooks::broken_nbvolumetet_volume_hdiv = [](int p) { return p; };
209 
210  if(continuity == CONTINUOUS)
211  CHKERR simpleInterface->addDomainField("U", space, base, 1);
212  else
213  CHKERR simpleInterface->addDomainBrokenField("U", space, base, 1);
214 
215  int order = 3;
216  CHKERR PetscOptionsGetInt(PETSC_NULL, "", "-order", &order, PETSC_NULL);
217  CHKERR simpleInterface->setFieldOrder("U", order);
218  CHKERR simpleInterface->setUp();
219 
220  auto bc_mng = mField.getInterface<BcManager>();
221  CHKERR bc_mng->removeSideDOFs(simpleInterface->getProblemName(), "ZERO_FLUX",
222  "U", SPACE_DIM, 0, 1, true);
223 
225 }
226 //! [Set up problem]
227 
228 //! [Set integration rule]
232 }
233 //! [Set integration rule]
234 
235 //! [Create common data]
237 //! [Create common data]
238 
239 //! [Boundary condition]
241 //! [Boundary condition]
242 
243 //! [Push operators to pipeline]
245 //! [Push operators to pipeline]
246 
247 //! [Solve]
248 MoFEMErrorCode Example::solveSystem() { return 0; }
249 
250 //! [Solve]
253 
254  PipelineManager *pipeline_mng = mField.getInterface<PipelineManager>();
255  pipeline_mng->getDomainLhsFE().reset();
256 
257  auto post_proc_fe = boost::make_shared<MyPostProc>(mField);
258  post_proc_fe->generateReferenceElementMesh();
259  pipeline_mng->getDomainRhsFE() = post_proc_fe;
260 
261  if (SPACE_DIM == 2) {
262  if (space == HCURL) {
263  auto jac_ptr = boost::make_shared<MatrixDouble>();
264  post_proc_fe->getOpPtrVector().push_back(
265  new OpCalculateHOJacForFace(jac_ptr));
266  post_proc_fe->getOpPtrVector().push_back(new OpMakeHdivFromHcurl());
267  post_proc_fe->getOpPtrVector().push_back(
269  }
270  }
271 
272  switch (space) {
273  case H1:
274  case L2:
275 
276  {
277 
279 
280  auto u_ptr = boost::make_shared<VectorDouble>();
281  post_proc_fe->getOpPtrVector().push_back(
282  new OpCalculateScalarFieldValues("U", u_ptr));
283  post_proc_fe->getOpPtrVector().push_back(
284 
285  new OpPPMap(
286 
287  post_proc_fe->getPostProcMesh(), post_proc_fe->getMapGaussPts(),
288 
289  {{"U", u_ptr}},
290 
291  {},
292 
293  {},
294 
295  {}
296 
297  )
298 
299  );
300  } break;
301  case HCURL:
302  case HDIV:
303 
304  {
306 
308  post_proc_fe->getOpPtrVector(), {space});
309  auto u_ptr = boost::make_shared<MatrixDouble>();
310  post_proc_fe->getOpPtrVector().push_back(
311  new OpCalculateHVecVectorField<3>("U", u_ptr));
312 
313  post_proc_fe->getOpPtrVector().push_back(
314 
315  new OpPPMap(
316 
317  post_proc_fe->getPostProcMesh(), post_proc_fe->getMapGaussPts(),
318 
319  {},
320 
321  {{"U", u_ptr}},
322 
323  {},
324 
325  {}
326 
327  )
328 
329  );
330  } break;
331  default:
332  break;
333  }
334 
335  auto scale_tag_val = [&]() {
337  auto &post_proc_mesh = post_proc_fe->getPostProcMesh();
338  Range nodes;
339  CHKERR post_proc_mesh.get_entities_by_type(0, MBVERTEX, nodes);
340  Tag th;
341  CHKERR post_proc_mesh.tag_get_handle("U", th);
342  int length;
343  CHKERR post_proc_mesh.tag_get_length(th, length);
344  std::vector<double> data(nodes.size() * length);
345  CHKERR post_proc_mesh.tag_get_data(th, nodes, &*data.begin());
346  double max_v = 0;
347  for (int i = 0; i != nodes.size(); ++i) {
348  double v = 0;
349  for (int d = 0; d != length; ++d)
350  v += pow(data[length * i + d], 2);
351  v = std::sqrt(v);
352  max_v = std::max(max_v, v);
353  }
354  for (auto &v : data)
355  v /= max_v;
356  CHKERR post_proc_mesh.tag_set_data(th, nodes, &*data.begin());
358  };
359 
360  auto prb_ptr = mField.get_problem(simpleInterface->getProblemName());
361 
362  size_t nb = 0;
363  auto dofs_ptr = prb_ptr->getNumeredRowDofsPtr();
364 
365  for (auto dof_ptr : (*dofs_ptr)) {
366  MOFEM_LOG("PLOTBASE", Sev::verbose) << *dof_ptr;
367  auto &val = const_cast<double &>(dof_ptr->getFieldData());
368  val = 1;
369  CHKERR pipeline_mng->loopFiniteElements();
370  CHKERR scale_tag_val();
371  CHKERR post_proc_fe->writeFile(
372  "out_base_dof_" + boost::lexical_cast<std::string>(nb) + ".h5m");
373  CHKERR post_proc_fe->getPostProcMesh().delete_mesh();
374  val = 0;
375  ++nb;
376  };
377 
379 }
380 //! [Postprocess results]
381 
382 //! [Check results]
384 //! [Check results]
385 
386 int main(int argc, char *argv[]) {
387 
388  // Initialisation of MoFEM/PETSc and MOAB data structures
389  MoFEM::Core::Initialize(&argc, &argv, (char *)0, help);
390 
391  try {
392 
393  //! [Register MoFEM discrete manager in PETSc]
394  DMType dm_name = "DMMOFEM";
395  CHKERR DMRegister_MoFEM(dm_name);
396  //! [Register MoFEM discrete manager in PETSc
397 
398  // Add logging channel for example
399  auto core_log = logging::core::get();
400  core_log->add_sink(
401  LogManager::createSink(LogManager::getStrmWorld(), "PLOTBASE"));
402  LogManager::setLog("PLOTBASE");
403  MOFEM_LOG_TAG("PLOTBASE", "plotbase");
404 
405  //! [Create MoAB]
406  moab::Core mb_instance; ///< mesh database
407  moab::Interface &moab = mb_instance; ///< mesh database interface
408  //! [Create MoAB]
409 
410  //! [Create MoFEM]
411  MoFEM::Core core(moab); ///< finite element database
412  MoFEM::Interface &m_field = core; ///< finite element database insterface
413  //! [Create MoFEM]
414 
415  //! [Example]
416  Example ex(m_field);
417  CHKERR ex.runProblem();
418  //! [Example]
419  }
420  CATCH_ERRORS;
421 
423 
424  return 0;
425 }
426 
429  moab::Core core_ref;
430  moab::Interface &moab_ref = core_ref;
431 
432  char ref_mesh_file_name[255];
433 
434  if (SPACE_DIM == 2)
435  strcpy(ref_mesh_file_name, "ref_mesh2d.h5m");
436  else if (SPACE_DIM == 3)
437  strcpy(ref_mesh_file_name, "ref_mesh3d.h5m");
438  else
439  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
440  "Dimension not implemented");
441 
442  CHKERR PetscOptionsGetString(PETSC_NULL, "", "-ref_file", ref_mesh_file_name,
443  255, PETSC_NULL);
444  CHKERR moab_ref.load_file(ref_mesh_file_name, 0, "");
445 
446  // Get elements
447  Range elems;
448  CHKERR moab_ref.get_entities_by_dimension(0, SPACE_DIM, elems);
449 
450  // Add mid-nodes on edges
451  EntityHandle meshset;
452  CHKERR moab_ref.create_meshset(MESHSET_SET, meshset);
453  CHKERR moab_ref.add_entities(meshset, elems);
454  CHKERR moab_ref.convert_entities(meshset, true, false, false);
455  CHKERR moab_ref.delete_entities(&meshset, 1);
456 
457  // Get nodes on the mesh
458  Range elem_nodes;
459  CHKERR moab_ref.get_connectivity(elems, elem_nodes, false);
460 
461  // Map node entity and Gauss pint number
462  std::map<EntityHandle, int> nodes_pts_map;
463 
464  // Set gauss points coordinates from the reference mesh
465  gaussPts.resize(SPACE_DIM + 1, elem_nodes.size(), false);
466  gaussPts.clear();
467  Range::iterator nit = elem_nodes.begin();
468  for (int gg = 0; nit != elem_nodes.end(); nit++, gg++) {
469  double coords[3];
470  CHKERR moab_ref.get_coords(&*nit, 1, coords);
471  for (auto d : {0, 1, 2})
472  gaussPts(d, gg) = coords[d];
473  nodes_pts_map[*nit] = gg;
474  }
475 
476  if (SPACE_DIM == 2) {
477  // Set size of adjacency matrix (note ho order nodes 3 nodes and 3 nodes on
478  // edges)
479  refEleMap.resize(elems.size(), 3 + 3);
480  } else if (SPACE_DIM == 3) {
481  refEleMap.resize(elems.size(), 4 + 6);
482  }
483 
484  // Set adjacency matrix
485  Range::iterator tit = elems.begin();
486  for (int tt = 0; tit != elems.end(); ++tit, ++tt) {
487  const EntityHandle *conn;
488  int num_nodes;
489  CHKERR moab_ref.get_connectivity(*tit, conn, num_nodes, false);
490  for (int nn = 0; nn != num_nodes; ++nn) {
491  refEleMap(tt, nn) = nodes_pts_map[conn[nn]];
492  }
493  }
494 
496 }
497 
500 
501  const int num_nodes = gaussPts.size2();
502 
503  // Calculate shape functions
504 
505  switch (numeredEntFiniteElementPtr->getEntType()) {
506  case MBTRI:
507  shapeFunctions.resize(num_nodes, 3);
508  CHKERR Tools::shapeFunMBTRI(&*shapeFunctions.data().begin(),
509  &gaussPts(0, 0), &gaussPts(1, 0), num_nodes);
510  break;
511  case MBQUAD: {
512  shapeFunctions.resize(num_nodes, 4);
513  for (int gg = 0; gg != num_nodes; gg++) {
514  double ksi = gaussPts(0, gg);
515  double eta = gaussPts(1, gg);
516  shapeFunctions(gg, 0) = N_MBQUAD0(ksi, eta);
517  shapeFunctions(gg, 1) = N_MBQUAD1(ksi, eta);
518  shapeFunctions(gg, 2) = N_MBQUAD2(ksi, eta);
519  shapeFunctions(gg, 3) = N_MBQUAD3(ksi, eta);
520  }
521  } break;
522  case MBTET: {
523  shapeFunctions.resize(num_nodes, 8);
524  CHKERR Tools::shapeFunMBTET(&*shapeFunctions.data().begin(),
525  &gaussPts(0, 0), &gaussPts(1, 0),
526  &gaussPts(2, 0), num_nodes);
527  } break;
528  case MBHEX: {
529  shapeFunctions.resize(num_nodes, 8);
530  for (int gg = 0; gg != num_nodes; gg++) {
531  double ksi = gaussPts(0, gg);
532  double eta = gaussPts(1, gg);
533  double zeta = gaussPts(2, gg);
534  shapeFunctions(gg, 0) = N_MBHEX0(ksi, eta, zeta);
535  shapeFunctions(gg, 1) = N_MBHEX1(ksi, eta, zeta);
536  shapeFunctions(gg, 2) = N_MBHEX2(ksi, eta, zeta);
537  shapeFunctions(gg, 3) = N_MBHEX3(ksi, eta, zeta);
538  shapeFunctions(gg, 4) = N_MBHEX4(ksi, eta, zeta);
539  shapeFunctions(gg, 5) = N_MBHEX5(ksi, eta, zeta);
540  shapeFunctions(gg, 6) = N_MBHEX6(ksi, eta, zeta);
541  shapeFunctions(gg, 7) = N_MBHEX7(ksi, eta, zeta);
542  }
543  } break;
544  default:
545  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
546  "Not implemented element type");
547  }
548 
549  // Create physical nodes
550 
551  // MoAB interface allowing for creating nodes and elements in the bulk
552  ReadUtilIface *iface;
553  CHKERR getPostProcMesh().query_interface(iface);
554 
555  std::vector<double *> arrays; /// pointers to memory allocated by MoAB for
556  /// storing X, Y, and Z coordinates
557  EntityHandle startv; // Starting handle for vertex
558  // Allocate memory for num_nodes, and return starting handle, and access to
559  // memort.
560  CHKERR iface->get_node_coords(3, num_nodes, 0, startv, arrays);
561 
562  mapGaussPts.resize(gaussPts.size2());
563  for (int gg = 0; gg != num_nodes; ++gg)
564  mapGaussPts[gg] = startv + gg;
565 
566  Tag th;
567  int def_in_the_loop = -1;
568  CHKERR getPostProcMesh().tag_get_handle("NB_IN_THE_LOOP", 1, MB_TYPE_INTEGER,
569  th, MB_TAG_CREAT | MB_TAG_SPARSE,
570  &def_in_the_loop);
571 
572  // Create physical elements
573 
574  const int num_el = refEleMap.size1();
575  const int num_nodes_on_ele = refEleMap.size2();
576 
577  EntityHandle starte; // Starting handle to first created element
578  EntityHandle *conn; // Access to MOAB memory with connectivity of elements
579 
580  // Create tris/tets in the bulk in MoAB database
581  if (SPACE_DIM == 2)
582  CHKERR iface->get_element_connect(num_el, num_nodes_on_ele, MBTRI, 0,
583  starte, conn);
584  else if (SPACE_DIM == 3)
585  CHKERR iface->get_element_connect(num_el, num_nodes_on_ele, MBTET, 0,
586  starte, conn);
587  else
588  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
589  "Dimension not implemented");
590 
591  // At this point elements (memory for elements) is allocated, at code bellow
592  // actual connectivity of elements is set.
593  for (unsigned int tt = 0; tt != refEleMap.size1(); ++tt) {
594  for (int nn = 0; nn != num_nodes_on_ele; ++nn)
595  conn[num_nodes_on_ele * tt + nn] = mapGaussPts[refEleMap(tt, nn)];
596  }
597 
598  // Finalise elements creation. At that point MOAB updates adjacency tables,
599  // and elements are ready to use.
600  CHKERR iface->update_adjacencies(starte, num_el, num_nodes_on_ele, conn);
601 
602  auto physical_elements = Range(starte, starte + num_el - 1);
603  CHKERR getPostProcMesh().tag_clear_data(th, physical_elements, &(nInTheLoop));
604 
605  EntityHandle fe_ent = numeredEntFiniteElementPtr->getEnt();
606  int fe_num_nodes;
607  {
608  const EntityHandle *conn;
609  mField.get_moab().get_connectivity(fe_ent, conn, fe_num_nodes, true);
610  coords.resize(3 * fe_num_nodes, false);
611  CHKERR mField.get_moab().get_coords(conn, fe_num_nodes, &coords[0]);
612  }
613 
614  // Set physical coordinates to physical nodes
617  &*shapeFunctions.data().begin());
618 
620  arrays[0], arrays[1], arrays[2]);
621  const double *t_coords_ele_x = &coords[0];
622  const double *t_coords_ele_y = &coords[1];
623  const double *t_coords_ele_z = &coords[2];
624  for (int gg = 0; gg != num_nodes; ++gg) {
626  t_coords_ele_x, t_coords_ele_y, t_coords_ele_z);
627  t_coords(i) = 0;
628  for (int nn = 0; nn != fe_num_nodes; ++nn) {
629  t_coords(i) += t_n * t_ele_coords(i);
630  for (auto ii : {0, 1, 2})
631  if (std::abs(t_coords(ii)) < std::numeric_limits<float>::epsilon())
632  t_coords(ii) = 0;
633  ++t_ele_coords;
634  ++t_n;
635  }
636  ++t_coords;
637  }
638 
640 }
641 
644  ParallelComm *pcomm_post_proc_mesh =
645  ParallelComm::get_pcomm(coreMeshPtr.get(), MYPCOMM_INDEX);
646  if (pcomm_post_proc_mesh != NULL)
647  delete pcomm_post_proc_mesh;
649 };
650 
653 
654  auto resolve_shared_ents = [&]() {
656 
657  ParallelComm *pcomm_post_proc_mesh =
658  ParallelComm::get_pcomm(&(getPostProcMesh()), MYPCOMM_INDEX);
659  if (pcomm_post_proc_mesh == NULL) {
660  // wrapRefMeshComm =
661  // boost::make_shared<WrapMPIComm>(T::mField.get_comm(), false);
662  pcomm_post_proc_mesh = new ParallelComm(
663  &(getPostProcMesh()),
664  PETSC_COMM_WORLD /*(T::wrapRefMeshComm)->get_comm()*/);
665  }
666 
667  CHKERR pcomm_post_proc_mesh->resolve_shared_ents(0);
668 
670  };
671 
672  CHKERR resolve_shared_ents();
673 
675 }
MoFEMFunctionReturnHot
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:460
help
static char help[]
Definition: plot_base.cpp:12
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
Example::checkResults
MoFEMErrorCode checkResults()
[Postprocess results]
Definition: dynamic_first_order_con_law.cpp:1205
EXECUTABLE_DIMENSION
#define EXECUTABLE_DIMENSION
Definition: plastic.cpp:13
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
Example::assembleSystem
MoFEMErrorCode assembleSystem()
[Push operators to pipeline]
Definition: dynamic_first_order_con_law.cpp:647
MoFEM::CoreTmp< 0 >
Core (interface) class.
Definition: Core.hpp:82
H1
@ H1
continuous field
Definition: definitions.h:85
MyPostProc
Definition: plot_base.cpp:34
FTensor::Tensor1
Definition: Tensor1_value.hpp:8
MyPostProc::setGaussPts
MoFEMErrorCode setGaussPts(int order)
Definition: plot_base.cpp:498
Example::Example
Example(MoFEM::Interface &m_field)
Definition: plot_base.cpp:50
EntityHandle
MoFEM::AinsworthOrderHooks::broken_nbvolumetet_volume_hdiv
static boost::function< int(int)> broken_nbvolumetet_volume_hdiv
Definition: Hdiv.hpp:29
MoFEM::OpSetContravariantPiolaTransformOnFace2D
OpSetContravariantPiolaTransformOnFace2DImpl< 2 > OpSetContravariantPiolaTransformOnFace2D
Definition: UserDataOperators.hpp:3076
N_MBHEX5
#define N_MBHEX5(x, y, z)
Definition: fem_tools.h:76
L2
@ L2
field with C-1 continuity
Definition: definitions.h:88
MoFEM::AinsworthOrderHooks::broken_nbvolumetet_face_hdiv
static boost::function< int(int)> broken_nbvolumetet_face_hdiv
Definition: Hdiv.hpp:28
MoFEM::PipelineManager::loopFiniteElements
MoFEMErrorCode loopFiniteElements(SmartPetscObj< DM > dm=nullptr)
Iterate finite elements.
Definition: PipelineManager.cpp:19
N_MBHEX0
#define N_MBHEX0(x, y, z)
Definition: fem_tools.h:71
main
int main(int argc, char *argv[])
[Check results]
Definition: plot_base.cpp:386
MoFEM::Exceptions::MoFEMErrorCode
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
Definition: Exceptions.hpp:56
MoFEM::Types::MatrixDouble
UBlasMatrix< double > MatrixDouble
Definition: Types.hpp:77
LASTCONTINUITY
@ LASTCONTINUITY
Definition: definitions.h:102
MoFEM::th
Tag th
Definition: Projection10NodeCoordsOnField.cpp:122
MoFEM::PipelineManager
PipelineManager interface.
Definition: PipelineManager.hpp:24
MoFEM.hpp
N_MBHEX7
#define N_MBHEX7(x, y, z)
Definition: fem_tools.h:78
MoFEM::CoreTmp< 0 >::Finalize
static MoFEMErrorCode Finalize()
Checks for options to be called at the conclusion of the program.
Definition: Core.cpp:112
MoFEM::Simple
Simple interface for fast problem set-up.
Definition: Simple.hpp:27
N_MBQUAD2
#define N_MBQUAD2(x, y)
quad shape function
Definition: fem_tools.h:59
Example::base
FieldApproximationBase base
Definition: plot_base.cpp:68
MOFEM_IMPOSSIBLE_CASE
@ MOFEM_IMPOSSIBLE_CASE
Definition: definitions.h:35
zeta
double zeta
Viscous hardening.
Definition: plastic.cpp:126
order
constexpr int order
Definition: dg_projection.cpp:18
MoFEM::DeprecatedCoreInterface
Deprecated interface functions.
Definition: DeprecatedCoreInterface.hpp:16
MoFEM::Interface
DeprecatedCoreInterface Interface
Definition: Interface.hpp:2010
MoFEM::AinsworthOrderHooks::broken_nbfacetri_edge_hdiv
static boost::function< int(int)> broken_nbfacetri_edge_hdiv
Definition: Hdiv.hpp:25
Example::readMesh
MoFEMErrorCode readMesh()
[Run problem]
Definition: dynamic_first_order_con_law.cpp:463
MoFEM::PostProcBrokenMeshInMoab
Definition: PostProcBrokenMeshInMoabBase.hpp:667
FieldSpace
FieldSpace
approximation spaces
Definition: definitions.h:82
N_MBHEX3
#define N_MBHEX3(x, y, z)
Definition: fem_tools.h:74
MoFEM::AinsworthOrderHooks::broken_nbvolumetet_edge_hdiv
static boost::function< int(int)> broken_nbvolumetet_edge_hdiv
Definition: Hdiv.hpp:27
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:548
MoFEM
implementation of Data Operators for Forces and Sources
Definition: Common.hpp:10
MoFEM::BcManager
Simple interface for fast problem set-up.
Definition: BcManager.hpp:25
FieldContinuity
FieldContinuity
Field continuity.
Definition: definitions.h:99
MoFEM::OpCalculateHVecVectorField
Get vector field for H-div approximation.
Definition: UserDataOperators.hpp:2115
Example
[Example]
Definition: plastic.cpp:177
OpPPMap
OpPostProcMapInMoab< SPACE_DIM, SPACE_DIM > OpPPMap
Definition: photon_diffusion.cpp:29
Example::solveSystem
MoFEMErrorCode solveSystem()
[Solve]
Definition: dynamic_first_order_con_law.cpp:893
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
eta
double eta
Definition: free_surface.cpp:170
N_MBHEX1
#define N_MBHEX1(x, y, z)
Definition: fem_tools.h:72
N_MBQUAD1
#define N_MBQUAD1(x, y)
quad shape function
Definition: fem_tools.h:58
MoFEM::FaceElementForcesAndSourcesCore
Face finite element.
Definition: FaceElementForcesAndSourcesCore.hpp:23
AINSWORTH_LOBATTO_BASE
@ AINSWORTH_LOBATTO_BASE
Definition: definitions.h:62
Example::setIntegrationRules
MoFEMErrorCode setIntegrationRules()
[Set up problem]
Definition: plot_base.cpp:229
MOFEM_LOG_TAG
#define MOFEM_LOG_TAG(channel, tag)
Tag channel.
Definition: LogManager.hpp:339
MyPostProc::generateReferenceElementMesh
MoFEMErrorCode generateReferenceElementMesh()
Definition: plot_base.cpp:427
MoFEM::AinsworthOrderHooks::broken_nbfacetri_face_hdiv
static boost::function< int(int)> broken_nbfacetri_face_hdiv
Definition: Hdiv.hpp:26
MoFEM::VolumeElementForcesAndSourcesCore
Volume finite element base.
Definition: VolumeElementForcesAndSourcesCore.hpp:26
i
FTensor::Index< 'i', SPACE_DIM > i
Definition: hcurl_divergence_operator_2d.cpp:27
Example::space
FieldSpace space
Definition: plot_base.cpp:69
FTensor::Index< 'i', 3 >
AINSWORTH_BERNSTEIN_BEZIER_BASE
@ AINSWORTH_BERNSTEIN_BEZIER_BASE
Definition: definitions.h:64
MoFEM::AddHOOps
Add operators pushing bases from local to physical configuration.
Definition: HODataOperators.hpp:413
SPACE_DIM
constexpr int SPACE_DIM
Definition: plot_base.cpp:26
v
const double v
phase velocity of light in medium (cm/ns)
Definition: initial_diffusion.cpp:40
ElementsAndOps
Definition: child_and_parent.cpp:18
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
MOFEM_LOG
#define MOFEM_LOG(channel, severity)
Log.
Definition: LogManager.hpp:308
MyPostProc::postProcess
MoFEMErrorCode postProcess()
Definition: plot_base.cpp:651
CATCH_ERRORS
#define CATCH_ERRORS
Catch errors.
Definition: definitions.h:385
DEMKOWICZ_JACOBI_BASE
@ DEMKOWICZ_JACOBI_BASE
Definition: definitions.h:66
FTensor::Tensor0
Definition: Tensor0.hpp:16
MyPostProc::refEleMap
ublas::matrix< int > refEleMap
Definition: plot_base.cpp:44
N_MBHEX2
#define N_MBHEX2(x, y, z)
Definition: fem_tools.h:73
MoFEM::Core
CoreTmp< 0 > Core
Definition: Core.hpp:1148
EntData
EntitiesFieldData::EntData EntData
Definition: plot_base.cpp:29
Example::setupProblem
MoFEMErrorCode setupProblem()
[Run problem]
Definition: plastic.cpp:233
UserDataOperator
ForcesAndSourcesCore::UserDataOperator UserDataOperator
Definition: HookeElement.hpp:75
Example::boundaryCondition
MoFEMErrorCode boundaryCondition()
[Set up problem]
Definition: dynamic_first_order_con_law.cpp:536
N_MBHEX4
#define N_MBHEX4(x, y, z)
Definition: fem_tools.h:75
DISCONTINUOUS
@ DISCONTINUOUS
Broken continuity (No effect on L2 space)
Definition: definitions.h:101
PostProcEle
PostProcBrokenMeshInMoab< DomainEle > PostProcEle
Definition: plot_base.cpp:32
AINSWORTH_LEGENDRE_BASE
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base .
Definition: definitions.h:60
HCURL
@ HCURL
field with continuous tangents
Definition: definitions.h:86
N_MBQUAD0
#define N_MBQUAD0(x, y)
quad shape function
Definition: fem_tools.h:57
MoFEM::PetscOptionsGetString
PetscErrorCode PetscOptionsGetString(PetscOptions *, const char pre[], const char name[], char str[], size_t size, PetscBool *set)
Definition: DeprecatedPetsc.hpp:172
MoFEM::PetscOptionsGetEList
PetscErrorCode PetscOptionsGetEList(PetscOptions *, const char pre[], const char name[], const char *const *list, PetscInt next, PetscInt *value, PetscBool *set)
Definition: DeprecatedPetsc.hpp:203
MoFEM::OpCalculateHOJacForFace
OpCalculateHOJacForFaceImpl< 2 > OpCalculateHOJacForFace
Definition: HODataOperators.hpp:264
sdf_hertz_2d_axisymm_plane.d
float d
Definition: sdf_hertz_2d_axisymm_plane.py:4
FieldApproximationBase
FieldApproximationBase
approximation base
Definition: definitions.h:58
Example::runProblem
MoFEMErrorCode runProblem()
[Run problem]
Definition: plastic.cpp:214
MyPostProc::preProcess
MoFEMErrorCode preProcess()
Definition: plot_base.cpp:642
N_MBQUAD3
#define N_MBQUAD3(x, y)
quad shape function
Definition: fem_tools.h:60
MoFEM::BitRefManager
Managing BitRefLevels.
Definition: BitRefManager.hpp:21
MoFEMFunctionBeginHot
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:453
MoFEM::OpMakeHdivFromHcurl
Make Hdiv space from Hcurl space in 2d.
Definition: UserDataOperators.hpp:2985
Example::createCommonData
MoFEMErrorCode createCommonData()
[Set up problem]
Definition: plastic.cpp:404
N_MBHEX6
#define N_MBHEX6(x, y, z)
Definition: fem_tools.h:77
DomainEle
ElementsAndOps< SPACE_DIM >::DomainEle DomainEle
Definition: child_and_parent.cpp:34
CONTINUOUS
@ CONTINUOUS
Regular field.
Definition: definitions.h:100
Example::outputResults
MoFEMErrorCode outputResults()
[Solve]
Definition: dynamic_first_order_con_law.cpp:1183
MyPostProc::shapeFunctions
MatrixDouble shapeFunctions
Definition: plot_base.cpp:45
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
HDIV
@ HDIV
field with continuous normal traction
Definition: definitions.h:87
MOFEM_NOT_IMPLEMENTED
@ MOFEM_NOT_IMPLEMENTED
Definition: definitions.h:32
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:698
MoFEM::PetscOptionsGetBool
PetscErrorCode PetscOptionsGetBool(PetscOptions *, const char pre[], const char name[], PetscBool *bval, PetscBool *set)
Definition: DeprecatedPetsc.hpp:182