7static char help[] =
"testing mesh refinement algorithm\n\n";
9int main(
int argc,
char *argv[]) {
15 PetscBool flg = PETSC_TRUE;
20 if (flg != PETSC_TRUE) {
22 "*** ERROR -my_file (MESH FILE NEEDED)");
25 moab::Core mb_instance;
26 moab::Interface &moab = mb_instance;
35 auto get_dim = [&]() {
39 0, 3, nb_ents_3d,
true);
45 0, 2, nb_ents_2d,
true);
61 if (
dim == 3 ||
dim == 2) {
66 "Dimension not handled by test");
72 auto refine_edges = [&](
auto bit0,
auto bit) {
80 Range edges_to_refine;
81 CHKERR moab.get_entities_by_type(*meshset_ptr, MBEDGE, edges_to_refine);
83 for (Range::iterator eit = edges_to_refine.begin();
84 eit != edges_to_refine.end(); eit++, ii++) {
87 CHKERR moab.add_entities(*meshset_ref_edges_ptr, &*eit, 1);
90 CHKERR refine->addVerticesInTheMiddleOfEdges(*meshset_ref_edges_ptr,
bit,
94 }
else if (
dim == 2) {
98 "Dimension not handled by test");
103 auto refine_ents_hanging_nodes = [&](
auto bit0,
auto bit) {
112 CHKERR moab.get_entities_by_dimension(*meshset_ptr,
dim, ents_dim);
114 for (Range::iterator eit = ents_dim.begin(); eit != ents_dim.end();
118 std::vector<EntityHandle> adj_ents;
119 CHKERR moab.get_adjacencies(&*eit, 1, 1,
false, adj_ents);
120 CHKERR moab.add_entities(*meshset_ref_edges_ptr, &*adj_ents.begin(),
125 CHKERR refine->addVerticesInTheMiddleOfEdges(*meshset_ref_edges_ptr,
bit,
126 false,
QUIET, 10000);
128 CHKERR refine->refineTetsHangingNodes(*meshset_ptr,
bit,
QUIET,
true);
129 }
else if (
dim == 2) {
130 CHKERR refine->refineTrisHangingNodes(*meshset_ptr,
bit,
QUIET,
true);
133 "Dimension not handled by test");
138 auto save_blessed_field = [&](
auto bit) {
141 std::ofstream myfile;
142 myfile.open(
"mesh_refine.txt");
148 CHKERR moab.get_entities_by_handle(*out_meshset_tet_ptr, tets);
151 for (Range::iterator tit = tets.begin(); tit != tets.end(); tit++) {
154 CHKERR moab.get_connectivity(*tit, conn, num_nodes,
true);
156 for (
int nn = 0; nn < num_nodes; nn++) {
158 myfile << conn[nn] <<
" ";
172 auto save_vtk = [&](
auto bit) {
177 CHKERR moab.write_file(
"out_mesh_refine.vtk",
"VTK",
"",
178 out_meshset_tet_ptr->get_ptr(), 1);
#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
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
#define CHKERR
Inline error check.
static MoFEMErrorCodeGeneric< PetscErrorCode > ierr
std::bitset< BITREFLEVEL_SIZE > BitRefLevel
Bit structure attached to each entity identifying to what mesh entity is attached.
implementation of Data Operators for Forces and Sources
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.
virtual moab::Interface & get_moab()=0
static MoFEMErrorCode Initialize(int *argc, char ***args, const char file[], const char help[])
Initializes the MoFEM database PETSc, MOAB and MPI.
static MoFEMErrorCode Finalize()
Checks for options to be called at the conclusion of the program.
Deprecated interface functions.
Mesh refinement interface.
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface refernce to pointer of interface.