Alpha_complex.h
1 /* This file is part of the Gudhi Library. The Gudhi library
2  * (Geometric Understanding in Higher Dimensions) is a generic C++
3  * library for computational topology.
4  *
5  * Author(s): Vincent Rouvreau
6  *
7  * Copyright (C) 2015 Inria
8  *
9  * This program is free software: you can redistribute it and/or modify
10  * it under the terms of the GNU General Public License as published by
11  * the Free Software Foundation, either version 3 of the License, or
12  * (at your option) any later version.
13  *
14  * This program is distributed in the hope that it will be useful,
15  * but WITHOUT ANY WARRANTY; without even the implied warranty of
16  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
17  * GNU General Public License for more details.
18  *
19  * You should have received a copy of the GNU General Public License
20  * along with this program. If not, see <http://www.gnu.org/licenses/>.
21  */
22 
23 #ifndef ALPHA_COMPLEX_H_
24 #define ALPHA_COMPLEX_H_
25 
26 #include <gudhi/Debug_utils.h>
27 // to construct Alpha_complex from a OFF file of points
28 #include <gudhi/Points_off_io.h>
29 
30 #include <stdlib.h>
31 #include <math.h> // isnan, fmax
32 
33 #include <CGAL/Delaunay_triangulation.h>
34 #include <CGAL/Epick_d.h>
35 #include <CGAL/Spatial_sort_traits_adapter_d.h>
36 #include <CGAL/property_map.h> // for CGAL::Identity_property_map
37 #include <CGAL/NT_converter.h>
38 
39 #include <iostream>
40 #include <vector>
41 #include <string>
42 #include <limits> // NaN
43 #include <map>
44 #include <utility> // std::pair
45 #include <stdexcept>
46 #include <numeric> // for std::iota
47 
48 namespace Gudhi {
49 
50 namespace alpha_complex {
51 
75 template<class Kernel = CGAL::Epick_d<CGAL::Dynamic_dimension_tag>>
77  public:
78  // Add an int in TDS to save point index in the structure
79  typedef CGAL::Triangulation_data_structure<typename Kernel::Dimension,
80  CGAL::Triangulation_vertex<Kernel, std::ptrdiff_t>,
81  CGAL::Triangulation_full_cell<Kernel> > TDS;
83  typedef CGAL::Delaunay_triangulation<Kernel, TDS> Delaunay_triangulation;
84 
86  typedef typename Kernel::Point_d Point_d;
89  typedef Kernel Geom_traits;
90 
91  private:
92  typedef typename Kernel::Compute_squared_radius_d Squared_Radius;
93  typedef typename Kernel::Side_of_bounded_sphere_d Is_Gabriel;
94  typedef typename Kernel::Point_dimension_d Point_Dimension;
95 
96  // Type required to compute squared radius, or side of bounded sphere on a vector of points.
97  typedef typename std::vector<Point_d> Vector_of_CGAL_points;
98 
99  // Vertex_iterator type from CGAL.
100  typedef typename Delaunay_triangulation::Vertex_iterator CGAL_vertex_iterator;
101 
102  // size_type type from CGAL.
103  typedef typename Delaunay_triangulation::size_type size_type;
104 
105  // Map type to switch from simplex tree vertex handle to CGAL vertex iterator.
106  typedef typename std::map< std::size_t, CGAL_vertex_iterator > Vector_vertex_iterator;
107 
108  private:
111  Vector_vertex_iterator vertex_handle_to_iterator_;
113  Delaunay_triangulation* triangulation_;
115  Kernel kernel_;
116 
117  public:
127  Alpha_complex(const std::string& off_file_name)
128  : triangulation_(nullptr) {
129  Gudhi::Points_off_reader<Point_d> off_reader(off_file_name);
130  if (!off_reader.is_valid()) {
131  std::cerr << "Alpha_complex - Unable to read file " << off_file_name << "\n";
132  exit(-1); // ----- >>
133  }
134 
135  init_from_range(off_reader.get_point_cloud());
136  }
137 
147  template<typename InputPointRange >
148  Alpha_complex(const InputPointRange& points)
149  : triangulation_(nullptr) {
150  init_from_range(points);
151  }
152 
156  delete triangulation_;
157  }
158 
159  // Forbid copy/move constructor/assignment operator
160  Alpha_complex(const Alpha_complex& other) = delete;
161  Alpha_complex& operator= (const Alpha_complex& other) = delete;
162  Alpha_complex (Alpha_complex&& other) = delete;
163  Alpha_complex& operator= (Alpha_complex&& other) = delete;
164 
171  const Point_d& get_point(std::size_t vertex) const {
172  return vertex_handle_to_iterator_.at(vertex)->point();
173  }
174 
179  std::size_t number_of_vertices() const {
180  return vertex_handle_to_iterator_.size();
181  }
182 
183  private:
184  template<typename InputPointRange >
185  void init_from_range(const InputPointRange& points) {
186  auto first = std::begin(points);
187  auto last = std::end(points);
188 
189  if (first != last) {
190  // point_dimension function initialization
191  Point_Dimension point_dimension = kernel_.point_dimension_d_object();
192 
193  // Delaunay triangulation is point dimension.
194  triangulation_ = new Delaunay_triangulation(point_dimension(*first));
195 
196  std::vector<Point_d> point_cloud(first, last);
197 
198  // Creates a vector {0, 1, ..., N-1}
199  std::vector<std::ptrdiff_t> indices(boost::counting_iterator<std::ptrdiff_t>(0),
200  boost::counting_iterator<std::ptrdiff_t>(point_cloud.size()));
201 
202  typedef boost::iterator_property_map<typename std::vector<Point_d>::iterator,
203  CGAL::Identity_property_map<std::ptrdiff_t>> Point_property_map;
204  typedef CGAL::Spatial_sort_traits_adapter_d<Kernel, Point_property_map> Search_traits_d;
205 
206  CGAL::spatial_sort(indices.begin(), indices.end(), Search_traits_d(std::begin(point_cloud)));
207 
208  typename Delaunay_triangulation::Full_cell_handle hint;
209  for (auto index : indices) {
210  typename Delaunay_triangulation::Vertex_handle pos = triangulation_->insert(point_cloud[index], hint);
211  // Save index value as data to retrieve it after insertion
212  pos->data() = index;
213  hint = pos->full_cell();
214  }
215  // --------------------------------------------------------------------------------------------
216  // double map to retrieve simplex tree vertex handles from CGAL vertex iterator and vice versa
217  // Loop on triangulation vertices list
218  for (CGAL_vertex_iterator vit = triangulation_->vertices_begin(); vit != triangulation_->vertices_end(); ++vit) {
219  if (!triangulation_->is_infinite(*vit)) {
220 #ifdef DEBUG_TRACES
221  std::cout << "Vertex insertion - " << vit->data() << " -> " << vit->point() << std::endl;
222 #endif // DEBUG_TRACES
223  vertex_handle_to_iterator_.emplace(vit->data(), vit);
224  }
225  }
226  // --------------------------------------------------------------------------------------------
227  }
228  }
229 
230  public:
231  template <typename SimplicialComplexForAlpha>
232  bool create_complex(SimplicialComplexForAlpha& complex) {
234  return create_complex(complex, std::numeric_limits<Filtration_value>::infinity());
235  }
236 
252  template <typename SimplicialComplexForAlpha, typename Filtration_value>
253  bool create_complex(SimplicialComplexForAlpha& complex, Filtration_value max_alpha_square) {
254  // From SimplicialComplexForAlpha type required to insert into a simplicial complex (with or without subfaces).
255  typedef typename SimplicialComplexForAlpha::Vertex_handle Vertex_handle;
256  typedef typename SimplicialComplexForAlpha::Simplex_handle Simplex_handle;
257  typedef std::vector<Vertex_handle> Vector_vertex;
258 
259  if (triangulation_ == nullptr) {
260  std::cerr << "Alpha_complex cannot create_complex from a NULL triangulation\n";
261  return false; // ----- >>
262  }
263  if (triangulation_->maximal_dimension() < 1) {
264  std::cerr << "Alpha_complex cannot create_complex from a zero-dimension triangulation\n";
265  return false; // ----- >>
266  }
267  if (complex.num_vertices() > 0) {
268  std::cerr << "Alpha_complex create_complex - complex is not empty\n";
269  return false; // ----- >>
270  }
271 
272  // --------------------------------------------------------------------------------------------
273  // Simplex_tree construction from loop on triangulation finite full cells list
274  if (triangulation_->number_of_vertices() > 0) {
275  for (auto cit = triangulation_->finite_full_cells_begin(); cit != triangulation_->finite_full_cells_end(); ++cit) {
276  Vector_vertex vertexVector;
277 #ifdef DEBUG_TRACES
278  std::cout << "Simplex_tree insertion ";
279 #endif // DEBUG_TRACES
280  for (auto vit = cit->vertices_begin(); vit != cit->vertices_end(); ++vit) {
281  if (*vit != nullptr) {
282 #ifdef DEBUG_TRACES
283  std::cout << " " << (*vit)->data();
284 #endif // DEBUG_TRACES
285  // Vector of vertex construction for simplex_tree structure
286  vertexVector.push_back((*vit)->data());
287  }
288  }
289 #ifdef DEBUG_TRACES
290  std::cout << std::endl;
291 #endif // DEBUG_TRACES
292  // Insert each simplex and its subfaces in the simplex tree - filtration is NaN
293  complex.insert_simplex_and_subfaces(vertexVector, std::numeric_limits<double>::quiet_NaN());
294  }
295  }
296  // --------------------------------------------------------------------------------------------
297 
298  // --------------------------------------------------------------------------------------------
299  // Will be re-used many times
300  Vector_of_CGAL_points pointVector;
301  // ### For i : d -> 0
302  for (int decr_dim = triangulation_->maximal_dimension(); decr_dim >= 0; decr_dim--) {
303  // ### Foreach Sigma of dim i
304  for (Simplex_handle f_simplex : complex.skeleton_simplex_range(decr_dim)) {
305  int f_simplex_dim = complex.dimension(f_simplex);
306  if (decr_dim == f_simplex_dim) {
307  pointVector.clear();
308 #ifdef DEBUG_TRACES
309  std::cout << "Sigma of dim " << decr_dim << " is";
310 #endif // DEBUG_TRACES
311  for (auto vertex : complex.simplex_vertex_range(f_simplex)) {
312  pointVector.push_back(get_point(vertex));
313 #ifdef DEBUG_TRACES
314  std::cout << " " << vertex;
315 #endif // DEBUG_TRACES
316  }
317 #ifdef DEBUG_TRACES
318  std::cout << std::endl;
319 #endif // DEBUG_TRACES
320  // ### If filt(Sigma) is NaN : filt(Sigma) = alpha(Sigma)
321  if (std::isnan(complex.filtration(f_simplex))) {
322  Filtration_value alpha_complex_filtration = 0.0;
323  // No need to compute squared_radius on a single point - alpha is 0.0
324  if (f_simplex_dim > 0) {
325  // squared_radius function initialization
326  Squared_Radius squared_radius = kernel_.compute_squared_radius_d_object();
327  CGAL::NT_converter<typename Geom_traits::FT, Filtration_value> cv;
328 
329  alpha_complex_filtration = cv(squared_radius(pointVector.begin(), pointVector.end()));
330  }
331  complex.assign_filtration(f_simplex, alpha_complex_filtration);
332 #ifdef DEBUG_TRACES
333  std::cout << "filt(Sigma) is NaN : filt(Sigma) =" << complex.filtration(f_simplex) << std::endl;
334 #endif // DEBUG_TRACES
335  }
336  propagate_alpha_filtration(complex, f_simplex, decr_dim);
337  }
338  }
339  }
340  // --------------------------------------------------------------------------------------------
341 
342  // --------------------------------------------------------------------------------------------
343  // As Alpha value is an approximation, we have to make filtration non decreasing while increasing the dimension
345  // Remove all simplices that have a filtration value greater than max_alpha_square
346  complex.prune_above_filtration(max_alpha_square);
347  // --------------------------------------------------------------------------------------------
348  return true;
349  }
350 
351  private:
352  template <typename SimplicialComplexForAlpha, typename Simplex_handle>
353  void propagate_alpha_filtration(SimplicialComplexForAlpha& complex, Simplex_handle f_simplex, int decr_dim) {
354  // From SimplicialComplexForAlpha type required to assign filtration values.
356 #ifdef DEBUG_TRACES
357  typedef typename SimplicialComplexForAlpha::Vertex_handle Vertex_handle;
358 #endif // DEBUG_TRACES
359 
360  // ### Foreach Tau face of Sigma
361  for (auto f_boundary : complex.boundary_simplex_range(f_simplex)) {
362 #ifdef DEBUG_TRACES
363  std::cout << " | --------------------------------------------------\n";
364  std::cout << " | Tau ";
365  for (auto vertex : complex.simplex_vertex_range(f_boundary)) {
366  std::cout << vertex << " ";
367  }
368  std::cout << "is a face of Sigma\n";
369  std::cout << " | isnan(complex.filtration(Tau)=" << std::isnan(complex.filtration(f_boundary)) << std::endl;
370 #endif // DEBUG_TRACES
371  // ### If filt(Tau) is not NaN
372  if (!std::isnan(complex.filtration(f_boundary))) {
373  // ### filt(Tau) = fmin(filt(Tau), filt(Sigma))
374  Filtration_value alpha_complex_filtration = fmin(complex.filtration(f_boundary),
375  complex.filtration(f_simplex));
376  complex.assign_filtration(f_boundary, alpha_complex_filtration);
377 #ifdef DEBUG_TRACES
378  std::cout << " | filt(Tau) = fmin(filt(Tau), filt(Sigma)) = " << complex.filtration(f_boundary) << std::endl;
379 #endif // DEBUG_TRACES
380  // ### Else
381  } else {
382  // No need to compute is_gabriel for dimension <= 2
383  // i.e. : Sigma = (3,1) => Tau = 1
384  if (decr_dim > 1) {
385  // insert the Tau points in a vector for is_gabriel function
386  Vector_of_CGAL_points pointVector;
387 #ifdef DEBUG_TRACES
388  Vertex_handle vertexForGabriel = Vertex_handle();
389 #endif // DEBUG_TRACES
390  for (auto vertex : complex.simplex_vertex_range(f_boundary)) {
391  pointVector.push_back(get_point(vertex));
392  }
393  // Retrieve the Sigma point that is not part of Tau - parameter for is_gabriel function
394  Point_d point_for_gabriel;
395  for (auto vertex : complex.simplex_vertex_range(f_simplex)) {
396  point_for_gabriel = get_point(vertex);
397  if (std::find(pointVector.begin(), pointVector.end(), point_for_gabriel) == pointVector.end()) {
398 #ifdef DEBUG_TRACES
399  // vertex is not found in Tau
400  vertexForGabriel = vertex;
401 #endif // DEBUG_TRACES
402  // No need to continue loop
403  break;
404  }
405  }
406  // is_gabriel function initialization
407  Is_Gabriel is_gabriel = kernel_.side_of_bounded_sphere_d_object();
408  bool is_gab = is_gabriel(pointVector.begin(), pointVector.end(), point_for_gabriel)
409  != CGAL::ON_BOUNDED_SIDE;
410 #ifdef DEBUG_TRACES
411  std::cout << " | Tau is_gabriel(Sigma)=" << is_gab << " - vertexForGabriel=" << vertexForGabriel << std::endl;
412 #endif // DEBUG_TRACES
413  // ### If Tau is not Gabriel of Sigma
414  if (false == is_gab) {
415  // ### filt(Tau) = filt(Sigma)
416  Filtration_value alpha_complex_filtration = complex.filtration(f_simplex);
417  complex.assign_filtration(f_boundary, alpha_complex_filtration);
418 #ifdef DEBUG_TRACES
419  std::cout << " | filt(Tau) = filt(Sigma) = " << complex.filtration(f_boundary) << std::endl;
420 #endif // DEBUG_TRACES
421  }
422  }
423  }
424  }
425  }
426 };
427 
428 } // namespace alpha_complex
429 
430 namespace alphacomplex = alpha_complex;
431 
432 } // namespace Gudhi
433 
434 #endif // ALPHA_COMPLEX_H_
OFF file reader implementation in order to read points from an OFF file.
Definition: Points_off_io.h:134
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...
Simplex_vertex_range simplex_vertex_range(Simplex_handle const &simplex)
Returns a range over vertices of a given simplex.
Definition: SimplicialComplexForAlpha.h:26
const std::vector< Point_d > & get_point_cloud() const
Point cloud getter.
Definition: Points_off_io.h:170
~Alpha_complex()
Alpha_complex destructor deletes the Delaunay triangulation.
Definition: Alpha_complex.h:155
Kernel::Point_d Point_d
A point in Euclidean space.
Definition: Alpha_complex.h:86
Alpha_complex(const InputPointRange &points)
Alpha_complex constructor from a list of points.
Definition: Alpha_complex.h:148
Alpha_complex(const std::string &off_file_name)
Alpha_complex constructor from an OFF file name.
Definition: Alpha_complex.h:127
unspecified Vertex_handle
Definition: SimplicialComplexForAlpha.h:37
Value type for a filtration function on a cell complex.
Definition: FiltrationValue.h:32
bool is_valid() const
Returns if the OFF file read operation was successful or not.
Definition: Points_off_io.h:162
The concept SimplicialComplexForAlpha describes the requirements for a type to implement a simplicial...
Definition: SimplicialComplexForAlpha.h:33
void prune_above_filtration(Filtration_value filtration)
Kernel Geom_traits
Geometric traits class that provides the geometric types and predicates needed by Delaunay triangulat...
Definition: Alpha_complex.h:89
Alpha complex data structure.
Definition: Alpha_complex.h:76
CGAL::Delaunay_triangulation< Kernel, TDS > Delaunay_triangulation
A Delaunay triangulation of a set of points in .
Definition: Alpha_complex.h:83
bool create_complex(SimplicialComplexForAlpha &complex, Filtration_value max_alpha_square)
Inserts all Delaunay triangulation into the simplicial complex. It also computes the filtration value...
Definition: Alpha_complex.h:253
int assign_filtration(Simplex_handle simplex, Filtration_value filtration)
unspecified Simplex_handle
Definition: SimplicialComplexForAlpha.h:35
Boundary_simplex_range boundary_simplex_range(Simplex_handle const &simplex)
Returns a range over boundaries of a given simplex.
std::size_t number_of_vertices() const
number_of_vertices returns the number of vertices (same as the number of points). ...
Definition: Alpha_complex.h:179
unspecified Filtration_value
Definition: SimplicialComplexForAlpha.h:39
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:171
GUDHI  Version 2.2.0  - C++ library for Topological Data Analysis (TDA) and Higher Dimensional Geometry Understanding.  - Copyright : GPL v3 Generated on Thu Jun 14 2018 15:00:54 for GUDHI by Doxygen 1.8.13