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 <boost/container/static_vector.hpp>
42 
43 #include <iostream>
44 #include <vector>
45 #include <unordered_map>
46 #include <stdexcept>
47 #include <cstddef> // for std::size_t
48 #include <memory> // for std::unique_ptr
49 #include <type_traits> // for std::conditional and std::enable_if
50 #include <limits> // for numeric_limits<>
51 
52 // Make compilation fail - required for external projects - https://github.com/GUDHI/gudhi-devel/issues/10
53 #if CGAL_VERSION_NR < 1041101000
54 # error Alpha_complex_3d is only available for CGAL >= 4.11
55 #endif
56 
57 namespace Gudhi {
58 
59 namespace alpha_complex {
60 
61 thread_local double RELATIVE_PRECISION_OF_TO_DOUBLE = 0.00001;
62 
63 // Value_from_iterator returns the filtration value from an iterator on alpha shapes values
64 //
65 // FAST SAFE EXACT
66 // CGAL::to_double(*iterator) CGAL::to_double(*iterator) CGAL::to_double(iterator->exact())
67 
68 template <complexity Complexity>
69 struct Value_from_iterator {
70  template <typename Iterator>
71  static double perform(Iterator it) {
72  // Default behaviour
73  return CGAL::to_double(*it);
74  }
75 };
76 
77 template <>
78 struct Value_from_iterator<complexity::EXACT> {
79  template <typename Iterator>
80  static double perform(Iterator it) {
81  return CGAL::to_double(it->exact());
82  }
83 };
84 
120 template <complexity Complexity = complexity::SAFE, bool Weighted = false, bool Periodic = false>
122  // Epick = Exact_predicates_inexact_constructions_kernel
123  // Epeck = Exact_predicates_exact_constructions_kernel
124  // Exact_alpha_comparison_tag = exact version of CGAL Alpha_shape_3 and of its objects (Alpha_shape_vertex_base_3 and
125  // Alpha_shape_cell_base_3). Not available if weighted or periodic.
126  // Can be CGAL::Tag_false or CGAL::Tag_true. Default is False.
127  // cf. https://doc.cgal.org/latest/Alpha_shapes_3/classCGAL_1_1Alpha__shape__3.html
128  //
129  // We could use Epick + CGAL::Tag_true for not weighted nor periodic, but during benchmark, we found a bug
130  // https://github.com/CGAL/cgal/issues/3460
131  // This is the reason we only use Epick + CGAL::Tag_false, or Epeck
132  //
133  // FAST SAFE EXACT
134  // Epick + CGAL::Tag_false Epeck Epeck
135  using Predicates = typename std::conditional<(Complexity == complexity::FAST),
136  CGAL::Exact_predicates_inexact_constructions_kernel,
137  CGAL::Exact_predicates_exact_constructions_kernel>::type;
138 
139  // The other way to do a conditional type. Here there are 3 possibilities
140  template <typename Predicates, bool Weighted_version, bool Periodic_version>
141  struct Kernel_3 {};
142 
143  template <typename Predicates, bool Is_periodic>
144  struct Kernel_3<Predicates, Is_periodic, false> {
145  using Kernel = Predicates;
146  };
147 
148  template <typename Predicates>
149  struct Kernel_3<Predicates, false, true> {
150  using Kernel = CGAL::Periodic_3_Delaunay_triangulation_traits_3<Predicates>;
151  };
152  template <typename Predicates>
153  struct Kernel_3<Predicates, true, true> {
154  using Kernel = CGAL::Periodic_3_regular_triangulation_traits_3<Predicates>;
155  };
156 
157  public:
158  using Kernel = typename Kernel_3<Predicates, Weighted, Periodic>::Kernel;
159 
160  private:
161  using TdsVb = typename std::conditional<Periodic, CGAL::Periodic_3_triangulation_ds_vertex_base_3<>,
162  CGAL::Triangulation_ds_vertex_base_3<>>::type;
163 
164  using Tvb = typename std::conditional<Weighted, CGAL::Regular_triangulation_vertex_base_3<Kernel, TdsVb>,
165  CGAL::Triangulation_vertex_base_3<Kernel, TdsVb>>::type;
166 
167  using Vb = CGAL::Alpha_shape_vertex_base_3<Kernel, Tvb>;
168 
169  using TdsCb = typename std::conditional<Periodic, CGAL::Periodic_3_triangulation_ds_cell_base_3<>,
170  CGAL::Triangulation_ds_cell_base_3<>>::type;
171 
172  using Tcb = typename std::conditional<Weighted, CGAL::Regular_triangulation_cell_base_3<Kernel, TdsCb>,
173  CGAL::Triangulation_cell_base_3<Kernel, TdsCb>>::type;
174 
175  using Cb = CGAL::Alpha_shape_cell_base_3<Kernel, Tcb>;
176  using Tds = CGAL::Triangulation_data_structure_3<Vb, Cb>;
177 
178  // The other way to do a conditional type. Here there 4 possibilities, cannot use std::conditional
179  template <typename Kernel, typename Tds, bool Weighted_version, bool Periodic_version>
180  struct Triangulation_3 {};
181 
182  template <typename Kernel, typename Tds>
183  struct Triangulation_3<Kernel, Tds, false, false> {
184  using Dt = CGAL::Delaunay_triangulation_3<Kernel, Tds>;
185  using Weighted_point_3 = void;
186  };
187  template <typename Kernel, typename Tds>
188  struct Triangulation_3<Kernel, Tds, true, false> {
189  using Dt = CGAL::Regular_triangulation_3<Kernel, Tds>;
190  using Weighted_point_3 = typename Dt::Weighted_point;
191  };
192  template <typename Kernel, typename Tds>
193  struct Triangulation_3<Kernel, Tds, false, true> {
194  using Dt = CGAL::Periodic_3_Delaunay_triangulation_3<Kernel, Tds>;
195  using Weighted_point_3 = void;
196  };
197  template <typename Kernel, typename Tds>
198  struct Triangulation_3<Kernel, Tds, true, true> {
199  using Dt = CGAL::Periodic_3_regular_triangulation_3<Kernel, Tds>;
200  using Weighted_point_3 = typename Dt::Weighted_point;
201  };
202 
210  using Dt = typename Triangulation_3<Kernel, Tds, Weighted, Periodic>::Dt;
211 
212  public:
218  using Alpha_shape_3 = CGAL::Alpha_shape_3<Dt>;
219 
222  using FT = typename Alpha_shape_3::FT;
223 
233  using Bare_point_3 = typename Kernel::Point_3;
234 
246  using Weighted_point_3 = typename Triangulation_3<Kernel, Tds, Weighted, Periodic>::Weighted_point_3;
247 
251  using Point_3 = typename Alpha_shape_3::Point;
252 
253  private:
254  using Dispatch =
255  CGAL::Dispatch_output_iterator<CGAL::cpp11::tuple<CGAL::Object, FT>,
256  CGAL::cpp11::tuple<std::back_insert_iterator<std::vector<CGAL::Object>>,
257  std::back_insert_iterator<std::vector<FT>>>>;
258 
259  using Cell_handle = typename Alpha_shape_3::Cell_handle;
260  using Facet = typename Alpha_shape_3::Facet;
261  using Edge = typename Alpha_shape_3::Edge;
262  using Alpha_vertex_handle = typename Alpha_shape_3::Vertex_handle;
263  using Vertex_list = boost::container::static_vector<Alpha_vertex_handle, 4>;
264 
265  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::make_unique<Alpha_shape_3>(
280  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  // FIXME: this test is only valid if we have a forward range
302  GUDHI_CHECK(boost::size(weights) == boost::size(points),
303  std::invalid_argument("Points number in range different from weights range number"));
304 
305  auto weighted_points_3 = boost::range::combine(points, weights)
306  | boost::adaptors::transformed([](auto const&t){return Weighted_point_3(boost::get<0>(t), boost::get<1>(t));});
307 
308  alpha_shape_3_ptr_ = std::make_unique<Alpha_shape_3>(
309  std::begin(weighted_points_3), std::end(weighted_points_3), 0, Alpha_shape_3::GENERAL);
310  }
311 
332  template <typename InputPointRange>
333  Alpha_complex_3d(const InputPointRange& points, FT x_min, FT y_min,
334  FT z_min, FT x_max, FT y_max, FT z_max) {
335  static_assert(Periodic, "This constructor is not available for non-periodic versions of Alpha_complex_3d");
336  // Checking if the cuboid is the same in x,y and z direction. If not, CGAL will not process it.
337  GUDHI_CHECK(
338  (x_max - x_min == y_max - y_min) && (x_max - x_min == z_max - z_min) && (z_max - z_min == y_max - y_min),
339  std::invalid_argument("The size of the cuboid in every directions is not the same."));
340 
341  // Define the periodic cube
342  Dt pdt(typename Kernel::Iso_cuboid_3(x_min, y_min, z_min, x_max, y_max, z_max));
343  // Heuristic for inserting large point sets (if pts is reasonably large)
344  pdt.insert(std::begin(points), std::end(points), true);
345  // As pdt won't be modified anymore switch to 1-sheeted cover if possible
346  if (!pdt.is_triangulation_in_1_sheet()) {
347  throw std::invalid_argument("Unable to construct a triangulation within a single periodic domain.");
348  }
349  pdt.convert_to_1_sheeted_covering();
350 
351  // alpha shape construction from points. CGAL has a strange behavior in REGULARIZED mode. This is the default mode
352  // Maybe need to set it to GENERAL mode
353  alpha_shape_3_ptr_ = std::make_unique<Alpha_shape_3>(pdt, 0, Alpha_shape_3::GENERAL);
354  }
355 
380  template <typename InputPointRange, typename WeightRange>
381  Alpha_complex_3d(const InputPointRange& points, WeightRange weights, FT x_min, FT y_min,
382  FT z_min, FT x_max, FT y_max, FT z_max) {
383  static_assert(Weighted, "This constructor is not available for non-weighted versions of Alpha_complex_3d");
384  static_assert(Periodic, "This constructor is not available for non-periodic versions of Alpha_complex_3d");
385  // FIXME: this test is only valid if we have a forward range
386  GUDHI_CHECK(boost::size(weights) == boost::size(points),
387  std::invalid_argument("Points number in range different from weights range number"));
388  // Checking if the cuboid is the same in x,y and z direction. If not, CGAL will not process it.
389  GUDHI_CHECK(
390  (x_max - x_min == y_max - y_min) && (x_max - x_min == z_max - z_min) && (z_max - z_min == y_max - y_min),
391  std::invalid_argument("The size of the cuboid in every directions is not the same."));
392 
393 #ifdef GUDHI_DEBUG
394  // Defined in GUDHI_DEBUG to avoid unused variable warning for GUDHI_CHECK
395  FT maximal_possible_weight = 0.015625 * (x_max - x_min) * (x_max - x_min);
396 #endif
397 
398  auto weighted_points_3 = boost::range::combine(points, weights)
399  | boost::adaptors::transformed([=](auto const&t){
400  auto w = boost::get<1>(t);
401  GUDHI_CHECK((w < maximal_possible_weight) && (w >= 0),
402  std::invalid_argument("Invalid weight " + std::to_string(w) +
403  ". Must be non-negative and less than maximal possible weight = 1/64*cuboid length squared."));
404  return Weighted_point_3(boost::get<0>(t), w);
405  });
406 
407  // Define the periodic cube
408  Dt pdt(typename Kernel::Iso_cuboid_3(x_min, y_min, z_min, x_max, y_max, z_max));
409  // Heuristic for inserting large point sets (if pts is reasonably large)
410  pdt.insert(std::begin(weighted_points_3), std::end(weighted_points_3), true);
411  // As pdt won't be modified anymore switch to 1-sheeted cover if possible
412  if (!pdt.is_triangulation_in_1_sheet()) {
413  throw std::invalid_argument("Unable to construct a triangulation within a single periodic domain.");
414  }
415  pdt.convert_to_1_sheeted_covering();
416 
417  // alpha shape construction from points. CGAL has a strange behavior in REGULARIZED mode. This is the default mode
418  // Maybe need to set it to GENERAL mode
419  alpha_shape_3_ptr_ = std::make_unique<Alpha_shape_3>(pdt, 0, Alpha_shape_3::GENERAL);
420  }
421 
436  template <typename SimplicialComplexForAlpha3d,
438  bool create_complex(SimplicialComplexForAlpha3d& complex,
439  Filtration_value max_alpha_square = std::numeric_limits<Filtration_value>::infinity()) {
440  if (complex.num_vertices() > 0) {
441  std::cerr << "Alpha_complex_3d create_complex - complex is not empty\n";
442  return false; // ----- >>
443  }
444 
445  using Complex_vertex_handle = typename SimplicialComplexForAlpha3d::Vertex_handle;
446  using Simplex_tree_vector_vertex = std::vector<Complex_vertex_handle>;
447 
448 #ifdef DEBUG_TRACES
449  std::size_t count_vertices = 0;
450  std::size_t count_edges = 0;
451  std::size_t count_facets = 0;
452  std::size_t count_cells = 0;
453 #endif // DEBUG_TRACES
454  std::vector<CGAL::Object> objects;
455  std::vector<FT> alpha_values;
456 
457  Dispatch dispatcher = CGAL::dispatch_output<CGAL::Object, FT>(std::back_inserter(objects),
458  std::back_inserter(alpha_values));
459 
460  alpha_shape_3_ptr_->filtration_with_alpha_values(dispatcher);
461 #ifdef DEBUG_TRACES
462  std::clog << "filtration_with_alpha_values returns : " << objects.size() << " objects" << std::endl;
463 #endif // DEBUG_TRACES
464  if (objects.size() == 0) {
465  std::cerr << "Alpha_complex_3d create_complex - no triangulation as points are on a 2d plane\n";
466  return false; // ----- >>
467  }
468 
469  using Alpha_value_iterator = typename std::vector<FT>::const_iterator;
470  Alpha_value_iterator alpha_value_iterator = alpha_values.begin();
471  for (auto object_iterator : objects) {
472  Vertex_list vertex_list;
473 
474  // Retrieve Alpha shape vertex list from object
475  if (const Cell_handle* cell = CGAL::object_cast<Cell_handle>(&object_iterator)) {
476  for (auto i = 0; i < 4; i++) {
477 #ifdef DEBUG_TRACES
478  std::clog << "from cell[" << i << "] - Point coordinates (" << (*cell)->vertex(i)->point() << ")"
479  << std::endl;
480 #endif // DEBUG_TRACES
481  vertex_list.push_back((*cell)->vertex(i));
482  }
483 #ifdef DEBUG_TRACES
484  count_cells++;
485 #endif // DEBUG_TRACES
486  } else if (const Facet* facet = CGAL::object_cast<Facet>(&object_iterator)) {
487  for (auto i = 0; i < 4; i++) {
488  if ((*facet).second != i) {
489 #ifdef DEBUG_TRACES
490  std::clog << "from facet=[" << i << "] - Point coordinates (" << (*facet).first->vertex(i)->point() << ")"
491  << std::endl;
492 #endif // DEBUG_TRACES
493  vertex_list.push_back((*facet).first->vertex(i));
494  }
495  }
496 #ifdef DEBUG_TRACES
497  count_facets++;
498 #endif // DEBUG_TRACES
499  } else if (const Edge* edge = CGAL::object_cast<Edge>(&object_iterator)) {
500  for (auto i : {(*edge).second, (*edge).third}) {
501 #ifdef DEBUG_TRACES
502  std::clog << "from edge[" << i << "] - Point coordinates (" << (*edge).first->vertex(i)->point() << ")"
503  << std::endl;
504 #endif // DEBUG_TRACES
505  vertex_list.push_back((*edge).first->vertex(i));
506  }
507 #ifdef DEBUG_TRACES
508  count_edges++;
509 #endif // DEBUG_TRACES
510  } else if (const Alpha_vertex_handle* vertex = CGAL::object_cast<Alpha_vertex_handle>(&object_iterator)) {
511 #ifdef DEBUG_TRACES
512  count_vertices++;
513  std::clog << "from vertex - Point coordinates (" << (*vertex)->point() << ")" << std::endl;
514 #endif // DEBUG_TRACES
515  vertex_list.push_back((*vertex));
516  }
517  // Construction of the vector of simplex_tree vertex from list of alpha_shapes vertex
518  Simplex_tree_vector_vertex the_simplex;
519  for (auto the_alpha_shape_vertex : vertex_list) {
520  auto the_map_iterator = map_cgal_simplex_tree.find(the_alpha_shape_vertex);
521  if (the_map_iterator == map_cgal_simplex_tree.end()) {
522  // alpha shape not found
523  Complex_vertex_handle vertex = map_cgal_simplex_tree.size();
524 #ifdef DEBUG_TRACES
525  std::clog << "Point (" << the_alpha_shape_vertex->point() << ") not found - insert new vertex id " << vertex
526  << std::endl;
527 #endif // DEBUG_TRACES
528  the_simplex.push_back(vertex);
529  map_cgal_simplex_tree.emplace(the_alpha_shape_vertex, vertex);
530  } else {
531  // alpha shape found
532  Complex_vertex_handle vertex = the_map_iterator->second;
533 #ifdef DEBUG_TRACES
534  std::clog << "Point (" << the_alpha_shape_vertex->point() << ") found as vertex id " << vertex << std::endl;
535 #endif // DEBUG_TRACES
536  the_simplex.push_back(vertex);
537  }
538  }
539  // Construction of the simplex_tree
540  Filtration_value filtr = Value_from_iterator<Complexity>::perform(alpha_value_iterator);
541 
542 #ifdef DEBUG_TRACES
543  std::clog << "filtration = " << filtr << std::endl;
544 #endif // DEBUG_TRACES
545  complex.insert_simplex(the_simplex, static_cast<Filtration_value>(filtr));
546  GUDHI_CHECK(alpha_value_iterator != alpha_values.end(), "CGAL provided more simplices than values");
547  ++alpha_value_iterator;
548  }
549 
550 #ifdef DEBUG_TRACES
551  std::clog << "vertices \t" << count_vertices << std::endl;
552  std::clog << "edges \t\t" << count_edges << std::endl;
553  std::clog << "facets \t\t" << count_facets << std::endl;
554  std::clog << "cells \t\t" << count_cells << std::endl;
555 #endif // DEBUG_TRACES
556  // --------------------------------------------------------------------------------------------
557  // As Alpha value is an approximation, we have to make filtration non decreasing while increasing the dimension
559  // Remove all simplices that have a filtration value greater than max_alpha_square
560  complex.prune_above_filtration(max_alpha_square);
561  // --------------------------------------------------------------------------------------------
562  return true;
563  }
564 
571  const Point_3& get_point(std::size_t vertex) {
572  if (map_cgal_simplex_tree.size() != cgal_vertex_iterator_vector.size()) {
573  cgal_vertex_iterator_vector.resize(map_cgal_simplex_tree.size());
574  for (auto map_iterator : map_cgal_simplex_tree) {
575  cgal_vertex_iterator_vector[map_iterator.second] = map_iterator.first;
576  }
577 
578  }
579  auto cgal_vertex_iterator = cgal_vertex_iterator_vector.at(vertex);
580  return cgal_vertex_iterator->point();
581  }
582 
583  private:
584  // use of a unique_ptr on cgal Alpha_shape_3, as copy and default constructor is not available - no need to be freed
585  std::unique_ptr<Alpha_shape_3> alpha_shape_3_ptr_;
586 
587  // Map type to switch from CGAL vertex iterator to simplex tree vertex handle.
588  std::unordered_map<Alpha_vertex_handle, std::size_t> map_cgal_simplex_tree;
589  // Vector type to switch from simplex tree vertex handle to CGAL vertex iterator.
590  std::vector<Alpha_vertex_handle> cgal_vertex_iterator_vector;
591 };
592 
593 } // namespace alpha_complex
594 
595 } // namespace Gudhi
596 
597 #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:438
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:571
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:333
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:381
Alpha complex data structure for 3d specific case.
Definition: Alpha_complex_3d.h:121
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:251
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:246
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:218
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:222
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:233
GUDHI  Version 3.4.1  - C++ library for Topological Data Analysis (TDA) and Higher Dimensional Geometry Understanding.  - Copyright : MIT Generated on Fri Jan 22 2021 09:41:15 for GUDHI by Doxygen 1.8.13