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
utilities.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_utilities_h
16#define dealii_cgal_utilities_h
17
18#include <deal.II/base/config.h>
19
21
23
24#ifdef DEAL_II_WITH_CGAL
26
28
29# include <deal.II/grid/tria.h>
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/Cartesian.h>
38# include <CGAL/Complex_2_in_triangulation_3.h>
39# include <CGAL/Exact_predicates_exact_constructions_kernel.h>
40# include <CGAL/Exact_predicates_inexact_constructions_kernel.h>
41# include <CGAL/Kernel_traits.h>
42# include <CGAL/Mesh_complex_3_in_triangulation_3.h>
43# include <CGAL/Mesh_criteria_3.h>
44# include <CGAL/Mesh_triangulation_3.h>
45// Disable a warning that we get with gcc-13 about a potential uninitialized
46// usage of an <anonymous> lambda function in this external CGAL header.
48# include <CGAL/Polygon_mesh_processing/corefinement.h>
50# include <CGAL/Polygon_mesh_processing/measure.h>
51# include <CGAL/Polygon_mesh_processing/remesh.h>
52# include <CGAL/Polygon_mesh_processing/triangulate_faces.h>
53# include <CGAL/Polyhedral_mesh_domain_with_features_3.h>
54# include <CGAL/Simple_cartesian.h>
55# include <CGAL/Surface_mesh.h>
56# include <CGAL/Triangulation_3.h>
57# include <CGAL/boost/graph/copy_face_graph.h>
58# include <CGAL/boost/graph/selection.h>
59# include <CGAL/convex_hull_3.h>
60# include <CGAL/make_mesh_3.h>
61# include <CGAL/make_surface_mesh.h>
62
63# include <fstream>
64# include <limits>
65# include <type_traits>
66
67
68
82namespace CGALWrappers
83{
84# ifdef CGAL_CONCURRENT_MESH_3
85 using ConcurrencyTag = CGAL::Parallel_tag;
86# else
87 using ConcurrencyTag = CGAL::Sequential_tag;
88# endif
89
100 {
109
117
124
131 };
132
133
134
145 template <typename C3t3>
146 void
148 CGAL::Surface_mesh<typename C3t3::Point::Point> &surface_mesh,
149 C3t3 &triangulation,
150 const AdditionalData<3> &data = AdditionalData<3>{});
151
195 template <typename CGALPointType>
196 void
198 const CGAL::Surface_mesh<CGALPointType> &surface_mesh_1,
199 const CGAL::Surface_mesh<CGALPointType> &surface_mesh_2,
200 const BooleanOperation &boolean_operation,
201 CGAL::Surface_mesh<CGALPointType> &output_surface_mesh);
202
214 template <typename CGALTriangulationType>
216 compute_quadrature(const CGALTriangulationType &tria,
217 const unsigned int degree);
218
234 template <int dim0, int dim1, int spacedim>
237 const typename ::Triangulation<dim0, spacedim>::cell_iterator &cell0,
238 const typename ::Triangulation<dim1, spacedim>::cell_iterator &cell1,
239 const unsigned int degree,
240 const BooleanOperation &bool_op,
241 const Mapping<dim0, spacedim> &mapping0 =
244 const Mapping<dim1, spacedim> &mapping1 =
247
263 template <int dim0, int dim1, int spacedim>
266 const typename ::Triangulation<dim0, spacedim>::cell_iterator &cell0,
267 const typename ::Triangulation<dim1, spacedim>::cell_iterator &cell1,
268 const unsigned int degree,
269 const Mapping<dim0, spacedim> &mapping0,
270 const Mapping<dim1, spacedim> &mapping1);
271
300 template <typename CGALPointType>
301 void
302 remesh_surface(CGAL::Surface_mesh<CGALPointType> &surface_mesh,
303 const AdditionalData<3> &data = AdditionalData<3>{});
304} // namespace CGALWrappers
305
306# ifndef DOXYGEN
307// Template implementations
308namespace CGALWrappers
309{
310 template <typename C3t3>
311 void
313 CGAL::Surface_mesh<typename C3t3::Point::Point> &surface_mesh,
314 C3t3 &triangulation,
315 const AdditionalData<3> &data)
316 {
317 using CGALPointType = typename C3t3::Point::Point;
318 Assert(CGAL::is_closed(surface_mesh),
319 ExcMessage("The surface mesh must be closed."));
320
321 using K = typename CGAL::Kernel_traits<CGALPointType>::Kernel;
322 using Mesh_domain = CGAL::Polyhedral_mesh_domain_with_features_3<
323 K,
324 CGAL::Surface_mesh<CGALPointType>>;
325 using Tr = typename CGAL::
326 Mesh_triangulation_3<Mesh_domain, CGAL::Default, ConcurrencyTag>::type;
327 using Mesh_criteria = CGAL::Mesh_criteria_3<Tr>;
328
329 CGAL::Polygon_mesh_processing::triangulate_faces(surface_mesh);
330 Mesh_domain domain(surface_mesh);
331 domain.detect_features();
332 Mesh_criteria criteria(CGAL::parameters::facet_size = data.facet_size,
333 CGAL::parameters::facet_angle = data.facet_angle,
334 CGAL::parameters::facet_distance =
335 data.facet_distance,
336 CGAL::parameters::cell_radius_edge_ratio =
337 data.cell_radius_edge_ratio,
338 CGAL::parameters::cell_size = data.cell_size);
339 // Mesh generation
340 triangulation = CGAL::make_mesh_3<C3t3>(domain,
341 criteria,
342 CGAL::parameters::no_perturb(),
343 CGAL::parameters::no_exude());
344 }
345
346
347
348 template <typename CGALPointType>
349 void
351 const CGAL::Surface_mesh<CGALPointType> &surface_mesh_1,
352 const CGAL::Surface_mesh<CGALPointType> &surface_mesh_2,
353 const BooleanOperation &boolean_operation,
354 CGAL::Surface_mesh<CGALPointType> &output_surface_mesh)
355 {
356 Assert(output_surface_mesh.is_empty(),
358 "output_surface_mesh must be empty upon calling this function"));
359 Assert(CGAL::is_closed(surface_mesh_1),
361 "The input surface_mesh_1 must be a closed surface mesh."));
362 Assert(CGAL::is_closed(surface_mesh_2),
364 "The input surface_mesh_2 must be a closed surface mesh."));
365 Assert(CGAL::is_triangle_mesh(surface_mesh_1),
366 ExcMessage("The first CGAL mesh must be triangulated."));
367 Assert(CGAL::is_triangle_mesh(surface_mesh_2),
368 ExcMessage("The second CGAL mesh must be triangulated."));
369
370 bool res = false;
371 auto surf_1 = surface_mesh_1;
372 auto surf_2 = surface_mesh_2;
373 namespace PMP = CGAL::Polygon_mesh_processing;
374 switch (boolean_operation)
375 {
377 res = PMP::corefine_and_compute_union(surf_1,
378 surf_2,
379 output_surface_mesh);
380 break;
382 res = PMP::corefine_and_compute_intersection(surf_1,
383 surf_2,
384 output_surface_mesh);
385 break;
387 res = PMP::corefine_and_compute_difference(surf_1,
388 surf_2,
389 output_surface_mesh);
390 break;
392 PMP::corefine(surf_1,
393 surf_2); // both surfaces are corefined
394 output_surface_mesh = std::move(surf_1);
395 res = true;
396 break;
397 default:
398 output_surface_mesh.clear();
399 break;
400 }
401 (void)res;
402 Assert(res,
403 ExcMessage("The boolean operation was not successfully computed."));
404 }
405
406
407
408 template <typename CGALTriangulationType>
409 ::Quadrature<CGALTriangulationType::Point::Ambient_dimension::value>
410 compute_quadrature(const CGALTriangulationType &tria,
411 const unsigned int degree)
412 {
413 Assert(tria.is_valid(), ExcMessage("The triangulation is not valid."));
414 Assert(CGALTriangulationType::Point::Ambient_dimension::value == 3,
416 Assert(degree > 0,
417 ExcMessage("The degree of the Quadrature formula is not positive."));
418
419 constexpr int spacedim =
420 CGALTriangulationType::Point::Ambient_dimension::value;
421 std::vector<std::array<::Point<spacedim>, spacedim + 1>>
422 vec_of_simplices; // tets
423
424 std::array<::Point<spacedim>, spacedim + 1> simplex;
425 for (const auto &f : tria.finite_cell_handles())
426 {
427 for (unsigned int i = 0; i < (spacedim + 1); ++i)
428 {
429 simplex[i] =
430 cgal_point_to_dealii_point<spacedim>(f->vertex(i)->point());
431 }
432
433 vec_of_simplices.push_back(simplex);
434 }
435
436 return QGaussSimplex<spacedim>(degree).mapped_quadrature(vec_of_simplices);
437 }
438
439
440
441 template <int dim0, int dim1, int spacedim>
442 ::Quadrature<spacedim>
444 const typename ::Triangulation<dim0, spacedim>::cell_iterator &cell0,
445 const typename ::Triangulation<dim1, spacedim>::cell_iterator &cell1,
446 const unsigned int degree,
447 const BooleanOperation &bool_op,
448 const Mapping<dim0, spacedim> &mapping0,
449 const Mapping<dim1, spacedim> &mapping1)
450 {
451 Assert(dim0 == 3 && dim1 == 3 && spacedim == 3,
452 ExcNotImplemented("2d geometries are not yet supported."));
453 if (dim1 > dim0)
454 {
456 cell1,
457 cell0,
458 degree,
459 bool_op,
460 mapping1,
461 mapping0); // This function works for dim1<=dim0, so swap them if this
462 // is not the case.
463 }
465 {
467 cell0, cell1, degree, mapping0, mapping1);
468 }
469 else
470 {
471 using K = CGAL::Exact_predicates_inexact_constructions_kernel;
472 using CGALPoint = CGAL::Point_3<K>;
473 using CGALTriangulation = CGAL::Triangulation_3<K>;
474 CGAL::Surface_mesh<CGALPoint> surface_1, surface_2, out_surface;
475 dealii_cell_to_cgal_surface_mesh(cell0, mapping0, surface_1);
476 dealii_cell_to_cgal_surface_mesh(cell1, mapping1, surface_2);
477 // They have to be triangle meshes
478 CGAL::Polygon_mesh_processing::triangulate_faces(surface_1);
479 CGAL::Polygon_mesh_processing::triangulate_faces(surface_2);
480 Assert(CGAL::is_triangle_mesh(surface_1) &&
481 CGAL::is_triangle_mesh(surface_2),
482 ExcMessage("The surface must be a triangle mesh."));
483 compute_boolean_operation(surface_1, surface_2, bool_op, out_surface);
484
485 CGALTriangulation tria;
486 tria.insert(out_surface.points().begin(), out_surface.points().end());
487 return compute_quadrature(tria, degree);
488 }
489 }
490
491
492
493 template <int dim0, int dim1, int spacedim>
494 ::Quadrature<spacedim>
496 const typename ::Triangulation<dim0, spacedim>::cell_iterator &cell0,
497 const typename ::Triangulation<dim1, spacedim>::cell_iterator &cell1,
498 const unsigned int degree,
499 const Mapping<dim0, spacedim> &mapping0,
500 const Mapping<dim1, spacedim> &mapping1)
501 {
502 Assert(dim0 == 3 && dim1 == 3 && spacedim == 3,
503 ExcNotImplemented("2d geometries are not yet supported."));
504 using K = CGAL::Exact_predicates_inexact_constructions_kernel;
505 using CGALPoint = CGAL::Point_3<K>;
506 using CGALTriangulation = CGAL::Triangulation_3<K>;
507
508 CGAL::Surface_mesh<CGALPoint> surface_1, surface_2, out_surface;
509 dealii_cell_to_cgal_surface_mesh(cell0, mapping0, surface_1);
510 dealii_cell_to_cgal_surface_mesh(cell1, mapping1, surface_2);
511 // They have to be triangle meshes
512 CGAL::Polygon_mesh_processing::triangulate_faces(surface_1);
513 CGAL::Polygon_mesh_processing::triangulate_faces(surface_2);
514
516 surface_2,
518 out_surface);
519 CGAL::Surface_mesh<CGALPoint> dummy;
520 CGALTriangulation tr;
521 CGAL::convex_hull_3(out_surface.points().begin(),
522 out_surface.points().end(),
523 dummy);
524 tr.insert(dummy.points().begin(), dummy.points().end());
525 return compute_quadrature(tr, degree);
526 }
527
528
529
530 template <typename CGALPointType>
531 void
532 remesh_surface(CGAL::Surface_mesh<CGALPointType> &cgal_mesh,
533 const AdditionalData<3> &data)
534 {
535 using K = CGAL::Exact_predicates_inexact_constructions_kernel;
536 using Mesh_domain = CGAL::Polyhedral_mesh_domain_with_features_3<K>;
537 // Polyhedron type
538 using Polyhedron = CGAL::Mesh_polyhedron_3<K>::type;
539 // Triangulation
540 using Tr = CGAL::Mesh_triangulation_3<Mesh_domain>::type;
541 using C3t3 =
542 CGAL::Mesh_complex_3_in_triangulation_3<Tr,
543 Mesh_domain::Corner_index,
544 Mesh_domain::Curve_index>;
545 // Criteria
546 using Mesh_criteria = CGAL::Mesh_criteria_3<Tr>;
547
548 Polyhedron poly;
549 CGAL::copy_face_graph(cgal_mesh, poly);
550 // Create a vector with only one element: the pointer to the
551 // polyhedron.
552 std::vector<Polyhedron *> poly_ptrs_vector(1, &poly);
553 // Create a polyhedral domain, with only one polyhedron,
554 // and no "bounding polyhedron", so the volumetric part of the
555 // domain will be empty.
556 Mesh_domain domain(poly_ptrs_vector.begin(), poly_ptrs_vector.end());
557 // Get sharp features
558 domain.detect_features(); // includes detection of borders
559
560 // Mesh criteria
561 const double default_value_edge_size = std::numeric_limits<double>::max();
562 if (data.edge_size > 0 &&
563 std::abs(data.edge_size - default_value_edge_size) > 1e-12)
564 {
565 Mesh_criteria criteria(CGAL::parameters::edge_size = data.edge_size,
566 CGAL::parameters::facet_size = data.facet_size,
567 CGAL::parameters::facet_angle = data.facet_angle,
568 CGAL::parameters::facet_distance =
569 data.facet_distance,
570 CGAL::parameters::cell_radius_edge_ratio =
571 data.cell_radius_edge_ratio,
572 CGAL::parameters::cell_size = data.cell_size);
573 // Mesh generation
574 C3t3 c3t3 = CGAL::make_mesh_3<C3t3>(domain,
575 criteria,
576 CGAL::parameters::no_perturb(),
577 CGAL::parameters::no_exude());
578 cgal_mesh.clear();
579 CGAL::facets_in_complex_3_to_triangle_mesh(c3t3, cgal_mesh);
580 }
581 else if (std::abs(data.edge_size - default_value_edge_size) <= 1e-12)
582 {
583 Mesh_criteria criteria(CGAL::parameters::facet_size = data.facet_size,
584 CGAL::parameters::facet_angle = data.facet_angle,
585 CGAL::parameters::facet_distance =
586 data.facet_distance,
587 CGAL::parameters::cell_radius_edge_ratio =
588 data.cell_radius_edge_ratio,
589 CGAL::parameters::cell_size = data.cell_size);
590 // Mesh generation
591 C3t3 c3t3 = CGAL::make_mesh_3<C3t3>(domain,
592 criteria,
593 CGAL::parameters::no_perturb(),
594 CGAL::parameters::no_exude());
595 cgal_mesh.clear();
596 CGAL::facets_in_complex_3_to_triangle_mesh(c3t3, cgal_mesh);
597 }
598 else
599 {
601 }
602 }
603
604
605
612 template <int spacedim>
613 void
614 resort_dealii_vertices_to_cgal_order(const unsigned int structdim,
615 std::vector<Point<spacedim>> &vertices)
616 {
617 // Mark the two arguments as "used" because some compilers complain about
618 // arguments used only within an 'if constexpr' block.
619 (void)structdim;
620 (void)vertices;
621
622 if constexpr (spacedim == 2)
623 if (ReferenceCell::n_vertices_to_type(structdim, vertices.size()) ==
625 std::swap(vertices[2], vertices[3]);
626 }
627
628
629
637 template <int dim, int spacedim>
638 std::vector<Point<spacedim>>
639 get_vertices_in_cgal_order(
640 const typename ::Triangulation<dim, spacedim>::cell_iterator &cell,
641 const Mapping<dim, spacedim> &mapping)
642 {
643 // Elements have to be rectangular or simplices
644 const unsigned int n_vertices = cell->n_vertices();
645 Assert((n_vertices == ReferenceCells::get_hypercube<dim>().n_vertices()) ||
646 (n_vertices == ReferenceCells::get_simplex<dim>().n_vertices()),
648
649 std::vector<Point<spacedim>> ordered_vertices(n_vertices);
650 std::copy_n(mapping.get_vertices(cell).begin(),
651 n_vertices,
652 ordered_vertices.begin());
653
654 resort_dealii_vertices_to_cgal_order(dim, ordered_vertices);
655
656 return ordered_vertices;
657 }
658
659
660
669 template <int dim, int spacedim>
670 std::vector<Point<spacedim>>
671 get_vertices_in_cgal_order(
672 const typename ::Triangulation<dim, spacedim>::cell_iterator &cell,
673 const unsigned int face_no,
674 const Mapping<dim, spacedim> &mapping)
675 {
676 // Elements have to be rectangular or simplices
677 const unsigned int n_vertices = cell->face(face_no)->n_vertices();
678 Assert(
679 (n_vertices == ReferenceCells::get_hypercube<dim - 1>().n_vertices()) ||
680 (n_vertices == ReferenceCells::get_simplex<dim - 1>().n_vertices()),
682
683 std::vector<Point<spacedim>> ordered_vertices(n_vertices);
684 std::copy_n(mapping.get_vertices(cell, face_no).begin(),
685 n_vertices,
686 ordered_vertices.begin());
687
688 resort_dealii_vertices_to_cgal_order(dim - 1, ordered_vertices);
689
690 return ordered_vertices;
691 }
692
693
694
695} // namespace CGALWrappers
696# endif
697
699
700#else
701
702// Make sure the scripts that create the C++20 module input files have
703// something to latch on if the preprocessor #ifdef above would
704// otherwise lead to an empty content of the file.
707
708#endif
709#endif
Abstract base class for mapping classes.
Definition mapping.h:320
virtual boost::container::small_vector< Point< spacedim >, ReferenceCells::max_n_vertices< dim >() > get_vertices(const typename Triangulation< dim, spacedim >::cell_iterator &cell) const
static ReferenceCell n_vertices_to_type(const int dim, const unsigned int n_vertices)
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:40
#define DEAL_II_DISABLE_EXTRA_DIAGNOSTICS
Definition config.h:603
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:41
#define DEAL_II_ENABLE_EXTRA_DIAGNOSTICS
Definition config.h:647
#define DEAL_II_ASSERT_UNREACHABLE()
static ::ExceptionBase & ExcNotImplemented()
#define Assert(cond, exc)
static ::ExceptionBase & ExcMessage(std::string arg1)
const Mapping< dim, spacedim > & get_default_linear_mapping(const Triangulation< dim, spacedim > &triangulation)
Definition mapping.cc:316
MappingQ< dim, spacedim > StaticMappingQ1< dim, spacedim >::mapping
Definition mapping_q1.h:104
void remesh_surface(CGAL::Surface_mesh< CGALPointType > &surface_mesh, const AdditionalData< 3 > &data=AdditionalData< 3 >{})
::Quadrature< CGALTriangulationType::Point::Ambient_dimension::value > compute_quadrature(const CGALTriangulationType &tria, const unsigned int degree)
std::vector< CGAL::Polygon_with_holes_2< KernelType > > compute_boolean_operation(const CGAL::Polygon_with_holes_2< KernelType > &polygon_1, const CGAL::Polygon_with_holes_2< KernelType > &polygon_2, const BooleanOperation &boolean_operation)
CGAL::Sequential_tag ConcurrencyTag
Definition utilities.h:87
CGAL::Exact_predicates_exact_constructions_kernel_with_sqrt K
::Quadrature< spacedim > compute_quadrature_on_boolean_operation(const typename ::Triangulation< dim0, spacedim >::cell_iterator &cell0, const typename ::Triangulation< dim1, spacedim >::cell_iterator &cell1, const unsigned int degree, const BooleanOperation &bool_op, const Mapping< dim0, spacedim > &mapping0=(ReferenceCells::get_hypercube< dim0 >() .template get_default_linear_mapping< dim0, spacedim >()), const Mapping< dim1, spacedim > &mapping1=(ReferenceCells::get_hypercube< dim1 >() .template get_default_linear_mapping< dim1, spacedim >()))
void cgal_surface_mesh_to_cgal_triangulation(CGAL::Surface_mesh< typename C3t3::Point::Point > &surface_mesh, C3t3 &triangulation, const AdditionalData< 3 > &data=AdditionalData< 3 >{})
void dealii_cell_to_cgal_surface_mesh(const typename ::Triangulation< dim, spacedim >::cell_iterator &cell, const ::Mapping< dim, spacedim > &mapping, CGAL::Surface_mesh< CGALPointType > &mesh)
::Quadrature< spacedim > compute_quadrature_on_intersection(const typename ::Triangulation< dim0, spacedim >::cell_iterator &cell0, const typename ::Triangulation< dim1, spacedim >::cell_iterator &cell1, const unsigned int degree, const Mapping< dim0, spacedim > &mapping0, const Mapping< dim1, spacedim > &mapping1)
void simplex(Triangulation< dim, dim > &tria, const std::vector< Point< dim > > &vertices)
constexpr const ReferenceCell & get_hypercube()
constexpr ReferenceCell Quadrilateral
constexpr const ReferenceCell & get_simplex()
::VectorizedArray< Number, width > abs(const ::VectorizedArray< Number, width > &)