v0.13.2
Loading...
Searching...
No Matches
PostProcBrokenMeshInMoabBase.cpp
Go to the documentation of this file.
1/**
2 * @file PostProc.cpp
3 * @brief Post processing elements and operators
4 *
5 * @copyright Copyright (c) 2022
6 *
7 */
8
9#include <MoFEM.hpp>
10
11namespace MoFEM {
12
14 : hoNodes(PETSC_TRUE), defMaxLevel(0), optPrefix(""), countEle(0),
15 countVertEle(0), nbVertices(0) {}
16
19
20 CHKERR PetscOptionsGetInt(optPrefix.c_str(), "-max_post_proc_ref_level",
21 &defMaxLevel, PETSC_NULL);
22 CHKERR PetscOptionsGetBool(optPrefix.c_str(), "-max_post_ho_nodes", &hoNodes,
23 PETSC_NULL);
24
25 if (defMaxLevel < 0)
26 SETERRQ(PETSC_COMM_WORLD, MOFEM_INVALID_DATA,
27 "Wrong parameter -max_post_proc_ref_level "
28 "should be positive number");
29
31};
32
35
36 CHKERR getOptions();
37
38 const int max_level = defMaxLevel;
39
40 moab::Core core_ref;
41 moab::Interface &moab_ref = core_ref;
42
43 auto create_reference_element = [&moab_ref]() {
45 constexpr double base_coords[] = {0, 0, 0, 1, 0, 0, 0, 1, 0, 0, 0, 1};
46 EntityHandle nodes[4];
47 for (int nn = 0; nn < 4; nn++) {
48 CHKERR
49 moab_ref.create_vertex(&base_coords[3 * nn], nodes[nn]);
50 }
51 EntityHandle tet;
52 CHKERR moab_ref.create_element(MBTET, nodes, 4, tet);
54 };
55
56 MoFEM::CoreTmp<-1> m_core_ref(moab_ref, PETSC_COMM_SELF, -2);
57 MoFEM::Interface &m_field_ref = m_core_ref;
58
59 auto refine_ref_tetrahedron = [this, &m_field_ref, max_level]() {
61 // seed ref mofem database by setting bit ref level to reference
62 // tetrahedron
63 CHKERR
64 m_field_ref.getInterface<BitRefManager>()->setBitRefLevelByDim(
65 0, 3, BitRefLevel().set(0));
66 for (int ll = 0; ll != max_level; ++ll) {
67 MOFEM_TAG_AND_LOG_C("WORLD", Sev::noisy, "PostProc", "Refine Level %d",
68 ll);
69 Range edges;
70 CHKERR m_field_ref.getInterface<BitRefManager>()
71 ->getEntitiesByTypeAndRefLevel(BitRefLevel().set(ll),
72 BitRefLevel().set(), MBEDGE, edges);
73 Range tets;
74 CHKERR m_field_ref.getInterface<BitRefManager>()
75 ->getEntitiesByTypeAndRefLevel(BitRefLevel().set(ll),
76 BitRefLevel(ll).set(), MBTET, tets);
77 // refine mesh
78 MeshRefinement *m_ref;
79 CHKERR m_field_ref.getInterface(m_ref);
80 CHKERR m_ref->addVerticesInTheMiddleOfEdges(edges,
81 BitRefLevel().set(ll + 1));
82 CHKERR m_ref->refineTets(tets, BitRefLevel().set(ll + 1));
83 }
85 };
86
87 auto get_ref_gauss_pts_and_shape_functions = [this, max_level, &moab_ref,
88 &m_field_ref]() {
90 for (int ll = 0; ll != max_level + 1; ++ll) {
91 Range tets;
92 CHKERR
93 m_field_ref.getInterface<BitRefManager>()->getEntitiesByTypeAndRefLevel(
94 BitRefLevel().set(ll), BitRefLevel().set(ll), MBTET, tets);
95 if (hoNodes) {
96 EntityHandle meshset;
97 CHKERR moab_ref.create_meshset(MESHSET_SET, meshset);
98 CHKERR moab_ref.add_entities(meshset, tets);
99 CHKERR moab_ref.convert_entities(meshset, true, false, false);
100 CHKERR moab_ref.delete_entities(&meshset, 1);
101 }
102 Range elem_nodes;
103 CHKERR moab_ref.get_connectivity(tets, elem_nodes, false);
104
105 auto &gauss_pts = levelGaussPtsOnRefMesh[ll];
106 gauss_pts.resize(elem_nodes.size(), 4, false);
107 std::map<EntityHandle, int> little_map;
108 Range::iterator nit = elem_nodes.begin();
109 for (int gg = 0; nit != elem_nodes.end(); nit++, gg++) {
110 CHKERR moab_ref.get_coords(&*nit, 1, &gauss_pts(gg, 0));
111 little_map[*nit] = gg;
112 }
113 gauss_pts = trans(gauss_pts);
114
115 auto &ref_tets = levelRef[ll];
116 Range::iterator tit = tets.begin();
117 for (int tt = 0; tit != tets.end(); ++tit, ++tt) {
118 const EntityHandle *conn;
119 int num_nodes;
120 CHKERR moab_ref.get_connectivity(*tit, conn, num_nodes, false);
121 if (tt == 0) {
122 ref_tets.resize(tets.size(), num_nodes);
123 }
124 for (int nn = 0; nn != num_nodes; ++nn) {
125 ref_tets(tt, nn) = little_map[conn[nn]];
126 }
127 }
128
129 auto &shape_functions = levelShapeFunctions[ll];
130 shape_functions.resize(elem_nodes.size(), 4);
131 CHKERR ShapeMBTET(&*shape_functions.data().begin(), &gauss_pts(0, 0),
132 &gauss_pts(1, 0), &gauss_pts(2, 0), elem_nodes.size());
133 }
135 };
136
137 levelRef.resize(max_level + 1);
138 levelGaussPtsOnRefMesh.resize(max_level + 1);
139 levelShapeFunctions.resize(max_level + 1);
140
141 CHKERR create_reference_element();
142 CHKERR refine_ref_tetrahedron();
143 CHKERR get_ref_gauss_pts_and_shape_functions();
144
146}
147
150
151 CHKERR getOptions();
152#ifndef NDEBUG
153 if (defMaxLevel > 0)
154 MOFEM_LOG("WORLD", Sev::warning)
155 << "Refinement for hexes is not implemented";
156#endif
157
158 moab::Core core_ref;
159 moab::Interface &moab_ref = core_ref;
160
161 auto create_reference_element = [&moab_ref]() {
163 constexpr double base_coords[] = {0, 0, 0, 1, 0, 0, 1, 1, 0, 0, 1, 0,
164 0, 0, 1, 1, 0, 1, 1, 1, 1, 0, 1, 1};
165 EntityHandle nodes[8];
166 for (int nn = 0; nn < 8; nn++)
167 CHKERR moab_ref.create_vertex(&base_coords[3 * nn], nodes[nn]);
168 EntityHandle hex;
169 CHKERR moab_ref.create_element(MBHEX, nodes, 8, hex);
171 };
172
173 auto add_ho_nodes = [&]() {
175 Range hexes;
176 CHKERR moab_ref.get_entities_by_type(0, MBHEX, hexes, true);
177 EntityHandle meshset;
178 CHKERR moab_ref.create_meshset(MESHSET_SET, meshset);
179 CHKERR moab_ref.add_entities(meshset, hexes);
180 CHKERR moab_ref.convert_entities(meshset, true, true, true);
181 CHKERR moab_ref.delete_entities(&meshset, 1);
183 };
184
185 auto set_gauss_pts = [&](std::map<EntityHandle, int> &little_map) {
187 Range hexes;
188 CHKERR moab_ref.get_entities_by_type(0, MBHEX, hexes, true);
189 Range hexes_nodes;
190 CHKERR moab_ref.get_connectivity(hexes, hexes_nodes, false);
191 auto &gauss_pts = levelGaussPtsOnRefMesh[0];
192 gauss_pts.resize(hexes_nodes.size(), 4, false);
193 size_t gg = 0;
194 for (auto node : hexes_nodes) {
195 CHKERR moab_ref.get_coords(&node, 1, &gauss_pts(gg, 0));
196 little_map[node] = gg;
197 ++gg;
198 }
199 gauss_pts = trans(gauss_pts);
201 };
202
203 auto set_ref_hexes = [&](std::map<EntityHandle, int> &little_map) {
205 Range hexes;
206 CHKERR moab_ref.get_entities_by_type(0, MBHEX, hexes, true);
207 size_t hh = 0;
208 auto &ref_hexes = levelRef[0];
209 for (auto hex : hexes) {
210 const EntityHandle *conn;
211 int num_nodes;
212 CHKERR moab_ref.get_connectivity(hex, conn, num_nodes, false);
213 if (ref_hexes.size2() != num_nodes) {
214 ref_hexes.resize(hexes.size(), num_nodes);
215 }
216 for (int nn = 0; nn != num_nodes; ++nn) {
217 ref_hexes(hh, nn) = little_map[conn[nn]];
218 }
219 ++hh;
220 }
222 };
223
224 auto set_shape_functions = [&]() {
226 auto &gauss_pts = levelGaussPtsOnRefMesh[0];
227 auto &shape_functions = levelShapeFunctions[0];
228 const auto nb_gauss_pts = gauss_pts.size2();
229 shape_functions.resize(nb_gauss_pts, 8);
230 for (int gg = 0; gg != nb_gauss_pts; ++gg) {
231 const double ksi = gauss_pts(0, gg);
232 const double zeta = gauss_pts(1, gg);
233 const double eta = gauss_pts(2, gg);
234 shape_functions(gg, 0) = N_MBHEX0(ksi, zeta, eta);
235 shape_functions(gg, 1) = N_MBHEX1(ksi, zeta, eta);
236 shape_functions(gg, 2) = N_MBHEX2(ksi, zeta, eta);
237 shape_functions(gg, 3) = N_MBHEX3(ksi, zeta, eta);
238 shape_functions(gg, 4) = N_MBHEX4(ksi, zeta, eta);
239 shape_functions(gg, 5) = N_MBHEX5(ksi, zeta, eta);
240 shape_functions(gg, 6) = N_MBHEX6(ksi, zeta, eta);
241 shape_functions(gg, 7) = N_MBHEX7(ksi, zeta, eta);
242 }
244 };
245
246 levelRef.resize(1);
247 levelGaussPtsOnRefMesh.resize(1);
248 levelShapeFunctions.resize(1);
249
250 CHKERR create_reference_element();
251 if (hoNodes)
252 CHKERR add_ho_nodes();
253 std::map<EntityHandle, int> little_map;
254 CHKERR set_gauss_pts(little_map);
255 CHKERR set_ref_hexes(little_map);
256 CHKERR set_shape_functions();
257
259}
260
263
264 CHKERR getOptions();
265 const int max_level = defMaxLevel;
266
267 moab::Core core_ref;
268 moab::Interface &moab_ref = core_ref;
269
270 auto create_reference_element = [&moab_ref]() {
272 constexpr double base_coords[] = {
273
274 0, 0,
275 0,
276
277 1, 0,
278 0,
279
280 0, 1,
281 0
282
283 };
284 EntityHandle nodes[3];
285 for (int nn = 0; nn != 3; ++nn)
286 CHKERR moab_ref.create_vertex(&base_coords[3 * nn], nodes[nn]);
287 EntityHandle tri;
288 CHKERR moab_ref.create_element(MBTRI, nodes, 3, tri);
289
290 Range edges;
291 CHKERR moab_ref.get_adjacencies(&tri, 1, 1, true, edges,
292 moab::Interface::UNION);
293
295 };
296
297 CHKERR create_reference_element();
298
299 MoFEM::CoreTmp<-1> m_core_ref(moab_ref, PETSC_COMM_SELF, -2);
300 MoFEM::Interface &m_field_ref = m_core_ref;
301
302 auto refine_ref_triangles = [this, &m_field_ref, max_level]() {
304 // seed ref mofem database by setting bit ref level to reference
305 // tetrahedron
306 CHKERR
307 m_field_ref.getInterface<BitRefManager>()->setBitRefLevelByDim(
308 0, 2, BitRefLevel().set(0));
309
310 for (int ll = 0; ll != max_level; ++ll) {
311 MOFEM_TAG_AND_LOG_C("WORLD", Sev::noisy, "PostProc", "Refine Level %d",
312 ll);
313 Range edges;
314 CHKERR m_field_ref.getInterface<BitRefManager>()
315 ->getEntitiesByTypeAndRefLevel(BitRefLevel().set(ll),
316 BitRefLevel().set(), MBEDGE, edges);
317 Range tris;
318 CHKERR m_field_ref.getInterface<BitRefManager>()
319 ->getEntitiesByTypeAndRefLevel(BitRefLevel().set(ll),
320 BitRefLevel(ll).set(), MBTRI, tris);
321 // refine mesh
322 auto m_ref = m_field_ref.getInterface<MeshRefinement>();
323 CHKERR m_ref->addVerticesInTheMiddleOfEdges(edges,
324 BitRefLevel().set(ll + 1));
325 CHKERR m_ref->refineTris(tris, BitRefLevel().set(ll + 1));
326 }
328 };
329
330 // auto set_gauss_pts = [&](std::map<EntityHandle, int> &little_map) {
331 // MoFEMFunctionBegin;
332 // Range faces;
333 // CHKERR moab_ref.get_entities_by_type(0, MBTRI, faces, true);
334 // Range faces_nodes;
335 // CHKERR moab_ref.get_connectivity(faces, faces_nodes, false);
336 // auto &gauss_pts = levelGaussPtsOnRefMesh[0];
337 // gauss_pts.resize(faces_nodes.size(), 4, false);
338 // size_t gg = 0;
339 // for (auto node : faces_nodes) {
340 // CHKERR moab_ref.get_coords(&node, 1, &gauss_pts(gg, 0));
341 // little_map[node] = gg;
342 // ++gg;
343 // }
344 // gauss_pts = trans(gauss_pts);
345 // MoFEMFunctionReturn(0);
346 // };
347
348 auto get_ref_gauss_pts_and_shape_functions = [this, max_level, &moab_ref,
349 &m_field_ref]() {
351 for (int ll = 0; ll != max_level + 1; ++ll) {
352
353 Range tris;
354 CHKERR
355 m_field_ref.getInterface<BitRefManager>()->getEntitiesByTypeAndRefLevel(
356 BitRefLevel().set(ll), BitRefLevel().set(ll), MBTRI, tris);
357
358 if (hoNodes) {
359 EntityHandle meshset;
360 CHKERR moab_ref.create_meshset(MESHSET_SET, meshset);
361 CHKERR moab_ref.add_entities(meshset, tris);
362 CHKERR moab_ref.convert_entities(meshset, true, false, false);
363 CHKERR moab_ref.delete_entities(&meshset, 1);
364 }
365
366 Range elem_nodes;
367 CHKERR moab_ref.get_connectivity(tris, elem_nodes, false);
368
369 auto &gauss_pts = levelGaussPtsOnRefMesh[ll];
370 gauss_pts.resize(elem_nodes.size(), 3, false);
371 std::map<EntityHandle, int> little_map;
372 Range::iterator nit = elem_nodes.begin();
373 for (int gg = 0; nit != elem_nodes.end(); nit++, gg++) {
374 CHKERR moab_ref.get_coords(&*nit, 1, &gauss_pts(gg, 0));
375 little_map[*nit] = gg;
376 }
377 gauss_pts = trans(gauss_pts);
378
379 auto &ref_tris = levelRef[ll];
380 Range::iterator tit = tris.begin();
381 for (int tt = 0; tit != tris.end(); ++tit, ++tt) {
382 const EntityHandle *conn;
383 int num_nodes;
384 CHKERR moab_ref.get_connectivity(*tit, conn, num_nodes, false);
385 if (tt == 0) {
386 ref_tris.resize(tris.size(), num_nodes);
387 }
388 for (int nn = 0; nn != num_nodes; ++nn) {
389 ref_tris(tt, nn) = little_map[conn[nn]];
390 }
391 }
392
393 auto &shape_functions = levelShapeFunctions[ll];
394 shape_functions.resize(elem_nodes.size(), 3);
395 CHKERR ShapeMBTRI(&*shape_functions.data().begin(), &gauss_pts(0, 0),
396 &gauss_pts(1, 0), elem_nodes.size());
397 }
399 };
400
401 levelRef.resize(max_level + 1);
402 levelGaussPtsOnRefMesh.resize(max_level + 1);
403 levelShapeFunctions.resize(max_level + 1);
404
405 CHKERR refine_ref_triangles();
406 CHKERR get_ref_gauss_pts_and_shape_functions();
407
409}
410
413
414 CHKERR getOptions();
415#ifndef NDEBUG
416 if (defMaxLevel > 0)
417 MOFEM_LOG("WORLD", Sev::warning)
418 << "Refinement for quad is not implemented";
419#endif
420
421 moab::Core core_ref;
422 moab::Interface &moab_ref = core_ref;
423
424 auto create_reference_element = [&moab_ref]() {
426 constexpr double base_coords[] = {
427
428 0, 0,
429 0,
430
431 1, 0,
432 0,
433
434 1, 1,
435 0,
436
437 0, 1,
438 0
439
440 };
441 EntityHandle nodes[4];
442 for (int nn = 0; nn < 4; nn++)
443 CHKERR moab_ref.create_vertex(&base_coords[3 * nn], nodes[nn]);
444 EntityHandle quad;
445 CHKERR moab_ref.create_element(MBQUAD, nodes, 4, quad);
447 };
448
449 auto add_ho_nodes = [&]() {
451 Range quads;
452 CHKERR moab_ref.get_entities_by_type(0, MBQUAD, quads, true);
453 EntityHandle meshset;
454 CHKERR moab_ref.create_meshset(MESHSET_SET, meshset);
455 CHKERR moab_ref.add_entities(meshset, quads);
456 CHKERR moab_ref.convert_entities(meshset, true, true, true);
457 CHKERR moab_ref.delete_entities(&meshset, 1);
459 };
460
461 auto set_gauss_pts = [&](std::map<EntityHandle, int> &little_map) {
463 Range quads;
464 CHKERR moab_ref.get_entities_by_type(0, MBQUAD, quads, true);
465 Range quads_nodes;
466 CHKERR moab_ref.get_connectivity(quads, quads_nodes, false);
467 auto &gauss_pts = levelGaussPtsOnRefMesh[0];
468 gauss_pts.resize(quads_nodes.size(), 4, false);
469 size_t gg = 0;
470 for (auto node : quads_nodes) {
471 CHKERR moab_ref.get_coords(&node, 1, &gauss_pts(gg, 0));
472 little_map[node] = gg;
473 ++gg;
474 }
475 gauss_pts = trans(gauss_pts);
477 };
478
479 auto set_ref_quads = [&](std::map<EntityHandle, int> &little_map) {
481 Range quads;
482 CHKERR moab_ref.get_entities_by_type(0, MBQUAD, quads, true);
483 size_t hh = 0;
484 auto &ref_quads = levelRef[0];
485 for (auto quad : quads) {
486 const EntityHandle *conn;
487 int num_nodes;
488 CHKERR moab_ref.get_connectivity(quad, conn, num_nodes, false);
489 if (ref_quads.size2() != num_nodes) {
490 ref_quads.resize(quads.size(), num_nodes);
491 }
492 for (int nn = 0; nn != num_nodes; ++nn) {
493 ref_quads(hh, nn) = little_map[conn[nn]];
494 }
495 ++hh;
496 }
498 };
499
500 auto set_shape_functions = [&]() {
502 auto &gauss_pts = levelGaussPtsOnRefMesh[0];
503 auto &shape_functions = levelShapeFunctions[0];
504 const auto nb_gauss_pts = gauss_pts.size2();
505 shape_functions.resize(nb_gauss_pts, 4);
506 for (int gg = 0; gg != nb_gauss_pts; ++gg) {
507 const double ksi = gauss_pts(0, gg);
508 const double zeta = gauss_pts(1, gg);
509 shape_functions(gg, 0) = N_MBQUAD0(ksi, zeta);
510 shape_functions(gg, 1) = N_MBQUAD1(ksi, zeta);
511 shape_functions(gg, 2) = N_MBQUAD2(ksi, zeta);
512 shape_functions(gg, 3) = N_MBQUAD3(ksi, zeta);
513 }
515 };
516
517 levelRef.resize(1);
518 levelGaussPtsOnRefMesh.resize(1);
519 levelShapeFunctions.resize(1);
520
521 CHKERR create_reference_element();
522 if (hoNodes)
523 CHKERR add_ho_nodes();
524 std::map<EntityHandle, int> little_map;
525 CHKERR set_gauss_pts(little_map);
526 CHKERR set_ref_quads(little_map);
527 CHKERR set_shape_functions();
528
530}
531
534
535 CHKERR getOptions();
536#ifndef NDEBUG
537 if (defMaxLevel > 0)
538 MOFEM_LOG("WORLD", Sev::warning)
539 << "Refinement for edges is not implemented";
540#endif
541
542 auto set_gauss_pts = [&](std::map<EntityHandle, int> &little_map) {
544
545 int nb_nodes = 2;
546 if (hoNodes)
547 nb_nodes = 3;
548
549 auto &gauss_pts = levelGaussPtsOnRefMesh[0];
550 gauss_pts.resize(2, nb_nodes, false);
551 gauss_pts.clear();
552
553 int nn = 0;
554 for (; nn != 2; ++nn) {
555 gauss_pts(0, nn) = static_cast<double>(nn);
556 little_map[nn] = nn;
557 }
558
559 if (nn < nb_nodes) {
560 gauss_pts(0, nn) = 0.5;
561 little_map[nn] = 2;
562 }
563
565 };
566
567 auto set_ref_edges = [&](std::map<EntityHandle, int> &little_map) {
569
570 int level = 0;
571 int nb_edges = level + 1;
572
573 int nb_nodes = 2;
574 if (hoNodes)
575 nb_nodes = 3;
576
577 auto &ref_edges = levelRef[level];
578 ref_edges.resize(nb_edges, nb_nodes, false);
579
580 for (int ee = 0; ee != nb_edges; ++ee) {
581 int nn = 0;
582 for (; nn != 2; ++nn) {
583 ref_edges(ee, nn) = nb_nodes * ee + nn;
584 }
585 if (nn < nb_nodes) {
586 ref_edges(ee, nn) = nb_nodes * ee + 2;
587 }
588 }
589
591 };
592
593 auto set_shape_functions = [&]() {
595 auto &gauss_pts = levelGaussPtsOnRefMesh[0];
596 auto &shape_functions = levelShapeFunctions[0];
597 const auto nb_gauss_pts = gauss_pts.size2();
598 shape_functions.resize(nb_gauss_pts, 2);
599 for (int gg = 0; gg != nb_gauss_pts; ++gg) {
600 const double ksi = gauss_pts(0, gg);
601 shape_functions(gg, 0) = N_MBEDGE0(ksi);
602 shape_functions(gg, 1) = N_MBEDGE1(ksi);
603 }
605 };
606
607 levelRef.resize(1);
608 levelGaussPtsOnRefMesh.resize(1);
609 levelShapeFunctions.resize(1);
610
611 std::map<EntityHandle, int> little_map;
612 CHKERR set_gauss_pts(little_map);
613 CHKERR set_ref_edges(little_map);
614 CHKERR set_shape_functions();
615
617}
618
619} // namespace MoFEM
#define MOFEM_TAG_AND_LOG_C(channel, severity, tag, format,...)
Tag and log in channel.
Definition: LogManager.hpp:370
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:346
@ MOFEM_INVALID_DATA
Definition: definitions.h:36
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:416
#define CHKERR
Inline error check.
Definition: definitions.h:535
#define N_MBQUAD3(x, y)
quad shape function
Definition: fem_tools.h:60
#define N_MBHEX7(x, y, z)
Definition: fem_tools.h:78
#define N_MBHEX3(x, y, z)
Definition: fem_tools.h:74
#define N_MBHEX5(x, y, z)
Definition: fem_tools.h:76
#define N_MBEDGE0(x)
edge shape function
Definition: fem_tools.h:105
#define N_MBHEX4(x, y, z)
Definition: fem_tools.h:75
PetscErrorCode ShapeMBTET(double *N, const double *G_X, const double *G_Y, const double *G_Z, int DIM)
calculate shape functions
Definition: fem_tools.c:306
#define N_MBHEX0(x, y, z)
Definition: fem_tools.h:71
#define N_MBHEX6(x, y, z)
Definition: fem_tools.h:77
#define N_MBHEX2(x, y, z)
Definition: fem_tools.h:73
#define N_MBQUAD0(x, y)
quad shape function
Definition: fem_tools.h:57
#define N_MBHEX1(x, y, z)
Definition: fem_tools.h:72
#define N_MBQUAD2(x, y)
quad shape function
Definition: fem_tools.h:59
#define N_MBQUAD1(x, y)
quad shape function
Definition: fem_tools.h:58
#define N_MBEDGE1(x)
edge shape function
Definition: fem_tools.h:106
PetscErrorCode ShapeMBTRI(double *N, const double *X, const double *Y, const int G_DIM)
calculate shape functions on triangle
Definition: fem_tools.c:182
constexpr double eta
#define MOFEM_LOG(channel, severity)
Log.
Definition: LogManager.hpp:308
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
Definition: Exceptions.hpp:56
std::bitset< BITREFLEVEL_SIZE > BitRefLevel
Bit structure attached to each entity identifying to what mesh entity is attached.
Definition: Types.hpp:40
implementation of Data Operators for Forces and Sources
Definition: Common.hpp:10
PetscErrorCode PetscOptionsGetInt(PetscOptions *, const char pre[], const char name[], PetscInt *ivalue, PetscBool *set)
PetscErrorCode PetscOptionsGetBool(PetscOptions *, const char pre[], const char name[], PetscBool *bval, PetscBool *set)
double zeta
Definition: plastic.cpp:107
Managing BitRefLevels.
Deprecated interface functions.
Mesh refinement interface.
PetscBool hoNodes
int defMaxLevel
PostProcGenerateRefMeshBase()
MoFEMErrorCode getOptions()
std::string optPrefix
Element for postprocessing. Uses MoAB to generate post-processing mesh.
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface refernce to pointer of interface.