Alpha_complex_3d.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) 2018 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_3D_H_
13 #define ALPHA_COMPLEX_3D_H_
14 
15 #include <boost/version.hpp>
16 #include <boost/variant.hpp>
17 
18 #include <gudhi/Debug_utils.h>
19 #include <gudhi/Alpha_complex_options.h>
20 
21 #include <CGAL/Exact_predicates_inexact_constructions_kernel.h>
22 #include <CGAL/Exact_predicates_exact_constructions_kernel.h>
23 #include <CGAL/Delaunay_triangulation_3.h>
24 #include <CGAL/Periodic_3_Delaunay_triangulation_traits_3.h>
25 #include <CGAL/Periodic_3_Delaunay_triangulation_3.h>
26 #include <CGAL/Periodic_3_regular_triangulation_traits_3.h>
27 #include <CGAL/Periodic_3_regular_triangulation_3.h>
28 #include <CGAL/Regular_triangulation_3.h>
29 #include <CGAL/Alpha_shape_3.h>
30 #include <CGAL/Alpha_shape_cell_base_3.h>
31 #include <CGAL/Alpha_shape_vertex_base_3.h>
32 
33 #include <CGAL/Object.h>
34 #include <CGAL/tuple.h>
35 #include <CGAL/iterator.h>
36 #include <CGAL/version.h> // for CGAL_VERSION_NR
37 
38 #include <Eigen/src/Core/util/Macros.h> // for EIGEN_VERSION_AT_LEAST
39 
40 #include <boost/container/static_vector.hpp>
41 
42 #include <iostream>
43 #include <vector>
44 #include <unordered_map>
45 #include <stdexcept>
46 #include <cstddef>
47 #include <memory> // for std::unique_ptr
48 #include <type_traits> // for std::conditional and std::enable_if
49 #include <limits> // for numeric_limits<>
50 
51 // Make compilation fail - required for external projects - https://github.com/GUDHI/gudhi-devel/issues/10
52 #if CGAL_VERSION_NR < 1041101000
53 # error Alpha_complex_3d is only available for CGAL >= 4.11
54 #endif
55 
56 #if !EIGEN_VERSION_AT_LEAST(3,1,0)
57 # error Alpha_complex_3d is only available for Eigen3 >= 3.1.0 installed with CGAL
58 #endif
59 
60 namespace Gudhi {
61 
62 namespace alpha_complex {
63 
64 #ifdef GUDHI_CAN_USE_CXX11_THREAD_LOCAL
65 thread_local
66 #endif // GUDHI_CAN_USE_CXX11_THREAD_LOCAL
67  double RELATIVE_PRECISION_OF_TO_DOUBLE = 0.00001;
68 
69 // Value_from_iterator returns the filtration value from an iterator on alpha shapes values
70 //
71 // FAST SAFE EXACT
72 // CGAL::to_double(*iterator) CGAL::to_double(*iterator) CGAL::to_double(iterator->exact())
73 
74 template <complexity Complexity>
75 struct Value_from_iterator {
76  template <typename Iterator>
77  static double perform(Iterator it) {
78  // Default behaviour
79  return CGAL::to_double(*it);
80  }
81 };
82 
83 template <>
84 struct Value_from_iterator<complexity::EXACT> {
85  template <typename Iterator>
86  static double perform(Iterator it) {
87  return CGAL::to_double(it->exact());
88  }
89 };
90 
126 template <complexity Complexity = complexity::SAFE, bool Weighted = false, bool Periodic = false>
128  // Epick = Exact_predicates_inexact_constructions_kernel
129  // Epeck = Exact_predicates_exact_constructions_kernel
130  // Exact_alpha_comparison_tag = exact version of CGAL Alpha_shape_3 and of its objects (Alpha_shape_vertex_base_3 and
131  // Alpha_shape_cell_base_3). Not available if weighted or periodic.
132  // Can be CGAL::Tag_false or CGAL::Tag_true. Default is False.
133  // cf. https://doc.cgal.org/latest/Alpha_shapes_3/classCGAL_1_1Alpha__shape__3.html
134  //
135  // We could use Epick + CGAL::Tag_true for not weighted nor periodic, but during benchmark, we found a bug
136  // https://github.com/CGAL/cgal/issues/3460
137  // This is the reason we only use Epick + CGAL::Tag_false, or Epeck
138  //
139  // FAST SAFE EXACT
140  // Epick + CGAL::Tag_false Epeck Epeck
141  using Predicates = typename std::conditional<(Complexity == complexity::FAST),
142  CGAL::Exact_predicates_inexact_constructions_kernel,
143  CGAL::Exact_predicates_exact_constructions_kernel>::type;
144 
145  // The other way to do a conditional type. Here there are 3 possibilities
146  template <typename Predicates, bool Weighted_version, bool Periodic_version>
147  struct Kernel_3 {};
148 
149  template <typename Predicates, bool Is_periodic>
150  struct Kernel_3<Predicates, Is_periodic, false> {
151  using Kernel = Predicates;
152  };
153 
154  template <typename Predicates>
155  struct Kernel_3<Predicates, false, true> {
156  using Kernel = CGAL::Periodic_3_Delaunay_triangulation_traits_3<Predicates>;
157  };
158  template <typename Predicates>
159  struct Kernel_3<Predicates, true, true> {
160  using Kernel = CGAL::Periodic_3_regular_triangulation_traits_3<Predicates>;
161  };
162 
163  using Kernel = typename Kernel_3<Predicates, Weighted, Periodic>::Kernel;
164 
165  using TdsVb = typename std::conditional<Periodic, CGAL::Periodic_3_triangulation_ds_vertex_base_3<>,
166  CGAL::Triangulation_ds_vertex_base_3<>>::type;
167 
168  using Tvb = typename std::conditional<Weighted, CGAL::Regular_triangulation_vertex_base_3<Kernel, TdsVb>,
169  CGAL::Triangulation_vertex_base_3<Kernel, TdsVb>>::type;
170 
171  using Vb = CGAL::Alpha_shape_vertex_base_3<Kernel, Tvb>;
172 
173  using TdsCb = typename std::conditional<Periodic, CGAL::Periodic_3_triangulation_ds_cell_base_3<>,
174  CGAL::Triangulation_ds_cell_base_3<>>::type;
175 
176  using Tcb = typename std::conditional<Weighted, CGAL::Regular_triangulation_cell_base_3<Kernel, TdsCb>,
177  CGAL::Triangulation_cell_base_3<Kernel, TdsCb>>::type;
178 
179  using Cb = CGAL::Alpha_shape_cell_base_3<Kernel, Tcb>;
180  using Tds = CGAL::Triangulation_data_structure_3<Vb, Cb>;
181 
182  // The other way to do a conditional type. Here there 4 possibilities, cannot use std::conditional
183  template <typename Kernel, typename Tds, bool Weighted_version, bool Periodic_version>
184  struct Triangulation_3 {};
185 
186  template <typename Kernel, typename Tds>
187  struct Triangulation_3<Kernel, Tds, false, false> {
188  using Dt = CGAL::Delaunay_triangulation_3<Kernel, Tds>;
189  using Weighted_point_3 = void;
190  };
191  template <typename Kernel, typename Tds>
192  struct Triangulation_3<Kernel, Tds, true, false> {
193  using Dt = CGAL::Regular_triangulation_3<Kernel, Tds>;
194  using Weighted_point_3 = typename Dt::Weighted_point;
195  };
196  template <typename Kernel, typename Tds>
197  struct Triangulation_3<Kernel, Tds, false, true> {
198  using Dt = CGAL::Periodic_3_Delaunay_triangulation_3<Kernel, Tds>;
199  using Weighted_point_3 = void;
200  };
201  template <typename Kernel, typename Tds>
202  struct Triangulation_3<Kernel, Tds, true, true> {
203  using Dt = CGAL::Periodic_3_regular_triangulation_3<Kernel, Tds>;
204  using Weighted_point_3 = typename Dt::Weighted_point;
205  };
206 
214  using Dt = typename Triangulation_3<Kernel, Tds, Weighted, Periodic>::Dt;
215 
216  public:
222  using Alpha_shape_3 = CGAL::Alpha_shape_3<Dt>;
223 
226  using FT = typename Alpha_shape_3::FT;
227 
236  using Point_3 = typename Kernel::Point_3;
237 
250  using Weighted_point_3 = typename Triangulation_3<Kernel, Tds, Weighted, Periodic>::Weighted_point_3;
251 
252  private:
253  using Dispatch =
254  CGAL::Dispatch_output_iterator<CGAL::cpp11::tuple<CGAL::Object, FT>,
255  CGAL::cpp11::tuple<std::back_insert_iterator<std::vector<CGAL::Object>>,
256  std::back_insert_iterator<std::vector<FT>>>>;
257 
258  using Cell_handle = typename Alpha_shape_3::Cell_handle;
259  using Facet = typename Alpha_shape_3::Facet;
260  using Edge = typename Alpha_shape_3::Edge;
261  using Alpha_vertex_handle = typename Alpha_shape_3::Vertex_handle;
262  using Vertex_list = boost::container::static_vector<Alpha_vertex_handle, 4>;
263 
264  public:
275  template <typename InputPointRange>
276  Alpha_complex_3d(const InputPointRange& points) {
277  static_assert(!Periodic, "This constructor is not available for periodic versions of Alpha_complex_3d");
278 
279  alpha_shape_3_ptr_ = std::unique_ptr<Alpha_shape_3>(
280  new Alpha_shape_3(std::begin(points), std::end(points), 0, Alpha_shape_3::GENERAL));
281  }
282 
297  template <typename InputPointRange, typename WeightRange>
298  Alpha_complex_3d(const InputPointRange& points, WeightRange weights) {
299  static_assert(Weighted, "This constructor is not available for non-weighted versions of Alpha_complex_3d");
300  static_assert(!Periodic, "This constructor is not available for periodic versions of Alpha_complex_3d");
301  GUDHI_CHECK((weights.size() == points.size()),
302  std::invalid_argument("Points number in range different from weights range number"));
303 
304  std::vector<Weighted_point_3> weighted_points_3;
305 
306  std::size_t index = 0;
307  weighted_points_3.reserve(points.size());
308  while ((index < weights.size()) && (index < points.size())) {
309  weighted_points_3.push_back(Weighted_point_3(points[index], weights[index]));
310  index++;
311  }
312 
313  alpha_shape_3_ptr_ = std::unique_ptr<Alpha_shape_3>(
314  new Alpha_shape_3(std::begin(weighted_points_3), std::end(weighted_points_3), 0, Alpha_shape_3::GENERAL));
315  }
316 
338  template <typename InputPointRange>
339  Alpha_complex_3d(const InputPointRange& points, FT x_min, FT y_min,
340  FT z_min, FT x_max, FT y_max, FT z_max) {
341  static_assert(Periodic, "This constructor is not available for non-periodic versions of Alpha_complex_3d");
342  // Checking if the cuboid is the same in x,y and z direction. If not, CGAL will not process it.
343  GUDHI_CHECK(
344  (x_max - x_min == y_max - y_min) && (x_max - x_min == z_max - z_min) && (z_max - z_min == y_max - y_min),
345  std::invalid_argument("The size of the cuboid in every directions is not the same."));
346 
347  // Define the periodic cube
348  Dt pdt(typename Kernel::Iso_cuboid_3(x_min, y_min, z_min, x_max, y_max, z_max));
349  // Heuristic for inserting large point sets (if pts is reasonably large)
350  pdt.insert(std::begin(points), std::end(points), true);
351  // As pdt won't be modified anymore switch to 1-sheeted cover if possible
352  if (!pdt.is_triangulation_in_1_sheet()) {
353  throw std::invalid_argument("Unable to construct a triangulation within a single periodic domain.");
354  }
355  pdt.convert_to_1_sheeted_covering();
356 
357  // alpha shape construction from points. CGAL has a strange behavior in REGULARIZED mode. This is the default mode
358  // Maybe need to set it to GENERAL mode
359  alpha_shape_3_ptr_ = std::unique_ptr<Alpha_shape_3>(new Alpha_shape_3(pdt, 0, Alpha_shape_3::GENERAL));
360  }
361 
386  template <typename InputPointRange, typename WeightRange>
387  Alpha_complex_3d(const InputPointRange& points, WeightRange weights, FT x_min, FT y_min,
388  FT z_min, FT x_max, FT y_max, FT z_max) {
389  static_assert(Weighted, "This constructor is not available for non-weighted versions of Alpha_complex_3d");
390  static_assert(Periodic, "This constructor is not available for non-periodic versions of Alpha_complex_3d");
391  GUDHI_CHECK((weights.size() == points.size()),
392  std::invalid_argument("Points number in range different from weights range number"));
393  // Checking if the cuboid is the same in x,y and z direction. If not, CGAL will not process it.
394  GUDHI_CHECK(
395  (x_max - x_min == y_max - y_min) && (x_max - x_min == z_max - z_min) && (z_max - z_min == y_max - y_min),
396  std::invalid_argument("The size of the cuboid in every directions is not the same."));
397 
398  std::vector<Weighted_point_3> weighted_points_3;
399 
400  std::size_t index = 0;
401  weighted_points_3.reserve(points.size());
402 
403 #ifdef GUDHI_DEBUG
404  // Defined in GUDHI_DEBUG to avoid unused variable warning for GUDHI_CHECK
405  FT maximal_possible_weight = 0.015625 * (x_max - x_min) * (x_max - x_min);
406 #endif
407 
408  while ((index < weights.size()) && (index < points.size())) {
409  GUDHI_CHECK((weights[index] < maximal_possible_weight) && (weights[index] >= 0),
410  std::invalid_argument("Invalid weight at index " + std::to_string(index + 1) +
411  ". Must be positive and less than maximal possible weight = 1/64*cuboid length "
412  "squared, which is not an acceptable input."));
413  weighted_points_3.push_back(Weighted_point_3(points[index], weights[index]));
414  index++;
415  }
416 
417  // Define the periodic cube
418  Dt pdt(typename Kernel::Iso_cuboid_3(x_min, y_min, z_min, x_max, y_max, z_max));
419  // Heuristic for inserting large point sets (if pts is reasonably large)
420  pdt.insert(std::begin(weighted_points_3), std::end(weighted_points_3), true);
421  // As pdt won't be modified anymore switch to 1-sheeted cover if possible
422  if (!pdt.is_triangulation_in_1_sheet()) {
423  throw std::invalid_argument("Unable to construct a triangulation within a single periodic domain.");
424  }
425  pdt.convert_to_1_sheeted_covering();
426 
427  // alpha shape construction from points. CGAL has a strange behavior in REGULARIZED mode. This is the default mode
428  // Maybe need to set it to GENERAL mode
429  alpha_shape_3_ptr_ = std::unique_ptr<Alpha_shape_3>(new Alpha_shape_3(pdt, 0, Alpha_shape_3::GENERAL));
430  }
431 
446  template <typename SimplicialComplexForAlpha3d,
448  bool create_complex(SimplicialComplexForAlpha3d& complex,
449  Filtration_value max_alpha_square = std::numeric_limits<Filtration_value>::infinity()) {
450  if (complex.num_vertices() > 0) {
451  std::cerr << "Alpha_complex_3d create_complex - complex is not empty\n";
452  return false; // ----- >>
453  }
454 
455  // using Filtration_value = typename SimplicialComplexForAlpha3d::Filtration_value;
456  using Complex_vertex_handle = typename SimplicialComplexForAlpha3d::Vertex_handle;
457  using Alpha_shape_simplex_tree_map = std::unordered_map<Alpha_vertex_handle, Complex_vertex_handle>;
458  using Simplex_tree_vector_vertex = std::vector<Complex_vertex_handle>;
459 
460 #ifdef DEBUG_TRACES
461  std::size_t count_vertices = 0;
462  std::size_t count_edges = 0;
463  std::size_t count_facets = 0;
464  std::size_t count_cells = 0;
465 #endif // DEBUG_TRACES
466  std::vector<CGAL::Object> objects;
467  std::vector<FT> alpha_values;
468 
469  Dispatch dispatcher = CGAL::dispatch_output<CGAL::Object, FT>(std::back_inserter(objects),
470  std::back_inserter(alpha_values));
471 
472  alpha_shape_3_ptr_->filtration_with_alpha_values(dispatcher);
473 #ifdef DEBUG_TRACES
474  std::cout << "filtration_with_alpha_values returns : " << objects.size() << " objects" << std::endl;
475 #endif // DEBUG_TRACES
476 
477  Alpha_shape_simplex_tree_map map_cgal_simplex_tree;
478  using Alpha_value_iterator = typename std::vector<FT>::const_iterator;
479  Alpha_value_iterator alpha_value_iterator = alpha_values.begin();
480  for (auto object_iterator : objects) {
481  Vertex_list vertex_list;
482 
483  // Retrieve Alpha shape vertex list from object
484  if (const Cell_handle* cell = CGAL::object_cast<Cell_handle>(&object_iterator)) {
485  for (auto i = 0; i < 4; i++) {
486 #ifdef DEBUG_TRACES
487  std::cout << "from cell[" << i << "]=" << (*cell)->vertex(i)->point() << std::endl;
488 #endif // DEBUG_TRACES
489  vertex_list.push_back((*cell)->vertex(i));
490  }
491 #ifdef DEBUG_TRACES
492  count_cells++;
493 #endif // DEBUG_TRACES
494  } else if (const Facet* facet = CGAL::object_cast<Facet>(&object_iterator)) {
495  for (auto i = 0; i < 4; i++) {
496  if ((*facet).second != i) {
497 #ifdef DEBUG_TRACES
498  std::cout << "from facet=[" << i << "]" << (*facet).first->vertex(i)->point() << std::endl;
499 #endif // DEBUG_TRACES
500  vertex_list.push_back((*facet).first->vertex(i));
501  }
502  }
503 #ifdef DEBUG_TRACES
504  count_facets++;
505 #endif // DEBUG_TRACES
506  } else if (const Edge* edge = CGAL::object_cast<Edge>(&object_iterator)) {
507  for (auto i : {(*edge).second, (*edge).third}) {
508 #ifdef DEBUG_TRACES
509  std::cout << "from edge[" << i << "]=" << (*edge).first->vertex(i)->point() << std::endl;
510 #endif // DEBUG_TRACES
511  vertex_list.push_back((*edge).first->vertex(i));
512  }
513 #ifdef DEBUG_TRACES
514  count_edges++;
515 #endif // DEBUG_TRACES
516  } else if (const Alpha_vertex_handle* vertex = CGAL::object_cast<Alpha_vertex_handle>(&object_iterator)) {
517 #ifdef DEBUG_TRACES
518  count_vertices++;
519  std::cout << "from vertex=" << (*vertex)->point() << std::endl;
520 #endif // DEBUG_TRACES
521  vertex_list.push_back((*vertex));
522  }
523  // Construction of the vector of simplex_tree vertex from list of alpha_shapes vertex
524  Simplex_tree_vector_vertex the_simplex;
525  for (auto the_alpha_shape_vertex : vertex_list) {
526  auto the_map_iterator = map_cgal_simplex_tree.find(the_alpha_shape_vertex);
527  if (the_map_iterator == map_cgal_simplex_tree.end()) {
528  // alpha shape not found
529  Complex_vertex_handle vertex = map_cgal_simplex_tree.size();
530 #ifdef DEBUG_TRACES
531  std::cout << "vertex [" << the_alpha_shape_vertex->point() << "] not found - insert " << vertex << std::endl;
532 #endif // DEBUG_TRACES
533  the_simplex.push_back(vertex);
534  map_cgal_simplex_tree.emplace(the_alpha_shape_vertex, vertex);
535  } else {
536  // alpha shape found
537  Complex_vertex_handle vertex = the_map_iterator->second;
538 #ifdef DEBUG_TRACES
539  std::cout << "vertex [" << the_alpha_shape_vertex->point() << "] found in " << vertex << std::endl;
540 #endif // DEBUG_TRACES
541  the_simplex.push_back(vertex);
542  }
543  }
544  // Construction of the simplex_tree
545  Filtration_value filtr = Value_from_iterator<Complexity>::perform(alpha_value_iterator);
546 
547 #ifdef DEBUG_TRACES
548  std::cout << "filtration = " << filtr << std::endl;
549 #endif // DEBUG_TRACES
550  complex.insert_simplex(the_simplex, static_cast<Filtration_value>(filtr));
551  GUDHI_CHECK(alpha_value_iterator != alpha_values.end(), "CGAL provided more simplices than values");
552  ++alpha_value_iterator;
553  }
554 
555 #ifdef DEBUG_TRACES
556  std::cout << "vertices \t" << count_vertices << std::endl;
557  std::cout << "edges \t\t" << count_edges << std::endl;
558  std::cout << "facets \t\t" << count_facets << std::endl;
559  std::cout << "cells \t\t" << count_cells << std::endl;
560 #endif // DEBUG_TRACES
561  // --------------------------------------------------------------------------------------------
562  // As Alpha value is an approximation, we have to make filtration non decreasing while increasing the dimension
564  // Remove all simplices that have a filtration value greater than max_alpha_square
565  complex.prune_above_filtration(max_alpha_square);
566  // --------------------------------------------------------------------------------------------
567  return true;
568  }
569 
570  private:
571  // use of a unique_ptr on cgal Alpha_shape_3, as copy and default constructor is not available - no need to be freed
572  std::unique_ptr<Alpha_shape_3> alpha_shape_3_ptr_;
573 };
574 
575 } // namespace alpha_complex
576 
577 } // namespace Gudhi
578 
579 #endif // ALPHA_COMPLEX_3D_H_
bool create_complex(SimplicialComplexForAlpha3d &complex, Filtration_value max_alpha_square=std::numeric_limits< Filtration_value >::infinity())
Inserts all Delaunay triangulation into the simplicial complex. It also computes the filtration value...
Definition: Alpha_complex_3d.h:448
Alpha_complex_3d(const InputPointRange &points, FT x_min, FT y_min, FT z_min, FT x_max, FT y_max, FT z_max)
Alpha_complex constructor from a list of points and an iso-cuboid coordinates.
Definition: Alpha_complex_3d.h:339
Alpha_complex_3d(const InputPointRange &points, WeightRange weights, FT x_min, FT y_min, FT z_min, FT x_max, FT y_max, FT z_max)
Alpha_complex constructor from a list of points, associated weights and an iso-cuboid coordinates...
Definition: Alpha_complex_3d.h:387
Alpha complex data structure for 3d specific case.
Definition: Alpha_complex_3d.h:127
Alpha_complex_3d(const InputPointRange &points, WeightRange weights)
Alpha_complex constructor from a list of points and associated weights.
Definition: Alpha_complex_3d.h:298
Definition: SimplicialComplexForAlpha.h:14
typename Triangulation_3< Kernel, Tds, Weighted, Periodic >::Weighted_point_3 Weighted_point_3
Gives public access to the Weighted_point_3 type. A Weighted point can be constructed as follows: ...
Definition: Alpha_complex_3d.h:250
Value type for a filtration function on a cell complex.
Definition: FiltrationValue.h:20
void prune_above_filtration(Filtration_value filtration)
The concept SimplicialComplexForAlpha3d describes the requirements for a type to implement a simplici...
Definition: SimplicialComplexForAlpha3d.h:21
CGAL::Alpha_shape_3< Dt > Alpha_shape_3
The CGAL 3D Alpha Shapes type.
Definition: Alpha_complex_3d.h:222
Alpha_complex_3d(const InputPointRange &points)
Alpha_complex constructor from a list of points.
Definition: Alpha_complex_3d.h:276
unspecified Filtration_value
Definition: SimplicialComplexForAlpha3d.h:25
typename Alpha_shape_3::FT FT
The alpha values type. Must be compatible with double.
Definition: Alpha_complex_3d.h:226
complexity
Alpha complex complexity template parameter possible values.
Definition: Alpha_complex_options.h:23
typename Kernel::Point_3 Point_3
Gives public access to the Point_3 type. Here is a Point_3 constructor example:
Definition: Alpha_complex_3d.h:236
void insert_simplex(std::vector< Vertex_handle > const &vertex_range, Filtration_value filtration)
Inserts a simplex from a given simplex (represented by a vector of Vertex_handle) in the simplicial c...
unspecified Vertex_handle
Definition: SimplicialComplexForAlpha3d.h:23
GUDHI  Version 3.0.0  - C++ library for Topological Data Analysis (TDA) and Higher Dimensional Geometry Understanding.  - Copyright : MIT Generated on Tue Sep 24 2019 09:57:50 for GUDHI by Doxygen 1.8.13