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 #include <boost/range/size.hpp>
18 #include <boost/range/combine.hpp>
19 #include <boost/range/adaptor/transformed.hpp>
20 
21 #include <gudhi/Debug_utils.h>
22 #include <gudhi/Alpha_complex_options.h>
23 
24 #include <CGAL/Exact_predicates_inexact_constructions_kernel.h>
25 #include <CGAL/Exact_predicates_exact_constructions_kernel.h>
26 #include <CGAL/Delaunay_triangulation_3.h>
27 #include <CGAL/Periodic_3_Delaunay_triangulation_traits_3.h>
28 #include <CGAL/Periodic_3_Delaunay_triangulation_3.h>
29 #include <CGAL/Periodic_3_regular_triangulation_traits_3.h>
30 #include <CGAL/Periodic_3_regular_triangulation_3.h>
31 #include <CGAL/Regular_triangulation_3.h>
32 #include <CGAL/Alpha_shape_3.h>
33 #include <CGAL/Alpha_shape_cell_base_3.h>
34 #include <CGAL/Alpha_shape_vertex_base_3.h>
35 
36 #include <CGAL/Object.h>
37 #include <CGAL/tuple.h>
38 #include <CGAL/iterator.h>
39 #include <CGAL/version.h> // for CGAL_VERSION_NR
40 
41 #include <Eigen/src/Core/util/Macros.h> // for EIGEN_VERSION_AT_LEAST
42 
43 #include <boost/container/static_vector.hpp>
44 
45 #include <iostream>
46 #include <vector>
47 #include <unordered_map>
48 #include <stdexcept>
49 #include <cstddef> // for std::size_t
50 #include <memory> // for std::unique_ptr
51 #include <type_traits> // for std::conditional and std::enable_if
52 #include <limits> // for numeric_limits<>
53 
54 // Make compilation fail - required for external projects - https://github.com/GUDHI/gudhi-devel/issues/10
55 #if CGAL_VERSION_NR < 1041101000
56 # error Alpha_complex_3d is only available for CGAL >= 4.11
57 #endif
58 
59 #if !EIGEN_VERSION_AT_LEAST(3,1,0)
60 # error Alpha_complex_3d is only available for Eigen3 >= 3.1.0 installed with CGAL
61 #endif
62 
63 namespace Gudhi {
64 
65 namespace alpha_complex {
66 
67 thread_local 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 
237  using Bare_point_3 = typename Kernel::Point_3;
238 
250  using Weighted_point_3 = typename Triangulation_3<Kernel, Tds, Weighted, Periodic>::Weighted_point_3;
251 
255  using Point_3 = typename Alpha_shape_3::Point;
256 
257  private:
258  using Dispatch =
259  CGAL::Dispatch_output_iterator<CGAL::cpp11::tuple<CGAL::Object, FT>,
260  CGAL::cpp11::tuple<std::back_insert_iterator<std::vector<CGAL::Object>>,
261  std::back_insert_iterator<std::vector<FT>>>>;
262 
263  using Cell_handle = typename Alpha_shape_3::Cell_handle;
264  using Facet = typename Alpha_shape_3::Facet;
265  using Edge = typename Alpha_shape_3::Edge;
266  using Alpha_vertex_handle = typename Alpha_shape_3::Vertex_handle;
267  using Vertex_list = boost::container::static_vector<Alpha_vertex_handle, 4>;
268 
269  public:
279  template <typename InputPointRange>
280  Alpha_complex_3d(const InputPointRange& points) {
281  static_assert(!Periodic, "This constructor is not available for periodic versions of Alpha_complex_3d");
282 
283  alpha_shape_3_ptr_ = std::make_unique<Alpha_shape_3>(
284  std::begin(points), std::end(points), 0, Alpha_shape_3::GENERAL);
285  }
286 
301  template <typename InputPointRange, typename WeightRange>
302  Alpha_complex_3d(const InputPointRange& points, WeightRange weights) {
303  static_assert(Weighted, "This constructor is not available for non-weighted versions of Alpha_complex_3d");
304  static_assert(!Periodic, "This constructor is not available for periodic versions of Alpha_complex_3d");
305  // FIXME: this test is only valid if we have a forward range
306  GUDHI_CHECK(boost::size(weights) == boost::size(points),
307  std::invalid_argument("Points number in range different from weights range number"));
308 
309  auto weighted_points_3 = boost::range::combine(points, weights)
310  | boost::adaptors::transformed([](auto const&t){return Weighted_point_3(boost::get<0>(t), boost::get<1>(t));});
311 
312  alpha_shape_3_ptr_ = std::make_unique<Alpha_shape_3>(
313  std::begin(weighted_points_3), std::end(weighted_points_3), 0, Alpha_shape_3::GENERAL);
314  }
315 
336  template <typename InputPointRange>
337  Alpha_complex_3d(const InputPointRange& points, FT x_min, FT y_min,
338  FT z_min, FT x_max, FT y_max, FT z_max) {
339  static_assert(Periodic, "This constructor is not available for non-periodic versions of Alpha_complex_3d");
340  // Checking if the cuboid is the same in x,y and z direction. If not, CGAL will not process it.
341  GUDHI_CHECK(
342  (x_max - x_min == y_max - y_min) && (x_max - x_min == z_max - z_min) && (z_max - z_min == y_max - y_min),
343  std::invalid_argument("The size of the cuboid in every directions is not the same."));
344 
345  // Define the periodic cube
346  Dt pdt(typename Kernel::Iso_cuboid_3(x_min, y_min, z_min, x_max, y_max, z_max));
347  // Heuristic for inserting large point sets (if pts is reasonably large)
348  pdt.insert(std::begin(points), std::end(points), true);
349  // As pdt won't be modified anymore switch to 1-sheeted cover if possible
350  if (!pdt.is_triangulation_in_1_sheet()) {
351  throw std::invalid_argument("Unable to construct a triangulation within a single periodic domain.");
352  }
353  pdt.convert_to_1_sheeted_covering();
354 
355  // alpha shape construction from points. CGAL has a strange behavior in REGULARIZED mode. This is the default mode
356  // Maybe need to set it to GENERAL mode
357  alpha_shape_3_ptr_ = std::make_unique<Alpha_shape_3>(pdt, 0, Alpha_shape_3::GENERAL);
358  }
359 
384  template <typename InputPointRange, typename WeightRange>
385  Alpha_complex_3d(const InputPointRange& points, WeightRange weights, FT x_min, FT y_min,
386  FT z_min, FT x_max, FT y_max, FT z_max) {
387  static_assert(Weighted, "This constructor is not available for non-weighted versions of Alpha_complex_3d");
388  static_assert(Periodic, "This constructor is not available for non-periodic versions of Alpha_complex_3d");
389  // FIXME: this test is only valid if we have a forward range
390  GUDHI_CHECK(boost::size(weights) == boost::size(points),
391  std::invalid_argument("Points number in range different from weights range number"));
392  // Checking if the cuboid is the same in x,y and z direction. If not, CGAL will not process it.
393  GUDHI_CHECK(
394  (x_max - x_min == y_max - y_min) && (x_max - x_min == z_max - z_min) && (z_max - z_min == y_max - y_min),
395  std::invalid_argument("The size of the cuboid in every directions is not the same."));
396 
397 #ifdef GUDHI_DEBUG
398  // Defined in GUDHI_DEBUG to avoid unused variable warning for GUDHI_CHECK
399  FT maximal_possible_weight = 0.015625 * (x_max - x_min) * (x_max - x_min);
400 #endif
401 
402  auto weighted_points_3 = boost::range::combine(points, weights)
403  | boost::adaptors::transformed([=](auto const&t){
404  auto w = boost::get<1>(t);
405  GUDHI_CHECK((w < maximal_possible_weight) && (w >= 0),
406  std::invalid_argument("Invalid weight " + std::to_string(w) +
407  ". Must be non-negative and less than maximal possible weight = 1/64*cuboid length squared."));
408  return Weighted_point_3(boost::get<0>(t), w);
409  });
410 
411  // Define the periodic cube
412  Dt pdt(typename Kernel::Iso_cuboid_3(x_min, y_min, z_min, x_max, y_max, z_max));
413  // Heuristic for inserting large point sets (if pts is reasonably large)
414  pdt.insert(std::begin(weighted_points_3), std::end(weighted_points_3), true);
415  // As pdt won't be modified anymore switch to 1-sheeted cover if possible
416  if (!pdt.is_triangulation_in_1_sheet()) {
417  throw std::invalid_argument("Unable to construct a triangulation within a single periodic domain.");
418  }
419  pdt.convert_to_1_sheeted_covering();
420 
421  // alpha shape construction from points. CGAL has a strange behavior in REGULARIZED mode. This is the default mode
422  // Maybe need to set it to GENERAL mode
423  alpha_shape_3_ptr_ = std::make_unique<Alpha_shape_3>(pdt, 0, Alpha_shape_3::GENERAL);
424  }
425 
440  template <typename SimplicialComplexForAlpha3d,
442  bool create_complex(SimplicialComplexForAlpha3d& complex,
443  Filtration_value max_alpha_square = std::numeric_limits<Filtration_value>::infinity()) {
444  if (complex.num_vertices() > 0) {
445  std::cerr << "Alpha_complex_3d create_complex - complex is not empty\n";
446  return false; // ----- >>
447  }
448 
449  using Complex_vertex_handle = typename SimplicialComplexForAlpha3d::Vertex_handle;
450  using Simplex_tree_vector_vertex = std::vector<Complex_vertex_handle>;
451 
452 #ifdef DEBUG_TRACES
453  std::size_t count_vertices = 0;
454  std::size_t count_edges = 0;
455  std::size_t count_facets = 0;
456  std::size_t count_cells = 0;
457 #endif // DEBUG_TRACES
458  std::vector<CGAL::Object> objects;
459  std::vector<FT> alpha_values;
460 
461  Dispatch dispatcher = CGAL::dispatch_output<CGAL::Object, FT>(std::back_inserter(objects),
462  std::back_inserter(alpha_values));
463 
464  alpha_shape_3_ptr_->filtration_with_alpha_values(dispatcher);
465 #ifdef DEBUG_TRACES
466  std::clog << "filtration_with_alpha_values returns : " << objects.size() << " objects" << std::endl;
467 #endif // DEBUG_TRACES
468  if (objects.size() == 0) {
469  std::cerr << "Alpha_complex_3d create_complex - no triangulation as points are on a 2d plane\n";
470  return false; // ----- >>
471  }
472 
473  using Alpha_value_iterator = typename std::vector<FT>::const_iterator;
474  Alpha_value_iterator alpha_value_iterator = alpha_values.begin();
475  for (auto object_iterator : objects) {
476  Vertex_list vertex_list;
477 
478  // Retrieve Alpha shape vertex list from object
479  if (const Cell_handle* cell = CGAL::object_cast<Cell_handle>(&object_iterator)) {
480  for (auto i = 0; i < 4; i++) {
481 #ifdef DEBUG_TRACES
482  std::clog << "from cell[" << i << "] - Point coordinates (" << (*cell)->vertex(i)->point() << ")"
483  << std::endl;
484 #endif // DEBUG_TRACES
485  vertex_list.push_back((*cell)->vertex(i));
486  }
487 #ifdef DEBUG_TRACES
488  count_cells++;
489 #endif // DEBUG_TRACES
490  } else if (const Facet* facet = CGAL::object_cast<Facet>(&object_iterator)) {
491  for (auto i = 0; i < 4; i++) {
492  if ((*facet).second != i) {
493 #ifdef DEBUG_TRACES
494  std::clog << "from facet=[" << i << "] - Point coordinates (" << (*facet).first->vertex(i)->point() << ")"
495  << std::endl;
496 #endif // DEBUG_TRACES
497  vertex_list.push_back((*facet).first->vertex(i));
498  }
499  }
500 #ifdef DEBUG_TRACES
501  count_facets++;
502 #endif // DEBUG_TRACES
503  } else if (const Edge* edge = CGAL::object_cast<Edge>(&object_iterator)) {
504  for (auto i : {(*edge).second, (*edge).third}) {
505 #ifdef DEBUG_TRACES
506  std::clog << "from edge[" << i << "] - Point coordinates (" << (*edge).first->vertex(i)->point() << ")"
507  << std::endl;
508 #endif // DEBUG_TRACES
509  vertex_list.push_back((*edge).first->vertex(i));
510  }
511 #ifdef DEBUG_TRACES
512  count_edges++;
513 #endif // DEBUG_TRACES
514  } else if (const Alpha_vertex_handle* vertex = CGAL::object_cast<Alpha_vertex_handle>(&object_iterator)) {
515 #ifdef DEBUG_TRACES
516  count_vertices++;
517  std::clog << "from vertex - Point coordinates (" << (*vertex)->point() << ")" << std::endl;
518 #endif // DEBUG_TRACES
519  vertex_list.push_back((*vertex));
520  }
521  // Construction of the vector of simplex_tree vertex from list of alpha_shapes vertex
522  Simplex_tree_vector_vertex the_simplex;
523  for (auto the_alpha_shape_vertex : vertex_list) {
524  auto the_map_iterator = map_cgal_simplex_tree.find(the_alpha_shape_vertex);
525  if (the_map_iterator == map_cgal_simplex_tree.end()) {
526  // alpha shape not found
527  Complex_vertex_handle vertex = map_cgal_simplex_tree.size();
528 #ifdef DEBUG_TRACES
529  std::clog << "Point (" << the_alpha_shape_vertex->point() << ") not found - insert new vertex id " << vertex
530  << std::endl;
531 #endif // DEBUG_TRACES
532  the_simplex.push_back(vertex);
533  map_cgal_simplex_tree.emplace(the_alpha_shape_vertex, vertex);
534  } else {
535  // alpha shape found
536  Complex_vertex_handle vertex = the_map_iterator->second;
537 #ifdef DEBUG_TRACES
538  std::clog << "Point (" << the_alpha_shape_vertex->point() << ") found as vertex id " << vertex << std::endl;
539 #endif // DEBUG_TRACES
540  the_simplex.push_back(vertex);
541  }
542  }
543  // Construction of the simplex_tree
544  Filtration_value filtr = Value_from_iterator<Complexity>::perform(alpha_value_iterator);
545 
546 #ifdef DEBUG_TRACES
547  std::clog << "filtration = " << filtr << std::endl;
548 #endif // DEBUG_TRACES
549  complex.insert_simplex(the_simplex, static_cast<Filtration_value>(filtr));
550  GUDHI_CHECK(alpha_value_iterator != alpha_values.end(), "CGAL provided more simplices than values");
551  ++alpha_value_iterator;
552  }
553 
554 #ifdef DEBUG_TRACES
555  std::clog << "vertices \t" << count_vertices << std::endl;
556  std::clog << "edges \t\t" << count_edges << std::endl;
557  std::clog << "facets \t\t" << count_facets << std::endl;
558  std::clog << "cells \t\t" << count_cells << std::endl;
559 #endif // DEBUG_TRACES
560  // --------------------------------------------------------------------------------------------
561  // As Alpha value is an approximation, we have to make filtration non decreasing while increasing the dimension
563  // Remove all simplices that have a filtration value greater than max_alpha_square
564  complex.prune_above_filtration(max_alpha_square);
565  // --------------------------------------------------------------------------------------------
566  return true;
567  }
568 
575  const Point_3& get_point(std::size_t vertex) {
576  if (map_cgal_simplex_tree.size() != cgal_vertex_iterator_vector.size()) {
577  cgal_vertex_iterator_vector.resize(map_cgal_simplex_tree.size());
578  for (auto map_iterator : map_cgal_simplex_tree) {
579  cgal_vertex_iterator_vector[map_iterator.second] = map_iterator.first;
580  }
581 
582  }
583  auto cgal_vertex_iterator = cgal_vertex_iterator_vector.at(vertex);
584  return cgal_vertex_iterator->point();
585  }
586 
587  private:
588  // use of a unique_ptr on cgal Alpha_shape_3, as copy and default constructor is not available - no need to be freed
589  std::unique_ptr<Alpha_shape_3> alpha_shape_3_ptr_;
590 
591  // Map type to switch from CGAL vertex iterator to simplex tree vertex handle.
592  std::unordered_map<Alpha_vertex_handle, std::size_t> map_cgal_simplex_tree;
593  // Vector type to switch from simplex tree vertex handle to CGAL vertex iterator.
594  std::vector<Alpha_vertex_handle> cgal_vertex_iterator_vector;
595 };
596 
597 } // namespace alpha_complex
598 
599 } // namespace Gudhi
600 
601 #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:442
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:575
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:337
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:385
Alpha complex data structure for 3d specific case.
Definition: Alpha_complex_3d.h:127
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:255
Alpha_complex_3d(const InputPointRange &points, WeightRange weights)
Alpha_complex constructor from a list of points and associated weights.
Definition: Alpha_complex_3d.h:302
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:280
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
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:237
GUDHI  Version 3.3.0  - C++ library for Topological Data Analysis (TDA) and Higher Dimensional Geometry Understanding.  - Copyright : MIT Generated on Tue Aug 11 2020 11:09:13 for GUDHI by Doxygen 1.8.13