libMesh
Public Member Functions | Protected Member Functions | Static Protected Member Functions | Protected Attributes | List of all members
libMesh::NetGenMeshInterface Class Reference

Class NetGenMeshInterface provides an interface for tetrahedralization of meshes using the NetGen library. More...

#include <mesh_netgen_interface.h>

Inheritance diagram for libMesh::NetGenMeshInterface:
[legend]

Public Member Functions

 NetGenMeshInterface (UnstructuredMesh &mesh)
 Constructor. More...
 
virtual ~NetGenMeshInterface () override=default
 Empty destructor. More...
 
virtual void triangulate () override
 Method invokes NetGen library to compute a tetrahedralization. More...
 
Realdesired_volume ()
 Sets and/or gets the desired tetrahedron volume. More...
 
bool & smooth_after_generating ()
 Sets/gets flag which tells whether to do two steps of Laplace mesh smoothing after generating the grid. More...
 
ElemTypeelem_type ()
 Sets and/or gets the desired element type. More...
 
void attach_hole_list (std::unique_ptr< std::vector< std::unique_ptr< UnstructuredMesh >>> holes)
 Attaches a vector of Mesh pointers defining holes which will be meshed around. More...
 

Protected Member Functions

unsigned check_hull_integrity ()
 This function checks the integrity of the current set of elements in the Mesh to see if they comprise a hull, that is: More...
 
void process_hull_integrity_result (unsigned result)
 This function prints an informative message and crashes based on the output of the check_hull_integrity() function. More...
 
void delete_2D_hull_elements ()
 Delete original convex hull elements from the Mesh after performing a Delaunay tetrahedralization. More...
 

Static Protected Member Functions

static BoundingBox volume_to_surface_mesh (UnstructuredMesh &mesh)
 Remove volume elements from the given mesh, after converting their outer boundary faces to surface elements. More...
 

Protected Attributes

MeshSerializer _serializer
 Tetgen only operates on serial meshes. More...
 
Real _desired_volume
 The desired volume for the elements in the resulting mesh. More...
 
bool _smooth_after_generating
 Flag which tells whether we should smooth the mesh after it is generated. More...
 
ElemType _elem_type
 The exact type of tetrahedra we intend to construct. More...
 
UnstructuredMesh_mesh
 Local reference to the mesh we are working with. More...
 
std::unique_ptr< std::vector< std::unique_ptr< UnstructuredMesh > > > _holes
 A pointer to a vector of meshes each defining a hole. More...
 

Detailed Description

Class NetGenMeshInterface provides an interface for tetrahedralization of meshes using the NetGen library.

For information about TetGen cf. NetGen github repository.

Author
Roy H. Stogner
Date
2024

Definition at line 57 of file mesh_netgen_interface.h.

Constructor & Destructor Documentation

◆ NetGenMeshInterface()

libMesh::NetGenMeshInterface::NetGenMeshInterface ( UnstructuredMesh mesh)
explicit

Constructor.

Takes a reference to the mesh containing the triangulated surface which is to be tetrahedralized.

Definition at line 69 of file mesh_netgen_interface.C.

69  :
72 {
73 }
MeshBase & mesh
MeshTetInterface(UnstructuredMesh &mesh)
Constructor.
MeshSerializer _serializer
Tetgen only operates on serial meshes.

◆ ~NetGenMeshInterface()

virtual libMesh::NetGenMeshInterface::~NetGenMeshInterface ( )
overridevirtualdefault

Empty destructor.

Member Function Documentation

◆ attach_hole_list()

void libMesh::MeshTetInterface::attach_hole_list ( std::unique_ptr< std::vector< std::unique_ptr< UnstructuredMesh >>>  holes)
inherited

Attaches a vector of Mesh pointers defining holes which will be meshed around.

We use unique_ptr here because we expect that we may need to modify these meshes internally.

Definition at line 125 of file mesh_tet_interface.C.

Referenced by MeshTetTest::testHole(), and MeshTetTest::testSphereShell().

126 {
127  _holes = std::move(holes);
128 }
std::unique_ptr< std::vector< std::unique_ptr< UnstructuredMesh > > > _holes
A pointer to a vector of meshes each defining a hole.

◆ check_hull_integrity()

unsigned libMesh::MeshTetInterface::check_hull_integrity ( )
protectedinherited

This function checks the integrity of the current set of elements in the Mesh to see if they comprise a hull, that is:

  • If they are all TRI3 elements
  • They all have non-nullptr neighbors
Returns
  • 0 if the mesh forms a valid hull
  • 1 if a non-TRI3 element is found
  • 2 if an element with a nullptr-neighbor is found

Definition at line 338 of file mesh_tet_interface.C.

References libMesh::MeshTetInterface::_mesh, libMesh::MeshBase::n_elem(), libMesh::Elem::neighbor_ptr_range(), libMesh::TRI3, and libMesh::Elem::type().

Referenced by triangulate(), and libMesh::TetGenMeshInterface::triangulate_conformingDelaunayMesh_carvehole().

339 {
340  // Check for easy return: if the Mesh is empty (i.e. if
341  // somebody called triangulate_conformingDelaunayMesh on
342  // a Mesh with no elements, then hull integrity check must
343  // fail...
344  if (_mesh.n_elem() == 0)
345  return 3;
346 
347  for (auto & elem : this->_mesh.element_ptr_range())
348  {
349  // Check for proper element type
350  if (elem->type() != TRI3)
351  {
352  //libmesh_error_msg("ERROR: Some of the elements in the original mesh were not TRI3!");
353  return 1;
354  }
355 
356  for (auto neigh : elem->neighbor_ptr_range())
357  {
358  if (neigh == nullptr)
359  {
360  // libmesh_error_msg("ERROR: Non-convex hull, cannot be tetrahedralized.");
361  return 2;
362  }
363  }
364  }
365 
366  // If we made it here, return success!
367  return 0;
368 }
virtual dof_id_type n_elem() const =0
UnstructuredMesh & _mesh
Local reference to the mesh we are working with.

◆ delete_2D_hull_elements()

void libMesh::MeshTetInterface::delete_2D_hull_elements ( )
protectedinherited

Delete original convex hull elements from the Mesh after performing a Delaunay tetrahedralization.

Definition at line 390 of file mesh_tet_interface.C.

References libMesh::MeshTetInterface::_mesh, libMesh::MeshBase::delete_elem(), libMesh::MeshBase::get_boundary_info(), libMesh::BoundaryInfo::regenerate_id_sets(), libMesh::TRI3, and libMesh::Elem::type().

Referenced by libMesh::TetGenMeshInterface::triangulate_conformingDelaunayMesh_carvehole().

391 {
392  for (auto & elem : this->_mesh.element_ptr_range())
393  {
394  // Check for proper element type. Yes, we legally delete elements while
395  // iterating over them because no entries from the underlying container
396  // are actually erased.
397  if (elem->type() == TRI3)
398  _mesh.delete_elem(elem);
399  }
400 
401  // We just removed any boundary info associated with hull element
402  // edges, so let's update the boundary id caches.
404 }
const BoundaryInfo & get_boundary_info() const
The information about boundary ids on the mesh.
Definition: mesh_base.h:165
virtual void delete_elem(Elem *e)=0
Removes element e from the mesh.
void regenerate_id_sets()
Clears and regenerates the cached sets of ids.
UnstructuredMesh & _mesh
Local reference to the mesh we are working with.

◆ desired_volume()

Real& libMesh::MeshTetInterface::desired_volume ( )
inlineinherited

Sets and/or gets the desired tetrahedron volume.

Set to zero to disable volume constraint.

Definition at line 68 of file mesh_tet_interface.h.

References libMesh::MeshTetInterface::_desired_volume.

Referenced by main(), and MeshTetTest::testSphereShell().

68 {return _desired_volume;}
Real _desired_volume
The desired volume for the elements in the resulting mesh.

◆ elem_type()

ElemType& libMesh::MeshTetInterface::elem_type ( )
inlineinherited

Sets and/or gets the desired element type.

This should be a Tet type.

Definition at line 81 of file mesh_tet_interface.h.

References libMesh::MeshTetInterface::_elem_type.

81 {return _elem_type;}
ElemType _elem_type
The exact type of tetrahedra we intend to construct.

◆ process_hull_integrity_result()

void libMesh::MeshTetInterface::process_hull_integrity_result ( unsigned  result)
protectedinherited

This function prints an informative message and crashes based on the output of the check_hull_integrity() function.

It is a separate function so that you can check hull integrity without crashing if you desire.

Definition at line 372 of file mesh_tet_interface.C.

References libMesh::err.

Referenced by libMesh::TetGenMeshInterface::triangulate_conformingDelaunayMesh_carvehole().

373 {
374  if (result != 0)
375  {
376  libMesh::err << "Error! Conforming Delaunay mesh tetrahedralization requires a convex hull." << std::endl;
377 
378  if (result==1)
379  {
380  libMesh::err << "Non-TRI3 elements were found in the input Mesh. ";
381  libMesh::err << "A constrained Delaunay triangulation requires a convex hull of TRI3 elements." << std::endl;
382  }
383 
384  libmesh_error_msg("Consider calling TetGenMeshInterface::pointset_convexhull() followed by Mesh::find_neighbors() first.");
385  }
386 }
OStreamProxy err

◆ smooth_after_generating()

bool& libMesh::MeshTetInterface::smooth_after_generating ( )
inlineinherited

Sets/gets flag which tells whether to do two steps of Laplace mesh smoothing after generating the grid.

False by default (for compatibility with old TetGenMeshInterface behavior).

Definition at line 75 of file mesh_tet_interface.h.

References libMesh::MeshTetInterface::_smooth_after_generating.

bool _smooth_after_generating
Flag which tells whether we should smooth the mesh after it is generated.

◆ triangulate()

void libMesh::NetGenMeshInterface::triangulate ( )
overridevirtual

Method invokes NetGen library to compute a tetrahedralization.

Implements libMesh::MeshTetInterface.

Definition at line 77 of file mesh_netgen_interface.C.

References libMesh::MeshTetInterface::_desired_volume, libMesh::MeshTetInterface::_holes, libMesh::MeshTetInterface::_mesh, libMesh::MeshBase::add_elem(), libMesh::MeshBase::add_point(), libMesh::MeshCommunication::broadcast(), libMesh::Elem::build_with_id(), libMesh::MeshTetInterface::check_hull_integrity(), libMesh::BoundingBox::contains(), libMesh::MeshBase::delete_elem(), libMesh::MeshBase::get_boundary_info(), libMesh::DofObject::id(), libMesh::libmesh_assert(), libMesh::libmesh_ignore(), libMesh::make_range(), libMesh::MeshTools::n_elem(), libMesh::MeshBase::n_nodes(), libMesh::MeshBase::node_ptr(), libMesh::Utility::pow(), libMesh::MeshBase::prepare_for_use(), libMesh::ParallelObject::processor_id(), libMesh::TET4, libMesh::TRI3, libMesh::TRI6, libMesh::TRI7, and libMesh::MeshTetInterface::volume_to_surface_mesh().

Referenced by main().

78 {
79  using namespace nglib;
80 
81  // We're hoping to do volume_to_surface_mesh in parallel at least,
82  // but then we'll need to serialize any hole meshes to rank 0 so it
83  // can use them in serial.
84 
85  const BoundingBox mesh_bb =
87 
88  std::vector<MeshSerializer> hole_serializers;
89  if (_holes)
90  for (std::unique_ptr<UnstructuredMesh> & hole : *_holes)
91  {
92  const BoundingBox hole_bb =
94 
95  libmesh_error_msg_if
96  (!mesh_bb.contains(hole_bb),
97  "Found hole with bounding box " << hole_bb <<
98  "\nextending outside of mesh bounding box " << mesh_bb);
99 
100  hole_serializers.emplace_back
101  (*hole, /* need_serial */ true,
102  /* serial_only_needed_on_proc_0 */ true);
103  }
104 
105  // If we're not rank 0, we're just going to wait for rank 0 to call
106  // Netgen, then receive its data afterward, we're not going to hope
107  // that Netgen does the exact same thing on every processor.
108  if (this->_mesh.processor_id() != 0)
109  {
110  // We don't need our holes anymore. Delete their serializers
111  // first to avoid dereferencing dangling pointers.
112  hole_serializers.clear();
113  if (_holes)
114  _holes->clear();
115 
116  // Receive the mesh data rank 0 will send later, then fix it up
117  // together
118  MeshCommunication().broadcast(this->_mesh);
119  this->_mesh.prepare_for_use();
120  return;
121  }
122 
123  this->check_hull_integrity();
124 
125  Ng_Meshing_Parameters params;
126 
127  // Override any default parameters we might need to, to avoid
128  // inserting nodes we don't want.
129  params.uselocalh = false;
130  params.minh = 0;
131  params.elementsperedge = 1;
132  params.elementspercurve = 1;
133  params.closeedgeenable = false;
134  params.closeedgefact = 0;
135  params.minedgelenenable = false;
136  params.minedgelen = 0;
137 
138  // Try to get a no-extra-nodes mesh if we're asked to, or try to
139  // translate our desired volume into NetGen terms otherwise.
140  //
141  // Spoiler alert: all we can do is try; NetGen uses a marching front
142  // algorithm that can insert extra nodes despite all my best
143  // efforts.
144  if (_desired_volume == 0) // shorthand for "no refinement"
145  {
146  params.maxh = std::numeric_limits<double>::max();
147  params.fineness = 0; // "coarse" in the docs
148  params.grading = 1; // "aggressive local grading" to avoid smoothing??
149 
150  // Turning off optimization steps avoids another opportunity for
151  // Netgen to try to add more nodes.
152  params.optsteps_3d = 0;
153  }
154  else
155  params.maxh = double(std::pow(_desired_volume, 1./3.));
156 
157  // Keep track of how NetGen copies of nodes map back to our original
158  // nodes, so we can connect new elements to nodes correctly.
159  std::unordered_map<int, dof_id_type> ng_to_libmesh_id;
160 
161  auto handle_ng_result = [](Ng_Result result) {
162  static const std::vector<std::string> result_types =
163  {"Netgen error", "Netgen success", "Netgen surface input error",
164  "Netgen volume failure", "Netgen STL input error",
165  "Netgen surface failure", "Netgen file not found"};
166 
167  if (result+1 >= 0 &&
168  std::size_t(result+1) < result_types.size())
169  libmesh_error_msg_if
170  (result, "Ng_GenerateVolumeMesh failed: " <<
171  result_types[result+1]);
172  else
173  libmesh_error_msg
174  ("Ng_GenerateVolumeMesh failed with an unknown error code");
175  };
176 
177  WrappedNgMesh ngmesh;
178 
179  // Create surface mesh in the WrappedNgMesh
180  {
181  // NetGen appears to use ONE-BASED numbering for its nodes, and
182  // since it doesn't return an id when adding nodes we'll have to
183  // track the numbering ourselves.
184  int ng_id = 1;
185 
186  auto create_surface_component =
187  [this, &ng_id, &ng_to_libmesh_id, &ngmesh]
188  (UnstructuredMesh & srcmesh, bool hole_mesh)
189  {
190  // Keep track of what nodes we've already added to the Netgen
191  // mesh vs what nodes we need to add. We'll keep track by id,
192  // not by point location. I don't know if Netgen can handle
193  // multiple nodes with the same point location, but if they can
194  // it's not going to be *us* who breaks that feature.
195  std::unordered_map<dof_id_type, int> libmesh_to_ng_id;
196 
197  // Keep track of what nodes we've already added to the main
198  // mesh from a hole mesh.
199  std::unordered_map<dof_id_type, dof_id_type> hole_to_main_mesh_id;
200 
201  // Use a separate array for passing points to NetGen, just in case
202  // we're not using double-precision ourselves.
203  std::array<double, 3> point_val;
204 
205  // And an array for element vertices
206  std::array<int, 3> elem_nodes;
207 
208  for (const auto * elem : srcmesh.element_ptr_range())
209  {
210  // If someone has triangles we can't triangulate, we have a
211  // problem
212  if (elem->type() == TRI6 ||
213  elem->type() == TRI7)
214  libmesh_not_implemented_msg
215  ("Netgen tetrahedralization currently only supports TRI3 boundaries");
216 
217  // If someone has non-triangles, let's just ignore them.
218  if (elem->type() != TRI3)
219  continue;
220 
221  for (int ni : make_range(3))
222  {
223  // Just using the "invert_trigs" option in NetGen params
224  // doesn't work for me, so we'll have to have properly
225  // oriented the tris earlier.
226  auto & elem_node = hole_mesh ? elem_nodes[2-ni] : elem_nodes[ni];
227 
228  const Node & n = elem->node_ref(ni);
229  auto n_id = n.id();
230  if (hole_mesh)
231  {
232  if (auto it = hole_to_main_mesh_id.find(n_id);
233  it != hole_to_main_mesh_id.end())
234  {
235  n_id = it->second;
236  }
237  else
238  {
239  Node * n_new = this->_mesh.add_point(n);
240  const dof_id_type n_new_id = n_new->id();
241  hole_to_main_mesh_id.emplace(n_id, n_new_id);
242  n_id = n_new_id;
243  }
244  }
245 
246  if (auto it = libmesh_to_ng_id.find(n_id);
247  it != libmesh_to_ng_id.end())
248  {
249  const int existing_ng_id = it->second;
250  elem_node = existing_ng_id;
251  }
252  else
253  {
254  for (auto i : make_range(3))
255  point_val[i] = double(n(i));
256 
257  Ng_AddPoint(ngmesh, point_val.data());
258 
259  ng_to_libmesh_id[ng_id] = n_id;
260  libmesh_to_ng_id[n_id] = ng_id;
261  elem_node = ng_id;
262  ++ng_id;
263  }
264  }
265 
266  Ng_AddSurfaceElement(ngmesh, NG_TRIG, elem_nodes.data());
267  }
268  };
269 
270  create_surface_component(this->_mesh, false);
271 
272  if (_holes)
273  for (const std::unique_ptr<UnstructuredMesh> & h : *_holes)
274  create_surface_component(*h, true);
275  }
276 
277  auto result = Ng_GenerateVolumeMesh(ngmesh, &params);
278  handle_ng_result(result);
279 
280  const int n_elem = Ng_GetNE(ngmesh);
281 
282  libmesh_error_msg_if (n_elem <= 0,
283  "NetGen failed to generate any tetrahedra");
284 
285  const dof_id_type n_points = Ng_GetNP(ngmesh);
286  const dof_id_type old_nodes = this->_mesh.n_nodes();
287 
288  // Netgen may have generated new interior nodes
289  if (n_points != old_nodes)
290  {
291  std::array<double, 3> point_val;
292 
293  // We should only be getting new nodes if we asked for them
294  if (!_desired_volume)
295  {
296  std::cout <<
297  "NetGen output " << n_points <<
298  " points when we gave it " <<
299  old_nodes << " and disabled refinement\n" <<
300  "If new interior points are acceptable in your mesh, please set\n" <<
301  "a non-zero desired_volume to indicate that. If new interior\n" <<
302  "points are not acceptable in your mesh, you may need a different\n" <<
303  "(non-advancing-front?) mesh generator." << std::endl;
304  libmesh_error();
305  }
306  else
307  for (auto i : make_range(old_nodes, n_points))
308  {
309  // i+1 since ng uses ONE-BASED numbering
310  Ng_GetPoint (ngmesh, i+1, point_val.data());
311  const Point p(point_val[0], point_val[1], point_val[2]);
312  Node * n_new = this->_mesh.add_point(p);
313  const dof_id_type n_new_id = n_new->id();
314  ng_to_libmesh_id[i+1] = n_new_id;
315  }
316  }
317 
318  for (auto * elem : this->_mesh.element_ptr_range())
319  this->_mesh.delete_elem(elem);
320 
321  BoundaryInfo * bi = & this->_mesh.get_boundary_info();
322 
323  for (auto i : make_range(n_elem))
324  {
325  // Enough data to return even a Tet10 without a segfault if nglib
326  // went nuts
327  int ngnodes[11];
328 
329  // i+1 since we must be 1-based with these ids too...
330  Ng_Volume_Element_Type ngtype =
331  Ng_GetVolumeElement(ngmesh, i+1, ngnodes);
332 
333  // But really nglib shouldn't go nuts
334  libmesh_assert(ngtype == NG_TET);
335  libmesh_ignore(ngtype);
336 
337  auto elem = this->_mesh.add_elem(Elem::build_with_id(TET4, i));
338  for (auto n : make_range(4))
339  {
340  const dof_id_type node_id =
341  libmesh_map_find(ng_to_libmesh_id, ngnodes[n]);
342  elem->set_node(n, this->_mesh.node_ptr(node_id));
343  }
344 
345  // NetGen and we disagree about node numbering orientation
346  elem->orient(bi);
347  }
348 
349  // We don't need our holes anymore. Delete their serializers
350  // first to avoid dereferencing dangling pointers.
351  hole_serializers.clear();
352  if (_holes)
353  _holes->clear();
354 
355  // Send our data to other ranks, then fix it up together
356  MeshCommunication().broadcast(this->_mesh);
357  this->_mesh.prepare_for_use();
358 }
std::unique_ptr< std::vector< std::unique_ptr< UnstructuredMesh > > > _holes
A pointer to a vector of meshes each defining a hole.
unsigned check_hull_integrity()
This function checks the integrity of the current set of elements in the Mesh to see if they comprise...
dof_id_type n_elem(const MeshBase::const_element_iterator &begin, const MeshBase::const_element_iterator &end)
Count up the number of elements of a specific type (as defined by an iterator range).
Definition: mesh_tools.C:969
void prepare_for_use(const bool skip_renumber_nodes_and_elements, const bool skip_find_neighbors)
Prepare a newly ecreated (or read) mesh for use.
Definition: mesh_base.C:759
const BoundaryInfo & get_boundary_info() const
The information about boundary ids on the mesh.
Definition: mesh_base.h:165
virtual Node * add_point(const Point &p, const dof_id_type id=DofObject::invalid_id, const processor_id_type proc_id=DofObject::invalid_processor_id)=0
Add a new Node at Point p to the end of the vertex array, with processor_id procid.
void libmesh_ignore(const Args &...)
T pow(const T &x)
Definition: utility.h:328
virtual void delete_elem(Elem *e)=0
Removes element e from the mesh.
virtual Elem * add_elem(Elem *e)=0
Add elem e to the end of the element array.
libmesh_assert(ctx)
static BoundingBox volume_to_surface_mesh(UnstructuredMesh &mesh)
Remove volume elements from the given mesh, after converting their outer boundary faces to surface el...
Real _desired_volume
The desired volume for the elements in the resulting mesh.
static std::unique_ptr< Elem > build_with_id(const ElemType type, dof_id_type id)
Calls the build() method above with a nullptr parent, and additionally sets the newly-created Elem&#39;s ...
Definition: elem.C:558
IntRange< T > make_range(T beg, T end)
The 2-parameter make_range() helper function returns an IntRange<T> when both input parameters are of...
Definition: int_range.h:140
virtual const Node * node_ptr(const dof_id_type i) const =0
processor_id_type processor_id() const
UnstructuredMesh & _mesh
Local reference to the mesh we are working with.
virtual dof_id_type n_nodes() const =0
uint8_t dof_id_type
Definition: id_types.h:67

◆ volume_to_surface_mesh()

BoundingBox libMesh::MeshTetInterface::volume_to_surface_mesh ( UnstructuredMesh mesh)
staticprotectedinherited

Remove volume elements from the given mesh, after converting their outer boundary faces to surface elements.

Returns the bounding box of the mesh; this is useful for detecting misplaced holes later.

Definition at line 131 of file mesh_tet_interface.C.

References libMesh::MeshBase::add_elem(), libMesh::MeshTools::Modification::all_tri(), libMesh::Elem::build_side_ptr(), libMesh::MeshBase::delete_elem(), libMesh::Elem::dim(), libMesh::Elem::flip(), libMesh::MeshBase::get_boundary_info(), libMesh::DofObject::id(), libMesh::MeshBase::is_prepared(), libMesh::MeshBase::is_replicated(), libMesh::Elem::level(), libMesh::libmesh_assert(), libMesh::Elem::low_order_key(), libMesh::make_range(), libMesh::MeshBase::max_elem_id(), mesh, libMesh::Elem::n_sides(), libMesh::Elem::neighbor_ptr(), libMesh::Elem::node_ref_range(), libMesh::MeshBase::parallel_max_unique_id(), libMesh::MeshBase::prepare_for_use(), libMesh::Real, libMesh::remote_elem, libMesh::DofObject::set_id(), libMesh::Elem::set_interior_parent(), libMesh::Elem::set_neighbor(), libMesh::DofObject::set_unique_id(), libMesh::Elem::side_index_range(), libMesh::Elem::side_ptr(), libMesh::BoundingBox::union_with(), and libMesh::MeshTools::valid_is_prepared().

Referenced by triangulate().

132 {
133  // If we've been handed an unprepared mesh then we need to be made
134  // aware of that and fix that; we're relying on neighbor pointers.
136 
137  if (!mesh.is_prepared())
139 
140  // We'll return a bounding box for use by subclasses in basic sanity checks.
141  BoundingBox surface_bb;
142 
143  // First convert all volume boundaries to surface elements; this
144  // gives us a manifold bounding the mesh, though it may not be a
145  // connected manifold even if the volume mesh was connected.
146  {
147  // Make sure ids are in sync and valid on a DistributedMesh
148  const dof_id_type max_orig_id = mesh.max_elem_id();
149 #ifdef LIBMESH_ENABLE_UNIQUE_ID
150  const unique_id_type max_unique_id = mesh.parallel_max_unique_id();
151 #endif
152 
153  // Change this if we add arbitrary polyhedra...
154  const dof_id_type max_sides = 6;
155 
156  std::unordered_set<Elem *> elems_to_delete;
157 
158  std::vector<std::unique_ptr<Elem>> elems_to_add;
159 
160  // Convert all faces to surface elements
161  for (auto * elem : mesh.active_element_ptr_range())
162  {
163  libmesh_error_msg_if (elem->dim() < 2,
164  "Cannot use meshes with 0D or 1D elements to define a volume");
165 
166  // If we've already got 2D elements then those are (part of)
167  // our surface.
168  if (elem->dim() == 2)
169  continue;
170 
171  // 3D elements will be removed after we've extracted their
172  // surface faces.
173  elems_to_delete.insert(elem);
174 
175  for (auto s : make_range(elem->n_sides()))
176  {
177  // If there's a neighbor on this side then there's not a
178  // boundary
179  if (elem->neighbor_ptr(s))
180  {
181  // We're not supporting AMR meshes here yet
182  if (elem->level() != elem->neighbor_ptr(s)->level())
183  libmesh_not_implemented_msg
184  ("Tetrahedralizaton of adapted meshes is not currently supported");
185  continue;
186  }
187 
188  elems_to_add.push_back(elem->build_side_ptr(s));
189  Elem * side_elem = elems_to_add.back().get();
190 
191  // Wipe the interior_parent before it can become a
192  // dangling pointer later
193  side_elem->set_interior_parent(nullptr);
194 
195  // If the mesh is replicated then its automatic id
196  // setting is fine. If not, then we need unambiguous ids
197  // independent of element traversal.
198  if (!mesh.is_replicated())
199  {
200  side_elem->set_id(max_orig_id + max_sides*elem->id() + s);
201 #ifdef LIBMESH_ENABLE_UNIQUE_ID
202  side_elem->set_unique_id(max_unique_id + max_sides*elem->id() + s);
203 #endif
204  }
205  }
206  }
207 
208  // If the mesh is replicated then its automatic neighbor finding
209  // is fine. If not, then we need to insert them ourselves, but
210  // it's easy because we can use the fact (from our implementation
211  // above) that our new elements have no parents or children, plus
212  // the fact (from the tiny fraction of homology I understand) that
213  // a manifold boundary is a manifold with no boundary.
214  //
215  // See UnstructuredMesh::find_neighbors() for more explanation of
216  // (a more complicated version of) the algorithm here.
217  if (!mesh.is_replicated())
218  {
219  typedef dof_id_type key_type;
220  typedef std::pair<Elem *, unsigned char> val_type;
221  typedef std::unordered_multimap<key_type, val_type> map_type;
222  map_type side_to_elem_map;
223 
224  std::unique_ptr<Elem> my_side, their_side;
225 
226  for (auto & elem : elems_to_add)
227  {
228  for (auto s : elem->side_index_range())
229  {
230  if (elem->neighbor_ptr(s))
231  continue;
232  const dof_id_type key = elem->low_order_key(s);
233  auto bounds = side_to_elem_map.equal_range(key);
234  if (bounds.first != bounds.second)
235  {
236  elem->side_ptr(my_side, s);
237  while (bounds.first != bounds.second)
238  {
239  Elem * potential_neighbor = bounds.first->second.first;
240  const unsigned int ns = bounds.first->second.second;
241  potential_neighbor->side_ptr(their_side, ns);
242  if (*my_side == *their_side)
243  {
244  elem->set_neighbor(s, potential_neighbor);
245  potential_neighbor->set_neighbor(ns, elem.get());
246  side_to_elem_map.erase (bounds.first);
247  break;
248  }
249  ++bounds.first;
250  }
251 
252  if (!elem->neighbor_ptr(s))
253  side_to_elem_map.emplace
254  (key, std::make_pair(elem.get(), cast_int<unsigned char>(s)));
255  }
256  }
257  }
258 
259  // At this point we *should* have a match for everything, so
260  // anything we don't have a match for is remote.
261  for (auto & elem : elems_to_add)
262  for (auto s : elem->side_index_range())
263  if (!elem->neighbor_ptr(s))
264  elem->set_neighbor(s, const_cast<RemoteElem*>(remote_elem));
265  }
266 
267  // Remove volume and edge elements
268  for (Elem * elem : elems_to_delete)
269  mesh.delete_elem(elem);
270 
271  // Add the new elements outside the loop so we don't risk
272  // invalidating iterators.
273  for (auto & elem : elems_to_add)
274  mesh.add_elem(std::move(elem));
275  }
276 
277  // Fix up neighbor pointers, element counts, etc.
279 
280  // We're making tets; we need to start with tris
282 
283  // Partition surface into connected components. At this point I'm
284  // finally going to give up and serialize, because at least we got
285  // from 3D down to 2D first, and because I don't want to have to
286  // turn flood_component into a while loop with a parallel sync in
287  // the middle, and because we do have to serialize *eventually*
288  // anyways unless we get a parallel tetrahedralizer backend someday.
289  MeshSerializer mesh_serializer(mesh);
290 
291  std::vector<std::unordered_set<Elem *>> components;
292  std::unordered_set<Elem *> in_component;
293 
294  for (auto * elem : mesh.element_ptr_range())
295  if (!in_component.count(elem))
296  components.emplace_back(flood_component(in_component, elem));
297 
298  const std::unordered_set<Elem *> * biggest_component = nullptr;
299  Real biggest_six_vol = 0;
300  for (const auto & component : components)
301  {
302  Real six_vol = six_times_signed_volume(component);
303  if (std::abs(six_vol) > std::abs(biggest_six_vol))
304  {
305  biggest_six_vol = six_vol;
306  biggest_component = &component;
307  }
308  }
309 
310  if (!biggest_component)
311  libmesh_error_msg("No non-zero-volume component found among " <<
312  components.size() << " boundary components");
313 
314  for (const auto & component : components)
315  if (&component != biggest_component)
316  {
317  for (Elem * elem: component)
318  mesh.delete_elem(elem);
319  }
320  else
321  {
322  for (Elem * elem: component)
323  {
324  if (biggest_six_vol < 0)
325  elem->flip(&mesh.get_boundary_info());
326 
327  for (auto & node : elem->node_ref_range())
328  surface_bb.union_with(node);
329  }
330  }
331 
333 
334  return surface_bb;
335 }
unique_id_type & set_unique_id()
Definition: dof_object.h:858
bool is_prepared() const
Definition: mesh_base.h:198
virtual unique_id_type parallel_max_unique_id() const =0
bool valid_is_prepared(const MeshBase &mesh)
A function for testing whether a mesh&#39;s cached is_prepared() setting is not a false positive...
Definition: mesh_tools.C:1182
void prepare_for_use(const bool skip_renumber_nodes_and_elements, const bool skip_find_neighbors)
Prepare a newly ecreated (or read) mesh for use.
Definition: mesh_base.C:759
This is the base class from which all geometric element types are derived.
Definition: elem.h:94
MeshBase & mesh
const BoundaryInfo & get_boundary_info() const
The information about boundary ids on the mesh.
Definition: mesh_base.h:165
void set_interior_parent(Elem *p)
Sets the pointer to the element&#39;s interior_parent.
Definition: elem.C:1248
dof_id_type & set_id()
Definition: dof_object.h:836
void all_tri(MeshBase &mesh)
Subdivides any non-simplex elements in a Mesh to produce simplex (triangular in 2D, tetrahedral in 3D) elements.
virtual void delete_elem(Elem *e)=0
Removes element e from the mesh.
virtual Elem * add_elem(Elem *e)=0
Add elem e to the end of the element array.
virtual dof_id_type max_elem_id() const =0
libmesh_assert(ctx)
void set_neighbor(const unsigned int i, Elem *n)
Assigns n as the neighbor.
Definition: elem.h:2618
Defines a Cartesian bounding box by the two corner extremum.
Definition: bounding_box.h:40
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
virtual std::unique_ptr< Elem > side_ptr(unsigned int i)=0
Temporarily serialize a DistributedMesh for non-distributed-mesh capable code paths.
virtual bool is_replicated() const
Definition: mesh_base.h:233
IntRange< T > make_range(T beg, T end)
The 2-parameter make_range() helper function returns an IntRange<T> when both input parameters are of...
Definition: int_range.h:140
void union_with(const Point &p)
Enlarges this bounding box to include the given point.
Definition: bounding_box.h:286
uint8_t unique_id_type
Definition: id_types.h:86
uint8_t dof_id_type
Definition: id_types.h:67
const RemoteElem * remote_elem
Definition: remote_elem.C:57

Member Data Documentation

◆ _desired_volume

Real libMesh::MeshTetInterface::_desired_volume
protectedinherited

The desired volume for the elements in the resulting mesh.

Unlimited (indicated by 0) by default

Definition at line 138 of file mesh_tet_interface.h.

Referenced by libMesh::MeshTetInterface::desired_volume(), triangulate(), and libMesh::TetGenMeshInterface::triangulate_pointset().

◆ _elem_type

ElemType libMesh::MeshTetInterface::_elem_type
protectedinherited

The exact type of tetrahedra we intend to construct.

Definition at line 149 of file mesh_tet_interface.h.

Referenced by libMesh::MeshTetInterface::elem_type().

◆ _holes

std::unique_ptr<std::vector<std::unique_ptr<UnstructuredMesh> > > libMesh::MeshTetInterface::_holes
protectedinherited

A pointer to a vector of meshes each defining a hole.

If this is nullptr, there are no holes!

Definition at line 160 of file mesh_tet_interface.h.

Referenced by triangulate().

◆ _mesh

UnstructuredMesh& libMesh::MeshTetInterface::_mesh
protectedinherited

◆ _serializer

MeshSerializer libMesh::NetGenMeshInterface::_serializer
protected

Tetgen only operates on serial meshes.

Definition at line 83 of file mesh_netgen_interface.h.

◆ _smooth_after_generating

bool libMesh::MeshTetInterface::_smooth_after_generating
protectedinherited

The documentation for this class was generated from the following files: