Skip to content

Commit

Permalink
Add CGAL Surface_mesh support
Browse files Browse the repository at this point in the history
  • Loading branch information
fdrmrc committed Apr 28, 2022
1 parent c02086a commit 591a201
Show file tree
Hide file tree
Showing 17 changed files with 672 additions and 0 deletions.
5 changes: 5 additions & 0 deletions cmake/config/template-arguments.in
Original file line number Diff line number Diff line change
Expand Up @@ -277,3 +277,8 @@ SYM_RANKS := { 2; 4 }
OUTPUT_FLAG_TYPES := { DXFlags; UcdFlags; GnuplotFlags; PovrayFlags; EpsFlags;
GmvFlags; TecplotFlags; VtkFlags; SvgFlags;
Deal_II_IntermediateFlags }

// CGAL Kernels
CGAL_KERNELS := {CGAL::Simple_cartesian<double>; CGAL::Exact_predicates_exact_constructions_kernel;
CGAL::Exact_predicates_exact_constructions_kernel_with_sqrt;
CGAL::Exact_predicates_inexact_constructions_kernel }
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) 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)
45 changes: 45 additions & 0 deletions cmake/modules/FindCGAL.cmake
Original file line number Diff line number Diff line change
@@ -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
)
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)
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
60 changes: 60 additions & 0 deletions include/deal.II/cgal/surface_mesh.h
Original file line number Diff line number Diff line change
@@ -0,0 +1,60 @@
// ---------------------------------------------------------------------
//
// 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_surface_mesh_h
#define dealii_cgal_surface_mesh_h

#include <deal.II/base/config.h>

#include <deal.II/fe/mapping.h>

#include <deal.II/grid/tria.h>

#include <deal.II/cgal/utilities.h>

#ifdef DEAL_II_WITH_CGAL
# include <CGAL/Surface_mesh.h>

DEAL_II_NAMESPACE_OPEN

namespace CGALWrappers
{
/**
* Build a CGAL::Surface_mesh from a deal.II cell. The class Surface_mesh is
* an implementation of a halfedge data structure and can be used to represent
* a polyhedral surface. It is an edge-centered data structure
* capable of maintaining incidence information of vertices, edges, and faces.
* Each edge is represented by two halfedges with opposite orientation. The
* orientation of a face is chosen so that the halfedges around a
* face are oriented counterclockwise.
*
* @param[in] cell The input deal.II cell iterator
* @param[in] mapping The mapping used to map the vertices of the cell
* @param[out] mesh The output CGAL::Surface_mesh
*/
template <typename CGALPointType, int dim, int spacedim>
void
to_cgal_mesh(
const typename dealii::Triangulation<dim, spacedim>::cell_iterator &cell,
const dealii::Mapping<dim, spacedim> & mapping,
CGAL::Surface_mesh<CGALPointType> & mesh);
} // namespace CGALWrappers



DEAL_II_NAMESPACE_CLOSE

#endif
#endif
114 changes: 114 additions & 0 deletions include/deal.II/cgal/utilities.h
Original file line number Diff line number Diff line change
@@ -0,0 +1,114 @@
// ---------------------------------------------------------------------
//
// 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/Exact_predicates_exact_constructions_kernel.h>
# include <CGAL/Exact_predicates_exact_constructions_kernel_with_sqrt.h>
# include <CGAL/Exact_predicates_inexact_constructions_kernel.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
1 change: 1 addition & 0 deletions source/CMakeLists.txt
Original file line number Diff line number Diff line change
Expand Up @@ -56,6 +56,7 @@ ADD_SUBDIRECTORY(fe)
ADD_SUBDIRECTORY(dofs)
ADD_SUBDIRECTORY(lac)
ADD_SUBDIRECTORY(base)
ADD_SUBDIRECTORY(cgal)
ADD_SUBDIRECTORY(gmsh)
ADD_SUBDIRECTORY(grid)
ADD_SUBDIRECTORY(hp)
Expand Down
35 changes: 35 additions & 0 deletions source/cgal/CMakeLists.txt
Original file line number Diff line number Diff line change
@@ -0,0 +1,35 @@
## ---------------------------------------------------------------------
##
## Copyright (C) 2012 - 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.
##
## ---------------------------------------------------------------------


INCLUDE_DIRECTORIES(BEFORE ${CMAKE_CURRENT_BINARY_DIR})


SET(_src
surface_mesh.cc
)



SET(_inst
surface_mesh.inst.in
)

FILE(GLOB _header
${CMAKE_SOURCE_DIR}/include/deal.II/cgal/*.h
)

DEAL_II_ADD_LIBRARY(obj_cgal OBJECT ${_src} ${_header} ${_inst})
EXPAND_INSTANTIATIONS(obj_cgal "${_inst}")

0 comments on commit 591a201

Please sign in to comment.