Loading...
Searching...
No Matches
Alpha_complex.h
1/* This file is part of the Gudhi Library - https://gudhi.inria.fr/ - which is released under MIT.
2 * See file LICENSE or go to https://gudhi.inria.fr/licensing/ for full license details.
3 * Author(s): Vincent Rouvreau
4 *
5 * Copyright (C) 2015 Inria
6 *
7 * Modification(s):
8 * - 2019/08 Vincent Rouvreau: Fix issue #10 for CGAL and Eigen3
9 * - YYYY/MM Author: Description of the modification
10 */
11
12#ifndef ALPHA_COMPLEX_H_
13#define ALPHA_COMPLEX_H_
14
15#include <gudhi/Alpha_complex/Alpha_kernel_d.h>
16#include <gudhi/Debug_utils.h>
17// to construct Alpha_complex from a OFF file of points
18#include <gudhi/Points_off_io.h>
19
20#include <cmath> // isnan, fmax
21#include <memory> // for std::unique_ptr
22#include <cstddef> // for std::size_t
23
24#include <CGAL/Delaunay_triangulation.h>
25#include <CGAL/Regular_triangulation.h> // aka. Weighted Delaunay triangulation
26#include <CGAL/Epeck_d.h> // For EXACT or SAFE version
27#include <CGAL/Epick_d.h> // For FAST version
28#include <CGAL/Spatial_sort_traits_adapter_d.h>
29#include <CGAL/property_map.h> // for CGAL::Identity_property_map
30#include <CGAL/version.h> // for CGAL_VERSION_NR
31#include <CGAL/NT_converter.h>
32
33#include <Eigen/src/Core/util/Macros.h> // for EIGEN_VERSION_AT_LEAST
34
35#include <boost/range/size.hpp>
36#include <boost/range/combine.hpp>
37#include <boost/range/adaptor/transformed.hpp>
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#include <algorithm> // for std::sort
48
49// Make compilation fail - required for external projects - https://github.com/GUDHI/gudhi-devel/issues/10
50#if CGAL_VERSION_NR < 1041101000
51# error Alpha_complex is only available for CGAL >= 4.11
52#endif
53
54#if !EIGEN_VERSION_AT_LEAST(3,1,0)
55# error Alpha_complex is only available for Eigen3 >= 3.1.0 installed with CGAL
56#endif
57
58namespace Gudhi {
59
60namespace alpha_complex {
61
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; };
64
102template<class Kernel = CGAL::Epeck_d<CGAL::Dynamic_dimension_tag>, bool Weighted = false>
104 private:
105 // Vertex_handle internal type (required by triangulation_ and vertices_).
106 using Internal_vertex_handle = std::ptrdiff_t;
107
108 public:
110 using Geom_traits = std::conditional_t<Weighted, CGAL::Regular_triangulation_traits_adapter<Kernel>, Kernel>;
111
112 // Add an int in TDS to save point index in the structure
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> >;
116
118 using Triangulation = std::conditional_t<Weighted, CGAL::Regular_triangulation<Kernel, TDS>,
119 CGAL::Delaunay_triangulation<Kernel, TDS>>;
120
123
124 // Numeric type of coordinates in the kernel
125 using FT = typename A_kernel_d::FT;
126
130 using Sphere = typename A_kernel_d::Sphere;
131
133 using Point_d = typename Geom_traits::Point_d;
134
135 private:
136 // Vertex_iterator type from CGAL.
137 using CGAL_vertex_iterator = typename Triangulation::Vertex_iterator;
138
139 // Structure to switch from simplex tree vertex handle to CGAL vertex iterator.
140 using Vector_vertex_iterator = std::vector< CGAL_vertex_iterator >;
141
142 private:
145 Vector_vertex_iterator vertex_handle_to_iterator_;
147 std::unique_ptr<Triangulation> triangulation_;
149 A_kernel_d kernel_;
153 std::vector<Internal_vertex_handle> vertices_;
154
156 std::vector<Sphere> cache_, old_cache_;
157
158 public:
168 Alpha_complex(const std::string& off_file_name) {
169 Gudhi::Points_off_reader<Point_d> off_reader(off_file_name);
170 if (!off_reader.is_valid()) {
171 std::cerr << "Alpha_complex - Unable to read file " << off_file_name << "\n";
172 exit(-1); // ----- >>
173 }
174
175 init_from_range(off_reader.get_point_cloud());
176 }
177
188 template<typename InputPointRange >
189 Alpha_complex(const InputPointRange& points) {
190 init_from_range(points);
191 }
192
205 template <typename InputPointRange, typename WeightRange>
206 Alpha_complex(const InputPointRange& points, WeightRange weights) {
207 static_assert(Weighted, "This constructor is not available for non-weighted versions of Alpha_complex");
208 // FIXME: this test is only valid if we have a forward range
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);
214 }
215
216 // Forbid copy/move constructor/assignment operator
217 Alpha_complex(const Alpha_complex& other) = delete;
218 Alpha_complex& operator= (const Alpha_complex& other) = delete;
219 Alpha_complex (Alpha_complex&& other) = delete;
220 Alpha_complex& operator= (Alpha_complex&& other) = delete;
221
224 std::size_t num_vertices() const {
225 if (triangulation_ == nullptr)
226 return 0;
227 else
228 return triangulation_->number_of_vertices();
229 }
230
237 const Point_d& get_point(std::size_t vertex) const {
238 auto it = vertex_handle_to_iterator_.at(vertex);
239 if (it == nullptr) throw std::out_of_range("This vertex is missing, maybe hidden by a duplicate or another heavier point.");
240 return it->point();
241 }
242
243 private:
244 template<typename InputPointRange >
245 void init_from_range(const InputPointRange& points) {
246 #if CGAL_VERSION_NR < 1050000000
247 if (Is_Epeck_D<Kernel>::value)
248 std::cerr << "It is strongly advised to use a CGAL version >= 5.0 with Epeck_d Kernel for performance reasons."
249 << std::endl;
250 #endif
251
252#if CGAL_VERSION_NR < 1050101000
253 // Make compilation fail if weighted and CGAL < 5.1
254 static_assert(!Weighted, "Weighted Alpha_complex is only available for CGAL >= 5.1");
255#endif
256
257 auto first = std::begin(points);
258 auto last = std::end(points);
259
260 if (first != last) {
261 // Delaunay triangulation init with point dimension.
262 triangulation_ = std::make_unique<Triangulation>(kernel_.get_dimension(*first));
263
264 std::vector<Point_d> point_cloud(first, last);
265
266 // Creates a vector {0, 1, ..., N-1}
267 std::vector<Internal_vertex_handle> indices(boost::counting_iterator<Internal_vertex_handle>(0),
268 boost::counting_iterator<Internal_vertex_handle>(point_cloud.size()));
269
270 using Point_property_map = boost::iterator_property_map<typename std::vector<Point_d>::iterator,
271 CGAL::Identity_property_map<Internal_vertex_handle>>;
272 using Search_traits_d = CGAL::Spatial_sort_traits_adapter_d<Geom_traits, Point_property_map>;
273
274 CGAL::spatial_sort(indices.begin(), indices.end(), Search_traits_d(std::begin(point_cloud)));
275
276 typename Triangulation::Full_cell_handle hint;
277 for (auto index : indices) {
278 typename Triangulation::Vertex_handle pos = triangulation_->insert(point_cloud[index], hint);
279 if (pos != nullptr) {
280 // Save index value as data to retrieve it after insertion
281 pos->data() = index;
282 hint = pos->full_cell();
283 }
284 }
285 // --------------------------------------------------------------------------------------------
286 // structure to retrieve CGAL points from vertex handle - one vertex handle per point.
287 // Needs to be constructed before as vertex handles arrives in no particular order.
288 vertex_handle_to_iterator_.resize(point_cloud.size());
289 // List of sorted unique vertices in the triangulation. We take advantage of the existing loop to construct it
290 // Vertices list avoids quadratic complexity with the Simplex_tree. We should not fill it up with Toplex_map e.g.
291 vertices_.reserve(triangulation_->number_of_vertices());
292 // Loop on triangulation vertices list
293 for (CGAL_vertex_iterator vit = triangulation_->vertices_begin(); vit != triangulation_->vertices_end(); ++vit) {
294 if (!triangulation_->is_infinite(*vit)) {
295#ifdef DEBUG_TRACES
296 std::clog << "Vertex insertion - " << vit->data() << " -> " << vit->point() << std::endl;
297#endif // DEBUG_TRACES
298 vertex_handle_to_iterator_[vit->data()] = vit;
299 vertices_.push_back(vit->data());
300 }
301 }
302 std::sort(vertices_.begin(), vertices_.end());
303 // --------------------------------------------------------------------------------------------
304 }
305 }
306
313 const Point_d& get_point_(std::size_t vertex) const {
314 return vertex_handle_to_iterator_[vertex]->point();
315 }
316
318 template<class SimplicialComplexForAlpha>
319 auto& get_cache(SimplicialComplexForAlpha& cplx, typename SimplicialComplexForAlpha::Simplex_handle s) {
320 auto k = cplx.key(s);
321 if(k==cplx.null_key()){
322 k = cache_.size();
323 cplx.assign_key(s, k);
324 // Using a transform_range is slower, currently.
325 thread_local std::vector<Point_d> v;
326 v.clear();
327 for (auto vertex : cplx.simplex_vertex_range(s))
328 v.push_back(get_point_(vertex));
329 cache_.emplace_back(kernel_.get_sphere(v.cbegin(), v.cend()));
330 }
331 return cache_[k];
332 }
333
335 template<class SimplicialComplexForAlpha>
336 auto radius(SimplicialComplexForAlpha& cplx, typename SimplicialComplexForAlpha::Simplex_handle s) {
337 auto k = cplx.key(s);
338 if(k!=cplx.null_key())
339 return kernel_.get_squared_radius(old_cache_[k]);
340 // Using a transform_range is slower, currently.
341 thread_local std::vector<Point_d> v;
342 v.clear();
343 for (auto vertex : cplx.simplex_vertex_range(s))
344 v.push_back(get_point_(vertex));
345 return kernel_.get_squared_radius(v.cbegin(), v.cend());
346 }
347
348 public:
372 template <typename SimplicialComplexForAlpha,
375 Filtration_value max_alpha_square = std::numeric_limits<Filtration_value>::infinity(),
376 bool exact = false,
377 bool default_filtration_value = false) {
378 // Filtration_value must be capable to represent the special value "Not-A-Number"
379 static_assert(std::numeric_limits<Filtration_value>::has_quiet_NaN);
380 // To support more general types for Filtration_value
381 using std::isnan;
382
383 // From SimplicialComplexForAlpha type required to insert into a simplicial complex (with or without subfaces).
385 using Simplex_handle = typename SimplicialComplexForAlpha::Simplex_handle;
386 using Vector_vertex = std::vector<Vertex_handle>;
387
388 if (triangulation_ == nullptr) {
389 std::cerr << "Alpha_complex cannot create_complex from a NULL triangulation\n";
390 return false; // ----- >>
391 }
392 if (triangulation_->maximal_dimension() < 1) {
393 std::cerr << "Alpha_complex cannot create_complex from a zero-dimension triangulation\n";
394 return false; // ----- >>
395 }
396 if (complex.num_vertices() > 0) {
397 std::cerr << "Alpha_complex create_complex - complex is not empty\n";
398 return false; // ----- >>
399 }
400
401 // --------------------------------------------------------------------------------------------
402 // Simplex_tree construction from loop on triangulation finite full cells list
403 if (num_vertices() > 0) {
404 std::vector<Vertex_handle> one_vertex(1);
405 for (auto vertex : vertices_) {
406#ifdef DEBUG_TRACES
407 std::clog << "SimplicialComplex insertion " << vertex << std::endl;
408#endif // DEBUG_TRACES
409 one_vertex[0] = vertex;
410 complex.insert_simplex_and_subfaces(one_vertex, std::numeric_limits<Filtration_value>::quiet_NaN());
411 }
412
413 for (auto cit = triangulation_->finite_full_cells_begin();
414 cit != triangulation_->finite_full_cells_end();
415 ++cit) {
416 Vector_vertex vertexVector;
417#ifdef DEBUG_TRACES
418 std::clog << "SimplicialComplex insertion ";
419#endif // DEBUG_TRACES
420 for (auto vit = cit->vertices_begin(); vit != cit->vertices_end(); ++vit) {
421 if (*vit != nullptr) {
422#ifdef DEBUG_TRACES
423 std::clog << " " << (*vit)->data();
424#endif // DEBUG_TRACES
425 // Vector of vertex construction for simplex_tree structure
426 vertexVector.push_back((*vit)->data());
427 }
428 }
429#ifdef DEBUG_TRACES
430 std::clog << std::endl;
431#endif // DEBUG_TRACES
432 // Insert each simplex and its subfaces in the simplex tree - filtration is NaN
433 complex.insert_simplex_and_subfaces(vertexVector, std::numeric_limits<Filtration_value>::quiet_NaN());
434 }
435 }
436 // --------------------------------------------------------------------------------------------
437
438 if (!default_filtration_value) {
439 CGAL::NT_converter<FT, Filtration_value> cgal_converter;
440 // --------------------------------------------------------------------------------------------
441 // ### For i : d -> 0
442 for (int decr_dim = triangulation_->maximal_dimension(); decr_dim >= 0; decr_dim--) {
443 // ### Foreach Sigma of dim i
444 for (Simplex_handle f_simplex : complex.skeleton_simplex_range(decr_dim)) {
445 int f_simplex_dim = complex.dimension(f_simplex);
446 if (decr_dim == f_simplex_dim) {
447 // ### If filt(Sigma) is NaN : filt(Sigma) = alpha(Sigma)
448 if (isnan(complex.filtration(f_simplex))) {
449 Filtration_value alpha_complex_filtration = 0.0;
450 // No need to compute squared_radius on a non-weighted single point - alpha is 0.0
451 if (Weighted || f_simplex_dim > 0) {
452 auto const& sqrad = radius(complex, f_simplex);
453#if CGAL_VERSION_NR >= 1050000000
454 if(exact) CGAL::exact(sqrad);
455#endif
456 alpha_complex_filtration = cgal_converter(sqrad);
457 }
458 complex.assign_filtration(f_simplex, alpha_complex_filtration);
459#ifdef DEBUG_TRACES
460 std::clog << "filt(Sigma) is NaN : filt(Sigma) =" << complex.filtration(f_simplex) << std::endl;
461#endif // DEBUG_TRACES
462 }
463 // No need to propagate further, unweighted points all have value 0
464 if (decr_dim > !Weighted)
465 propagate_alpha_filtration(complex, f_simplex);
466 }
467 }
468 old_cache_ = std::move(cache_);
469 cache_.clear();
470 }
471 // --------------------------------------------------------------------------------------------
472
473 // --------------------------------------------------------------------------------------------
474 if (!exact)
475 // As Alpha value is an approximation, we have to make filtration non decreasing while increasing the dimension
476 // Only in not exact version, cf. https://github.com/GUDHI/gudhi-devel/issues/57
478 // Remove all simplices that have a filtration value greater than max_alpha_square
479 complex.prune_above_filtration(max_alpha_square);
480 // --------------------------------------------------------------------------------------------
481 }
482 return true;
483 }
484
485 private:
486 template <typename SimplicialComplexForAlpha, typename Simplex_handle>
487 void propagate_alpha_filtration(SimplicialComplexForAlpha& complex, Simplex_handle f_simplex) {
488 // From SimplicialComplexForAlpha type required to assign filtration values.
490 // To support more general types for Filtration_value
491 using std::isnan;
492
493 // ### Foreach Tau face of Sigma
494 for (auto face_opposite_vertex : complex.boundary_opposite_vertex_simplex_range(f_simplex)) {
495 auto f_boundary = face_opposite_vertex.first;
496#ifdef DEBUG_TRACES
497 std::clog << " | --------------------------------------------------\n";
498 std::clog << " | Tau ";
499 for (auto vertex : complex.simplex_vertex_range(f_boundary)) {
500 std::clog << vertex << " ";
501 }
502 std::clog << "is a face of Sigma\n";
503 std::clog << " | isnan(complex.filtration(Tau)=" << isnan(complex.filtration(f_boundary)) << std::endl;
504#endif // DEBUG_TRACES
505 // ### If filt(Tau) is not NaN
506 if (!isnan(complex.filtration(f_boundary))) {
507 // ### filt(Tau) = fmin(filt(Tau), filt(Sigma))
508 Filtration_value alpha_complex_filtration = fmin(complex.filtration(f_boundary),
509 complex.filtration(f_simplex));
510 complex.assign_filtration(f_boundary, alpha_complex_filtration);
511#ifdef DEBUG_TRACES
512 std::clog << " | filt(Tau) = fmin(filt(Tau), filt(Sigma)) = " << complex.filtration(f_boundary) << std::endl;
513#endif // DEBUG_TRACES
514 // ### Else
515 } else {
516 auto const& cache=get_cache(complex, f_boundary);
517 bool is_gab = kernel_.is_gabriel(cache, get_point_(face_opposite_vertex.second));
518#ifdef DEBUG_TRACES
519 std::clog << " | Tau is_gabriel(Sigma)=" << is_gab << " - vertexForGabriel=" << face_opposite_vertex.second << std::endl;
520#endif // DEBUG_TRACES
521 // ### If Tau is not Gabriel of Sigma
522 if (false == is_gab) {
523 // ### filt(Tau) = filt(Sigma)
524 Filtration_value alpha_complex_filtration = complex.filtration(f_simplex);
525 complex.assign_filtration(f_boundary, alpha_complex_filtration);
526#ifdef DEBUG_TRACES
527 std::clog << " | filt(Tau) = filt(Sigma) = " << complex.filtration(f_boundary) << std::endl;
528#endif // DEBUG_TRACES
529 }
530 }
531 }
532 }
533};
534
535} // namespace alpha_complex
536
537namespace alphacomplex = alpha_complex;
538
539} // namespace Gudhi
540
541#endif // ALPHA_COMPLEX_H_
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:374
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)
void prune_above_filtration(Filtration_value filtration)
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