deal.II version 9.7.0
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Loading...
Searching...
No Matches
triangulation.h
Go to the documentation of this file.
1// ------------------------------------------------------------------------
2//
3// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 2022 - 2025 by the deal.II authors
5//
6// This file is part of the deal.II library.
7//
8// Part of the source code is dual licensed under Apache-2.0 WITH
9// LLVM-exception OR LGPL-2.1-or-later. Detailed license information
10// governing the source code and code contributions can be found in
11// LICENSE.md and CONTRIBUTING.md at the top level directory of deal.II.
12//
13// ------------------------------------------------------------------------
14
15#ifndef dealii_cgal_triangulation_h
16#define dealii_cgal_triangulation_h
17
18#include <deal.II/base/config.h>
19
22
24
25#include <deal.II/grid/tria.h>
27
28#ifdef DEAL_II_WITH_CGAL
30
31# include <boost/hana.hpp>
32
33# include <CGAL/version.h>
34# if CGAL_VERSION_MAJOR >= 6
35# include <CGAL/Installation/internal/disable_deprecation_warnings_and_errors.h>
36# endif
37# include <CGAL/Complex_2_in_triangulation_3.h>
38# include <CGAL/IO/facets_in_complex_2_to_triangle_mesh.h>
39# include <CGAL/Implicit_surface_3.h>
40# include <CGAL/Labeled_mesh_domain_3.h>
41# include <CGAL/Mesh_complex_3_in_triangulation_3.h>
42# include <CGAL/Mesh_criteria_3.h>
43# include <CGAL/Mesh_triangulation_3.h>
44# include <CGAL/Polyhedron_3.h>
45# include <CGAL/Surface_mesh.h>
46# include <CGAL/Surface_mesh_default_triangulation_3.h>
47# include <CGAL/Triangulation_2.h>
48# include <CGAL/Triangulation_3.h>
49# include <CGAL/make_mesh_3.h>
50# include <CGAL/make_surface_mesh.h>
51
53
54namespace CGALWrappers
55{
86 template <int spacedim, typename CGALTriangulation>
87 void
88 add_points_to_cgal_triangulation(const std::vector<Point<spacedim>> &points,
89 CGALTriangulation &triangulation);
90
122 template <typename CGALTriangulation, int dim, int spacedim>
123 void
124 cgal_triangulation_to_dealii_triangulation(
125 const CGALTriangulation &cgal_triangulation,
126 Triangulation<dim, spacedim> &dealii_triangulation);
127
151 template <typename CGALTriangulationType,
152 typename CornerIndexType,
153 typename CurveIndexType>
154 void
155 cgal_triangulation_to_dealii_triangulation(
156 const CGAL::Mesh_complex_3_in_triangulation_3<CGALTriangulationType,
157 CornerIndexType,
158 CurveIndexType>
159 &cgal_triangulation,
160 Triangulation<3> &dealii_triangulation);
161
183 template <typename CGAL_MeshType>
184 void
185 cgal_surface_mesh_to_dealii_triangulation(const CGAL_MeshType &cgal_mesh,
186 Triangulation<2, 3> &triangulation);
187
188
189
190# ifndef DOXYGEN
191 // Template implementation
192 template <int spacedim, typename CGALTriangulation>
193 void
194 add_points_to_cgal_triangulation(const std::vector<Point<spacedim>> &points,
195 CGALTriangulation &triangulation)
196 {
197 Assert(triangulation.is_valid(),
199 "The triangulation you pass to this function should be a valid "
200 "CGAL triangulation."));
201
202 std::vector<typename CGALTriangulation::Point> cgal_points(points.size());
203 std::transform(points.begin(),
204 points.end(),
205 cgal_points.begin(),
206 [](const auto &p) {
207 return CGALWrappers::dealii_point_to_cgal_point<
208 typename CGALTriangulation::Point>(p);
209 });
210
211 triangulation.insert(cgal_points.begin(), cgal_points.end());
212 Assert(triangulation.is_valid(),
214 "The Triangulation is no longer valid after inserting the points. "
215 "Bailing out."));
216 }
217
218
219
220 template <typename CGALTriangulationType,
221 typename CornerIndexType,
222 typename CurveIndexType>
223 void
224 cgal_triangulation_to_dealii_triangulation(
225 const CGAL::Mesh_complex_3_in_triangulation_3<CGALTriangulationType,
226 CornerIndexType,
227 CurveIndexType>
228 &cgal_triangulation,
229 Triangulation<3> &dealii_triangulation)
230 {
231 using C3T3 = CGAL::Mesh_complex_3_in_triangulation_3<CGALTriangulationType,
232 CornerIndexType,
233 CurveIndexType>;
234
235 // Extract all vertices first
236 std::vector<Point<3>> dealii_vertices;
237 std::map<typename C3T3::Vertex_handle, unsigned int>
238 cgal_to_dealii_vertex_map;
239
240 std::size_t inum = 0;
241 for (auto vit = cgal_triangulation.triangulation().finite_vertices_begin();
242 vit != cgal_triangulation.triangulation().finite_vertices_end();
243 ++vit)
244 {
245 if (vit->in_dimension() <= -1)
246 continue;
247 cgal_to_dealii_vertex_map[vit] = inum++;
248 dealii_vertices.emplace_back(
249 CGALWrappers::cgal_point_to_dealii_point<3>(vit->point()));
250 }
251
252 // Now build cell connectivity
253 std::vector<CellData<3>> cells;
254 for (auto cgal_cell = cgal_triangulation.cells_in_complex_begin();
255 cgal_cell != cgal_triangulation.cells_in_complex_end();
256 ++cgal_cell)
257 {
258 CellData<3> cell(ReferenceCells::Tetrahedron.n_vertices());
259 for (unsigned int i = 0; i < 4; ++i)
260 cell.vertices[i] = cgal_to_dealii_vertex_map[cgal_cell->vertex(i)];
261 cell.manifold_id = cgal_triangulation.subdomain_index(cgal_cell);
262 cells.push_back(cell);
263 }
264
265 // Do the same with surface patches, if possible
266 SubCellData subcell_data;
267 if constexpr (std::is_integral_v<typename C3T3::Surface_patch_index>)
268 {
269 for (auto face = cgal_triangulation.facets_in_complex_begin();
270 face != cgal_triangulation.facets_in_complex_end();
271 ++face)
272 {
273 const auto &[cgal_cell, cgal_vertex_face_index] = *face;
274 CellData<2> dealii_face(ReferenceCells::Triangle.n_vertices());
275 // A face is identified by a cell and by the index within the cell
276 // of the opposite vertex. Loop over vertices, and retain only those
277 // that belong to this face.
278 int j = 0;
279 for (int i = 0; i < 4; ++i)
280 if (i != cgal_vertex_face_index)
281 dealii_face.vertices[j++] =
282 cgal_to_dealii_vertex_map[cgal_cell->vertex(i)];
283 dealii_face.manifold_id =
284 cgal_triangulation.surface_patch_index(cgal_cell,
285 cgal_vertex_face_index);
286 subcell_data.boundary_quads.emplace_back(dealii_face);
287 }
288 }
289 // and curves
290 if constexpr (std::is_integral_v<typename C3T3::Curve_index>)
291 {
292 for (auto edge = cgal_triangulation.edges_in_complex_begin();
293 edge != cgal_triangulation.edges_in_complex_end();
294 ++edge)
295 {
296 const auto &[cgal_cell, v1, v2] = *edge;
297 CellData<1> dealii_edge(ReferenceCells::Line.n_vertices());
298 dealii_edge.vertices[0] =
299 cgal_to_dealii_vertex_map[cgal_cell->vertex(v1)];
300 dealii_edge.vertices[1] =
301 cgal_to_dealii_vertex_map[cgal_cell->vertex(v2)];
302 dealii_edge.manifold_id =
303 cgal_triangulation.curve_index(cgal_cell->vertex(v1),
304 cgal_cell->vertex(v2));
305 subcell_data.boundary_lines.emplace_back(dealii_edge);
306 }
307 }
308
309 dealii_triangulation.create_triangulation(dealii_vertices,
310 cells,
311 subcell_data);
312 }
313
314
315
316 template <typename CGALTriangulation, int dim, int spacedim>
317 void
318 cgal_triangulation_to_dealii_triangulation(
319 const CGALTriangulation &cgal_triangulation,
320 Triangulation<dim, spacedim> &dealii_triangulation)
321 {
322 AssertThrow(cgal_triangulation.dimension() == dim,
323 ExcMessage("The dimension of the input CGAL triangulation (" +
324 std::to_string(cgal_triangulation.dimension()) +
325 ") does not match the dimension of the output "
326 "deal.II triangulation (" +
327 std::to_string(dim) + ")."));
328
329 Assert(dealii_triangulation.n_cells() == 0,
330 ExcMessage("The output triangulation object needs to be empty."));
331
332 // deal.II storage data structures
333 std::vector<Point<spacedim>> vertices(
334 cgal_triangulation.number_of_vertices());
335 std::vector<CellData<dim>> cells;
336 SubCellData subcell_data;
337
338 // CGAL storage data structures
339 std::map<typename CGALTriangulation::Vertex_handle, unsigned int>
340 vertex_map;
341 {
342 unsigned int i = 0;
343 for (const auto &v : cgal_triangulation.finite_vertex_handles())
344 {
345 vertices[i] =
346 CGALWrappers::cgal_point_to_dealii_point<spacedim>(v->point());
347 vertex_map[v] = i++;
348 }
349 }
350
351 const auto has_faces = boost::hana::is_valid(
352 [](auto &&obj) -> decltype(obj.finite_face_handles()) {});
353
354 const auto has_cells = boost::hana::is_valid(
355 [](auto &&obj) -> decltype(obj.finite_cell_handles()) {});
356
357 // Different loops for Triangulation_2 and Triangulation_3 types.
358 if constexpr (decltype(has_faces(cgal_triangulation)){})
359 {
360 // This is a non-degenerate Triangulation_2
361 if (cgal_triangulation.dimension() == 2)
362 for (const auto &f : cgal_triangulation.finite_face_handles())
363 {
364 CellData<dim> cell(ReferenceCells::Triangle.n_vertices());
365 for (unsigned int i = 0;
367 ++i)
368 cell.vertices[i] = vertex_map[f->vertex(i)];
369 cells.push_back(cell);
370 }
371 else if (cgal_triangulation.dimension() == 1)
372 // This is a degenerate Triangulation_2, made of edges
373 for (const auto &e : cgal_triangulation.finite_edges())
374 {
375 // An edge is identified by a face and a vertex index in the
376 // face
377 const auto &f = e.first;
378 const auto &i = e.second;
379 CellData<dim> cell(ReferenceCells::Line.n_vertices());
380 unsigned int id = 0;
381 // Since an edge is identified by a face (a triangle) and the
382 // index of the opposite vertex to the edge, we can use this
383 // logic to infer the indices of the vertices of the edge: loop
384 // over all vertices, and keep only those that are not the
385 // opposite vertex of the edge.
386 for (unsigned int j = 0;
388 ++j)
389 if (j != i)
390 cell.vertices[id++] = vertex_map[f->vertex(j)];
391 cells.push_back(cell);
392 }
393 else
394 {
396 }
397 }
398 else if constexpr (decltype(has_cells(cgal_triangulation)){})
399 {
400 // This is a non-degenerate Triangulation_3
401 if (cgal_triangulation.dimension() == 3)
402 for (const auto &c : cgal_triangulation.finite_cell_handles())
403 {
404 CellData<dim> cell(ReferenceCells::Tetrahedron.n_vertices());
405 for (unsigned int i = 0;
407 ++i)
408 cell.vertices[i] = vertex_map[c->vertex(i)];
409 cells.push_back(cell);
410 }
411 else if (cgal_triangulation.dimension() == 2)
412 // This is a degenerate Triangulation_3, made of triangles
413 for (const auto &facet : cgal_triangulation.finite_facets())
414 {
415 // A facet is identified by a cell and the opposite vertex index
416 // in the face
417 const auto &c = facet.first;
418 const auto &i = facet.second;
419 CellData<dim> cell(ReferenceCells::Triangle.n_vertices());
420 unsigned int id = 0;
421 // Since a face is identified by a cell (a tetrahedron) and the
422 // index of the opposite vertex to the face, we can use this
423 // logic to infer the indices of the vertices of the face: loop
424 // over all vertices, and keep only those that are not the
425 // opposite vertex of the face.
426 for (unsigned int j = 0;
428 ++j)
429 if (j != i)
430 cell.vertices[id++] = vertex_map[c->vertex(j)];
431 cells.push_back(cell);
432 }
433 else if (cgal_triangulation.dimension() == 1)
434 // This is a degenerate Triangulation_3, made of edges
435 for (const auto &edge : cgal_triangulation.finite_edges())
436 {
437 // An edge is identified by a cell and its two vertices
438 const auto &[c, i, j] = edge;
439 CellData<dim> cell(ReferenceCells::Line.n_vertices());
440 cell.vertices[0] = vertex_map[c->vertex(i)];
441 cell.vertices[1] = vertex_map[c->vertex(j)];
442 cells.push_back(cell);
443 }
444 else
445 {
447 }
448 }
449 dealii_triangulation.create_triangulation(vertices, cells, subcell_data);
450 }
451
452
453
454 template <typename CGAL_MeshType>
455 void
456 cgal_surface_mesh_to_dealii_triangulation(const CGAL_MeshType &cgal_mesh,
457 Triangulation<2, 3> &triangulation)
458 {
459 Assert(triangulation.n_cells() == 0,
461 "Triangulation must be empty upon calling this function."));
462
463 const auto is_surface_mesh =
464 boost::hana::is_valid([](auto &&obj) -> decltype(obj.faces()) {});
465
466 const auto is_polyhedral =
467 boost::hana::is_valid([](auto &&obj) -> decltype(obj.facets_begin()) {});
468
469 // Collect Vertices and cells
470 std::vector<::Point<3>> vertices;
471 std::vector<CellData<2>> cells;
472 SubCellData subcell_data;
473
474 // Different loops for Polyhedron or Surface_mesh types
475 if constexpr (decltype(is_surface_mesh(cgal_mesh)){})
476 {
477 AssertThrow(cgal_mesh.num_vertices() > 0,
478 ExcMessage("CGAL surface mesh is empty."));
479 vertices.reserve(cgal_mesh.num_vertices());
480 std::map<typename CGAL_MeshType::Vertex_index, unsigned int> vertex_map;
481 {
482 unsigned int i = 0;
483 for (const auto &v : cgal_mesh.vertices())
484 {
485 vertices.emplace_back(CGALWrappers::cgal_point_to_dealii_point<3>(
486 cgal_mesh.point(v)));
487 vertex_map[v] = i++;
488 }
489 }
490
491 // Collect CellData
492 for (const auto &face : cgal_mesh.faces())
493 {
494 const auto face_vertices =
495 CGAL::vertices_around_face(cgal_mesh.halfedge(face), cgal_mesh);
496
497 AssertThrow(face_vertices.size() == 3 || face_vertices.size() == 4,
498 ExcMessage("Only triangle or quadrilateral surface "
499 "meshes are supported in deal.II"));
500
501 CellData<2> c(face_vertices.size());
502 auto it_vertex = c.vertices.begin();
503 for (const auto &v : face_vertices)
504 *(it_vertex++) = vertex_map[v];
505
506 // Make sure the numberfing is consistent with the one in deal.II
507 if (face_vertices.size() == 4)
508 std::swap(c.vertices[3], c.vertices[2]);
509
510 cells.emplace_back(c);
511 }
512 }
513 else if constexpr (decltype(is_polyhedral(cgal_mesh)){})
514 {
515 AssertThrow(cgal_mesh.size_of_vertices() > 0,
516 ExcMessage("CGAL surface mesh is empty."));
517 vertices.reserve(cgal_mesh.size_of_vertices());
518 std::map<decltype(cgal_mesh.vertices_begin()), unsigned int> vertex_map;
519 {
520 unsigned int i = 0;
521 for (auto it = cgal_mesh.vertices_begin();
522 it != cgal_mesh.vertices_end();
523 ++it)
524 {
525 vertices.emplace_back(
526 CGALWrappers::cgal_point_to_dealii_point<3>(it->point()));
527 vertex_map[it] = i++;
528 }
529 }
530
531 // Loop over faces of Polyhedron, fill CellData
532 for (auto face = cgal_mesh.facets_begin();
533 face != cgal_mesh.facets_end();
534 ++face)
535 {
536 auto j = face->facet_begin();
537 const unsigned int vertices_per_face = CGAL::circulator_size(j);
539 ExcMessage("Only triangle or quadrilateral surface "
540 "meshes are supported in deal.II. You "
541 "tried to read a mesh where a face has " +
542 std::to_string(vertices_per_face) +
543 " vertices per face."));
544
545 CellData<2> c(vertices_per_face);
546 auto it = c.vertices.begin();
547 for (unsigned int i = 0; i < vertices_per_face; ++i)
548 {
549 *(it++) = vertex_map[j->vertex()];
550 ++j;
551 }
552
553 if (vertices_per_face == 4)
554 std::swap(c.vertices[3], c.vertices[2]);
555
556 cells.emplace_back(c);
557 }
558 }
559 else
560 {
561 AssertThrow(false,
563 "Unsupported CGAL surface triangulation type."));
564 }
565 triangulation.create_triangulation(vertices, cells, subcell_data);
566 }
567} // namespace CGALWrappers
568# endif // doxygen
569
571
572#else
573
574// Make sure the scripts that create the C++20 module input files have
575// something to latch on if the preprocessor #ifdef above would
576// otherwise lead to an empty content of the file.
579
580#endif
581#endif
unsigned int n_vertices() const
virtual void create_triangulation(const std::vector< Point< spacedim > > &vertices, const std::vector< CellData< dim > > &cells, const SubCellData &subcelldata)
unsigned int n_cells() const
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:40
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:41
#define DEAL_II_ASSERT_UNREACHABLE()
constexpr unsigned int GeometryInfo< dim >::vertices_per_face
#define Assert(cond, exc)
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcMessage(std::string arg1)
#define AssertThrow(cond, exc)
SymmetricTensor< 2, dim, Number > e(const Tensor< 2, dim, Number > &F)
constexpr ReferenceCell Triangle
constexpr ReferenceCell Tetrahedron
constexpr ReferenceCell Line
std::vector< CellData< 2 > > boundary_quads
Definition cell_data.h:248
std::vector< CellData< 1 > > boundary_lines
Definition cell_data.h:232