v0.14.0
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), countEle(0), countVertEle(0),
15 nbVertices(0) {}
16
19
20 CHKERR PetscOptionsGetInt(prefix.size() ? prefix.c_str() : PETSC_NULL,
21 "-max_post_proc_ref_level", &defMaxLevel,
22 PETSC_NULL);
23 CHKERR PetscOptionsGetBool(prefix.size() ? prefix.c_str() : PETSC_NULL,
24 "-max_post_ho_nodes", &hoNodes, PETSC_NULL);
25
26 if (defMaxLevel < 0)
27 SETERRQ(PETSC_COMM_WORLD, MOFEM_INVALID_DATA,
28 "Wrong parameter -max_post_proc_ref_level "
29 "should be positive number");
30
32};
33
36
37 const int max_level = defMaxLevel;
38
39 moab::Core core_ref;
40 moab::Interface &moab_ref = core_ref;
41
42 auto create_reference_element = [&moab_ref]() {
44 constexpr double base_coords[] = {0, 0, 0, 1, 0, 0, 0, 1, 0, 0, 0, 1};
45 EntityHandle nodes[4];
46 for (int nn = 0; nn < 4; nn++) {
47 CHKERR
48 moab_ref.create_vertex(&base_coords[3 * nn], nodes[nn]);
49 }
50 EntityHandle tet;
51 CHKERR moab_ref.create_element(MBTET, nodes, 4, tet);
53 };
54
55 MoFEM::CoreTmp<-1> m_core_ref(moab_ref, PETSC_COMM_SELF, -2);
56 MoFEM::Interface &m_field_ref = m_core_ref;
57
58 auto refine_ref_tetrahedron = [this, &m_field_ref, max_level]() {
60 // seed ref mofem database by setting bit ref level to reference
61 // tetrahedron
62 CHKERR
63 m_field_ref.getInterface<BitRefManager>()->setBitRefLevelByDim(
64 0, 3, BitRefLevel().set(0));
65 for (int ll = 0; ll != max_level; ++ll) {
66 MOFEM_TAG_AND_LOG_C("WORLD", Sev::noisy, "PostProc", "Refine Level %d",
67 ll);
68 Range edges;
69 CHKERR m_field_ref.getInterface<BitRefManager>()
70 ->getEntitiesByTypeAndRefLevel(BitRefLevel().set(ll),
71 BitRefLevel().set(), MBEDGE, edges);
72 Range tets;
73 CHKERR m_field_ref.getInterface<BitRefManager>()
74 ->getEntitiesByTypeAndRefLevel(BitRefLevel().set(ll),
75 BitRefLevel(ll).set(), MBTET, tets);
76 // refine mesh
77 MeshRefinement *m_ref;
78 CHKERR m_field_ref.getInterface(m_ref);
79 CHKERR m_ref->addVerticesInTheMiddleOfEdges(edges,
80 BitRefLevel().set(ll + 1));
81 CHKERR m_ref->refineTets(tets, BitRefLevel().set(ll + 1));
82 }
84 };
85
86 auto get_ref_gauss_pts_and_shape_functions = [this, max_level, &moab_ref,
87 &m_field_ref]() {
89 for (int ll = 0; ll != max_level + 1; ++ll) {
90 Range tets;
91 CHKERR
92 m_field_ref.getInterface<BitRefManager>()->getEntitiesByTypeAndRefLevel(
93 BitRefLevel().set(ll), BitRefLevel().set(ll), MBTET, tets);
94 if (hoNodes) {
95 EntityHandle meshset;
96 CHKERR moab_ref.create_meshset(MESHSET_SET, meshset);
97 CHKERR moab_ref.add_entities(meshset, tets);
98 CHKERR moab_ref.convert_entities(meshset, true, false, false);
99 CHKERR moab_ref.delete_entities(&meshset, 1);
100 }
101 Range elem_nodes;
102 CHKERR moab_ref.get_connectivity(tets, elem_nodes, false);
103
104 auto &gauss_pts = levelGaussPtsOnRefMesh[ll];
105 gauss_pts.resize(elem_nodes.size(), 4, false);
106 std::map<EntityHandle, int> little_map;
107 Range::iterator nit = elem_nodes.begin();
108 for (int gg = 0; nit != elem_nodes.end(); nit++, gg++) {
109 CHKERR moab_ref.get_coords(&*nit, 1, &gauss_pts(gg, 0));
110 little_map[*nit] = gg;
111 }
112 gauss_pts = trans(gauss_pts);
113
114 auto &ref_tets = levelRef[ll];
115 Range::iterator tit = tets.begin();
116 for (int tt = 0; tit != tets.end(); ++tit, ++tt) {
117 const EntityHandle *conn;
118 int num_nodes;
119 CHKERR moab_ref.get_connectivity(*tit, conn, num_nodes, false);
120 if (tt == 0) {
121 ref_tets.resize(tets.size(), num_nodes);
122 }
123 for (int nn = 0; nn != num_nodes; ++nn) {
124 ref_tets(tt, nn) = little_map[conn[nn]];
125 }
126 }
127
128 auto &shape_functions = levelShapeFunctions[ll];
129 shape_functions.resize(elem_nodes.size(), 4);
130 CHKERR ShapeMBTET(&*shape_functions.data().begin(), &gauss_pts(0, 0),
131 &gauss_pts(1, 0), &gauss_pts(2, 0), elem_nodes.size());
132 }
134 };
135
136 levelRef.resize(max_level + 1);
137 levelGaussPtsOnRefMesh.resize(max_level + 1);
138 levelShapeFunctions.resize(max_level + 1);
139
140 CHKERR create_reference_element();
141 CHKERR refine_ref_tetrahedron();
142 CHKERR get_ref_gauss_pts_and_shape_functions();
143
145}
146
149
150#ifndef NDEBUG
151 if (defMaxLevel > 0)
152 MOFEM_LOG("WORLD", Sev::warning)
153 << "Refinement for hexes is not implemented";
154#endif
155
156 moab::Core core_ref;
157 moab::Interface &moab_ref = core_ref;
158
159 auto create_reference_element = [&moab_ref]() {
161 constexpr double base_coords[] = {0, 0, 0, 1, 0, 0, 1, 1, 0, 0, 1, 0,
162 0, 0, 1, 1, 0, 1, 1, 1, 1, 0, 1, 1};
163 EntityHandle nodes[8];
164 for (int nn = 0; nn < 8; nn++)
165 CHKERR moab_ref.create_vertex(&base_coords[3 * nn], nodes[nn]);
166 EntityHandle hex;
167 CHKERR moab_ref.create_element(MBHEX, nodes, 8, hex);
169 };
170
171 auto add_ho_nodes = [&]() {
173 Range hexes;
174 CHKERR moab_ref.get_entities_by_type(0, MBHEX, hexes, true);
175 EntityHandle meshset;
176 CHKERR moab_ref.create_meshset(MESHSET_SET, meshset);
177 CHKERR moab_ref.add_entities(meshset, hexes);
178 CHKERR moab_ref.convert_entities(meshset, true, true, true);
179 CHKERR moab_ref.delete_entities(&meshset, 1);
181 };
182
183 auto set_gauss_pts = [&](std::map<EntityHandle, int> &little_map) {
185 Range hexes;
186 CHKERR moab_ref.get_entities_by_type(0, MBHEX, hexes, true);
187 Range hexes_nodes;
188 CHKERR moab_ref.get_connectivity(hexes, hexes_nodes, false);
189 auto &gauss_pts = levelGaussPtsOnRefMesh[0];
190 gauss_pts.resize(hexes_nodes.size(), 4, false);
191 size_t gg = 0;
192 for (auto node : hexes_nodes) {
193 CHKERR moab_ref.get_coords(&node, 1, &gauss_pts(gg, 0));
194 little_map[node] = gg;
195 ++gg;
196 }
197 gauss_pts = trans(gauss_pts);
199 };
200
201 auto set_ref_hexes = [&](std::map<EntityHandle, int> &little_map) {
203 Range hexes;
204 CHKERR moab_ref.get_entities_by_type(0, MBHEX, hexes, true);
205 size_t hh = 0;
206 auto &ref_hexes = levelRef[0];
207 for (auto hex : hexes) {
208 const EntityHandle *conn;
209 int num_nodes;
210 CHKERR moab_ref.get_connectivity(hex, conn, num_nodes, false);
211 if (ref_hexes.size2() != num_nodes) {
212 ref_hexes.resize(hexes.size(), num_nodes);
213 }
214 for (int nn = 0; nn != num_nodes; ++nn) {
215 ref_hexes(hh, nn) = little_map[conn[nn]];
216 }
217 ++hh;
218 }
220 };
221
222 auto set_shape_functions = [&]() {
224 auto &gauss_pts = levelGaussPtsOnRefMesh[0];
225 auto &shape_functions = levelShapeFunctions[0];
226 const auto nb_gauss_pts = gauss_pts.size2();
227 shape_functions.resize(nb_gauss_pts, 8);
228 for (int gg = 0; gg != nb_gauss_pts; ++gg) {
229 const double ksi = gauss_pts(0, gg);
230 const double zeta = gauss_pts(1, gg);
231 const double eta = gauss_pts(2, gg);
232 shape_functions(gg, 0) = N_MBHEX0(ksi, zeta, eta);
233 shape_functions(gg, 1) = N_MBHEX1(ksi, zeta, eta);
234 shape_functions(gg, 2) = N_MBHEX2(ksi, zeta, eta);
235 shape_functions(gg, 3) = N_MBHEX3(ksi, zeta, eta);
236 shape_functions(gg, 4) = N_MBHEX4(ksi, zeta, eta);
237 shape_functions(gg, 5) = N_MBHEX5(ksi, zeta, eta);
238 shape_functions(gg, 6) = N_MBHEX6(ksi, zeta, eta);
239 shape_functions(gg, 7) = N_MBHEX7(ksi, zeta, eta);
240 }
242 };
243
244 levelRef.resize(1);
245 levelGaussPtsOnRefMesh.resize(1);
246 levelShapeFunctions.resize(1);
247
248 CHKERR create_reference_element();
249 if (hoNodes)
250 CHKERR add_ho_nodes();
251 std::map<EntityHandle, int> little_map;
252 CHKERR set_gauss_pts(little_map);
253 CHKERR set_ref_hexes(little_map);
254 CHKERR set_shape_functions();
255
257}
258
261
262 const int max_level = defMaxLevel;
263
264 moab::Core core_ref;
265 moab::Interface &moab_ref = core_ref;
266
267 auto create_reference_element = [&moab_ref]() {
269 constexpr double base_coords[] = {
270
271 0, 0,
272 0,
273
274 1, 0,
275 0,
276
277 0, 1,
278 0
279
280 };
281 EntityHandle nodes[3];
282 for (int nn = 0; nn != 3; ++nn)
283 CHKERR moab_ref.create_vertex(&base_coords[3 * nn], nodes[nn]);
284 EntityHandle tri;
285 CHKERR moab_ref.create_element(MBTRI, nodes, 3, tri);
286
287 Range edges;
288 CHKERR moab_ref.get_adjacencies(&tri, 1, 1, true, edges,
289 moab::Interface::UNION);
290
292 };
293
294 CHKERR create_reference_element();
295
296 MoFEM::CoreTmp<-1> m_core_ref(moab_ref, PETSC_COMM_SELF, -2);
297 MoFEM::Interface &m_field_ref = m_core_ref;
298
299 auto refine_ref_triangles = [this, &m_field_ref, max_level]() {
301 // seed ref mofem database by setting bit ref level to reference
302 // tetrahedron
303 CHKERR
304 m_field_ref.getInterface<BitRefManager>()->setBitRefLevelByDim(
305 0, 2, BitRefLevel().set(0));
306
307 for (int ll = 0; ll != max_level; ++ll) {
308 MOFEM_TAG_AND_LOG_C("WORLD", Sev::noisy, "PostProc", "Refine Level %d",
309 ll);
310 Range edges;
311 CHKERR m_field_ref.getInterface<BitRefManager>()
312 ->getEntitiesByTypeAndRefLevel(BitRefLevel().set(ll),
313 BitRefLevel().set(), MBEDGE, edges);
314 Range tris;
315 CHKERR m_field_ref.getInterface<BitRefManager>()
316 ->getEntitiesByTypeAndRefLevel(BitRefLevel().set(ll),
317 BitRefLevel(ll).set(), MBTRI, tris);
318 // refine mesh
319 auto m_ref = m_field_ref.getInterface<MeshRefinement>();
320 CHKERR m_ref->addVerticesInTheMiddleOfEdges(edges,
321 BitRefLevel().set(ll + 1));
322 CHKERR m_ref->refineTris(tris, BitRefLevel().set(ll + 1));
323 }
325 };
326
327 // auto set_gauss_pts = [&](std::map<EntityHandle, int> &little_map) {
328 // MoFEMFunctionBegin;
329 // Range faces;
330 // CHKERR moab_ref.get_entities_by_type(0, MBTRI, faces, true);
331 // Range faces_nodes;
332 // CHKERR moab_ref.get_connectivity(faces, faces_nodes, false);
333 // auto &gauss_pts = levelGaussPtsOnRefMesh[0];
334 // gauss_pts.resize(faces_nodes.size(), 4, false);
335 // size_t gg = 0;
336 // for (auto node : faces_nodes) {
337 // CHKERR moab_ref.get_coords(&node, 1, &gauss_pts(gg, 0));
338 // little_map[node] = gg;
339 // ++gg;
340 // }
341 // gauss_pts = trans(gauss_pts);
342 // MoFEMFunctionReturn(0);
343 // };
344
345 auto get_ref_gauss_pts_and_shape_functions = [this, max_level, &moab_ref,
346 &m_field_ref]() {
348 for (int ll = 0; ll != max_level + 1; ++ll) {
349
350 Range tris;
351 CHKERR
352 m_field_ref.getInterface<BitRefManager>()->getEntitiesByTypeAndRefLevel(
353 BitRefLevel().set(ll), BitRefLevel().set(ll), MBTRI, tris);
354
355 if (hoNodes) {
356 EntityHandle meshset;
357 CHKERR moab_ref.create_meshset(MESHSET_SET, meshset);
358 CHKERR moab_ref.add_entities(meshset, tris);
359 CHKERR moab_ref.convert_entities(meshset, true, false, false);
360 CHKERR moab_ref.delete_entities(&meshset, 1);
361 }
362
363 Range elem_nodes;
364 CHKERR moab_ref.get_connectivity(tris, elem_nodes, false);
365
366 auto &gauss_pts = levelGaussPtsOnRefMesh[ll];
367 gauss_pts.resize(elem_nodes.size(), 3, false);
368 std::map<EntityHandle, int> little_map;
369 Range::iterator nit = elem_nodes.begin();
370 for (int gg = 0; nit != elem_nodes.end(); nit++, gg++) {
371 CHKERR moab_ref.get_coords(&*nit, 1, &gauss_pts(gg, 0));
372 little_map[*nit] = gg;
373 }
374 gauss_pts = trans(gauss_pts);
375
376 auto &ref_tris = levelRef[ll];
377 Range::iterator tit = tris.begin();
378 for (int tt = 0; tit != tris.end(); ++tit, ++tt) {
379 const EntityHandle *conn;
380 int num_nodes;
381 CHKERR moab_ref.get_connectivity(*tit, conn, num_nodes, false);
382 if (tt == 0) {
383 ref_tris.resize(tris.size(), num_nodes);
384 }
385 for (int nn = 0; nn != num_nodes; ++nn) {
386 ref_tris(tt, nn) = little_map[conn[nn]];
387 }
388 }
389
390 auto &shape_functions = levelShapeFunctions[ll];
391 shape_functions.resize(elem_nodes.size(), 3);
392 CHKERR ShapeMBTRI(&*shape_functions.data().begin(), &gauss_pts(0, 0),
393 &gauss_pts(1, 0), elem_nodes.size());
394 }
396 };
397
398 levelRef.resize(max_level + 1);
399 levelGaussPtsOnRefMesh.resize(max_level + 1);
400 levelShapeFunctions.resize(max_level + 1);
401
402 CHKERR refine_ref_triangles();
403 CHKERR get_ref_gauss_pts_and_shape_functions();
404
406}
407
410
411#ifndef NDEBUG
412 if (defMaxLevel > 0)
413 MOFEM_LOG("WORLD", Sev::warning)
414 << "Refinement for quad is not implemented";
415#endif
416
417 moab::Core core_ref;
418 moab::Interface &moab_ref = core_ref;
419
420 auto create_reference_element = [&moab_ref]() {
422 constexpr double base_coords[] = {
423
424 0, 0,
425 0,
426
427 1, 0,
428 0,
429
430 1, 1,
431 0,
432
433 0, 1,
434 0
435
436 };
437 EntityHandle nodes[4];
438 for (int nn = 0; nn < 4; nn++)
439 CHKERR moab_ref.create_vertex(&base_coords[3 * nn], nodes[nn]);
440 EntityHandle quad;
441 CHKERR moab_ref.create_element(MBQUAD, nodes, 4, quad);
443 };
444
445 auto add_ho_nodes = [&]() {
447 Range quads;
448 CHKERR moab_ref.get_entities_by_type(0, MBQUAD, quads, true);
449 EntityHandle meshset;
450 CHKERR moab_ref.create_meshset(MESHSET_SET, meshset);
451 CHKERR moab_ref.add_entities(meshset, quads);
452 CHKERR moab_ref.convert_entities(meshset, true, true, true);
453 CHKERR moab_ref.delete_entities(&meshset, 1);
455 };
456
457 auto set_gauss_pts = [&](std::map<EntityHandle, int> &little_map) {
459 Range quads;
460 CHKERR moab_ref.get_entities_by_type(0, MBQUAD, quads, true);
461 Range quads_nodes;
462 CHKERR moab_ref.get_connectivity(quads, quads_nodes, false);
463 auto &gauss_pts = levelGaussPtsOnRefMesh[0];
464 gauss_pts.resize(quads_nodes.size(), 4, false);
465 size_t gg = 0;
466 for (auto node : quads_nodes) {
467 CHKERR moab_ref.get_coords(&node, 1, &gauss_pts(gg, 0));
468 little_map[node] = gg;
469 ++gg;
470 }
471 gauss_pts = trans(gauss_pts);
473 };
474
475 auto set_ref_quads = [&](std::map<EntityHandle, int> &little_map) {
477 Range quads;
478 CHKERR moab_ref.get_entities_by_type(0, MBQUAD, quads, true);
479 size_t hh = 0;
480 auto &ref_quads = levelRef[0];
481 for (auto quad : quads) {
482 const EntityHandle *conn;
483 int num_nodes;
484 CHKERR moab_ref.get_connectivity(quad, conn, num_nodes, false);
485 if (ref_quads.size2() != num_nodes) {
486 ref_quads.resize(quads.size(), num_nodes);
487 }
488 for (int nn = 0; nn != num_nodes; ++nn) {
489 ref_quads(hh, nn) = little_map[conn[nn]];
490 }
491 ++hh;
492 }
494 };
495
496 auto set_shape_functions = [&]() {
498 auto &gauss_pts = levelGaussPtsOnRefMesh[0];
499 auto &shape_functions = levelShapeFunctions[0];
500 const auto nb_gauss_pts = gauss_pts.size2();
501 shape_functions.resize(nb_gauss_pts, 4);
502 for (int gg = 0; gg != nb_gauss_pts; ++gg) {
503 const double ksi = gauss_pts(0, gg);
504 const double zeta = gauss_pts(1, gg);
505 shape_functions(gg, 0) = N_MBQUAD0(ksi, zeta);
506 shape_functions(gg, 1) = N_MBQUAD1(ksi, zeta);
507 shape_functions(gg, 2) = N_MBQUAD2(ksi, zeta);
508 shape_functions(gg, 3) = N_MBQUAD3(ksi, zeta);
509 }
511 };
512
513 levelRef.resize(1);
514 levelGaussPtsOnRefMesh.resize(1);
515 levelShapeFunctions.resize(1);
516
517 CHKERR create_reference_element();
518 if (hoNodes)
519 CHKERR add_ho_nodes();
520 std::map<EntityHandle, int> little_map;
521 CHKERR set_gauss_pts(little_map);
522 CHKERR set_ref_quads(little_map);
523 CHKERR set_shape_functions();
524
526}
527
530
531#ifndef NDEBUG
532 if (defMaxLevel > 0)
533 MOFEM_LOG("WORLD", Sev::warning)
534 << "Refinement for edges is not implemented";
535#endif
536
537 auto set_gauss_pts = [&](std::map<EntityHandle, int> &little_map) {
539
540 int nb_nodes = 2;
541 if (hoNodes)
542 nb_nodes = 3;
543
544 auto &gauss_pts = levelGaussPtsOnRefMesh[0];
545 gauss_pts.resize(2, nb_nodes, false);
546 gauss_pts.clear();
547
548 int nn = 0;
549 for (; nn != 2; ++nn) {
550 gauss_pts(0, nn) = static_cast<double>(nn);
551 little_map[nn] = nn;
552 }
553
554 if (nn < nb_nodes) {
555 gauss_pts(0, nn) = 0.5;
556 little_map[nn] = 2;
557 }
558
560 };
561
562 auto set_ref_edges = [&](std::map<EntityHandle, int> &little_map) {
564
565 int level = 0;
566 int nb_edges = level + 1;
567
568 int nb_nodes = 2;
569 if (hoNodes)
570 nb_nodes = 3;
571
572 auto &ref_edges = levelRef[level];
573 ref_edges.resize(nb_edges, nb_nodes, false);
574
575 for (int ee = 0; ee != nb_edges; ++ee) {
576 int nn = 0;
577 for (; nn != 2; ++nn) {
578 ref_edges(ee, nn) = nb_nodes * ee + nn;
579 }
580 if (nn < nb_nodes) {
581 ref_edges(ee, nn) = nb_nodes * ee + 2;
582 }
583 }
584
586 };
587
588 auto set_shape_functions = [&]() {
590 auto &gauss_pts = levelGaussPtsOnRefMesh[0];
591 auto &shape_functions = levelShapeFunctions[0];
592 const auto nb_gauss_pts = gauss_pts.size2();
593 shape_functions.resize(nb_gauss_pts, 2);
594 for (int gg = 0; gg != nb_gauss_pts; ++gg) {
595 const double ksi = gauss_pts(0, gg);
596 shape_functions(gg, 0) = N_MBEDGE0(ksi);
597 shape_functions(gg, 1) = N_MBEDGE1(ksi);
598 }
600 };
601
602 levelRef.resize(1);
603 levelGaussPtsOnRefMesh.resize(1);
604 levelShapeFunctions.resize(1);
605
606 std::map<EntityHandle, int> little_map;
607 CHKERR set_gauss_pts(little_map);
608 CHKERR set_ref_edges(little_map);
609 CHKERR set_shape_functions();
610
612}
613
614} // 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
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
Viscous hardening.
Definition: plastic.cpp:173
Managing BitRefLevels.
Deprecated interface functions.
Mesh refinement interface.
PetscBool hoNodes
int defMaxLevel
virtual MoFEMErrorCode getOptions(std::string prefix)
PostProcGenerateRefMeshBase()
Element for postprocessing. Uses MoAB to generate post-processing mesh.
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface refernce to pointer of interface.