]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Add header file for point conversion 13929/head
authorMarco Feder <marco.feder@sissa.it>
Wed, 8 Jun 2022 12:45:21 +0000 (14:45 +0200)
committerMarco Feder <marco.feder@sissa.it>
Wed, 8 Jun 2022 22:17:32 +0000 (00:17 +0200)
include/deal.II/cgal/point_conversion.h [new file with mode: 0644]
include/deal.II/cgal/surface_mesh.h
include/deal.II/cgal/utilities.h

diff --git a/include/deal.II/cgal/point_conversion.h b/include/deal.II/cgal/point_conversion.h
new file mode 100644 (file)
index 0000000..9ccf20f
--- /dev/null
@@ -0,0 +1,105 @@
+// ---------------------------------------------------------------------
+//
+// 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.
+//
+// ---------------------------------------------------------------------
+
+#ifndef dealii_cgal_point_conversion_h
+#define dealii_cgal_point_conversion_h
+
+#include <deal.II/base/config.h>
+
+#include <deal.II/base/point.h>
+
+#ifdef DEAL_II_WITH_CGAL
+
+#  include <CGAL/Cartesian.h>
+#  include <CGAL/Exact_predicates_exact_constructions_kernel.h>
+#  include <CGAL/Exact_predicates_inexact_constructions_kernel.h>
+#  include <CGAL/Simple_cartesian.h>
+
+
+DEAL_II_NAMESPACE_OPEN
+namespace CGALWrappers
+{
+  /**
+   * Convert from a deal.II Point to any compatible CGAL point.
+   *
+   * @tparam CGALPointType Any of the CGAL point types
+   * @tparam dim Dimension of the point
+   * @param [in] p An input deal.II Point<dim>
+   * @return CGALPointType A CGAL point
+   */
+  template <typename CGALPointType, int dim>
+  inline CGALPointType
+  dealii_point_to_cgal_point(const dealii::Point<dim> &p);
+
+  /**
+   * Convert from various CGAL point types to deal.II Point.
+   *
+   * @tparam dim Dimension of the point
+   * @tparam CGALPointType Any of the CGAL point types
+   * @param p An input CGAL point type
+   * @return dealii::Point<dim> The corresponding deal.II point.
+   */
+  template <int dim, typename CGALPointType>
+  inline dealii::Point<dim>
+  cgal_point_to_dealii_point(const CGALPointType &p);
+
+
+#  ifndef DOXYGEN
+  // Template implementations
+
+  template <typename CGALPointType, int dim>
+  inline CGALPointType
+  dealii_point_to_cgal_point(const dealii::Point<dim> &p)
+  {
+    constexpr int cdim = CGALPointType::Ambient_dimension::value;
+    static_assert(dim <= cdim, "Only dim <= cdim supported");
+    if constexpr (cdim == 1)
+      return CGALPointType(p[0]);
+    else if constexpr (cdim == 2)
+      return CGALPointType(p[0], dim > 1 ? p[1] : 0);
+    else if constexpr (cdim == 3)
+      return CGALPointType(p[0], dim > 1 ? p[1] : 0, dim > 2 ? p[2] : 0);
+    else
+      Assert(false, dealii::ExcNotImplemented());
+    return CGALPointType();
+  }
+
+
+
+  template <int dim, typename CGALPointType>
+  inline dealii::Point<dim>
+  cgal_point_to_dealii_point(const CGALPointType &p)
+  {
+    constexpr int cdim = CGALPointType::Ambient_dimension::value;
+    if constexpr (dim == 1)
+      return dealii::Point<dim>(CGAL::to_double(p.x()));
+    else if constexpr (dim == 2)
+      return dealii::Point<dim>(CGAL::to_double(p.x()),
+                                cdim > 1 ? CGAL::to_double(p.y()) : 0);
+    else if constexpr (dim == 3)
+      return dealii::Point<dim>(CGAL::to_double(p.x()),
+                                cdim > 1 ? CGAL::to_double(p.y()) : 0,
+                                cdim > 2 ? CGAL::to_double(p.z()) : 0);
+    else
+      Assert(false, dealii::ExcNotImplemented());
+  }
+} // namespace CGALWrappers
+
+#  endif
+
+DEAL_II_NAMESPACE_CLOSE
+
+#endif
+#endif
index 818be5c06997fad7a933d721750315a2a16c1778..87f20fd5a44d6dc8ce2e6d7d4f6e347e534ca81b 100644 (file)
 
 #include <deal.II/grid/tria.h>
 
-#include <deal.II/cgal/utilities.h>
 
 #ifdef DEAL_II_WITH_CGAL
 #  include <CGAL/Polygon_mesh_processing/stitch_borders.h>
 #  include <CGAL/Surface_mesh.h>
+#  include <deal.II/cgal/point_conversion.h>
 
 
 DEAL_II_NAMESPACE_OPEN
index d94bc384cfcdc1593a13a7b2dd1a5a11f6395fe9..719a3e7d7df9581bb95bed6c01c4c5ba7a0b3f1b 100644 (file)
@@ -209,30 +209,6 @@ namespace CGALWrappers
 
 
 
-  /**
-   * Convert from a deal.II Point to any compatible CGAL point.
-   *
-   * @tparam CGALPointType Any of the CGAL point types
-   * @tparam dim Dimension of the point
-   * @param [in] p An input deal.II Point<dim>
-   * @return CGALPointType A CGAL point
-   */
-  template <typename CGALPointType, int dim>
-  inline CGALPointType
-  dealii_point_to_cgal_point(const dealii::Point<dim> &p);
-
-  /**
-   * Convert from various CGAL point types to deal.II Point.
-   *
-   * @tparam dim Dimension of the point
-   * @tparam CGALPointType Any of the CGAL point types
-   * @param p An input CGAL point type
-   * @return dealii::Point<dim> The corresponding deal.II point.
-   */
-  template <int dim, typename CGALPointType>
-  inline dealii::Point<dim>
-  cgal_point_to_dealii_point(const CGALPointType &p);
-
   /**
    * Given a closed CGAL::Surface_mesh, this function fills the
    * region bounded by the surface with tets.
@@ -375,45 +351,6 @@ namespace CGALWrappers
 // Template implementations
 namespace CGALWrappers
 {
-  template <typename CGALPointType, int dim>
-  inline CGALPointType
-  dealii_point_to_cgal_point(const dealii::Point<dim> &p)
-  {
-    constexpr int cdim = CGALPointType::Ambient_dimension::value;
-    static_assert(dim <= cdim, "Only dim <= cdim supported");
-    if constexpr (cdim == 1)
-      return CGALPointType(p[0]);
-    else if constexpr (cdim == 2)
-      return CGALPointType(p[0], dim > 1 ? p[1] : 0);
-    else if constexpr (cdim == 3)
-      return CGALPointType(p[0], dim > 1 ? p[1] : 0, dim > 2 ? p[2] : 0);
-    else
-      Assert(false, dealii::ExcNotImplemented());
-    return CGALPointType();
-  }
-
-
-
-  template <int dim, typename CGALPointType>
-  inline dealii::Point<dim>
-  cgal_point_to_dealii_point(const CGALPointType &p)
-  {
-    constexpr int cdim = CGALPointType::Ambient_dimension::value;
-    if constexpr (dim == 1)
-      return dealii::Point<dim>(CGAL::to_double(p.x()));
-    else if constexpr (dim == 2)
-      return dealii::Point<dim>(CGAL::to_double(p.x()),
-                                cdim > 1 ? CGAL::to_double(p.y()) : 0);
-    else if constexpr (dim == 3)
-      return dealii::Point<dim>(CGAL::to_double(p.x()),
-                                cdim > 1 ? CGAL::to_double(p.y()) : 0,
-                                cdim > 2 ? CGAL::to_double(p.z()) : 0);
-    else
-      Assert(false, dealii::ExcNotImplemented());
-  }
-
-
-
   template <typename C3t3>
   void
   cgal_surface_mesh_to_cgal_triangulation(

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.