v0.14.0
node_merge.cpp

test node merging

Date
2022-12-11
/**
* @file node_merge.cpp
* @example node_merge.cpp
* @brief test node merging
* @date 2022-12-11
*
* @copyright Copyright (c) 2022
*
*/
#include <MoFEM.hpp>
using namespace MoFEM;
static char help[] = "...\n\n";
static int debug = 1;
int main(int argc, char *argv[]) {
MoFEM::Core::Initialize(&argc, &argv, (char *)0, help);
try {
moab::Core mb_instance;
moab::Interface &moab = mb_instance;
int rank;
MPI_Comm_rank(PETSC_COMM_WORLD, &rank);
// Read parameters from line command
PetscBool flg = PETSC_TRUE;
char mesh_file_name[255];
#if PETSC_VERSION_GE(3, 6, 4)
CHKERR PetscOptionsGetString(PETSC_NULL, "", "-my_file", mesh_file_name,
255, &flg);
#else
CHKERR PetscOptionsGetString(PETSC_NULL, PETSC_NULL, "-my_file",
mesh_file_name, 255, &flg);
#endif
if (flg != PETSC_TRUE) {
SETERRQ(PETSC_COMM_SELF, 1, "*** ERROR -my_file (MESH FILE NEEDED)");
}
// Read mesh to MOAB
const char *option;
option = ""; //"PARALLEL=BCAST;";//;DEBUG_IO";
CHKERR moab.load_file(mesh_file_name, 0, option);
// Create MoFEM (Joseph) databas
MoFEM::Core core(moab);
MoFEM::Interface &m_field = core;
BitRefLevel bit_level0;
bit_level0.set(0);
CHKERR m_field.getInterface<BitRefManager>()->setBitRefLevelByDim(
0, 3, bit_level0);
NodeMergerInterface *node_merger_iface;
CHKERR m_field.getInterface(node_merger_iface);
int ii = 1;
for (; ii < 2; ii++) {
BitRefLevel bit_level1;
bit_level1.set(ii);
Range edges;
// CHKERR moab.get_entities_by_type(0,MBEDGE,edges,false);
->getEntitiesByTypeAndRefLevel(bit_level0, BitRefLevel().set(),
MBEDGE, edges);
Range::iterator eit = edges.begin();
const EntityHandle *conn;
int num_nodes;
CHKERR moab.get_connectivity(*eit, conn, num_nodes, true);
CHKERR node_merger_iface->mergeNodes(conn[0], conn[1], bit_level1,
bit_level0);
EntityHandle meshset_level1;
CHKERR moab.create_meshset(MESHSET_SET, meshset_level1);
->getEntitiesByTypeAndRefLevel(bit_level1, BitRefLevel().set(), MBTET,
meshset_level1);
std::ostringstream ss;
ss << "node_merger_" << ii << ".vtk";
if (debug)
CHKERR moab.write_file(ss.str().c_str(), "VTK", "", &meshset_level1, 1);
bit_level0 = bit_level1;
}
Range tets;
CHKERR m_field.getInterface<BitRefManager>()->getEntitiesByTypeAndRefLevel(
BitRefLevel().set(ii - 1), BitRefLevel().set(), MBTET, tets);
std::cout << tets << std::endl;
if (tets.size() != 10) {
SETERRQ1(PETSC_COMM_SELF, MOFEM_ATOM_TEST_INVALID,
"diffrent number of tets than expected = %u", tets.size());
}
}
return 0;
}
MoFEM::UnknownInterface::getInterface
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.
Definition: UnknownInterface.hpp:93
MoFEM::CoreTmp< 0 >
Core (interface) class.
Definition: Core.hpp:82
EntityHandle
MoFEM::NodeMergerInterface::mergeNodes
MoFEMErrorCode mergeNodes(EntityHandle father, EntityHandle mother, Range &out_tets, Range *tets_ptr=NULL, const bool only_if_improve_quality=false, const double move=0, const int line_search=0, Tag th=NULL, const int verb=0)
merge nodes which sharing edge
Definition: NodeMerger.cpp:37
MoFEM.hpp
MoFEM::CoreTmp< 0 >::Finalize
static MoFEMErrorCode Finalize()
Checks for options to be called at the conclusion of the program.
Definition: Core.cpp:112
help
static char help[]
Definition: node_merge.cpp:15
MoFEM::NodeMergerInterface
Merge node by collapsing edge between them.
Definition: NodeMerger.hpp:18
MoFEM::DeprecatedCoreInterface
Deprecated interface functions.
Definition: DeprecatedCoreInterface.hpp:16
MoFEM::Interface
DeprecatedCoreInterface Interface
Definition: Interface.hpp:2010
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:548
MoFEM
implementation of Data Operators for Forces and Sources
Definition: Common.hpp:10
mesh_file_name
char mesh_file_name[255]
Definition: mesh_smoothing.cpp:22
Range
MoFEM::CoreTmp< 0 >::Initialize
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
CATCH_ERRORS
#define CATCH_ERRORS
Catch errors.
Definition: definitions.h:385
MoFEM::Core
CoreTmp< 0 > Core
Definition: Core.hpp:1148
MoFEM::PetscOptionsGetString
PetscErrorCode PetscOptionsGetString(PetscOptions *, const char pre[], const char name[], char str[], size_t size, PetscBool *set)
Definition: DeprecatedPetsc.hpp:172
MoFEM::BitRefManager
Managing BitRefLevels.
Definition: BitRefManager.hpp:21
MoFEM::Types::BitRefLevel
std::bitset< BITREFLEVEL_SIZE > BitRefLevel
Bit structure attached to each entity identifying to what mesh entity is attached.
Definition: Types.hpp:40
MOFEM_ATOM_TEST_INVALID
@ MOFEM_ATOM_TEST_INVALID
Definition: definitions.h:40
main
int main(int argc, char *argv[])
Definition: node_merge.cpp:18
debug
static int debug
Definition: node_merge.cpp:16