]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Added support for CGAL library.
authorMarco Feder <marco.feder@sissa.it>
Tue, 26 Apr 2022 16:32:00 +0000 (19:32 +0300)
committerLuca Heltai <luca.heltai@sissa.it>
Thu, 28 Apr 2022 16:35:48 +0000 (16:35 +0000)
Co-authored-by: Luca Heltai <luca.heltai@gmail.com>
cmake/configure/configure_cgal.cmake [new file with mode: 0644]
cmake/modules/FindCGAL.cmake [new file with mode: 0644]
doc/doxygen/options.dox.in
doc/news/changes/major/20220426FederHeltai [new file with mode: 0644]
include/deal.II/base/config.h.in
include/deal.II/cgal/utilities.h [new file with mode: 0644]
tests/cgal/CMakeLists.txt [new file with mode: 0644]
tests/cgal/cgal_point_01.cc [new file with mode: 0644]
tests/cgal/cgal_point_01.output [new file with mode: 0644]

diff --git a/cmake/configure/configure_cgal.cmake b/cmake/configure/configure_cgal.cmake
new file mode 100644 (file)
index 0000000..d13ac06
--- /dev/null
@@ -0,0 +1,20 @@
+## ---------------------------------------------------------------------
+##
+## Copyright (C) 2017 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.
+##
+## ---------------------------------------------------------------------
+
+#
+# Configuration for the CGAL library:
+#
+
+CONFIGURE_FEATURE(CGAL)
diff --git a/cmake/modules/FindCGAL.cmake b/cmake/modules/FindCGAL.cmake
new file mode 100644 (file)
index 0000000..34a7716
--- /dev/null
@@ -0,0 +1,45 @@
+## ---------------------------------------------------------------------
+##
+## Copyright (C) 2017 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.
+##
+## ---------------------------------------------------------------------
+
+#
+# Try to find the CGAL libraries
+#
+# This module exports
+#   
+#   CGAL_INCLUDE_DIRS
+#
+
+SET(CGAL_DIR "" CACHE PATH "An optional hint to a CGAL installation")
+SET_IF_EMPTY(CGAL_DIR "$ENV{CGAL_DIR}")
+
+IF(NOT "${CGAL_DIR}" STREQUAL "")
+  SET(CGAL_DIR ${CGAL_DIR})
+ENDIF()
+
+IF(DEAL_II_HAVE_CXX17)
+  # temporarily disable ${CMAKE_SOURCE_DIR}/cmake/modules for module lookup
+  LIST(REMOVE_ITEM CMAKE_MODULE_PATH ${CMAKE_SOURCE_DIR}/cmake/modules/)
+  FIND_PACKAGE(CGAL)
+  LIST(APPEND CMAKE_MODULE_PATH ${CMAKE_SOURCE_DIR}/cmake/modules/)
+ELSE(DEAL_II_HAVE_CXX17)
+  SET(CGAL_FOUND FALSE)
+  SET(CGAL_INCLUDE_DIRS "-NOTFOUND")
+  MESSAGE("-- WARNING: CGAL wrappers require cxx17. Disabling CGAL Support.")
+ENDIF(DEAL_II_HAVE_CXX17)
+
+DEAL_II_PACKAGE_HANDLE(CGAL
+  INCLUDE_DIRS REQUIRED CGAL_INCLUDE_DIRS
+  CLEAR CGAL_INCLUDE_DIRS
+  )
index 25a3c15df49f9504355a340ef4fb879a093bd569..308f28dd3c1c24cd0c371477a58e19db0533aa14 100644 (file)
@@ -194,6 +194,7 @@ PREDEFINED             = DOXYGEN=1 \
                          DEAL_II_ARPACK_WITH_PARPACK=1 \
                          DEAL_II_WITH_ASSIMP=1 \
                          DEAL_II_WITH_BOOST=1 \
+                         DEAL_II_WITH_CGAL=1 \
                          DEAL_II_WITH_TASKFLOW=1 \
                          DEAL_II_WITH_COMPLEX_VALUES=1 \
                          DEAL_II_WITH_CUDA=1 \
diff --git a/doc/news/changes/major/20220426FederHeltai b/doc/news/changes/major/20220426FederHeltai
new file mode 100644 (file)
index 0000000..e4d7d88
--- /dev/null
@@ -0,0 +1,3 @@
+New: Added support for the CGAL library (www.cgal.org).
+<br>
+(Marco Feder, Luca Heltai, 2022/04/26)
index 7d0e646b393f0c36e00d82176682021212ab3296..aff1915bcbaeddf3c97b0747bff05731a4b5d1cf 100644 (file)
@@ -39,6 +39,7 @@
 #cmakedefine DEAL_II_WITH_ARPACK
 #cmakedefine DEAL_II_WITH_ARBORX
 #cmakedefine DEAL_II_WITH_ASSIMP
+#cmakedefine DEAL_II_WITH_CGAL
 #cmakedefine DEAL_II_WITH_COMPLEX_VALUES
 #cmakedefine DEAL_II_WITH_CUDA
 #cmakedefine DEAL_II_WITH_GINKGO
diff --git a/include/deal.II/cgal/utilities.h b/include/deal.II/cgal/utilities.h
new file mode 100644 (file)
index 0000000..deba433
--- /dev/null
@@ -0,0 +1,110 @@
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2020 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_utilities_h
+#define dealii_cgal_utilities_h
+
+#include <deal.II/base/config.h>
+
+#ifdef DEAL_II_WITH_CGAL
+#  include <CGAL/Cartesian.h>
+#  include <CGAL/Simple_cartesian.h>
+
+DEAL_II_NAMESPACE_OPEN
+/**
+ * Interface to the Computational Geometry Algorithm Library (CGAL).
+ *
+ * CGAL is a software project that provides easy access to efficient and
+ * reliable geometric algorithms. The library offers data structures and
+ * algorithms like triangulations, Voronoi diagrams, Boolean operations on
+ * polygons and polyhedra, point set processing, arrangements of curves, surface
+ * and volume mesh generation, geometry processing,  alpha shapes, convex hull
+ * algorithms, shape reconstruction, AABB and KD trees...
+ *
+ * You can learn more about the CGAL library at https://www.cgal.org/
+ */
+namespace CGALWrappers
+{
+  /**
+   * Convert from 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
+  to_cgal(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
+   * @return dealii::Point<dim>
+   */
+  template <int dim, typename CGALPointType>
+  inline dealii::Point<dim>
+  to_dealii(const CGALPointType &p);
+} // namespace CGALWrappers
+
+#  ifndef DOXYGEN
+// Template implementations
+namespace CGALWrappers
+{
+  template <typename CGALPointType, int dim>
+  inline CGALPointType
+  to_cgal(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>
+  to_dealii(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
diff --git a/tests/cgal/CMakeLists.txt b/tests/cgal/CMakeLists.txt
new file mode 100644 (file)
index 0000000..9c78c28
--- /dev/null
@@ -0,0 +1,6 @@
+CMAKE_MINIMUM_REQUIRED(VERSION 3.1.0)
+INCLUDE(../setup_testsubproject.cmake)
+PROJECT(testsuite CXX)
+IF(DEAL_II_WITH_CGAL)
+  DEAL_II_PICKUP_TESTS()
+ENDIF()
diff --git a/tests/cgal/cgal_point_01.cc b/tests/cgal/cgal_point_01.cc
new file mode 100644 (file)
index 0000000..26ca14c
--- /dev/null
@@ -0,0 +1,46 @@
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2020 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.
+//
+// ---------------------------------------------------------------------
+
+// Test conversion from CGAL to deal.II Point and viceversa
+
+#include <deal.II/base/point.h>
+
+#include <CGAL/IO/io.h>
+#include <CGAL/Simple_cartesian.h>
+#include <deal.II/cgal/utilities.h>
+
+#include "../tests.h"
+
+using namespace CGALWrappers;
+
+int
+main()
+{
+  initlog();
+  using CGALPoint = CGAL::Point_3<CGAL::Simple_cartesian<double>>;
+  // Test conversion from deal.II Point to CGAL Point
+  {
+    const Point<3> p(1.0, 2.0, 3.0);
+    const auto     cgal_point = to_cgal<CGALPoint>(p);
+    deallog << "CGAL Point: " << cgal_point << std::endl;
+  }
+
+  // Test conversion from CGAL Point to deal.II Point
+  {
+    const CGALPoint cgal_point(1.0, 2.0, 3.0);
+    const auto      deal_ii_point = to_dealii<3>(cgal_point);
+    deallog << "deal.II Point: " << deal_ii_point << std::endl;
+  }
+}
diff --git a/tests/cgal/cgal_point_01.output b/tests/cgal/cgal_point_01.output
new file mode 100644 (file)
index 0000000..0da2203
--- /dev/null
@@ -0,0 +1,3 @@
+
+DEAL::CGAL Point: 1.00000 2.00000 3.00000
+DEAL::deal.II Point: 1.00000 2.00000 3.00000

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.