v0.14.0
Loading...
Searching...
No Matches
Public Member Functions | List of all members
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)
 
virtual MoFEMErrorCode generateReferenceElementMesh ()=0
 

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 34 of file PostProcBrokenMeshInMoabBase.cpp.

34 {
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}
#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
#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
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
PetscBool hoNodes
int defMaxLevel
std::vector< MatrixDouble > levelShapeFunctions
std::vector< MatrixDouble > levelGaussPtsOnRefMesh
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface refernce to pointer of interface.

◆ PostProcGenerateRefMeshBase()

MoFEM::PostProcGenerateRefMeshBase::PostProcGenerateRefMeshBase ( )

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