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
19
20 // <h1>Miscellaneous Example 6 - Meshing with LibMesh's TetGen and Triangle Interfaces</h1>
21 // \author John W. Peterson
22 // \date 2011
23 //
24 // LibMesh provides interfaces to both Triangle and TetGen for generating
25 // Delaunay triangulations and tetrahedralizations in two and three dimensions
26 // (respectively).
27
28 // Local header files
29 #include "libmesh/elem.h"
30 #include "libmesh/face_tri3.h"
31 #include "libmesh/mesh.h"
32 #include "libmesh/mesh_generation.h"
33 #include "libmesh/mesh_tetgen_interface.h"
34 #include "libmesh/mesh_triangle_holes.h"
35 #include "libmesh/mesh_triangle_interface.h"
36 #include "libmesh/node.h"
37 #include "libmesh/replicated_mesh.h"
38 #include "libmesh/utility.h" // libmesh_map_find()
39
40 // Bring in everything from the libMesh namespace
41 using namespace libMesh;
42
43 // Major functions called by main
44 void triangulate_domain(const Parallel::Communicator & comm);
45 void tetrahedralize_domain(const Parallel::Communicator & comm);
46
47 // Helper routine for tetrahedralize_domain(). Adds the points and elements
48 // of a convex hull generated by TetGen to the input mesh
49 void add_cube_convex_hull_to_mesh(MeshBase & mesh, Point lower_limit, Point upper_limit);
50
51
52
53
54 // Begin the main program.
main(int argc,char ** argv)55 int main (int argc, char ** argv)
56 {
57 // Initialize libMesh and any dependent libraries, like in example 2.
58 LibMeshInit init (argc, argv);
59
60 libmesh_example_requires(2 <= LIBMESH_DIM, "2D support");
61
62 libMesh::out << "Triangulating an L-shaped domain with holes" << std::endl;
63
64 // 1.) 2D triangulation of L-shaped domain with three holes of different shape
65 triangulate_domain(init.comm());
66
67 libmesh_example_requires(3 <= LIBMESH_DIM, "3D support");
68
69 libMesh::out << "Tetrahedralizing a prismatic domain with a hole" << std::endl;
70
71 // 2.) 3D tetrahedralization of rectangular domain with hole.
72 tetrahedralize_domain(init.comm());
73
74 return 0;
75 }
76
77
78
79
triangulate_domain(const Parallel::Communicator & comm)80 void triangulate_domain(const Parallel::Communicator & comm)
81 {
82 #ifdef LIBMESH_HAVE_TRIANGLE
83 // Use typedefs for slightly less typing.
84 typedef TriangleInterface::Hole Hole;
85 typedef TriangleInterface::PolygonHole PolygonHole;
86 typedef TriangleInterface::ArbitraryHole ArbitraryHole;
87
88 // Libmesh mesh that will eventually be created.
89 Mesh mesh(comm, 2);
90
91 // The points which make up the L-shape:
92 mesh.add_point(Point(0., 0.));
93 mesh.add_point(Point(0., -1.));
94 mesh.add_point(Point(-1., -1.));
95 mesh.add_point(Point(-1., 1.));
96 mesh.add_point(Point(1., 1.));
97 mesh.add_point(Point(1., 0.));
98
99 // Declare the TriangleInterface object. This is where
100 // we can set parameters of the triangulation and where the
101 // actual triangulate function lives.
102 TriangleInterface t(mesh);
103
104 // Customize the variables for the triangulation
105 t.desired_area() = .01;
106
107 // A Planar Straight Line Graph (PSLG) is essentially a list
108 // of segments which have to exist in the final triangulation.
109 // For an L-shaped domain, Triangle will compute the convex
110 // hull of boundary points if we do not specify the PSLG.
111 // The PSLG algorithm is also required for triangulating domains
112 // containing holes
113 t.triangulation_type() = TriangleInterface::PSLG;
114
115 // Turn on/off Laplacian mesh smoothing after generation.
116 // By default this is on.
117 t.smooth_after_generating() = true;
118
119 // Define holes...
120
121 // hole_1 is a circle (discretized by 50 points)
122 PolygonHole hole_1(Point(-0.5, 0.5), // center
123 0.25, // radius
124 50); // n. points
125
126 // hole_2 is itself a triangle
127 PolygonHole hole_2(Point(0.5, 0.5), // center
128 0.1, // radius
129 3); // n. points
130
131 // hole_3 is an ellipse of 100 points which we define here
132 Point ellipse_center(-0.5, -0.5);
133 const unsigned int n_ellipse_points=100;
134 std::vector<Point> ellipse_points(n_ellipse_points);
135 const Real
136 dtheta = 2*libMesh::pi / static_cast<Real>(n_ellipse_points),
137 a = .1,
138 b = .2;
139
140 for (unsigned int i=0; i<n_ellipse_points; ++i)
141 ellipse_points[i]= Point(ellipse_center(0)+a*cos(i*dtheta),
142 ellipse_center(1)+b*sin(i*dtheta));
143
144 ArbitraryHole hole_3(ellipse_center, ellipse_points);
145
146 // Create the vector of Hole*'s ...
147 std::vector<Hole *> holes;
148 holes.push_back(&hole_1);
149 holes.push_back(&hole_2);
150 holes.push_back(&hole_3);
151
152 // ... and attach it to the triangulator object
153 t.attach_hole_list(&holes);
154
155 // Triangulate!
156 t.triangulate();
157
158 // Write the result to file
159 mesh.write("delaunay_l_shaped_hole.e");
160
161 #else
162 // Avoid compiler warnings
163 libmesh_ignore(comm);
164 #endif // LIBMESH_HAVE_TRIANGLE
165 }
166
167
168
tetrahedralize_domain(const Parallel::Communicator & comm)169 void tetrahedralize_domain(const Parallel::Communicator & comm)
170 {
171 #ifdef LIBMESH_HAVE_TETGEN
172 // The algorithm is broken up into several steps:
173 // 1.) A convex hull is constructed for a rectangular hole.
174 // 2.) A convex hull is constructed for the domain exterior.
175 // 3.) Neighbor information is updated so TetGen knows there is a convex hull
176 // 4.) A vector of hole points is created.
177 // 5.) The domain is tetrahedralized, the mesh is written out, etc.
178
179 // The mesh we will eventually generate
180 ReplicatedMesh mesh(comm, 3);
181
182 // Lower and Upper bounding box limits for a rectangular hole within the unit cube.
183 Point hole_lower_limit(0.2, 0.2, 0.4);
184 Point hole_upper_limit(0.8, 0.8, 0.6);
185
186 // 1.) Construct a convex hull for the hole
187 add_cube_convex_hull_to_mesh(mesh, hole_lower_limit, hole_upper_limit);
188
189 // 2.) Generate elements comprising the outer boundary of the domain.
190 add_cube_convex_hull_to_mesh(mesh,
191 Point(0., 0., 0.),
192 Point(1., 1., 1.));
193
194 // 3.) Update neighbor information so that TetGen can verify there is a convex hull.
195 mesh.find_neighbors();
196
197 // 4.) Set up vector of hole points
198 std::vector<Point> hole(1);
199 hole[0] = Point(0.5*(hole_lower_limit + hole_upper_limit));
200
201 // 5.) Set parameters and tetrahedralize the domain
202
203 // 0 means "use TetGen default value"
204 double quality_constraint = 2.0;
205
206 // The volume constraint determines the max-allowed tetrahedral
207 // volume in the Mesh. TetGen will split cells which are larger than
208 // this size
209 double volume_constraint = 0.001;
210
211 // Construct the Delaunay tetrahedralization
212 TetGenMeshInterface t(mesh);
213
214 // In debug mode, set the tetgen switches
215 // -V (verbose) and
216 // -CC (check consistency and constrained Delaunay)
217 // In optimized mode, only switch Q (quiet) is set.
218 // For more options, see tetgen website:
219 // (http://wias-berlin.de/software/tetgen/1.5/doc/manual/manual005.html#cmd-m)
220 #ifdef DEBUG
221 t.set_switches("VCC");
222 #endif
223 t.triangulate_conformingDelaunayMesh_carvehole(hole,
224 quality_constraint,
225 volume_constraint);
226
227 // Find neighbors, etc in preparation for writing out the Mesh
228 mesh.prepare_for_use();
229
230 // Finally, write out the result
231 mesh.write("hole_3D.e");
232 #else
233 // Avoid compiler warnings
234 libmesh_ignore(comm);
235 #endif // LIBMESH_HAVE_TETGEN
236 }
237
238
239
240
241
242
243
244
245
246
247
248
249
add_cube_convex_hull_to_mesh(MeshBase & mesh,Point lower_limit,Point upper_limit)250 void add_cube_convex_hull_to_mesh(MeshBase & mesh,
251 Point lower_limit,
252 Point upper_limit)
253 {
254 #ifdef LIBMESH_HAVE_TETGEN
255 ReplicatedMesh cube_mesh(mesh.comm(), 3);
256
257 unsigned n_elem = 1;
258
259 MeshTools::Generation::build_cube(cube_mesh,
260 n_elem, n_elem, n_elem, // n. elements in each direction
261 lower_limit(0), upper_limit(0),
262 lower_limit(1), upper_limit(1),
263 lower_limit(2), upper_limit(2),
264 HEX8);
265
266 // The pointset_convexhull() algorithm will ignore the Hex8s
267 // in the Mesh, and just construct the triangulation
268 // of the convex hull.
269 TetGenMeshInterface t(cube_mesh);
270 t.pointset_convexhull();
271
272 // Now add all nodes from the boundary of the cube_mesh to the input mesh.
273
274 // Map from "node id in cube_mesh" -> "node id in mesh". Initially inserted
275 // with a dummy value, later to be assigned a value by the input mesh.
276 std::map<unsigned, unsigned> node_id_map;
277 typedef std::map<unsigned, unsigned>::iterator iterator;
278
279 for (auto & elem : cube_mesh.element_ptr_range())
280 for (auto s : elem->side_index_range())
281 if (elem->neighbor_ptr(s) == nullptr)
282 {
283 // Add the node IDs of this side to the set
284 std::unique_ptr<Elem> side = elem->side_ptr(s);
285
286 for (auto n : side->node_index_range())
287 node_id_map.emplace(side->node_id(n), /*dummy_value=*/0);
288 }
289
290 // For each node in the map, insert it into the input mesh and keep
291 // track of the ID assigned.
292 for (iterator it=node_id_map.begin(); it != node_id_map.end(); ++it)
293 {
294 // Id of the node in the cube mesh
295 unsigned id = (*it).first;
296
297 // Pointer to node in the cube mesh
298 Node & old_node = cube_mesh.node_ref(id);
299
300 // Add geometric point to input mesh
301 Node * new_node = mesh.add_point (old_node);
302
303 // Track ID value of new_node in map
304 (*it).second = new_node->id();
305 }
306
307 // With the points added and the map data structure in place, we are
308 // ready to add each TRI3 element of the cube_mesh to the input Mesh
309 // with proper node assignments
310 for (auto & old_elem : cube_mesh.element_ptr_range())
311 if (old_elem->type() == TRI3)
312 {
313 Elem * new_elem = mesh.add_elem(new Tri3);
314
315 // Assign nodes in new elements. Since this is an example,
316 // we'll do it in several steps.
317 for (auto i : old_elem->node_index_range())
318 {
319 // Mapping to node ID in input mesh
320 unsigned int new_node_id = libmesh_map_find(node_id_map, old_elem->node_id(i));
321
322 // Node pointer assigned from input mesh
323 new_elem->set_node(i) = mesh.node_ptr(new_node_id);
324 }
325 }
326 #else
327 // Avoid compiler warnings
328 libmesh_ignore(mesh, lower_limit, upper_limit);
329 #endif // LIBMESH_HAVE_TETGEN
330 }
331