v0.15.0
Loading...
Searching...
No Matches
Functions | Variables
mesh_refine_atom.cpp File Reference

Test mesh refinement by splitting edges. More...

#include <MoFEM.hpp>

Go to the source code of this file.

Functions

int main (int argc, char *argv[])
 

Variables

static char help [] = "testing mesh refinement algorithm\n\n"
 

Detailed Description

Test mesh refinement by splitting edges.

`~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~p

Author
Anonymous author under MoFEM contributors agreement (MiT license)

Definition in file mesh_refine_atom.cpp.

Function Documentation

◆ main()

int main ( int  argc,
char *  argv[] 
)

Definition at line 17 of file mesh_refine_atom.cpp.

17 {
18
19 MoFEM::Core::Initialize(&argc, &argv, (char *)0, help);
20
21 try {
22
23 PetscBool flg = PETSC_TRUE;
24 char mesh_file_name[255];
25 ierr = PetscOptionsGetString(PETSC_NULLPTR, "", "-my_file", mesh_file_name,
26 255, &flg);
28 if (flg != PETSC_TRUE) {
29 SETERRQ(PETSC_COMM_SELF, MOFEM_INVALID_DATA,
30 "*** ERROR -my_file (MESH FILE NEEDED)");
31 }
32
33 moab::Core mb_instance;
34 moab::Interface &moab = mb_instance;
35
36 const char *option;
37 option = "";
38 CHKERR moab.load_file(mesh_file_name, 0, option);
39
40 MoFEM::Core core(moab);
41 MoFEM::Interface &m_field = core;
42
43 auto get_dim = [&]() {
44 int dim;
45 int nb_ents_3d;
46 CHKERR m_field.get_moab().get_number_entities_by_dimension(
47 0, 3, nb_ents_3d, true);
48 if (nb_ents_3d > 0) {
49 dim = 3;
50 } else {
51 int nb_ents_2d;
52 CHKERR m_field.get_moab().get_number_entities_by_dimension(
53 0, 2, nb_ents_2d, true);
54 if (nb_ents_2d > 0) {
55 dim = 2;
56 } else {
57 dim = 1;
58 }
59 }
60 return dim;
61 };
62
63 auto dim = get_dim();
64
65 auto refine = m_field.getInterface<MeshRefinement>();
66
67 BitRefLevel bit_level0;
68 bit_level0.set(0);
69 if (dim == 3 || dim == 2) {
70 CHKERR m_field.getInterface<BitRefManager>()->setBitRefLevelByDim(
71 0, dim, bit_level0);
72 } else {
73 SETERRQ(PETSC_COMM_WORLD, MOFEM_ATOM_TEST_INVALID,
74 "Dimension not handled by test");
75 }
76
77 BitRefLevel bit_level1;
78 bit_level1.set(1);
79
80 auto refine_edges = [&](auto bit0, auto bit) {
82 auto meshset_ptr = get_temp_meshset_ptr(moab);
83 CHKERR m_field.getInterface<BitRefManager>()->getEntitiesByRefLevel(
84 bit0, BitRefLevel().set(), *meshset_ptr);
85
86 // random mesh refinement
87 auto meshset_ref_edges_ptr = get_temp_meshset_ptr(moab);
88 Range edges_to_refine;
89 CHKERR moab.get_entities_by_type(*meshset_ptr, MBEDGE, edges_to_refine);
90 int ii = 0;
91 for (Range::iterator eit = edges_to_refine.begin();
92 eit != edges_to_refine.end(); eit++, ii++) {
93 int numb = ii % 2;
94 if (numb == 0) {
95 CHKERR moab.add_entities(*meshset_ref_edges_ptr, &*eit, 1);
96 }
97 }
98 CHKERR refine->addVerticesInTheMiddleOfEdges(*meshset_ref_edges_ptr, bit,
99 false, QUIET, 10000);
100 if (dim == 3) {
101 CHKERR refine->refineTets(*meshset_ptr, bit, QUIET, true);
102 } else if (dim == 2) {
103 CHKERR refine->refineTris(*meshset_ptr, bit, QUIET, true);
104 } else {
105 SETERRQ(PETSC_COMM_WORLD, MOFEM_ATOM_TEST_INVALID,
106 "Dimension not handled by test");
107 }
109 };
110
111 auto refine_ents_hanging_nodes = [&](auto bit0, auto bit) {
113 auto meshset_ptr = get_temp_meshset_ptr(moab);
114 CHKERR m_field.getInterface<BitRefManager>()->getEntitiesByRefLevel(
115 bit0, BitRefLevel().set(), *meshset_ptr);
116
117 // random mesh refinement
118 auto meshset_ref_edges_ptr = get_temp_meshset_ptr(moab);
119 Range ents_dim;
120 CHKERR moab.get_entities_by_dimension(*meshset_ptr, dim, ents_dim);
121 int ii = 0;
122 for (Range::iterator eit = ents_dim.begin(); eit != ents_dim.end();
123 eit++, ii++) {
124 int numb = ii % 2;
125 if (numb == 0) {
126 std::vector<EntityHandle> adj_ents;
127 CHKERR moab.get_adjacencies(&*eit, 1, 1, false, adj_ents);
128 CHKERR moab.add_entities(*meshset_ref_edges_ptr, &*adj_ents.begin(),
129 adj_ents.size());
130 }
131 }
132
133 CHKERR refine->addVerticesInTheMiddleOfEdges(*meshset_ref_edges_ptr, bit,
134 false, QUIET, 10000);
135 if (dim == 3) {
136 CHKERR refine->refineTetsHangingNodes(*meshset_ptr, bit, QUIET, true);
137 } else if (dim == 2) {
138 CHKERR refine->refineTrisHangingNodes(*meshset_ptr, bit, QUIET, true);
139 } else {
140 SETERRQ(PETSC_COMM_WORLD, MOFEM_ATOM_TEST_INVALID,
141 "Dimension not handled by test");
142 }
144 };
145
146 auto save_blessed_field = [&](auto bit) {
148
149 std::ofstream myfile;
150 myfile.open("mesh_refine.txt");
151
152 auto out_meshset_tet_ptr = get_temp_meshset_ptr(moab);
153 CHKERR m_field.getInterface<BitRefManager>()->getEntitiesByDimAndRefLevel(
154 bit, BitRefLevel().set(), dim, *out_meshset_tet_ptr);
155 Range tets;
156 CHKERR moab.get_entities_by_handle(*out_meshset_tet_ptr, tets);
157 {
158 int ii = 0;
159 for (Range::iterator tit = tets.begin(); tit != tets.end(); tit++) {
160 int num_nodes;
161 const EntityHandle *conn;
162 CHKERR moab.get_connectivity(*tit, conn, num_nodes, true);
163
164 for (int nn = 0; nn < num_nodes; nn++) {
165 // cout << conn[nn] << " ";
166 myfile << conn[nn] << " ";
167 }
168 // cout << std::endl;
169 myfile << std::endl;
170 if (ii > 25)
171 break;
172 }
173 }
174
175 myfile.close();
176
178 };
179
180 auto save_vtk = [&](auto bit) {
182 auto out_meshset_tet_ptr = get_temp_meshset_ptr(moab);
183 CHKERR m_field.getInterface<BitRefManager>()->getEntitiesByDimAndRefLevel(
184 bit, BitRefLevel().set(), dim, *out_meshset_tet_ptr);
185 CHKERR moab.write_file("out_mesh_refine.vtk", "VTK", "",
186 out_meshset_tet_ptr->get_ptr(), 1);
188 };
189
190 CHKERR refine_edges(BitRefLevel().set(0), BitRefLevel().set(1));
191 CHKERR refine_ents_hanging_nodes(BitRefLevel().set(1),
192 BitRefLevel().set(2));
193 CHKERR save_blessed_field(BitRefLevel().set(2));
194 CHKERR save_vtk(BitRefLevel().set(2));
195 }
197
199
200 return 0;
201}
@ QUIET
#define CATCH_ERRORS
Catch errors.
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
#define CHKERRG(n)
Check error code of MoFEM/MOAB/PETSc function.
@ MOFEM_ATOM_TEST_INVALID
Definition definitions.h:40
@ MOFEM_INVALID_DATA
Definition definitions.h:36
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
#define CHKERR
Inline error check.
auto bit
set bit
static char help[]
char mesh_file_name[255]
auto get_dim(const Range &ents)
Definition BcManager.cpp:10
static MoFEMErrorCodeGeneric< PetscErrorCode > ierr
std::bitset< BITREFLEVEL_SIZE > BitRefLevel
Bit structure attached to each entity identifying to what mesh entity is attached.
Definition Types.hpp:40
PetscErrorCode PetscOptionsGetString(PetscOptions *, const char pre[], const char name[], char str[], size_t size, PetscBool *set)
auto get_temp_meshset_ptr(moab::Interface &moab)
Create smart pointer to temporary meshset.
Managing BitRefLevels.
virtual moab::Interface & get_moab()=0
Core (interface) class.
Definition Core.hpp:82
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
static MoFEMErrorCode Finalize()
Checks for options to be called at the conclusion of the program.
Definition Core.cpp:118
Deprecated interface functions.
Mesh refinement interface.
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.

Variable Documentation

◆ help

char help[] = "testing mesh refinement algorithm\n\n"
static

Definition at line 15 of file mesh_refine_atom.cpp.