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> // for std::size_t
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 thread_local double RELATIVE_PRECISION_OF_TO_DOUBLE = 0.00001;
65 
66 // Value_from_iterator returns the filtration value from an iterator on alpha shapes values
67 //
68 // FAST SAFE EXACT
69 // CGAL::to_double(*iterator) CGAL::to_double(*iterator) CGAL::to_double(iterator->exact())
70 
71 template <complexity Complexity>
72 struct Value_from_iterator {
73  template <typename Iterator>
74  static double perform(Iterator it) {
75  // Default behaviour
76  return CGAL::to_double(*it);
77  }
78 };
79 
80 template <>
81 struct Value_from_iterator<complexity::EXACT> {
82  template <typename Iterator>
83  static double perform(Iterator it) {
84  return CGAL::to_double(it->exact());
85  }
86 };
87 
123 template <complexity Complexity = complexity::SAFE, bool Weighted = false, bool Periodic = false>
125  // Epick = Exact_predicates_inexact_constructions_kernel
126  // Epeck = Exact_predicates_exact_constructions_kernel
127  // Exact_alpha_comparison_tag = exact version of CGAL Alpha_shape_3 and of its objects (Alpha_shape_vertex_base_3 and
128  // Alpha_shape_cell_base_3). Not available if weighted or periodic.
129  // Can be CGAL::Tag_false or CGAL::Tag_true. Default is False.
130  // cf. https://doc.cgal.org/latest/Alpha_shapes_3/classCGAL_1_1Alpha__shape__3.html
131  //
132  // We could use Epick + CGAL::Tag_true for not weighted nor periodic, but during benchmark, we found a bug
133  // https://github.com/CGAL/cgal/issues/3460
134  // This is the reason we only use Epick + CGAL::Tag_false, or Epeck
135  //
136  // FAST SAFE EXACT
137  // Epick + CGAL::Tag_false Epeck Epeck
138  using Predicates = typename std::conditional<(Complexity == complexity::FAST),
139  CGAL::Exact_predicates_inexact_constructions_kernel,
140  CGAL::Exact_predicates_exact_constructions_kernel>::type;
141 
142  // The other way to do a conditional type. Here there are 3 possibilities
143  template <typename Predicates, bool Weighted_version, bool Periodic_version>
144  struct Kernel_3 {};
145 
146  template <typename Predicates, bool Is_periodic>
147  struct Kernel_3<Predicates, Is_periodic, false> {
148  using Kernel = Predicates;
149  };
150 
151  template <typename Predicates>
152  struct Kernel_3<Predicates, false, true> {
153  using Kernel = CGAL::Periodic_3_Delaunay_triangulation_traits_3<Predicates>;
154  };
155  template <typename Predicates>
156  struct Kernel_3<Predicates, true, true> {
157  using Kernel = CGAL::Periodic_3_regular_triangulation_traits_3<Predicates>;
158  };
159 
160  using Kernel = typename Kernel_3<Predicates, Weighted, Periodic>::Kernel;
161 
162  using TdsVb = typename std::conditional<Periodic, CGAL::Periodic_3_triangulation_ds_vertex_base_3<>,
163  CGAL::Triangulation_ds_vertex_base_3<>>::type;
164 
165  using Tvb = typename std::conditional<Weighted, CGAL::Regular_triangulation_vertex_base_3<Kernel, TdsVb>,
166  CGAL::Triangulation_vertex_base_3<Kernel, TdsVb>>::type;
167 
168  using Vb = CGAL::Alpha_shape_vertex_base_3<Kernel, Tvb>;
169 
170  using TdsCb = typename std::conditional<Periodic, CGAL::Periodic_3_triangulation_ds_cell_base_3<>,
171  CGAL::Triangulation_ds_cell_base_3<>>::type;
172 
173  using Tcb = typename std::conditional<Weighted, CGAL::Regular_triangulation_cell_base_3<Kernel, TdsCb>,
174  CGAL::Triangulation_cell_base_3<Kernel, TdsCb>>::type;
175 
176  using Cb = CGAL::Alpha_shape_cell_base_3<Kernel, Tcb>;
177  using Tds = CGAL::Triangulation_data_structure_3<Vb, Cb>;
178 
179  // The other way to do a conditional type. Here there 4 possibilities, cannot use std::conditional
180  template <typename Kernel, typename Tds, bool Weighted_version, bool Periodic_version>
181  struct Triangulation_3 {};
182 
183  template <typename Kernel, typename Tds>
184  struct Triangulation_3<Kernel, Tds, false, false> {
185  using Dt = CGAL::Delaunay_triangulation_3<Kernel, Tds>;
186  using Weighted_point_3 = void;
187  };
188  template <typename Kernel, typename Tds>
189  struct Triangulation_3<Kernel, Tds, true, false> {
190  using Dt = CGAL::Regular_triangulation_3<Kernel, Tds>;
191  using Weighted_point_3 = typename Dt::Weighted_point;
192  };
193  template <typename Kernel, typename Tds>
194  struct Triangulation_3<Kernel, Tds, false, true> {
195  using Dt = CGAL::Periodic_3_Delaunay_triangulation_3<Kernel, Tds>;
196  using Weighted_point_3 = void;
197  };
198  template <typename Kernel, typename Tds>
199  struct Triangulation_3<Kernel, Tds, true, true> {
200  using Dt = CGAL::Periodic_3_regular_triangulation_3<Kernel, Tds>;
201  using Weighted_point_3 = typename Dt::Weighted_point;
202  };
203 
211  using Dt = typename Triangulation_3<Kernel, Tds, Weighted, Periodic>::Dt;
212 
213  public:
219  using Alpha_shape_3 = CGAL::Alpha_shape_3<Dt>;
220 
223  using FT = typename Alpha_shape_3::FT;
224 
234  using Bare_point_3 = typename Kernel::Point_3;
235 
247  using Weighted_point_3 = typename Triangulation_3<Kernel, Tds, Weighted, Periodic>::Weighted_point_3;
248 
252  using Point_3 = typename Alpha_shape_3::Point;
253 
254  private:
255  using Dispatch =
256  CGAL::Dispatch_output_iterator<CGAL::cpp11::tuple<CGAL::Object, FT>,
257  CGAL::cpp11::tuple<std::back_insert_iterator<std::vector<CGAL::Object>>,
258  std::back_insert_iterator<std::vector<FT>>>>;
259 
260  using Cell_handle = typename Alpha_shape_3::Cell_handle;
261  using Facet = typename Alpha_shape_3::Facet;
262  using Edge = typename Alpha_shape_3::Edge;
263  using Alpha_vertex_handle = typename Alpha_shape_3::Vertex_handle;
264  using Vertex_list = boost::container::static_vector<Alpha_vertex_handle, 4>;
265 
266  public:
276  template <typename InputPointRange>
277  Alpha_complex_3d(const InputPointRange& points) {
278  static_assert(!Periodic, "This constructor is not available for periodic versions of Alpha_complex_3d");
279 
280  alpha_shape_3_ptr_ = std::unique_ptr<Alpha_shape_3>(
281  new Alpha_shape_3(std::begin(points), std::end(points), 0, Alpha_shape_3::GENERAL));
282  }
283 
298  template <typename InputPointRange, typename WeightRange>
299  Alpha_complex_3d(const InputPointRange& points, WeightRange weights) {
300  static_assert(Weighted, "This constructor is not available for non-weighted versions of Alpha_complex_3d");
301  static_assert(!Periodic, "This constructor is not available for periodic versions of Alpha_complex_3d");
302  GUDHI_CHECK((weights.size() == points.size()),
303  std::invalid_argument("Points number in range different from weights range number"));
304 
305  std::vector<Weighted_point_3> weighted_points_3;
306 
307  std::size_t index = 0;
308  weighted_points_3.reserve(points.size());
309  while ((index < weights.size()) && (index < points.size())) {
310  weighted_points_3.push_back(Weighted_point_3(points[index], weights[index]));
311  index++;
312  }
313 
314  alpha_shape_3_ptr_ = std::unique_ptr<Alpha_shape_3>(
315  new Alpha_shape_3(std::begin(weighted_points_3), std::end(weighted_points_3), 0, Alpha_shape_3::GENERAL));
316  }
317 
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 Complex_vertex_handle = typename SimplicialComplexForAlpha3d::Vertex_handle;
456  using Simplex_tree_vector_vertex = std::vector<Complex_vertex_handle>;
457 
458 #ifdef DEBUG_TRACES
459  std::size_t count_vertices = 0;
460  std::size_t count_edges = 0;
461  std::size_t count_facets = 0;
462  std::size_t count_cells = 0;
463 #endif // DEBUG_TRACES
464  std::vector<CGAL::Object> objects;
465  std::vector<FT> alpha_values;
466 
467  Dispatch dispatcher = CGAL::dispatch_output<CGAL::Object, FT>(std::back_inserter(objects),
468  std::back_inserter(alpha_values));
469 
470  alpha_shape_3_ptr_->filtration_with_alpha_values(dispatcher);
471 #ifdef DEBUG_TRACES
472  std::clog << "filtration_with_alpha_values returns : " << objects.size() << " objects" << std::endl;
473 #endif // DEBUG_TRACES
474 
475  using Alpha_value_iterator = typename std::vector<FT>::const_iterator;
476  Alpha_value_iterator alpha_value_iterator = alpha_values.begin();
477  for (auto object_iterator : objects) {
478  Vertex_list vertex_list;
479 
480  // Retrieve Alpha shape vertex list from object
481  if (const Cell_handle* cell = CGAL::object_cast<Cell_handle>(&object_iterator)) {
482  for (auto i = 0; i < 4; i++) {
483 #ifdef DEBUG_TRACES
484  std::clog << "from cell[" << i << "] - Point coordinates (" << (*cell)->vertex(i)->point() << ")"
485  << std::endl;
486 #endif // DEBUG_TRACES
487  vertex_list.push_back((*cell)->vertex(i));
488  }
489 #ifdef DEBUG_TRACES
490  count_cells++;
491 #endif // DEBUG_TRACES
492  } else if (const Facet* facet = CGAL::object_cast<Facet>(&object_iterator)) {
493  for (auto i = 0; i < 4; i++) {
494  if ((*facet).second != i) {
495 #ifdef DEBUG_TRACES
496  std::clog << "from facet=[" << i << "] - Point coordinates (" << (*facet).first->vertex(i)->point() << ")"
497  << std::endl;
498 #endif // DEBUG_TRACES
499  vertex_list.push_back((*facet).first->vertex(i));
500  }
501  }
502 #ifdef DEBUG_TRACES
503  count_facets++;
504 #endif // DEBUG_TRACES
505  } else if (const Edge* edge = CGAL::object_cast<Edge>(&object_iterator)) {
506  for (auto i : {(*edge).second, (*edge).third}) {
507 #ifdef DEBUG_TRACES
508  std::clog << "from edge[" << i << "] - Point coordinates (" << (*edge).first->vertex(i)->point() << ")"
509  << 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::clog << "from vertex - Point coordinates (" << (*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::clog << "Point (" << the_alpha_shape_vertex->point() << ") not found - insert new vertex id " << vertex
532  << std::endl;
533 #endif // DEBUG_TRACES
534  the_simplex.push_back(vertex);
535  map_cgal_simplex_tree.emplace(the_alpha_shape_vertex, vertex);
536  } else {
537  // alpha shape found
538  Complex_vertex_handle vertex = the_map_iterator->second;
539 #ifdef DEBUG_TRACES
540  std::clog << "Point (" << the_alpha_shape_vertex->point() << ") found as vertex id " << vertex << std::endl;
541 #endif // DEBUG_TRACES
542  the_simplex.push_back(vertex);
543  }
544  }
545  // Construction of the simplex_tree
546  Filtration_value filtr = Value_from_iterator<Complexity>::perform(alpha_value_iterator);
547 
548 #ifdef DEBUG_TRACES
549  std::clog << "filtration = " << filtr << std::endl;
550 #endif // DEBUG_TRACES
551  complex.insert_simplex(the_simplex, static_cast<Filtration_value>(filtr));
552  GUDHI_CHECK(alpha_value_iterator != alpha_values.end(), "CGAL provided more simplices than values");
553  ++alpha_value_iterator;
554  }
555 
556 #ifdef DEBUG_TRACES
557  std::clog << "vertices \t" << count_vertices << std::endl;
558  std::clog << "edges \t\t" << count_edges << std::endl;
559  std::clog << "facets \t\t" << count_facets << std::endl;
560  std::clog << "cells \t\t" << count_cells << std::endl;
561 #endif // DEBUG_TRACES
562  // --------------------------------------------------------------------------------------------
563  // As Alpha value is an approximation, we have to make filtration non decreasing while increasing the dimension
565  // Remove all simplices that have a filtration value greater than max_alpha_square
566  complex.prune_above_filtration(max_alpha_square);
567  // --------------------------------------------------------------------------------------------
568  return true;
569  }
570 
577  const Point_3& get_point(std::size_t vertex) {
578  if (map_cgal_simplex_tree.size() != cgal_vertex_iterator_vector.size()) {
579  cgal_vertex_iterator_vector.resize(map_cgal_simplex_tree.size());
580  for (auto map_iterator : map_cgal_simplex_tree) {
581  cgal_vertex_iterator_vector[map_iterator.second] = map_iterator.first;
582  }
583 
584  }
585  auto cgal_vertex_iterator = cgal_vertex_iterator_vector.at(vertex);
586  return cgal_vertex_iterator->point();
587  }
588 
589  private:
590  // use of a unique_ptr on cgal Alpha_shape_3, as copy and default constructor is not available - no need to be freed
591  std::unique_ptr<Alpha_shape_3> alpha_shape_3_ptr_;
592 
593  // Map type to switch from CGAL vertex iterator to simplex tree vertex handle.
594  std::unordered_map<Alpha_vertex_handle, std::size_t> map_cgal_simplex_tree;
595  // Vector type to switch from simplex tree vertex handle to CGAL vertex iterator.
596  std::vector<Alpha_vertex_handle> cgal_vertex_iterator_vector;
597 };
598 
599 } // namespace alpha_complex
600 
601 } // namespace Gudhi
602 
603 #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
const Point_3 & get_point(std::size_t vertex)
get_point returns the point corresponding to the vertex given as parameter.
Definition: Alpha_complex_3d.h:577
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:124
typename Alpha_shape_3::Point Point_3
Alpha_complex_3d::Point_3 type is either a Alpha_complex_3d::Bare_point_3 (Weighted = false) or a Alp...
Definition: Alpha_complex_3d.h:252
Alpha_complex_3d(const InputPointRange &points, WeightRange weights)
Alpha_complex constructor from a list of points and associated weights.
Definition: Alpha_complex_3d.h:299
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:247
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:219
Alpha_complex_3d(const InputPointRange &points)
Alpha_complex constructor from a list of points.
Definition: Alpha_complex_3d.h:277
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:223
complexity
Alpha complex complexity template parameter possible values.
Definition: Alpha_complex_options.h:23
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
typename Kernel::Point_3 Bare_point_3
Gives public access to the Bare_point_3 (bare aka. unweighed) type. Here is a Bare_point_3 constructo...
Definition: Alpha_complex_3d.h:234
GUDHI  Version 3.2.0  - C++ library for Topological Data Analysis (TDA) and Higher Dimensional Geometry Understanding.  - Copyright : MIT Generated on Wed May 20 2020 07:53:46 for GUDHI by Doxygen 1.8.13