v0.14.0
uniform_mesh_refinement.cpp
Go to the documentation of this file.
1 /** \file uniform_mesh_refinement.cpp
2 
3  \brief Uniformly refine mesh
4 
5 */
6 
7 #include <MoFEM.hpp>
8 using namespace MoFEM;
9 
10 static char help[] = "Uniform mesh refinement\n\n"
11 
12  "Usage example:\n"
13 
14  "$ ./uniform_mesh_refinement -my_file mesh.h5m "
15  "-output_file refined_mesh.h5m"
16 
17  "\n\n";
18 
19 int main(int argc, char *argv[]) {
20 
21  MoFEM::Core::Initialize(&argc, &argv, (char *)0, help);
22 
23  try {
24 
25  // global variables
26  char mesh_file_name[255];
27  char mesh_out_file[255] = "out.h5m";
28  int dim = 3;
29  int nb_levels = 1;
30  PetscBool shift = PETSC_TRUE;
31  PetscBool flg_file = PETSC_FALSE;
32  PetscBool debug = PETSC_FALSE;
33 
34  ierr = PetscOptionsBegin(PETSC_COMM_WORLD, "", "none", "none");
35  CHKERRQ(ierr);
36 
37  CHKERR PetscOptionsString("-my_file", "mesh file name", "", "mesh.h5m",
38  mesh_file_name, 255, &flg_file);
39  if (flg_file != PETSC_TRUE)
40  CHKERR PetscOptionsString("-file_name", "mesh file name", "", "mesh.h5m",
41  mesh_file_name, 255, &flg_file);
42  CHKERR PetscOptionsString("-output_file", "output mesh file name", "",
43  "mesh.h5m", mesh_out_file, 255, PETSC_NULL);
44  CHKERR PetscOptionsInt("-dim", "entities dim", "", dim, &dim, PETSC_NULL);
45  CHKERR PetscOptionsInt("-nb_levels", "number of refinement levels", "",
46  nb_levels, &nb_levels, PETSC_NULL);
47  CHKERR PetscOptionsBool("-shift",
48  "shift bits, squash entities of refined elements",
49  "", shift, &shift, PETSC_NULL);
50  CHKERR PetscOptionsBool("-debug",
51  "write additional files with bit ref levels", "",
52  debug, &debug, PETSC_NULL);
53 
54  ierr = PetscOptionsEnd();
55  CHKERRQ(ierr);
56 
57  moab::Core mb_instance;
58  moab::Interface &moab = mb_instance;
59  ParallelComm *pcomm = ParallelComm::get_pcomm(&moab, MYPCOMM_INDEX);
60  if (pcomm == NULL)
61  pcomm = new ParallelComm(&moab, PETSC_COMM_WORLD);
62 
63  const char *option;
64  option = "";
65  CHKERR moab.load_file(mesh_file_name, 0, option);
66 
67  // Create MoFEM database
68  MoFEM::Core core(moab);
69  MoFEM::Interface &m_field = core;
70 
71  if (flg_file != PETSC_TRUE) {
72  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
73  "*** ERROR -my_file (-file_name) (mesh file needed)");
74  }
75 
76  BitRefManager *bit_ref_manager;
77  CHKERR m_field.getInterface(bit_ref_manager);
78 
79  auto bit = [](auto l) { return BitRefLevel().set(l); };
80 
81  CHKERR bit_ref_manager->setBitRefLevelByDim(0, dim, bit(0));
82 
83  for (auto l = 0; l != nb_levels; ++l) {
84  Range ents;
85  CHKERR bit_ref_manager->getEntitiesByDimAndRefLevel(
86  bit(l), BitRefLevel().set(), dim, ents);
87 
88  Range edges;
89  CHKERR moab.get_adjacencies(ents, 1, true, edges, moab::Interface::UNION);
90 
91  EntityHandle meshset_ref_edges;
92  CHKERR moab.create_meshset(MESHSET_SET, meshset_ref_edges);
93  CHKERR moab.add_entities(meshset_ref_edges, edges);
94 
95  MeshRefinement *refine = m_field.getInterface<MeshRefinement>();
96 
97  CHKERR refine->addVerticesInTheMiddleOfEdges(meshset_ref_edges,
98  bit(l + 1));
99  if (dim == 3) {
100  CHKERR refine->refineTets(ents, bit(l + 1));
101  } else if (dim == 2) {
102  CHKERR refine->refineTris(ents, bit(l + 1));
103  } else {
104  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
105  "Refinement implemented only for three and two dimensions");
106  }
107 
108  auto update_meshsets = [&]() {
110  // update cubit meshsets
111  for (_IT_CUBITMESHSETS_FOR_LOOP_(m_field, ciit)) {
112  EntityHandle cubit_meshset = ciit->meshset;
113  CHKERR bit_ref_manager->updateMeshsetByEntitiesChildren(
114  cubit_meshset, bit(l + 1), cubit_meshset, MBMAXTYPE, true);
115  }
117  };
118 
119  auto update_partition_sets = [&]() {
121 
122  ParallelComm *pcomm = ParallelComm::get_pcomm(
123  &m_field.get_moab(), m_field.get_basic_entity_data_ptr()->pcommID);
124  Tag part_tag = pcomm->part_tag();
125 
126  Range r_tagged_sets;
127  CHKERR m_field.get_moab().get_entities_by_type_and_tag(
128  0, MBENTITYSET, &part_tag, NULL, 1, r_tagged_sets,
129  moab::Interface::UNION);
130 
131  std::vector<EntityHandle> tagged_sets(r_tagged_sets.size());
132  std::copy(r_tagged_sets.begin(), r_tagged_sets.end(),
133  tagged_sets.begin());
134 
135  auto order_tagged_sets = [&]() {
137  std::vector<int> parts(tagged_sets.size());
138  CHKERR m_field.get_moab().tag_get_data(
139  part_tag, &*tagged_sets.begin(), tagged_sets.size(),
140  &*parts.begin());
141  map<int, EntityHandle> m_tagged_sets;
142  for (int p = 0; p != tagged_sets.size(); ++p) {
143  m_tagged_sets[parts[p]] = tagged_sets[p];
144  }
145  for (int p = 0; p != tagged_sets.size(); ++p) {
146  tagged_sets[p] = m_tagged_sets.at(p);
147  }
149  };
150 
151  auto add_children = [&]() {
153  std::vector<Range> part_ents(tagged_sets.size());
154 
155  for (int p = 0; p != tagged_sets.size(); ++p) {
156  Range ents;
157  CHKERR moab.get_entities_by_dimension(tagged_sets[p], dim, ents,
158  true);
159  CHKERR bit_ref_manager->filterEntitiesByRefLevel(
160  bit(l), BitRefLevel().set(), ents);
161 
162  Range children;
163  CHKERR bit_ref_manager->updateRangeByChildren(ents, children);
164  children = children.subset_by_dimension(dim);
165  CHKERR bit_ref_manager->filterEntitiesByRefLevel(
166  bit(l + 1), BitRefLevel().set(), children);
167 
168  Range adj;
169  for (auto d = 1; d != dim; ++d) {
170  CHKERR moab.get_adjacencies(children.subset_by_dimension(dim), d,
171  false, adj, moab::Interface::UNION);
172  }
173 
174  part_ents[p].merge(children);
175  part_ents[p].merge(adj);
176  }
177 
178  for (int p = 1; p != tagged_sets.size(); ++p) {
179  for (int pp = 0; pp != p; pp++) {
180  part_ents[p] = subtract(part_ents[p], part_ents[pp]);
181  }
182  }
183 
184  for (int p = 0; p != tagged_sets.size(); ++p) {
185  CHKERR moab.add_entities(tagged_sets[p], part_ents[p]);
186  CHKERR moab.tag_clear_data(part_tag, part_ents[p], &p);
187  }
188 
189  if (debug) {
190 
191  auto save_range = [&](const std::string name, const Range &r) {
193  auto meshset_ptr = get_temp_meshset_ptr(m_field.get_moab());
194  CHKERR m_field.get_moab().add_entities(*meshset_ptr, r);
195  CHKERR m_field.get_moab().write_file(name.c_str(), "VTK", "",
196  meshset_ptr->get_ptr(), 1);
198  };
199 
200  for (int p = 0; p != tagged_sets.size(); ++p) {
201  MOFEM_LOG("WORLD", Sev::inform)
202  << "Write part " << p << " level " << l;
203  Range ents;
204  CHKERR m_field.get_moab().get_entities_by_handle(tagged_sets[p],
205  ents, true);
206  CHKERR bit_ref_manager->filterEntitiesByRefLevel(
207  bit(l + 1), BitRefLevel().set(), ents);
208  CHKERR save_range("part" + boost::lexical_cast<std::string>(p) +
209  "_" + boost::lexical_cast<std::string>(l) +
210  ".vtk",
211  ents);
212  }
213  }
214 
216  };
217 
218  CHKERR order_tagged_sets();
219  CHKERR add_children();
220 
222  };
223 
224  CHKERR update_partition_sets();
225  CHKERR update_meshsets();
226 
227  CHKERR moab.delete_entities(&meshset_ref_edges, 1);
228  }
229 
230  if (debug) {
231  for (int l = 0; l <= nb_levels; ++l) {
232  MOFEM_LOG("WORLD", Sev::inform) << "Write level " << l;
233  CHKERR bit_ref_manager->writeBitLevel(
234  bit(l), BitRefLevel().set(),
235  ("level" + boost::lexical_cast<std::string>(l) + ".vtk").c_str(),
236  "VTK", "");
237  }
238  }
239 
240  if (shift == PETSC_TRUE) {
241  MOFEM_LOG("WORLD", Sev::inform) << "Shift bits";
242  CHKERR core.getInterface<BitRefManager>()->shiftRightBitRef(
243  nb_levels, BitRefLevel().set(), VERBOSE);
244  }
245 
246  CHKERR moab.write_file(mesh_out_file);
247  }
248  CATCH_ERRORS;
249 
251  CHKERRQ(ierr);
252 
253  return 0;
254 }
MoFEM::UnknownInterface::getInterface
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.
Definition: UnknownInterface.hpp:93
MYPCOMM_INDEX
#define MYPCOMM_INDEX
default communicator number PCOMM
Definition: definitions.h:228
MoFEM::MeshRefinement::refineTris
MoFEMErrorCode refineTris(const EntityHandle meshset, const BitRefLevel &bit, int verb=QUIET, const bool debug=false)
refine triangles in the meshset
Definition: MeshRefinement.cpp:921
MoFEM::CoreTmp< 0 >
Core (interface) class.
Definition: Core.hpp:82
EntityHandle
MoFEM::MeshRefinement::addVerticesInTheMiddleOfEdges
MoFEMErrorCode addVerticesInTheMiddleOfEdges(const EntityHandle meshset, const BitRefLevel &bit, const bool recursive=false, int verb=QUIET, EntityHandle start_v=0)
make vertices in the middle of edges in meshset and add them to refinement levels defined by bit
Definition: MeshRefinement.cpp:42
nb_levels
constexpr int nb_levels
Definition: level_set.cpp:58
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
MoFEM::BitRefManager::writeBitLevel
MoFEMErrorCode writeBitLevel(const BitRefLevel bit, const BitRefLevel mask, const char *file_name, const char *file_type, const char *options, const bool check_for_empty=true) const
Write bit ref level to file.
Definition: BitRefManager.cpp:634
sdf.r
int r
Definition: sdf.py:8
MoFEM::BitRefManager::filterEntitiesByRefLevel
MoFEMErrorCode filterEntitiesByRefLevel(const BitRefLevel bit, const BitRefLevel mask, Range &ents, int verb=QUIET) const
filter entities by bit ref level
Definition: BitRefManager.cpp:746
MoFEM::DeprecatedCoreInterface
Deprecated interface functions.
Definition: DeprecatedCoreInterface.hpp:16
MoFEM::Interface
DeprecatedCoreInterface Interface
Definition: Interface.hpp:2010
MoFEM::MeshRefinement::refineTets
MoFEMErrorCode refineTets(const EntityHandle meshset, const BitRefLevel &bit, int verb=QUIET, const bool debug=false)
refine TET in the meshset
Definition: MeshRefinement.cpp:197
VERBOSE
@ VERBOSE
Definition: definitions.h:222
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:548
MoFEM::CoreInterface::get_moab
virtual moab::Interface & get_moab()=0
MoFEM
implementation of Data Operators for Forces and Sources
Definition: Common.hpp:10
MoFEM::BitRefManager::setBitRefLevelByDim
MoFEMErrorCode setBitRefLevelByDim(const EntityHandle meshset, const int dim, const BitRefLevel bit, int verb=QUIET) const
Set the Bit Ref Level By Dim object.
Definition: BitRefManager.cpp:499
main
int main(int argc, char *argv[])
Definition: uniform_mesh_refinement.cpp:19
MoFEM::BitRefManager::getEntitiesByDimAndRefLevel
MoFEMErrorCode getEntitiesByDimAndRefLevel(const BitRefLevel bit, const BitRefLevel mask, const int dim, const EntityHandle meshset, int verb=0) const
add all ents from ref level given by bit to meshset
Definition: BitRefManager.cpp:822
MoFEM::MeshRefinement
Mesh refinement interface.
Definition: MeshRefinement.hpp:26
MoFEM::get_temp_meshset_ptr
auto get_temp_meshset_ptr(moab::Interface &moab)
Create smart pointer to temporary meshset.
Definition: Templates.hpp:1886
bit
auto bit
set bit
Definition: hanging_node_approx.cpp:75
MoFEM::CoreInterface::get_basic_entity_data_ptr
virtual boost::shared_ptr< BasicEntityData > & get_basic_entity_data_ptr()=0
Get pointer to basic entity data.
MoFEM::BitRefManager::updateMeshsetByEntitiesChildren
MoFEMErrorCode updateMeshsetByEntitiesChildren(const EntityHandle parent, const BitRefLevel &parent_bit, const BitRefLevel &parent_mask, const BitRefLevel &child_bit, const BitRefLevel &child_mask, const EntityHandle child, EntityType child_type, const bool recursive=false, int verb=0)
Get child entities form meshset containing parent entities.
Definition: BitRefManager.cpp:1029
help
static char help[]
Definition: uniform_mesh_refinement.cpp:10
debug
static const bool debug
Definition: dm_create_subdm.cpp:12
mesh_file_name
char mesh_file_name[255]
Definition: mesh_smoothing.cpp:22
MoFEM::BitRefManager::updateRangeByChildren
MoFEMErrorCode updateRangeByChildren(const Range &parent, Range &child, MoFEMTypes bh=MF_ZERO)
Update range by childrens.
Definition: BitRefManager.cpp:1144
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
save_range
auto save_range
Definition: thermo_elastic.cpp:144
MOFEM_LOG
#define MOFEM_LOG(channel, severity)
Log.
Definition: LogManager.hpp:308
CATCH_ERRORS
#define CATCH_ERRORS
Catch errors.
Definition: definitions.h:385
MoFEM::Core
CoreTmp< 0 > Core
Definition: Core.hpp:1148
MoFEM::Exceptions::ierr
static MoFEMErrorCodeGeneric< PetscErrorCode > ierr
Definition: Exceptions.hpp:76
MOFEM_DATA_INCONSISTENCY
@ MOFEM_DATA_INCONSISTENCY
Definition: definitions.h:31
sdf_hertz_2d_axisymm_plane.d
float d
Definition: sdf_hertz_2d_axisymm_plane.py:4
_IT_CUBITMESHSETS_FOR_LOOP_
#define _IT_CUBITMESHSETS_FOR_LOOP_(MESHSET_MANAGER, IT)
Iterator that loops over all the Cubit MeshSets in a moFEM field.
Definition: MeshsetsManager.hpp:34
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
MoFEMFunctionReturn
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:429
MOFEM_NOT_IMPLEMENTED
@ MOFEM_NOT_IMPLEMENTED
Definition: definitions.h:32
MoFEMFunctionBegin
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:359
l
FTensor::Index< 'l', 3 > l
Definition: matrix_function.cpp:21