v0.14.0
Functions | Variables
node_merge.cpp File Reference
#include <MoFEM.hpp>

Go to the source code of this file.

Functions

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

Variables

static char help [] = "...\n\n"
 
static int debug = 1
 

Function Documentation

◆ main()

int main ( int  argc,
char *  argv[] 
)
Examples
node_merge.cpp.

Definition at line 18 of file node_merge.cpp.

18  {
19 
20  MoFEM::Core::Initialize(&argc, &argv, (char *)0, help);
21 
22  try {
23 
24  moab::Core mb_instance;
25  moab::Interface &moab = mb_instance;
26  int rank;
27  MPI_Comm_rank(PETSC_COMM_WORLD, &rank);
28 
29  // Read parameters from line command
30  PetscBool flg = PETSC_TRUE;
31  char mesh_file_name[255];
32 #if PETSC_VERSION_GE(3, 6, 4)
33  CHKERR PetscOptionsGetString(PETSC_NULL, "", "-my_file", mesh_file_name,
34  255, &flg);
35 #else
36  CHKERR PetscOptionsGetString(PETSC_NULL, PETSC_NULL, "-my_file",
37  mesh_file_name, 255, &flg);
38 #endif
39  if (flg != PETSC_TRUE) {
40  SETERRQ(PETSC_COMM_SELF, 1, "*** ERROR -my_file (MESH FILE NEEDED)");
41  }
42 
43  // Read mesh to MOAB
44  const char *option;
45  option = ""; //"PARALLEL=BCAST;";//;DEBUG_IO";
46  CHKERR moab.load_file(mesh_file_name, 0, option);
47 
48  // Create MoFEM (Joseph) databas
49  MoFEM::Core core(moab);
50  MoFEM::Interface &m_field = core;
51 
52  BitRefLevel bit_level0;
53  bit_level0.set(0);
54  CHKERR m_field.getInterface<BitRefManager>()->setBitRefLevelByDim(
55  0, 3, bit_level0);
56 
57  NodeMergerInterface *node_merger_iface;
58  CHKERR m_field.getInterface(node_merger_iface);
59 
60  int ii = 1;
61  for (; ii < 2; ii++) {
62 
63  BitRefLevel bit_level1;
64  bit_level1.set(ii);
65 
66  Range edges;
67  // CHKERR moab.get_entities_by_type(0,MBEDGE,edges,false);
69  ->getEntitiesByTypeAndRefLevel(bit_level0, BitRefLevel().set(),
70  MBEDGE, edges);
71  Range::iterator eit = edges.begin();
72 
73  const EntityHandle *conn;
74  int num_nodes;
75  CHKERR moab.get_connectivity(*eit, conn, num_nodes, true);
76  CHKERR node_merger_iface->mergeNodes(conn[0], conn[1], bit_level1,
77  bit_level0);
78 
79  EntityHandle meshset_level1;
80  CHKERR moab.create_meshset(MESHSET_SET, meshset_level1);
82  ->getEntitiesByTypeAndRefLevel(bit_level1, BitRefLevel().set(), MBTET,
83  meshset_level1);
84 
85  std::ostringstream ss;
86  ss << "node_merger_" << ii << ".vtk";
87 
88  if (debug)
89  CHKERR moab.write_file(ss.str().c_str(), "VTK", "", &meshset_level1, 1);
90  bit_level0 = bit_level1;
91  }
92 
93  Range tets;
94  CHKERR m_field.getInterface<BitRefManager>()->getEntitiesByTypeAndRefLevel(
95  BitRefLevel().set(ii - 1), BitRefLevel().set(), MBTET, tets);
96 
97  std::cout << tets << std::endl;
98  if (tets.size() != 10) {
99  SETERRQ1(PETSC_COMM_SELF, MOFEM_ATOM_TEST_INVALID,
100  "diffrent number of tets than expected = %u", tets.size());
101  }
102  }
103  CATCH_ERRORS;
104 
106  return 0;
107 }

Variable Documentation

◆ debug

int debug = 1
static
Examples
node_merge.cpp.

Definition at line 16 of file node_merge.cpp.

◆ help

char help[] = "...\n\n"
static
Examples
node_merge.cpp.

Definition at line 15 of file node_merge.cpp.

MoFEM::UnknownInterface::getInterface
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface refernce 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::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:1975
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:535
mesh_file_name
char mesh_file_name[255]
Definition: mesh_smoothing.cpp:23
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:372
MoFEM::Core
CoreTmp< 0 > Core
Definition: Core.hpp:1094
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
debug
static int debug
Definition: node_merge.cpp:16