12#ifndef ALPHA_COMPLEX_H_ 
   13#define ALPHA_COMPLEX_H_ 
   15#include <gudhi/Alpha_complex/Alpha_kernel_d.h> 
   16#include <gudhi/Debug_utils.h> 
   18#include <gudhi/Points_off_io.h> 
   24#include <CGAL/Delaunay_triangulation.h> 
   25#include <CGAL/Regular_triangulation.h>   
   26#include <CGAL/Epeck_d.h>   
   27#include <CGAL/Epick_d.h>   
   28#include <CGAL/Spatial_sort_traits_adapter_d.h> 
   29#include <CGAL/property_map.h>   
   30#include <CGAL/version.h>   
   31#include <CGAL/NT_converter.h> 
   33#include <Eigen/src/Core/util/Macros.h>   
   35#include <boost/range/size.hpp> 
   36#include <boost/range/combine.hpp> 
   37#include <boost/range/adaptor/transformed.hpp> 
   50#if CGAL_VERSION_NR < 1041101000 
   51# error Alpha_complex is only available for CGAL >= 4.11 
   54#if !EIGEN_VERSION_AT_LEAST(3,1,0) 
   55# error Alpha_complex is only available for Eigen3 >= 3.1.0 installed with CGAL 
   60namespace alpha_complex {
 
   62template<
typename D> 
struct Is_Epeck_D { 
static const bool value = 
false; };
 
   63template<
typename D> 
struct Is_Epeck_D<CGAL::Epeck_d<D>> { 
static const bool value = 
true; };
 
  102template<
class Kernel = CGAL::Epeck_d<CGAL::Dynamic_dimension_tag>, 
bool Weighted = false>
 
  106  using Internal_vertex_handle = std::ptrdiff_t;
 
  110  using Geom_traits = std::conditional_t<Weighted, CGAL::Regular_triangulation_traits_adapter<Kernel>, Kernel>;
 
  113  using TDS = CGAL::Triangulation_data_structure<
typename Geom_traits::Dimension,
 
  114                                                 CGAL::Triangulation_vertex<Geom_traits, Internal_vertex_handle>,
 
  115                                                 CGAL::Triangulation_full_cell<Geom_traits> >;
 
  118  using Triangulation = std::conditional_t<Weighted, CGAL::Regular_triangulation<Kernel, TDS>,
 
  119                                                     CGAL::Delaunay_triangulation<Kernel, TDS>>;
 
  125  using FT = 
typename A_kernel_d::FT;
 
  130  using Sphere = 
typename A_kernel_d::Sphere;
 
  133  using Point_d = 
typename Geom_traits::Point_d;
 
  137  using CGAL_vertex_iterator = 
typename Triangulation::Vertex_iterator;
 
  140  using Vector_vertex_iterator = std::vector< CGAL_vertex_iterator >;
 
  145  Vector_vertex_iterator vertex_handle_to_iterator_;
 
  147  std::unique_ptr<Triangulation> triangulation_;
 
  153  std::vector<Internal_vertex_handle> vertices_;
 
  156  std::vector<Sphere> cache_, old_cache_;
 
  171      std::cerr << 
"Alpha_complex - Unable to read file " << off_file_name << 
"\n";
 
  188  template<
typename InputPo
intRange >
 
  190    init_from_range(points);
 
  205  template <
typename InputPo
intRange, 
typename WeightRange>
 
  207    static_assert(Weighted, 
"This constructor is not available for non-weighted versions of Alpha_complex");
 
  209    GUDHI_CHECK(boost::size(weights) == boost::size(points),
 
  210                std::invalid_argument(
"Points number in range different from weights range number"));
 
  211    auto weighted_points = boost::range::combine(points, weights)
 
  212      | boost::adaptors::transformed([](
auto const&t){
return Point_d(boost::get<0>(t), boost::get<1>(t));});
 
  213    init_from_range(weighted_points);
 
  225    if (triangulation_ == 
nullptr)
 
  228      return triangulation_->number_of_vertices();
 
  238    return vertex_handle_to_iterator_.at(vertex)->point();
 
  242  template<
typename InputPo
intRange >
 
  243  void init_from_range(
const InputPointRange& points) {
 
  244    #if CGAL_VERSION_NR < 1050000000 
  245    if (Is_Epeck_D<Kernel>::value)
 
  246      std::cerr << 
"It is strongly advised to use a CGAL version >= 5.0 with Epeck_d Kernel for performance reasons." 
  250#if CGAL_VERSION_NR < 1050101000 
  252    static_assert(!Weighted, 
"Weighted Alpha_complex is only available for CGAL >= 5.1");
 
  255    auto first = std::begin(points);
 
  256    auto last = std::end(points);
 
  260      triangulation_ = std::make_unique<Triangulation>(kernel_.get_dimension(*first));
 
  262      std::vector<Point_d> point_cloud(first, last);
 
  265      std::vector<Internal_vertex_handle> indices(boost::counting_iterator<Internal_vertex_handle>(0),
 
  266                                                  boost::counting_iterator<Internal_vertex_handle>(point_cloud.size()));
 
  268      using Point_property_map = boost::iterator_property_map<typename std::vector<Point_d>::iterator,
 
  269                                                              CGAL::Identity_property_map<Internal_vertex_handle>>;
 
  270      using Search_traits_d = CGAL::Spatial_sort_traits_adapter_d<Geom_traits, Point_property_map>;
 
  272      CGAL::spatial_sort(indices.begin(), indices.end(), Search_traits_d(std::begin(point_cloud)));
 
  274      typename Triangulation::Full_cell_handle hint;
 
  275      for (
auto index : indices) {
 
  276        typename Triangulation::Vertex_handle pos = triangulation_->insert(point_cloud[index], hint);
 
  277        if (pos != 
nullptr) {
 
  280          hint = pos->full_cell();
 
  286      vertex_handle_to_iterator_.resize(point_cloud.size());
 
  289      vertices_.reserve(triangulation_->number_of_vertices());
 
  291      for (CGAL_vertex_iterator vit = triangulation_->vertices_begin(); vit != triangulation_->vertices_end(); ++vit) {
 
  292        if (!triangulation_->is_infinite(*vit)) {
 
  294          std::clog << 
"Vertex insertion - " << vit->data() << 
" -> " << vit->point() << std::endl;
 
  296          vertex_handle_to_iterator_[vit->data()] = vit;
 
  297          vertices_.push_back(vit->data());
 
  300      std::sort(vertices_.begin(), vertices_.end());
 
  311  const Point_d& get_point_(std::size_t vertex)
 const {
 
  312    return vertex_handle_to_iterator_[vertex]->point();
 
  316  template<
class SimplicialComplexForAlpha>
 
  318    auto k = cplx.key(s);
 
  319    if(k==cplx.null_key()){
 
  321      cplx.assign_key(s, k);
 
  323      thread_local std::vector<Point_d> v;
 
  325      for (
auto vertex : cplx.simplex_vertex_range(s))
 
  326        v.push_back(get_point_(vertex));
 
  327      cache_.emplace_back(kernel_.get_sphere(v.cbegin(), v.cend()));
 
  333  template<
class SimplicialComplexForAlpha>
 
  335    auto k = cplx.key(s);
 
  336    if(k!=cplx.null_key())
 
  337      return kernel_.get_squared_radius(old_cache_[k]);
 
  339    thread_local std::vector<Point_d> v;
 
  341    for (
auto vertex : cplx.simplex_vertex_range(s))
 
  342      v.push_back(get_point_(vertex));
 
  343    return kernel_.get_squared_radius(v.cbegin(), v.cend());
 
  370  template <
typename SimplicialComplexForAlpha,
 
  373                      Filtration_value max_alpha_square = std::numeric_limits<Filtration_value>::infinity(),
 
  375                      bool default_filtration_value = 
false) {
 
  379    using Vector_vertex = std::vector<Vertex_handle>;
 
  381    if (triangulation_ == 
nullptr) {
 
  382      std::cerr << 
"Alpha_complex cannot create_complex from a NULL triangulation\n";
 
  385    if (triangulation_->maximal_dimension() < 1) {
 
  386      std::cerr << 
"Alpha_complex cannot create_complex from a zero-dimension triangulation\n";
 
  390      std::cerr << 
"Alpha_complex create_complex - complex is not empty\n";
 
  397      std::vector<Vertex_handle> one_vertex(1);
 
  398      for (
auto vertex : vertices_) {
 
  400        std::clog << 
"SimplicialComplex insertion " << vertex << std::endl;
 
  402        one_vertex[0] = vertex;
 
  406      for (
auto cit = triangulation_->finite_full_cells_begin();
 
  407           cit != triangulation_->finite_full_cells_end();
 
  409        Vector_vertex vertexVector;
 
  411        std::clog << 
"SimplicialComplex insertion ";
 
  413        for (
auto vit = cit->vertices_begin(); vit != cit->vertices_end(); ++vit) {
 
  414          if (*vit != 
nullptr) {
 
  416            std::clog << 
" " << (*vit)->data();
 
  419            vertexVector.push_back((*vit)->data());
 
  423        std::clog << std::endl;
 
  431    if (!default_filtration_value) {
 
  432      CGAL::NT_converter<FT, Filtration_value> cgal_converter;
 
  435      for (
int decr_dim = triangulation_->maximal_dimension(); decr_dim >= 0; decr_dim--) {
 
  438          int f_simplex_dim = complex.
dimension(f_simplex);
 
  439          if (decr_dim == f_simplex_dim) {
 
  441            if (std::isnan(complex.filtration(f_simplex))) {
 
  444              if (Weighted || f_simplex_dim > 0) {
 
  445                auto const& sqrad = radius(complex, f_simplex);
 
  446#if CGAL_VERSION_NR >= 1050000000 
  447                if(exact) CGAL::exact(sqrad);
 
  449                alpha_complex_filtration = cgal_converter(sqrad);
 
  453              std::clog << 
"filt(Sigma) is NaN : filt(Sigma) =" << complex.filtration(f_simplex) << std::endl;
 
  457            if (decr_dim > !Weighted)
 
  458              propagate_alpha_filtration(complex, f_simplex);
 
  461        old_cache_ = std::move(cache_);
 
  479  template <
typename SimplicialComplexForAlpha, 
typename Simplex_handle>
 
  485    for (
auto face_opposite_vertex : complex.boundary_opposite_vertex_simplex_range(f_simplex)) {
 
  486      auto f_boundary = face_opposite_vertex.first;
 
  488      std::clog << 
" | --------------------------------------------------\n";
 
  489      std::clog << 
" | Tau ";
 
  491        std::clog << vertex << 
" ";
 
  493      std::clog << 
"is a face of Sigma\n";
 
  494      std::clog << 
" | isnan(complex.filtration(Tau)=" << std::isnan(complex.filtration(f_boundary)) << std::endl;
 
  497      if (!std::isnan(complex.filtration(f_boundary))) {
 
  499        Filtration_value alpha_complex_filtration = fmin(complex.filtration(f_boundary),
 
  500                                                                             complex.filtration(f_simplex));
 
  503        std::clog << 
" | filt(Tau) = fmin(filt(Tau), filt(Sigma)) = " << complex.filtration(f_boundary) << std::endl;
 
  507        auto const& cache=get_cache(complex, f_boundary);
 
  508        bool is_gab = kernel_.is_gabriel(cache, get_point_(face_opposite_vertex.second));
 
  510        std::clog << 
" | Tau is_gabriel(Sigma)=" << is_gab << 
" - vertexForGabriel=" << face_opposite_vertex.second << std::endl;
 
  513        if (
false == is_gab) {
 
  518          std::clog << 
" | filt(Tau) = filt(Sigma) = " << complex.filtration(f_boundary) << std::endl;
 
  528namespace alphacomplex = alpha_complex;
 
OFF file reader implementation in order to read points from an OFF file.
Definition: Points_off_io.h:122
const std::vector< Point_d > & get_point_cloud() const
Point cloud getter.
Definition: Points_off_io.h:158
bool is_valid() const
Returns if the OFF file read operation was successful or not.
Definition: Points_off_io.h:150
Alpha complex data structure.
Definition: Alpha_complex.h:103
bool create_complex(SimplicialComplexForAlpha &complex, Filtration_value max_alpha_square=std::numeric_limits< Filtration_value >::infinity(), bool exact=false, bool default_filtration_value=false)
Inserts all Delaunay triangulation into the simplicial complex. It also computes the filtration value...
Definition: Alpha_complex.h:372
std::conditional_t< Weighted, CGAL::Regular_triangulation< Kernel, TDS >, CGAL::Delaunay_triangulation< Kernel, TDS > > Triangulation
A (Weighted or not) Delaunay triangulation of a set of points in .
Definition: Alpha_complex.h:119
typename A_kernel_d::Sphere Sphere
Sphere is a std::pair<Kernel::Point_d, Kernel::FT> (aka. circurmcenter and squared radius)....
Definition: Alpha_complex.h:130
Alpha_complex(const InputPointRange &points, WeightRange weights)
Alpha_complex constructor from a list of points and weights.
Definition: Alpha_complex.h:206
std::conditional_t< Weighted, CGAL::Regular_triangulation_traits_adapter< Kernel >, Kernel > Geom_traits
Geometric traits class that provides the geometric types and predicates needed by the triangulations.
Definition: Alpha_complex.h:110
const Point_d & get_point(std::size_t vertex) const
get_point returns the point corresponding to the vertex given as parameter.
Definition: Alpha_complex.h:237
Alpha_complex(const std::string &off_file_name)
Alpha_complex constructor from an OFF file name.
Definition: Alpha_complex.h:168
typename Geom_traits::Point_d Point_d
A point, or a weighted point in Euclidean space.
Definition: Alpha_complex.h:133
std::size_t num_vertices() const
Returns the number of finite vertices in the triangulation.
Definition: Alpha_complex.h:224
Alpha_complex(const InputPointRange &points)
Alpha_complex constructor from a list of points.
Definition: Alpha_complex.h:189
Alpha complex kernel container.
Definition: Alpha_kernel_d.h:42
Value type for a filtration function on a cell complex.
Definition: FiltrationValue.h:20
The concept SimplicialComplexForAlpha describes the requirements for a type to implement a simplicial...
Definition: SimplicialComplexForAlpha.h:21
Skeleton_simplex_range skeleton_simplex_range
Returns a range over the simplices of the skeleton of the simplicial complex, for a given dimension.
Definition: SimplicialComplexForAlpha.h:70
int assign_filtration(Simplex_handle simplex, Filtration_value filtration)
std::size_t num_vertices()
int dimension(Simplex_handle simplex)
void prune_above_filtration(Filtration_value filtration)
void make_filtration_non_decreasing()
Simplex_vertex_range simplex_vertex_range(Simplex_handle const &simplex)
Returns a range over vertices of a given simplex.
void insert_simplex_and_subfaces(std::vector< Vertex_handle > const &vertex_range, Filtration_value filtration)
Inserts a simplex with vertices from a given simplex (represented by a vector of Vertex_handle) in th...
unspecified Simplex_handle
Definition: SimplicialComplexForAlpha.h:23
unspecified Vertex_handle
Definition: SimplicialComplexForAlpha.h:25
unspecified Filtration_value
Definition: SimplicialComplexForAlpha.h:27
Handle type for the vertices of a cell complex.
Definition: VertexHandle.h:15