1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2020 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/unstructured_mesh.h"
32 #include "libmesh/utility.h" // binary_find
33 #include "libmesh/mesh_tetgen_wrapper.h"
34
35 namespace libMesh
36 {
37
38 //----------------------------------------------------------------------
39 // TetGenMeshInterface class members
TetGenMeshInterface(UnstructuredMesh & mesh)40 TetGenMeshInterface::TetGenMeshInterface (UnstructuredMesh & mesh) :
41 _mesh(mesh),
42 _serializer(_mesh),
43 _switches("Q")
44 {
45 }
46
set_switches(const std::string & switches)47 void TetGenMeshInterface::set_switches(const std::string & switches)
48 {
49 // set the tetgen switch manually:
50 // p = tetrahedralizes a piecewise linear complex (see definition in user manual)
51 // Q = quiet, no terminal output
52 // q = specify a minimum radius/edge ratio
53 // a = tetrahedron volume constraint
54 // V = verbose output
55 // for full list of options and their meaning: see the tetgen manual
56 // (http://wias-berlin.de/software/tetgen/1.5/doc/manual/manual005.html)
57 _switches = switches;
58 }
59
60
triangulate_pointset()61 void TetGenMeshInterface::triangulate_pointset ()
62 {
63 // class tetgen_wrapper allows library access on a basic level:
64 TetGenWrapper tetgen_wrapper;
65
66 // fill input structure with point set data:
67 this->fill_pointlist(tetgen_wrapper);
68
69 // Run TetGen triangulation method:
70 // Q = quiet, no terminal output
71 // V = verbose, more terminal output
72 // Note: if no switch is used, the input must be a list of 3D points
73 // (.node file) and the Delaunay tetrahedralization of this point set
74 // will be generated.
75
76 // Can we apply quality and volume constraints in
77 // triangulate_pointset()?. On at least one test problem,
78 // specifying any quality or volume constraints here causes tetgen
79 // to segfault down in the insphere method: a nullptr is passed
80 // to the routine.
81 std::ostringstream oss;
82 oss << _switches;
83 // oss << "V"; // verbose operation
84 //oss << "q" << std::fixed << 2.0; // quality constraint
85 //oss << "a" << std::fixed << 100.; // volume constraint
86 tetgen_wrapper.set_switches(oss.str());
87
88 // Run tetgen
89 tetgen_wrapper.run_tetgen();
90
91 // save elements to mesh structure, nodes will not be changed:
92 const unsigned int num_elements = tetgen_wrapper.get_numberoftetrahedra();
93
94 // Vector that temporarily holds the node labels defining element.
95 unsigned int node_labels[4];
96
97 for (unsigned int i=0; i<num_elements; ++i)
98 {
99 auto elem = Elem::build(TET4);
100
101 // Get the nodes associated with this element
102 for (auto j : elem->node_index_range())
103 node_labels[j] = tetgen_wrapper.get_element_node(i,j);
104
105 // Associate the nodes with this element
106 this->assign_nodes_to_elem(node_labels, elem.get());
107
108 // Finally, add this element to the mesh.
109 this->_mesh.add_elem(std::move(elem));
110 }
111 }
112
113
114
pointset_convexhull()115 void TetGenMeshInterface::pointset_convexhull ()
116 {
117 // class tetgen_wrapper allows library access on a basic level
118 TetGenWrapper tetgen_wrapper;
119
120 // Copy Mesh's node points into TetGen data structure
121 this->fill_pointlist(tetgen_wrapper);
122
123 // Run TetGen triangulation method:
124 // Q = quiet, no terminal output
125 // Note: if no switch is used, the input must be a list of 3D points
126 // (.node file) and the Delaunay tetrahedralization of this point set
127 // will be generated. In this particular function, we are throwing
128 // away the tetrahedra generated by TetGen, and keeping only the
129 // convex hull...
130 tetgen_wrapper.set_switches(_switches);
131 tetgen_wrapper.run_tetgen();
132 unsigned int num_elements = tetgen_wrapper.get_numberoftrifaces();
133
134 // Delete *all* old elements. Yes, we legally delete elements while
135 // iterating over them because no entries from the underlying container
136 // are actually erased.
137 for (auto & elem : this->_mesh.element_ptr_range())
138 this->_mesh.delete_elem (elem);
139
140 // We just removed any boundary info associated with element faces
141 // or edges, so let's update the boundary id caches.
142 this->_mesh.get_boundary_info().regenerate_id_sets();
143
144 // Add the 2D elements which comprise the convex hull back to the mesh.
145 // Vector that temporarily holds the node labels defining element.
146 unsigned int node_labels[3];
147
148 for (unsigned int i=0; i<num_elements; ++i)
149 {
150 auto elem = Elem::build(TRI3);
151
152 // Get node labels associated with this element
153 for (auto j : elem->node_index_range())
154 node_labels[j] = tetgen_wrapper.get_triface_node(i,j);
155
156 this->assign_nodes_to_elem(node_labels, elem.get());
157
158 // Finally, add this element to the mesh.
159 this->_mesh.add_elem(std::move(elem));
160 }
161 }
162
163
164
165
166
triangulate_conformingDelaunayMesh(double quality_constraint,double volume_constraint)167 void TetGenMeshInterface::triangulate_conformingDelaunayMesh (double quality_constraint,
168 double volume_constraint)
169 {
170 // start triangulation method with empty holes list:
171 std::vector<Point> noholes;
172 triangulate_conformingDelaunayMesh_carvehole(noholes, quality_constraint, volume_constraint);
173 }
174
175
176
triangulate_conformingDelaunayMesh_carvehole(const std::vector<Point> & holes,double quality_constraint,double volume_constraint)177 void TetGenMeshInterface::triangulate_conformingDelaunayMesh_carvehole (const std::vector<Point> & holes,
178 double quality_constraint,
179 double volume_constraint)
180 {
181 // Before calling this function, the Mesh must contain a convex hull
182 // of TRI3 elements which define the boundary.
183 unsigned hull_integrity_check = check_hull_integrity();
184
185 // Possibly die if hull integrity check failed
186 this->process_hull_integrity_result(hull_integrity_check);
187
188 // class tetgen_wrapper allows library access on a basic level
189 TetGenWrapper tetgen_wrapper;
190
191 // Copy Mesh's node points into TetGen data structure
192 this->fill_pointlist(tetgen_wrapper);
193
194 // >>> fill input structure "tetgenio" with facet data:
195 int facet_num = this->_mesh.n_elem();
196
197 // allocate memory in "tetgenio" structure:
198 tetgen_wrapper.allocate_facetlist
199 (facet_num, cast_int<int>(holes.size()));
200
201
202 // Set up tetgen data structures with existing facet information
203 // from the convex hull.
204 {
205 int insertnum = 0;
206 for (auto & elem : this->_mesh.element_ptr_range())
207 {
208 tetgen_wrapper.allocate_facet_polygonlist(insertnum, 1);
209 tetgen_wrapper.allocate_polygon_vertexlist(insertnum, 0, 3);
210
211 for (auto j : elem->node_index_range())
212 {
213 // We need to get the sequential index of elem->node_ptr(j), but
214 // it should already be stored in _sequential_to_libmesh_node_map...
215 unsigned libmesh_node_id = elem->node_id(j);
216
217 // The libmesh node IDs may not be sequential, but can we assume
218 // they are at least in order??? We will do so here.
219 std::vector<unsigned>::iterator node_iter =
220 Utility::binary_find(_sequential_to_libmesh_node_map.begin(),
221 _sequential_to_libmesh_node_map.end(),
222 libmesh_node_id);
223
224 // Check to see if not found: this could also indicate the sequential
225 // node map is not sorted...
226 libmesh_error_msg_if(node_iter == _sequential_to_libmesh_node_map.end(),
227 "Global node " << libmesh_node_id << " not found in sequential node map!");
228
229 int sequential_index = cast_int<int>
230 (std::distance(_sequential_to_libmesh_node_map.begin(),
231 node_iter));
232
233 // Debugging:
234 // libMesh::out << "libmesh_node_id=" << libmesh_node_id
235 // << ", sequential_index=" << sequential_index
236 // << std::endl;
237
238 tetgen_wrapper.set_vertex(insertnum, // facet number
239 0, // polygon (always 0)
240 j, // local vertex index in tetgen input
241 sequential_index);
242 }
243
244 // Go to next facet in polygonlist
245 insertnum++;
246 }
247 }
248
249
250
251 // fill hole list (if there are holes):
252 if (holes.size() > 0)
253 {
254 unsigned hole_index = 0;
255 for (Point ihole : holes)
256 tetgen_wrapper.set_hole(hole_index++,
257 REAL(ihole(0)),
258 REAL(ihole(1)),
259 REAL(ihole(2)));
260 }
261
262
263 // Run TetGen triangulation method
264
265 // Assemble switches: we append the user's switches (if any) to
266 // - 'p' tetrahedralize a piecewise linear complex
267 // - 'C' check consistency of mesh (avoid inverted elements)
268 // (see definition and further options in user manual
269 // http://wias-berlin.de/software/tetgen/1.5/doc/manual/manual005.html )
270 std::ostringstream oss;
271 oss << "pC";
272 oss << _switches;
273
274 if (quality_constraint != 0)
275 oss << "q" << std::fixed << quality_constraint;
276
277 if (volume_constraint != 0)
278 oss << "a" << std::fixed << volume_constraint;
279
280 std::string params = oss.str();
281
282 tetgen_wrapper.set_switches(params); // TetGen switches: Piecewise linear complex, Quiet mode
283 tetgen_wrapper.run_tetgen();
284
285 // => nodes:
286 unsigned int old_nodesnum = this->_mesh.n_nodes();
287 REAL x=0., y=0., z=0.;
288 const unsigned int num_nodes = tetgen_wrapper.get_numberofpoints();
289
290 // Debugging:
291 // libMesh::out << "Original mesh had " << old_nodesnum << " nodes." << std::endl;
292 // libMesh::out << "Reserving space for " << num_nodes << " total nodes." << std::endl;
293
294 // Reserve space for additional nodes in the node map
295 _sequential_to_libmesh_node_map.reserve(num_nodes);
296
297 // Add additional nodes to the Mesh.
298 // Original code had i<=num_nodes here (Note: the indexing is:
299 // foo[3*i], [3*i+1], [3*i+2]) But according to the TetGen docs, "In
300 // all cases, the first item in any array is stored starting at
301 // index [0]."
302 for (unsigned int i=old_nodesnum; i<num_nodes; i++)
303 {
304 // Fill in x, y, z values
305 tetgen_wrapper.get_output_node(i, x,y,z);
306
307 // Catch the node returned by add_point()... this will tell us the ID
308 // assigned by the Mesh.
309 Node * new_node = this->_mesh.add_point ( Point(x,y,z) );
310
311 // Store this new ID in our sequential-to-libmesh node mapping array
312 _sequential_to_libmesh_node_map.push_back( new_node->id() );
313 }
314
315 // Debugging:
316 // std::copy(_sequential_to_libmesh_node_map.begin(),
317 // _sequential_to_libmesh_node_map.end(),
318 // std::ostream_iterator<unsigned>(std::cout, " "));
319 // std::cout << std::endl;
320
321
322 // => tetrahedra:
323 const unsigned int num_elements = tetgen_wrapper.get_numberoftetrahedra();
324
325 // Vector that temporarily holds the node labels defining element connectivity.
326 unsigned int node_labels[4];
327
328 for (unsigned int i=0; i<num_elements; i++)
329 {
330 // TetGen only supports Tet4 elements.
331 auto elem = Elem::build(TET4);
332
333 // Fill up the the node_labels vector
334 for (auto j : elem->node_index_range())
335 node_labels[j] = tetgen_wrapper.get_element_node(i,j);
336
337 // Associate nodes with this element
338 this->assign_nodes_to_elem(node_labels, elem.get());
339
340 // Finally, add this element to the mesh
341 this->_mesh.add_elem(std::move(elem));
342 }
343
344 // Delete original convex hull elements. Is there ever a case where
345 // we should not do this?
346 this->delete_2D_hull_elements();
347 }
348
349
350
351
352
fill_pointlist(TetGenWrapper & wrapper)353 void TetGenMeshInterface::fill_pointlist(TetGenWrapper & wrapper)
354 {
355 // fill input structure with point set data:
356 wrapper.allocate_pointlist( this->_mesh.n_nodes() );
357
358 // Make enough space to store a mapping between the implied sequential
359 // node numbering used in tetgen and libmesh's (possibly) non-sequential
360 // numbering scheme.
361 _sequential_to_libmesh_node_map.clear();
362 _sequential_to_libmesh_node_map.resize( this->_mesh.n_nodes() );
363
364 {
365 unsigned index = 0;
366 for (auto & node : this->_mesh.node_ptr_range())
367 {
368 _sequential_to_libmesh_node_map[index] = node->id();
369 wrapper.set_node(index++,
370 REAL((*node)(0)),
371 REAL((*node)(1)),
372 REAL((*node)(2)));
373 }
374 }
375 }
376
377
378
379
380
assign_nodes_to_elem(unsigned * node_labels,Elem * elem)381 void TetGenMeshInterface::assign_nodes_to_elem(unsigned * node_labels, Elem * elem)
382 {
383 for (auto j : elem->node_index_range())
384 {
385 // Get the mapped node index to ask the Mesh for
386 unsigned mapped_node_id = _sequential_to_libmesh_node_map[ node_labels[j] ];
387
388 Node * current_node = this->_mesh.node_ptr( mapped_node_id );
389
390 elem->set_node(j) = current_node;
391 }
392 }
393
394
395
396
397
check_hull_integrity()398 unsigned TetGenMeshInterface::check_hull_integrity()
399 {
400 // Check for easy return: if the Mesh is empty (i.e. if
401 // somebody called triangulate_conformingDelaunayMesh on
402 // a Mesh with no elements, then hull integrity check must
403 // fail...
404 if (_mesh.n_elem() == 0)
405 return 3;
406
407 for (auto & elem : this->_mesh.element_ptr_range())
408 {
409 // Check for proper element type
410 if (elem->type() != TRI3)
411 {
412 //libmesh_error_msg("ERROR: Some of the elements in the original mesh were not TRI3!");
413 return 1;
414 }
415
416 for (auto neigh : elem->neighbor_ptr_range())
417 {
418 if (neigh == nullptr)
419 {
420 // libmesh_error_msg("ERROR: Non-convex hull, cannot be tetrahedralized.");
421 return 2;
422 }
423 }
424 }
425
426 // If we made it here, return success!
427 return 0;
428 }
429
430
431
432
433
process_hull_integrity_result(unsigned result)434 void TetGenMeshInterface::process_hull_integrity_result(unsigned result)
435 {
436 if (result != 0)
437 {
438 libMesh::err << "Error! Conforming Delaunay mesh tetrahedralization requires a convex hull." << std::endl;
439
440 if (result==1)
441 {
442 libMesh::err << "Non-TRI3 elements were found in the input Mesh. ";
443 libMesh::err << "A constrained Delaunay triangulation requires a convex hull of TRI3 elements." << std::endl;
444 }
445
446 libmesh_error_msg("Consider calling TetGenMeshInterface::pointset_convexhull() followed by Mesh::find_neighbors() first.");
447 }
448 }
449
450
451
452
delete_2D_hull_elements()453 void TetGenMeshInterface::delete_2D_hull_elements()
454 {
455 for (auto & elem : this->_mesh.element_ptr_range())
456 {
457 // Check for proper element type. Yes, we legally delete elements while
458 // iterating over them because no entries from the underlying container
459 // are actually erased.
460 if (elem->type() == TRI3)
461 _mesh.delete_elem(elem);
462 }
463
464 // We just removed any boundary info associated with hull element
465 // edges, so let's update the boundary id caches.
466 this->_mesh.get_boundary_info().regenerate_id_sets();
467 }
468
469
470
471 } // namespace libMesh
472
473
474 #endif // #ifdef LIBMESH_HAVE_TETGEN
475