v0.15.0
Loading...
Searching...
No Matches
MoFEM::PostProcGenerateRefMesh< MBTET > Struct Reference

#include "src/post_proc/PostProcBrokenMeshInMoabBase.hpp"

Inheritance diagram for MoFEM::PostProcGenerateRefMesh< MBTET >:
[legend]
Collaboration diagram for MoFEM::PostProcGenerateRefMesh< MBTET >:
[legend]

Public Member Functions

MoFEMErrorCode generateReferenceElementMesh ()
 
 PostProcGenerateRefMeshBase ()
 
- Public Member Functions inherited from MoFEM::PostProcGenerateRefMeshBase
 PostProcGenerateRefMeshBase ()
 
virtual ~PostProcGenerateRefMeshBase ()=default
 
virtual MoFEMErrorCode getOptions (std::string prefix)
 

Additional Inherited Members

- Public Attributes inherited from MoFEM::PostProcGenerateRefMeshBase
std::vector< MatrixDoublelevelShapeFunctions
 
std::vector< MatrixDoublelevelGaussPtsOnRefMesh
 
std::vector< ublas::matrix< int > > levelRef
 
EntityHandle startingVertEleHandle
 
std::vector< double * > verticesOnEleArrays
 
EntityHandle startingEleHandle
 
EntityHandleeleConn
 
int countEle
 
int countVertEle
 
int nbVertices
 
int nbEles
 
PetscBool hoNodes
 
int defMaxLevel
 

Detailed Description

Definition at line 65 of file PostProcBrokenMeshInMoabBase.hpp.

Member Function Documentation

◆ generateReferenceElementMesh()

MoFEMErrorCode MoFEM::PostProcGenerateRefMesh< MBTET >::generateReferenceElementMesh ( )
virtual

Implements MoFEM::PostProcGenerateRefMeshBase.

Definition at line 36 of file PostProcBrokenMeshInMoabBase.cpp.

36 {
38
39 const int max_level = defMaxLevel;
40
41 moab::Core core_ref;
42 moab::Interface &moab_ref = core_ref;
43
44 auto create_reference_element = [&moab_ref]() {
46 constexpr double base_coords[] = {0, 0, 0, 1, 0, 0, 0, 1, 0, 0, 0, 1};
47 EntityHandle nodes[4];
48 for (int nn = 0; nn < 4; nn++) {
49 CHKERR
50 moab_ref.create_vertex(&base_coords[3 * nn], nodes[nn]);
51 }
52 EntityHandle tet;
53 CHKERR moab_ref.create_element(MBTET, nodes, 4, tet);
55 };
56
57 MoFEM::CoreTmp<-1> m_core_ref(moab_ref, PETSC_COMM_SELF, -2);
58 MoFEM::Interface &m_field_ref = m_core_ref;
59
60 auto refine_ref_tetrahedron = [this, &m_field_ref, max_level]() {
62 // seed ref mofem database by setting bit ref level to reference
63 // tetrahedron
64 CHKERR
65 m_field_ref.getInterface<BitRefManager>()->setBitRefLevelByDim(
66 0, 3, BitRefLevel().set(0));
67 for (int ll = 0; ll != max_level; ++ll) {
68 MOFEM_TAG_AND_LOG_C("WORLD", Sev::noisy, "PostProc", "Refine Level %d",
69 ll);
70 Range edges;
71 CHKERR m_field_ref.getInterface<BitRefManager>()
72 ->getEntitiesByTypeAndRefLevel(BitRefLevel().set(ll),
73 BitRefLevel().set(), MBEDGE, edges);
74 Range tets;
75 CHKERR m_field_ref.getInterface<BitRefManager>()
76 ->getEntitiesByTypeAndRefLevel(BitRefLevel().set(ll),
77 BitRefLevel(ll).set(), MBTET, tets);
78 // refine mesh
79 MeshRefinement *m_ref;
80 CHKERR m_field_ref.getInterface(m_ref);
81 CHKERR m_ref->addVerticesInTheMiddleOfEdges(edges,
82 BitRefLevel().set(ll + 1));
83 CHKERR m_ref->refineTets(tets, BitRefLevel().set(ll + 1));
84 }
86 };
87
88 auto get_ref_gauss_pts_and_shape_functions = [this, max_level, &moab_ref,
89 &m_field_ref]() {
91 for (int ll = 0; ll != max_level + 1; ++ll) {
92 Range tets;
93 CHKERR
94 m_field_ref.getInterface<BitRefManager>()->getEntitiesByTypeAndRefLevel(
95 BitRefLevel().set(ll), BitRefLevel().set(ll), MBTET, tets);
96 if (hoNodes) {
97 EntityHandle meshset;
98 CHKERR moab_ref.create_meshset(MESHSET_SET, meshset);
99 CHKERR moab_ref.add_entities(meshset, tets);
100 CHKERR moab_ref.convert_entities(meshset, true, false, false);
101 CHKERR moab_ref.delete_entities(&meshset, 1);
102 }
103 Range elem_nodes;
104 CHKERR moab_ref.get_connectivity(tets, elem_nodes, false);
105
106 auto &gauss_pts = levelGaussPtsOnRefMesh[ll];
107 gauss_pts.resize(elem_nodes.size(), 4, false);
108 std::map<EntityHandle, int> little_map;
109 Range::iterator nit = elem_nodes.begin();
110 for (int gg = 0; nit != elem_nodes.end(); nit++, gg++) {
111 CHKERR moab_ref.get_coords(&*nit, 1, &gauss_pts(gg, 0));
112 little_map[*nit] = gg;
113 }
114 gauss_pts = trans(gauss_pts);
115
116 auto &ref_tets = levelRef[ll];
117 Range::iterator tit = tets.begin();
118 for (int tt = 0; tit != tets.end(); ++tit, ++tt) {
119 const EntityHandle *conn;
120 int num_nodes;
121 CHKERR moab_ref.get_connectivity(*tit, conn, num_nodes, false);
122 if (tt == 0) {
123 ref_tets.resize(tets.size(), num_nodes);
124 }
125 for (int nn = 0; nn != num_nodes; ++nn) {
126 ref_tets(tt, nn) = little_map[conn[nn]];
127 }
128 }
129
130 auto &shape_functions = levelShapeFunctions[ll];
131 shape_functions.resize(elem_nodes.size(), 4);
132 CHKERR ShapeMBTET(&*shape_functions.data().begin(), &gauss_pts(0, 0),
133 &gauss_pts(1, 0), &gauss_pts(2, 0), elem_nodes.size());
134 }
136 };
137
138 levelRef.resize(max_level + 1);
139 levelGaussPtsOnRefMesh.resize(max_level + 1);
140 levelShapeFunctions.resize(max_level + 1);
141
142 CHKERR create_reference_element();
143 CHKERR refine_ref_tetrahedron();
144 CHKERR get_ref_gauss_pts_and_shape_functions();
145
147}
#define MOFEM_TAG_AND_LOG_C(channel, severity, tag, format,...)
Tag and log in channel.
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
#define CHKERR
Inline error check.
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
std::bitset< BITREFLEVEL_SIZE > BitRefLevel
Bit structure attached to each entity identifying to what mesh entity is attached.
Definition Types.hpp:40
Deprecated interface functions.
std::vector< ublas::matrix< int > > levelRef
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.

◆ PostProcGenerateRefMeshBase()

MoFEM::PostProcGenerateRefMeshBase::PostProcGenerateRefMeshBase ( )

The documentation for this struct was generated from the following files: