]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Surface deal.II tria to volumetric deal.II tria 13767/head
authorMarco Feder <marco.feder@sissa.it>
Wed, 8 Jun 2022 08:12:25 +0000 (10:12 +0200)
committerMarco Feder <marco.feder@sissa.it>
Sat, 11 Jun 2022 09:45:45 +0000 (11:45 +0200)
include/deal.II/grid/grid_generator.h
source/grid/grid_generator_cgal.cc
tests/cgal/cgal_triangulation_07.cc [new file with mode: 0644]
tests/cgal/cgal_triangulation_07.output [new file with mode: 0644]

index 054e28f5b7a6442e493383d2290b2941c5288d9b..a0f20e4c7a368ef9ceef4a7bc14d1891af1b13bd 100644 (file)
@@ -1807,6 +1807,26 @@ namespace GridGenerator
                       CGALWrappers::AdditionalData<dim>{},
                     const Point<3> &interior_point    = Point<3>(),
                     const double &  outer_ball_radius = 1.0);
+
+  /**
+   * Create a deal.II Triangulation<3> out of a deal.II Triangulation<2,3>
+   * by filling it with tetrahedra.
+   *
+   * The last optional argument @p data can be used to pass additional
+   * arguments to the CGAL::Mesh_criteria_3 class (see
+   * https://doc.cgal.org/latest/Mesh_3/index.html for more information).
+   *
+   *
+   * @param [in] surface_tria The input deal.II Triangulation<2,3>.
+   * @param [out] vol_tria The output deal.II Triangulation<3>.
+   * @param[in] data Additional parameters to pass to the CGAL::make_mesh_3
+   * function.
+   */
+  void
+  surface_mesh_to_volumetric_mesh(const Triangulation<2, 3> &surface_tria,
+                                  Triangulation<3> &         vol_tria,
+                                  const CGALWrappers::AdditionalData<3> &data =
+                                    CGALWrappers::AdditionalData<3>{});
   ///@}
 
   /**
index ed9e7481aeb7d3758f8fafa87981743de77acb4b..a5a7fafa694010c73b440b01e46ded239a56f0f0 100644 (file)
@@ -143,6 +143,72 @@ namespace GridGenerator
     (void)outer_ball_radius;
     AssertThrow(false, ExcMessage("This function needs CGAL to be installed."));
 
+#  endif
+  }
+
+
+
+  void
+  surface_mesh_to_volumetric_mesh(const Triangulation<2, 3> &surface_tria,
+                                  Triangulation<3> &         vol_tria,
+                                  const CGALWrappers::AdditionalData<3> &data)
+  {
+#  ifdef DEAL_II_WITH_CGAL
+    Assert(
+      surface_tria.n_cells() > 0,
+      ExcMessage(
+        "The input triangulation cannot be empty when calling this function."));
+    Assert(
+      vol_tria.n_cells() == 0,
+      ExcMessage(
+        "The output triangulation must be empty when calling this function."));
+    using K       = CGAL::Exact_predicates_inexact_constructions_kernel;
+    using Point_3 = K::Point_3;
+
+    using Mesh_domain =
+      CGAL::Polyhedral_mesh_domain_with_features_3<K,
+                                                   CGAL::Surface_mesh<Point_3>>;
+    using Tr            = CGAL::Mesh_triangulation_3<Mesh_domain,
+                                          CGAL::Default,
+                                          CGALWrappers::ConcurrencyTag>::type;
+    using Mesh_criteria = CGAL::Mesh_criteria_3<Tr>;
+    using C3t3 =
+      CGAL::Mesh_complex_3_in_triangulation_3<Tr,
+                                              Mesh_domain::Corner_index,
+                                              Mesh_domain::Curve_index>;
+
+    CGAL::Surface_mesh<Point_3> mesh;
+    // This function "fills" the missing arrow of the following diagram.
+    //  Tria<2,3>                           Tria<3>
+    //      |                                 ^
+    //      |                                 |
+    //      |                                 |
+    //      |                                 |
+    //      V                                 |
+    // CGAL::Surface_mesh -----------> CGAL::C3t3
+    CGALWrappers::dealii_tria_to_cgal_surface_mesh(surface_tria, mesh);
+    CGAL::Polygon_mesh_processing::triangulate_faces(mesh);
+    CGAL::Polygon_mesh_processing::stitch_borders(mesh);
+    Mesh_domain domain(mesh);
+    domain.detect_features();
+    Mesh_criteria criteria(CGAL::parameters::facet_size  = data.facet_size,
+                           CGAL::parameters::facet_angle = data.facet_angle,
+                           CGAL::parameters::facet_distance =
+                             data.facet_distance,
+                           CGAL::parameters::cell_radius_edge_ratio =
+                             data.cell_radius_edge_ratio,
+                           CGAL::parameters::cell_size = data.cell_size);
+    const auto cgal_triangulation = CGAL::make_mesh_3<C3t3>(domain, criteria);
+    CGALWrappers::cgal_triangulation_to_dealii_triangulation(cgal_triangulation,
+                                                             vol_tria);
+
+#  else
+
+    (void)surface_tria;
+    (void)vol_tria;
+    (void)data;
+    AssertThrow(false, ExcMessage("This function needs CGAL to be installed."));
+
 #  endif
   }
 } // namespace GridGenerator
diff --git a/tests/cgal/cgal_triangulation_07.cc b/tests/cgal/cgal_triangulation_07.cc
new file mode 100644 (file)
index 0000000..76f45e6
--- /dev/null
@@ -0,0 +1,52 @@
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2022 by the deal.II authors
+//
+// This file is part of the deal.II library.
+//
+// The deal.II library is free software; you can use it, redistribute
+// it, and/or modify it under the terms of the GNU Lesser General
+// Public License as published by the Free Software Foundation; either
+// version 2.1 of the License, or (at your option) any later version.
+// The full text of the license can be found in the file LICENSE.md at
+// the top level directory of deal.II.
+//
+// ---------------------------------------------------------------------
+
+// Create a Tria<3> out of a Tria<2,3>.
+
+#include <deal.II/base/config.h>
+
+#include <deal.II/grid/grid_generator.h>
+#include <deal.II/grid/grid_out.h>
+#include <deal.II/grid/grid_tools.h>
+#include <deal.II/grid/tria.h>
+
+#include <deal.II/cgal/triangulation.h>
+
+#include "../tests.h"
+
+using namespace CGALWrappers;
+int
+main()
+{
+  initlog();
+  Triangulation<2, 3> surface_tria;
+  GridGenerator::hyper_sphere(surface_tria, {0., 1., 0.}, 1.);
+  surface_tria.refine_global(3);
+
+  Triangulation<3>  out_tria;
+  AdditionalData<3> data;
+  data.facet_size             = 0.1;
+  data.facet_distance         = 0.;
+  data.cell_radius_edge_ratio = 2.;
+  data.cell_size              = 0.1;
+  GridGenerator::surface_mesh_to_volumetric_mesh(surface_tria, out_tria, data);
+
+  GridOut       go;
+  std::ofstream out_tria_name("tria_surface_to_volumetric.vtk");
+  go.write_vtk(out_tria, out_tria_name);
+
+  remove("tria_surface_to_volumetric.vtk");
+  deallog << "OK" << std::endl;
+}
diff --git a/tests/cgal/cgal_triangulation_07.output b/tests/cgal/cgal_triangulation_07.output
new file mode 100644 (file)
index 0000000..0fd8fc1
--- /dev/null
@@ -0,0 +1,2 @@
+
+DEAL::OK

In the beginning the Universe was created. This has made a lot of people very angry and has been widely regarded as a bad move.

Douglas Adams


Typeset in Trocchi and Trocchi Bold Sans Serif.