]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Address some comments
authorMarco Feder <marco.feder@sissa.it>
Fri, 29 Apr 2022 22:58:03 +0000 (00:58 +0200)
committerLuca Heltai <luca.heltai@sissa.it>
Mon, 2 May 2022 10:30:32 +0000 (13:30 +0300)
doc/news/changes/minor/20220428MarcoFederLucaHeltai [deleted file]
include/deal.II/cgal/surface_mesh.h
source/cgal/surface_mesh.cc
source/cgal/surface_mesh.inst.in
tests/cgal/cgal_surface_mesh_01.cc

diff --git a/doc/news/changes/minor/20220428MarcoFederLucaHeltai b/doc/news/changes/minor/20220428MarcoFederLucaHeltai
deleted file mode 100644 (file)
index 1b6bdd0..0000000
+++ /dev/null
@@ -1,4 +0,0 @@
-New: Added function CGALWrappers::to_cgal_mesh() to convert a deal.II cell
-into a CGAL::Surface_mesh.
-<br>
-(Marco Feder, Luca Heltai, 2022/04/28)
index c4752ac2b607a407c1417d0ac8fcc7f1373ca0e0..bdafe40bc53b7ae64a754c6324aae1a165ee4077 100644 (file)
@@ -44,7 +44,7 @@ namespace CGALWrappers
    * More information on this class is available at
    * https://doc.cgal.org/latest/Surface_mesh/index.html
    *
-   * The function will throw an exception in dimension one. In dimension two it
+   * The function will throw an exception in dimension one. In dimension two, it
    * generates a surface mesh of the quadrilateral cell or of the triangle cell,
    * while in dimension three it will generate the surface mesh of the cell,
    * i.e., a polyhedral mesh containing the faces of the input cell.
@@ -62,7 +62,7 @@ namespace CGALWrappers
    */
   template <typename CGALPointType, int dim, int spacedim>
   void
-  to_cgal_mesh(
+  convert_to_cgal_surface_mesh(
     const typename dealii::Triangulation<dim, spacedim>::cell_iterator &cell,
     const dealii::Mapping<dim, spacedim> &                              mapping,
     CGAL::Surface_mesh<CGALPointType> &                                 mesh);
index 8a47384b37dd1b46bb0e9dbeb6f01ec6b71882b1..5ec2681d3c7a835d54b8329ff4f542bf4ca9d981 100644 (file)
@@ -35,7 +35,6 @@ namespace
     const unsigned                                nv = face->n_vertices();
     std::vector<typename CGAL_Mesh::Vertex_index> indices;
 
-
     switch (nv)
       {
         case 2:
@@ -80,11 +79,16 @@ namespace CGALWrappers
 {
   template <typename CGALPointType, int dim, int spacedim>
   void
-  to_cgal_mesh(const typename Triangulation<dim, spacedim>::cell_iterator &cell,
-               const Mapping<dim, spacedim> &     mapping,
-               CGAL::Surface_mesh<CGALPointType> &mesh)
+  convert_to_cgal_surface_mesh(
+    const typename Triangulation<dim, spacedim>::cell_iterator &cell,
+    const Mapping<dim, spacedim> &                              mapping,
+    CGAL::Surface_mesh<CGALPointType> &                         mesh)
   {
     Assert(dim > 1, ExcImpossibleInDim(dim));
+    Assert(
+      mesh.is_empty(),
+      ExcMessage(
+        "The CGAL::Surface_mesh object must be empty upon calling this function."));
     using Mesh           = CGAL::Surface_mesh<CGALPointType>;
     using Vertex_index   = typename Mesh::Vertex_index;
     const auto &vertices = mapping.get_vertices(cell);
@@ -108,6 +112,7 @@ namespace CGALWrappers
                   mesh,
                   (f % 2 == 0 || cell->n_vertices() != 8));
   }
+
   // explicit instantiations
 #    include "surface_mesh.inst"
 
index 9910460881a4b675f02ce161c456f118d543b3a8..3fe79dc6d4c95e889f2307ce378762d358a967ac 100644 (file)
@@ -18,7 +18,8 @@
 for (dim : DIMENSIONS; spacedim : SPACE_DIMENSIONS; cgal_kernel : CGAL_KERNELS)
   {
 #if dim <= spacedim
-    template void to_cgal_mesh<typename cgal_kernel::Point_3, dim, spacedim>(
+    template void
+    convert_to_cgal_surface_mesh<typename cgal_kernel::Point_3, dim, spacedim>(
       const typename Triangulation<dim, spacedim>::cell_iterator &cell,
       const Mapping<dim, spacedim> &                              mapping,
       CGAL::Surface_mesh<typename cgal_kernel::Point_3> &         mesh);
index 267945ff855e2bb1770e853350d181176cea7d93..5ce3152b66ed435e7a7496acc48aad2ea62a22ba 100644 (file)
@@ -13,7 +13,7 @@
 //
 // ---------------------------------------------------------------------
 
-// Convert a deal.II cell to a cgal Surface_mesh
+// Convert a deal.II cell to a cgal Surface_mesh.
 
 #include <deal.II/base/config.h>
 
@@ -38,21 +38,26 @@ void
 test()
 {
   deallog << "dim= " << dim << ",\t spacedim= " << spacedim << std::endl;
-  std::vector<std::vector<unsigned int>> d2t = {{}, {2}, {3, 4}, {4, 5, 6, 8}};
-  for (const auto nv : d2t[dim])
+  using namespace ReferenceCells;
+  std::vector<std::vector<ReferenceCell>> ref_cells = {
+    {},
+    {Line},
+    {Triangle, Quadrilateral},
+    {Tetrahedron, Pyramid, Wedge, Hexahedron}};
+  for (const auto r_cell : ref_cells[dim])
     {
       Triangulation<dim, spacedim>  tria;
       CGAL::Surface_mesh<CGALPoint> mesh;
 
-      const auto ref     = ReferenceCell::n_vertices_to_type(dim, nv);
-      const auto mapping = ref.template get_default_mapping<dim, spacedim>(1);
-      GridGenerator::reference_cell(tria, ref);
+      const auto mapping =
+        r_cell.template get_default_mapping<dim, spacedim>(1);
+      GridGenerator::reference_cell(tria, r_cell);
 
       const auto cell = tria.begin_active();
-      to_cgal_mesh(cell, *mapping, mesh);
+      convert_to_cgal_surface_mesh(cell, *mapping, mesh);
 
       Assert(mesh.is_valid(), dealii::ExcMessage("The CGAL mesh is not valid"));
-      deallog << "deal vertices: " << nv << ", cgal vertices "
+      deallog << "deal vertices: " << cell->n_vertices() << ", cgal vertices "
               << mesh.num_vertices() << std::endl;
       deallog << "deal faces: " << cell->n_faces() << ", cgal faces "
               << mesh.num_faces() << std::endl;

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.