Loading [MathJax]/extensions/tex2jax.js
libMesh
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
mesh_tetgen_interface.C
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2025 Benjamin S. Kirk, John W. Peterson, Roy H. Stogner
3 
4 // This library is free software; you can redistribute it and/or
5 // modify it under the terms of the GNU Lesser General Public
6 // License as published by the Free Software Foundation; either
7 // version 2.1 of the License, or (at your option) any later version.
8 
9 // This library is distributed in the hope that it will be useful,
10 // but WITHOUT ANY WARRANTY; without even the implied warranty of
11 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
12 // Lesser General Public License for more details.
13 
14 // You should have received a copy of the GNU Lesser General Public
15 // License along with this library; if not, write to the Free Software
16 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
17 
18 #include "libmesh/libmesh_config.h"
19 #ifdef LIBMESH_HAVE_TETGEN
20 
21 
22 // C++ includes
23 #include <sstream>
24 
25 // Local includes
26 #include "libmesh/mesh_tetgen_interface.h"
27 
28 #include "libmesh/boundary_info.h"
29 #include "libmesh/cell_tet4.h"
30 #include "libmesh/face_tri3.h"
31 #include "libmesh/mesh_smoother_laplace.h"
32 #include "libmesh/unstructured_mesh.h"
33 #include "libmesh/utility.h" // binary_find
34 #include "libmesh/mesh_tetgen_wrapper.h"
35 
36 namespace libMesh
37 {
38 
39 //----------------------------------------------------------------------
40 // TetGenMeshInterface class members
43  _serializer(_mesh),
44  _switches("Q")
45 {
46 }
47 
48 void TetGenMeshInterface::set_switches(std::string switches)
49 {
50  // set the tetgen switch manually:
51  // p = tetrahedralizes a piecewise linear complex (see definition in user manual)
52  // Q = quiet, no terminal output
53  // q = specify a minimum radius/edge ratio
54  // a = tetrahedron volume constraint
55  // V = verbose output
56  // for full list of options and their meaning: see the tetgen manual
57  // (http://wias-berlin.de/software/tetgen/1.5/doc/manual/manual005.html)
58  _switches = std::move(switches);
59 }
60 
61 
63 {
64  this->triangulate_pointset();
65 }
66 
67 
69 {
70  // class tetgen_wrapper allows library access on a basic level:
71  TetGenWrapper tetgen_wrapper;
72 
73  // fill input structure with point set data:
74  this->fill_pointlist(tetgen_wrapper);
75 
76  // Run TetGen triangulation method:
77  // Q = quiet, no terminal output
78  // V = verbose, more terminal output
79  // Note: if no switch is used, the input must be a list of 3D points
80  // (.node file) and the Delaunay tetrahedralization of this point set
81  // will be generated.
82 
83  // Can we apply quality and volume constraints in
84  // triangulate_pointset()?. On at least one test problem,
85  // specifying any quality or volume constraints here causes tetgen
86  // to segfault down in the insphere method: a nullptr is passed
87  // to the routine.
88  std::ostringstream oss;
89  oss << _switches;
90  // oss << "V"; // verbose operation
91  //oss << "q" << std::fixed << 2.0; // quality constraint
92  //oss << "a" << std::fixed << 100.; // volume constraint
93 
94  // But if the user wants refinement, let's do our best.
95  if (_desired_volume)
96  oss << "a" << std::fixed << _desired_volume; // volume constraint
97 
98  tetgen_wrapper.set_switches(oss.str());
99 
100  // Run tetgen
101  tetgen_wrapper.run_tetgen();
102 
103  // save elements to mesh structure, nodes will not be changed:
104  const unsigned int num_elements = tetgen_wrapper.get_numberoftetrahedra();
105 
106  // Vector that temporarily holds the node labels defining element.
107  unsigned int node_labels[4];
108 
109  for (unsigned int i=0; i<num_elements; ++i)
110  {
111  auto elem = Elem::build(TET4);
112 
113  // Get the nodes associated with this element
114  for (auto j : elem->node_index_range())
115  node_labels[j] = tetgen_wrapper.get_element_node(i,j);
116 
117  // Associate the nodes with this element
118  this->assign_nodes_to_elem(node_labels, elem.get());
119 
120  // Finally, add this element to the mesh.
121  this->_mesh.add_elem(std::move(elem));
122  }
123 
124  // To the naked eye, a few smoothing iterations usually looks better.
125  // We don't do this by default.
126  if (this->_smooth_after_generating)
127  LaplaceMeshSmoother(this->_mesh).smooth(2);
128 }
129 
130 
131 
133 {
134  // class tetgen_wrapper allows library access on a basic level
135  TetGenWrapper tetgen_wrapper;
136 
137  // Copy Mesh's node points into TetGen data structure
138  this->fill_pointlist(tetgen_wrapper);
139 
140  // Run TetGen triangulation method:
141  // Q = quiet, no terminal output
142  // Note: if no switch is used, the input must be a list of 3D points
143  // (.node file) and the Delaunay tetrahedralization of this point set
144  // will be generated. In this particular function, we are throwing
145  // away the tetrahedra generated by TetGen, and keeping only the
146  // convex hull...
147  tetgen_wrapper.set_switches(_switches);
148  tetgen_wrapper.run_tetgen();
149  unsigned int num_elements = tetgen_wrapper.get_numberoftrifaces();
150 
151  // Delete *all* old elements. Yes, we legally delete elements while
152  // iterating over them because no entries from the underlying container
153  // are actually erased.
154  for (auto & elem : this->_mesh.element_ptr_range())
155  this->_mesh.delete_elem (elem);
156 
157  // We just removed any boundary info associated with element faces
158  // or edges, so let's update the boundary id caches.
160 
161  // Add the 2D elements which comprise the convex hull back to the mesh.
162  // Vector that temporarily holds the node labels defining element.
163  unsigned int node_labels[3];
164 
165  for (unsigned int i=0; i<num_elements; ++i)
166  {
167  auto elem = Elem::build(TRI3);
168 
169  // Get node labels associated with this element
170  for (auto j : elem->node_index_range())
171  node_labels[j] = tetgen_wrapper.get_triface_node(i,j);
172 
173  this->assign_nodes_to_elem(node_labels, elem.get());
174 
175  // Finally, add this element to the mesh.
176  this->_mesh.add_elem(std::move(elem));
177  }
178 
179  // To the naked eye, a few smoothing iterations usually looks better.
180  // We don't do this by default.
181  if (this->_smooth_after_generating)
182  LaplaceMeshSmoother(this->_mesh).smooth(2);
183 }
184 
185 
186 
187 
188 
190  double volume_constraint)
191 {
192  // start triangulation method with empty holes list:
193  std::vector<Point> noholes;
194  triangulate_conformingDelaunayMesh_carvehole(noholes, quality_constraint, volume_constraint);
195 }
196 
197 
198 
200  double quality_constraint,
201  double volume_constraint)
202 {
203  // Before calling this function, the Mesh must contain a convex hull
204  // of TRI3 elements which define the boundary.
205  unsigned hull_integrity_check = check_hull_integrity();
206 
207  // Possibly die if hull integrity check failed
208  this->process_hull_integrity_result(hull_integrity_check);
209 
210  // class tetgen_wrapper allows library access on a basic level
211  TetGenWrapper tetgen_wrapper;
212 
213  // Copy Mesh's node points into TetGen data structure
214  this->fill_pointlist(tetgen_wrapper);
215 
216  // >>> fill input structure "tetgenio" with facet data:
217  int facet_num = this->_mesh.n_elem();
218 
219  // allocate memory in "tetgenio" structure:
220  tetgen_wrapper.allocate_facetlist
221  (facet_num, cast_int<int>(holes.size()));
222 
223 
224  // Set up tetgen data structures with existing facet information
225  // from the convex hull.
226  {
227  int insertnum = 0;
228  for (auto & elem : this->_mesh.element_ptr_range())
229  {
230  tetgen_wrapper.allocate_facet_polygonlist(insertnum, 1);
231  tetgen_wrapper.allocate_polygon_vertexlist(insertnum, 0, 3);
232 
233  for (auto j : elem->node_index_range())
234  {
235  // We need to get the sequential index of elem->node_ptr(j), but
236  // it should already be stored in _sequential_to_libmesh_node_map...
237  unsigned libmesh_node_id = elem->node_id(j);
238 
239  // The libmesh node IDs may not be sequential, but can we assume
240  // they are at least in order??? We will do so here.
241  std::vector<unsigned>::iterator node_iter =
244  libmesh_node_id);
245 
246  // Check to see if not found: this could also indicate the sequential
247  // node map is not sorted...
248  libmesh_error_msg_if(node_iter == _sequential_to_libmesh_node_map.end(),
249  "Global node " << libmesh_node_id << " not found in sequential node map!");
250 
251  int sequential_index = cast_int<int>
253  node_iter));
254 
255  // Debugging:
256  // libMesh::out << "libmesh_node_id=" << libmesh_node_id
257  // << ", sequential_index=" << sequential_index
258  // << std::endl;
259 
260  tetgen_wrapper.set_vertex(insertnum, // facet number
261  0, // polygon (always 0)
262  j, // local vertex index in tetgen input
263  sequential_index);
264  }
265 
266  // Go to next facet in polygonlist
267  insertnum++;
268  }
269  }
270 
271 
272 
273  // fill hole list (if there are holes):
274  if (holes.size() > 0)
275  {
276  unsigned hole_index = 0;
277  for (Point ihole : holes)
278  tetgen_wrapper.set_hole(hole_index++,
279  REAL(ihole(0)),
280  REAL(ihole(1)),
281  REAL(ihole(2)));
282  }
283 
284 
285  // Run TetGen triangulation method
286 
287  // Assemble switches: we append the user's switches (if any) to
288  // - 'p' tetrahedralize a piecewise linear complex
289  // - 'C' check consistency of mesh (avoid inverted elements)
290  // (see definition and further options in user manual
291  // http://wias-berlin.de/software/tetgen/1.5/doc/manual/manual005.html )
292  std::ostringstream oss;
293  oss << "pC";
294  oss << _switches;
295 
296  if (quality_constraint != 0)
297  oss << "q" << std::fixed << quality_constraint;
298 
299  if (volume_constraint != 0)
300  oss << "a" << std::fixed << volume_constraint;
301 
302  std::string params = oss.str();
303 
304  tetgen_wrapper.set_switches(params); // TetGen switches: Piecewise linear complex, Quiet mode
305  tetgen_wrapper.run_tetgen();
306 
307  // => nodes:
308  unsigned int old_nodesnum = this->_mesh.n_nodes();
309  REAL x=0., y=0., z=0.;
310  const unsigned int num_nodes = tetgen_wrapper.get_numberofpoints();
311 
312  // Debugging:
313  // libMesh::out << "Original mesh had " << old_nodesnum << " nodes." << std::endl;
314  // libMesh::out << "Reserving space for " << num_nodes << " total nodes." << std::endl;
315 
316  // Reserve space for additional nodes in the node map
317  _sequential_to_libmesh_node_map.reserve(num_nodes);
318 
319  // Add additional nodes to the Mesh.
320  // Original code had i<=num_nodes here (Note: the indexing is:
321  // foo[3*i], [3*i+1], [3*i+2]) But according to the TetGen docs, "In
322  // all cases, the first item in any array is stored starting at
323  // index [0]."
324  for (unsigned int i=old_nodesnum; i<num_nodes; i++)
325  {
326  // Fill in x, y, z values
327  tetgen_wrapper.get_output_node(i, x,y,z);
328 
329  // Catch the node returned by add_point()... this will tell us the ID
330  // assigned by the Mesh.
331  Node * new_node = this->_mesh.add_point ( Point(x,y,z) );
332 
333  // Store this new ID in our sequential-to-libmesh node mapping array
334  _sequential_to_libmesh_node_map.push_back( new_node->id() );
335  }
336 
337  // Debugging:
338  // std::copy(_sequential_to_libmesh_node_map.begin(),
339  // _sequential_to_libmesh_node_map.end(),
340  // std::ostream_iterator<unsigned>(std::cout, " "));
341  // std::cout << std::endl;
342 
343 
344  // => tetrahedra:
345  const unsigned int num_elements = tetgen_wrapper.get_numberoftetrahedra();
346 
347  // Vector that temporarily holds the node labels defining element connectivity.
348  unsigned int node_labels[4];
349 
350  for (unsigned int i=0; i<num_elements; i++)
351  {
352  // TetGen only supports Tet4 elements.
353  auto elem = Elem::build(TET4);
354 
355  // Fill up the the node_labels vector
356  for (auto j : elem->node_index_range())
357  node_labels[j] = tetgen_wrapper.get_element_node(i,j);
358 
359  // Associate nodes with this element
360  this->assign_nodes_to_elem(node_labels, elem.get());
361 
362  // Finally, add this element to the mesh
363  this->_mesh.add_elem(std::move(elem));
364  }
365 
366  // Delete original convex hull elements. Is there ever a case where
367  // we should not do this?
368  this->delete_2D_hull_elements();
369 
370  // To the naked eye, a few smoothing iterations usually looks better.
371  // We don't do this by default.
372  if (this->_smooth_after_generating)
374 }
375 
376 
377 
378 
379 
381 {
382  // fill input structure with point set data:
383  wrapper.allocate_pointlist( this->_mesh.n_nodes() );
384 
385  // Make enough space to store a mapping between the implied sequential
386  // node numbering used in tetgen and libmesh's (possibly) non-sequential
387  // numbering scheme.
390 
391  {
392  unsigned index = 0;
393  for (auto & node : this->_mesh.node_ptr_range())
394  {
395  _sequential_to_libmesh_node_map[index] = node->id();
396  wrapper.set_node(index++,
397  REAL((*node)(0)),
398  REAL((*node)(1)),
399  REAL((*node)(2)));
400  }
401  }
402 }
403 
404 
405 
406 
407 
408 void TetGenMeshInterface::assign_nodes_to_elem(unsigned * node_labels, Elem * elem)
409 {
410  for (auto j : elem->node_index_range())
411  {
412  // Get the mapped node index to ask the Mesh for
413  unsigned mapped_node_id = _sequential_to_libmesh_node_map[ node_labels[j] ];
414 
415  Node * current_node = this->_mesh.node_ptr( mapped_node_id );
416 
417  elem->set_node(j) = current_node;
418  }
419 }
420 
421 } // namespace libMesh
422 
423 
424 #endif // #ifdef LIBMESH_HAVE_TETGEN
virtual Node *& set_node(const unsigned int i)
Definition: elem.h:2494
A Node is like a Point, but with more information.
Definition: node.h:52
virtual void smooth() override
Redefinition of the smooth function from the base class.
void set_switches(std::string new_switches)
Method to set switches to tetgen, allowing for different behaviours.
int get_element_node(unsigned i, unsigned j)
unsigned check_hull_integrity()
This function checks the integrity of the current set of elements in the Mesh to see if they comprise...
void get_output_node(unsigned i, REAL &x, REAL &y, REAL &z)
int get_triface_node(unsigned i, unsigned j)
std::string _switches
Parameter controlling the behaviour of tetgen.
void triangulate_conformingDelaunayMesh(double quality_constraint=0., double volume_constraint=0.)
Method invokes TetGen library to compute a Delaunay tetrahedralization from the nodes point set...
This is the base class from which all geometric element types are derived.
Definition: elem.h:94
std::vector< unsigned > _sequential_to_libmesh_node_map
We should not assume libmesh nodes are numbered sequentially...
MeshBase & mesh
void set_vertex(unsigned i, unsigned j, unsigned k, int nodeindex)
Sets index of ith facet, jth polygon, kth vertex in the TetGen input.
TetGenMeshInterface(UnstructuredMesh &mesh)
Constructor.
The libMesh namespace provides an interface to certain functionality in the library.
void allocate_polygon_vertexlist(unsigned i, unsigned j, int numofvertices)
Allocates memory, sets number of vertices for polygon j, facet i in the TetGen input.
void allocate_pointlist(int numofpoints)
Allocates memory, sets number of nodes in the TetGen input.
const BoundaryInfo & get_boundary_info() const
The information about boundary ids on the mesh.
Definition: mesh_base.h:160
Real distance(const Point &p)
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.
This class defines the data structures necessary for Laplace smoothing.
void triangulate_pointset()
Method invokes TetGen library to compute a Delaunay tetrahedralization from the nodes point set...
The TetGenWrapper provides an interface for basic access to TetGen data structures and methods...
virtual void delete_elem(Elem *e)=0
Removes element e from the mesh.
void set_node(unsigned i, REAL x, REAL y, REAL z)
Sets coordinates of point i in the TetGen input.
dof_id_type id() const
Definition: dof_object.h:828
The UnstructuredMesh class is derived from the MeshBase class.
void triangulate_conformingDelaunayMesh_carvehole(const std::vector< Point > &holes, double quality_constraint=0., double volume_constraint=0.)
Method invokes TetGen library to compute a Delaunay tetrahedralization from the nodes point set...
virtual Elem * add_elem(Elem *e)=0
Add elem e to the end of the element array.
static std::unique_ptr< Elem > build(const ElemType type, Elem *p=nullptr)
Definition: elem.C:321
void assign_nodes_to_elem(unsigned *node_labels, Elem *elem)
Assigns the node IDs contained in the &#39;node_labels&#39; array to &#39;elem&#39;.
void run_tetgen()
Starts the triangulation.
void regenerate_id_sets()
Clears and regenerates the cached sets of ids.
void allocate_facet_polygonlist(unsigned i, int numofpolygons)
Allocates memory, sets number of polygons for facet i in the TetGen input.
void set_switches(std::string_view s)
Method to set TetGen commandline switches -p Tetrahedralizes a piecewise linear complex (...
void set_hole(unsigned i, REAL x, REAL y, REAL z)
Sets coordinates of hole i in the TetGen input.
Real _desired_volume
The desired volume for the elements in the resulting mesh.
void process_hull_integrity_result(unsigned result)
This function prints an informative message and crashes based on the output of the check_hull_integri...
ForwardIterator binary_find(ForwardIterator first, ForwardIterator last, const T &value)
The STL provides std::binary_search() which returns true or false depending on whether the searched-f...
Definition: utility.h:265
IntRange< unsigned short > node_index_range() const
Definition: elem.h:2605
void allocate_facetlist(int numoffacets, int numofholes)
Allocates memory, sets number of facets, holes in the TetGen input.
void fill_pointlist(TetGenWrapper &wrapper)
This function copies nodes from the _mesh into TetGen&#39;s pointlist.
void pointset_convexhull()
Method invokes TetGen library to compute a Delaunay tetrahedralization from the nodes point set...
virtual dof_id_type n_elem() const =0
virtual const Node * node_ptr(const dof_id_type i) const =0
virtual void triangulate() override
Method invokes TetGen library to compute a Delaunay tetrahedralization.
UnstructuredMesh & _mesh
Local reference to the mesh we are working with.
A Point defines a location in LIBMESH_DIM dimensional Real space.
Definition: point.h:39
void delete_2D_hull_elements()
Delete original convex hull elements from the Mesh after performing a Delaunay tetrahedralization.
virtual dof_id_type n_nodes() const =0
bool _smooth_after_generating
Flag which tells whether we should smooth the mesh after it is generated.
Class MeshTetInterface provides an abstract interface for tetrahedralization of meshes by subclasses...