Skip to content

Commit

Permalink
Merge pull request #13644 from luca-heltai/cgal-support
Browse files Browse the repository at this point in the history
Added support for CGAL library.
  • Loading branch information
peterrum committed May 2, 2022
2 parents 10e9568 + 92aa866 commit 82cf0e0
Show file tree
Hide file tree
Showing 12 changed files with 265 additions and 0 deletions.
1 change: 1 addition & 0 deletions .github/workflows/linux.yml
Original file line number Diff line number Diff line change
Expand Up @@ -65,6 +65,7 @@ jobs:
libp4est-dev \
trilinos-all-dev \
petsc-dev \
libcgal-dev \
libmetis-dev \
libhdf5-mpi-dev
- name: info
Expand Down
20 changes: 20 additions & 0 deletions cmake/configure/configure_cgal.cmake
Original file line number Diff line number Diff line change
@@ -0,0 +1,20 @@
## ---------------------------------------------------------------------
##
## 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.
##
## ---------------------------------------------------------------------

#
# Configuration for the CGAL library:
#

CONFIGURE_FEATURE(CGAL)
46 changes: 46 additions & 0 deletions cmake/modules/FindCGAL.cmake
Original file line number Diff line number Diff line change
@@ -0,0 +1,46 @@
## ---------------------------------------------------------------------
##
## 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.
##
## ---------------------------------------------------------------------

#
# 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 C++17. Disabling CGAL Support.")
ENDIF(DEAL_II_HAVE_CXX17)

DEAL_II_PACKAGE_HANDLE(CGAL
INCLUDE_DIRS REQUIRED CGAL_INCLUDE_DIRS
USER_INCLUDE_DIRS REQUIRED CGAL_INCLUDE_DIRS
CLEAR CGAL_INCLUDE_DIRS
)
1 change: 1 addition & 0 deletions doc/doxygen/options.dox.in
Original file line number Diff line number Diff line change
Expand Up @@ -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 \
Expand Down
3 changes: 3 additions & 0 deletions doc/news/changes/major/20220426FederHeltai
Original file line number Diff line number Diff line change
@@ -0,0 +1,3 @@
New: Added support for the CGAL library (www.cgal.org).
<br>
(Marco Feder, Luca Heltai, 2022/04/26)
26 changes: 26 additions & 0 deletions doc/readme.html
Original file line number Diff line number Diff line change
Expand Up @@ -455,6 +455,32 @@ <h4>Optional interfaces to other software packages</h4>
</p>
</dd>

<dt><a name="CGAL" />
<a href="https://www.cgal.org/" target="_top">CGAL</a>
</dt>
<dd>
<p>
<a href="https://www.cgal.org/" target="_top">CGAL</a>
is a software project that provides easy access to
efficient and reliable geometric algorithms in the form of a C++
library. CGAL is used in various areas needing geometric
computation, such as geographic information systems, computer
aided design, molecular biology, medical imaging, computer
graphics, and robotics.

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.

The CGAL Wrappers use c++17 features. In order to use CGAL, add
<code>-DCGAL_DIR=/path/to/cgal</code> to the deal.II CMake call,
and make sure that <code>DEAL_II_HAVE_CXX17</code> is enabled.
</p>
</dd>


<dt><a name="CUDA"></a><a href="https://developer.nvidia.com/cuda-zone">CUDA</a></dt>
<dd>
Expand Down
1 change: 1 addition & 0 deletions doc/users/cmake_dealii.html
Original file line number Diff line number Diff line change
Expand Up @@ -456,6 +456,7 @@ <h3>Feature configuration</h3>
DEAL_II_WITH_ARPACK
DEAL_II_WITH_ASSIMP
DEAL_II_WITH_BOOST
DEAL_II_WITH_CGAL
DEAL_II_WITH_COMPLEX_VALUES
DEAL_II_WITH_CUDA
DEAL_II_WITH_GINKGO
Expand Down
1 change: 1 addition & 0 deletions include/deal.II/base/config.h.in
Original file line number Diff line number Diff line change
Expand Up @@ -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
Expand Down
111 changes: 111 additions & 0 deletions include/deal.II/cgal/utilities.h
Original file line number Diff line number Diff line change
@@ -0,0 +1,111 @@
// ---------------------------------------------------------------------
//
// 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_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
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);
} // namespace CGALWrappers

# ifndef DOXYGEN
// 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());
}
} // namespace CGALWrappers
# endif

DEAL_II_NAMESPACE_CLOSE

#endif
#endif
6 changes: 6 additions & 0 deletions tests/cgal/CMakeLists.txt
Original file line number Diff line number Diff line change
@@ -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()
46 changes: 46 additions & 0 deletions tests/cgal/cgal_point_01.cc
Original file line number Diff line number Diff line change
@@ -0,0 +1,46 @@
// ---------------------------------------------------------------------
//
// 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.
//
// ---------------------------------------------------------------------

// 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> dealii_point(1.0, 2.0, 3.0);
const auto cgal_point = dealii_point_to_cgal_point<CGALPoint>(dealii_point);
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 dealii_point = cgal_point_to_dealii_point<3>(cgal_point);
deallog << "deal.II Point: " << dealii_point << std::endl;
}
}
3 changes: 3 additions & 0 deletions tests/cgal/cgal_point_01.output
Original file line number Diff line number Diff line change
@@ -0,0 +1,3 @@

DEAL::CGAL Point: 1.00000 2.00000 3.00000
DEAL::deal.II Point: 1.00000 2.00000 3.00000

0 comments on commit 82cf0e0

Please sign in to comment.