12 #ifndef ALPHA_COMPLEX_3D_H_ 13 #define ALPHA_COMPLEX_3D_H_ 15 #include <boost/version.hpp> 16 #include <boost/variant.hpp> 18 #include <gudhi/Debug_utils.h> 19 #include <gudhi/Alpha_complex_options.h> 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> 33 #include <CGAL/Object.h> 34 #include <CGAL/tuple.h> 35 #include <CGAL/iterator.h> 36 #include <CGAL/version.h> 38 #include <Eigen/src/Core/util/Macros.h> 40 #include <boost/container/static_vector.hpp> 44 #include <unordered_map> 48 #include <type_traits> 52 #if CGAL_VERSION_NR < 1041101000 53 # error Alpha_complex_3d is only available for CGAL >= 4.11 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 62 namespace alpha_complex {
64 thread_local
double RELATIVE_PRECISION_OF_TO_DOUBLE = 0.00001;
71 template <complexity Complexity>
72 struct Value_from_iterator {
73 template <
typename Iterator>
74 static double perform(Iterator it) {
76 return CGAL::to_double(*it);
82 template <
typename Iterator>
83 static double perform(Iterator it) {
84 return CGAL::to_double(it->exact());
123 template <complexity Complexity = complexity::SAFE,
bool Weighted = false,
bool Periodic = false>
138 using Predicates =
typename std::conditional<(Complexity ==
complexity::FAST),
139 CGAL::Exact_predicates_inexact_constructions_kernel,
140 CGAL::Exact_predicates_exact_constructions_kernel>::type;
143 template <
typename Predicates,
bool Weighted_version,
bool Periodic_version>
146 template <
typename Predicates,
bool Is_periodic>
147 struct Kernel_3<Predicates, Is_periodic, false> {
148 using Kernel = Predicates;
151 template <
typename Predicates>
152 struct Kernel_3<Predicates, false, true> {
153 using Kernel = CGAL::Periodic_3_Delaunay_triangulation_traits_3<Predicates>;
155 template <
typename Predicates>
156 struct Kernel_3<Predicates, true, true> {
157 using Kernel = CGAL::Periodic_3_regular_triangulation_traits_3<Predicates>;
160 using Kernel =
typename Kernel_3<Predicates, Weighted, Periodic>::Kernel;
162 using TdsVb =
typename std::conditional<Periodic, CGAL::Periodic_3_triangulation_ds_vertex_base_3<>,
163 CGAL::Triangulation_ds_vertex_base_3<>>::type;
165 using Tvb =
typename std::conditional<Weighted, CGAL::Regular_triangulation_vertex_base_3<Kernel, TdsVb>,
166 CGAL::Triangulation_vertex_base_3<Kernel, TdsVb>>::type;
168 using Vb = CGAL::Alpha_shape_vertex_base_3<Kernel, Tvb>;
170 using TdsCb =
typename std::conditional<Periodic, CGAL::Periodic_3_triangulation_ds_cell_base_3<>,
171 CGAL::Triangulation_ds_cell_base_3<>>::type;
173 using Tcb =
typename std::conditional<Weighted, CGAL::Regular_triangulation_cell_base_3<Kernel, TdsCb>,
174 CGAL::Triangulation_cell_base_3<Kernel, TdsCb>>::type;
176 using Cb = CGAL::Alpha_shape_cell_base_3<Kernel, Tcb>;
177 using Tds = CGAL::Triangulation_data_structure_3<Vb, Cb>;
180 template <
typename Kernel,
typename Tds,
bool Weighted_version,
bool Periodic_version>
181 struct Triangulation_3 {};
183 template <
typename Kernel,
typename Tds>
184 struct Triangulation_3<Kernel, Tds, false, false> {
185 using Dt = CGAL::Delaunay_triangulation_3<Kernel, Tds>;
188 template <
typename Kernel,
typename Tds>
189 struct Triangulation_3<Kernel, Tds, true, false> {
190 using Dt = CGAL::Regular_triangulation_3<Kernel, Tds>;
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>;
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>;
211 using Dt =
typename Triangulation_3<Kernel, Tds, Weighted, Periodic>::Dt;
223 using FT =
typename Alpha_shape_3::FT;
247 using Weighted_point_3 =
typename Triangulation_3<Kernel, Tds, Weighted, Periodic>::Weighted_point_3;
252 using Point_3 =
typename Alpha_shape_3::Point;
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>>>>;
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>;
276 template <
typename InputPo
intRange>
278 static_assert(!Periodic,
"This constructor is not available for periodic versions of Alpha_complex_3d");
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));
298 template <
typename InputPo
intRange,
typename WeightRange>
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"));
305 std::vector<Weighted_point_3> weighted_points_3;
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]));
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));
338 template <
typename InputPo
intRange>
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");
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."));
348 Dt pdt(
typename Kernel::Iso_cuboid_3(x_min, y_min, z_min, x_max, y_max, z_max));
350 pdt.insert(std::begin(points), std::end(points),
true);
352 if (!pdt.is_triangulation_in_1_sheet()) {
353 throw std::invalid_argument(
"Unable to construct a triangulation within a single periodic domain.");
355 pdt.convert_to_1_sheeted_covering();
359 alpha_shape_3_ptr_ = std::unique_ptr<Alpha_shape_3>(
new Alpha_shape_3(pdt, 0, Alpha_shape_3::GENERAL));
386 template <
typename InputPo
intRange,
typename WeightRange>
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"));
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."));
398 std::vector<Weighted_point_3> weighted_points_3;
400 std::size_t index = 0;
401 weighted_points_3.reserve(points.size());
405 FT maximal_possible_weight = 0.015625 * (x_max - x_min) * (x_max - x_min);
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]));
418 Dt pdt(
typename Kernel::Iso_cuboid_3(x_min, y_min, z_min, x_max, y_max, z_max));
420 pdt.insert(std::begin(weighted_points_3), std::end(weighted_points_3),
true);
422 if (!pdt.is_triangulation_in_1_sheet()) {
423 throw std::invalid_argument(
"Unable to construct a triangulation within a single periodic domain.");
425 pdt.convert_to_1_sheeted_covering();
429 alpha_shape_3_ptr_ = std::unique_ptr<Alpha_shape_3>(
new Alpha_shape_3(pdt, 0, Alpha_shape_3::GENERAL));
449 Filtration_value max_alpha_square = std::numeric_limits<Filtration_value>::infinity()) {
451 std::cerr <<
"Alpha_complex_3d create_complex - complex is not empty\n";
456 using Simplex_tree_vector_vertex = std::vector<Complex_vertex_handle>;
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;
467 Dispatch dispatcher = CGAL::dispatch_output<CGAL::Object, FT>(std::back_inserter(objects),
468 std::back_inserter(alpha_values));
470 alpha_shape_3_ptr_->filtration_with_alpha_values(dispatcher);
472 std::clog <<
"filtration_with_alpha_values returns : " << objects.size() <<
" objects" << std::endl;
473 #endif // DEBUG_TRACES 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;
481 if (
const Cell_handle* cell = CGAL::object_cast<Cell_handle>(&object_iterator)) {
482 for (
auto i = 0; i < 4; i++) {
484 std::clog <<
"from cell[" << i <<
"] - Point coordinates (" << (*cell)->vertex(i)->point() <<
")" 486 #endif // DEBUG_TRACES 487 vertex_list.push_back((*cell)->vertex(i));
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) {
496 std::clog <<
"from facet=[" << i <<
"] - Point coordinates (" << (*facet).first->vertex(i)->point() <<
")" 498 #endif // DEBUG_TRACES 499 vertex_list.push_back((*facet).first->vertex(i));
504 #endif // DEBUG_TRACES 505 }
else if (
const Edge* edge = CGAL::object_cast<Edge>(&object_iterator)) {
506 for (
auto i : {(*edge).second, (*edge).third}) {
508 std::clog <<
"from edge[" << i <<
"] - Point coordinates (" << (*edge).first->vertex(i)->point() <<
")" 510 #endif // DEBUG_TRACES 511 vertex_list.push_back((*edge).first->vertex(i));
515 #endif // DEBUG_TRACES 516 }
else if (
const Alpha_vertex_handle* vertex = CGAL::object_cast<Alpha_vertex_handle>(&object_iterator)) {
519 std::clog <<
"from vertex - Point coordinates (" << (*vertex)->point() <<
")" << std::endl;
520 #endif // DEBUG_TRACES 521 vertex_list.push_back((*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()) {
529 Complex_vertex_handle vertex = map_cgal_simplex_tree.size();
531 std::clog <<
"Point (" << the_alpha_shape_vertex->point() <<
") not found - insert new vertex id " << vertex
533 #endif // DEBUG_TRACES 534 the_simplex.push_back(vertex);
535 map_cgal_simplex_tree.emplace(the_alpha_shape_vertex, vertex);
538 Complex_vertex_handle vertex = the_map_iterator->second;
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);
546 Filtration_value filtr = Value_from_iterator<Complexity>::perform(alpha_value_iterator);
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;
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 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;
585 auto cgal_vertex_iterator = cgal_vertex_iterator_vector.at(vertex);
586 return cgal_vertex_iterator->point();
591 std::unique_ptr<Alpha_shape_3> alpha_shape_3_ptr_;
594 std::unordered_map<Alpha_vertex_handle, std::size_t> map_cgal_simplex_tree;
596 std::vector<Alpha_vertex_handle> cgal_vertex_iterator_vector;
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
std::size_t num_vertices()
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
void make_filtration_non_decreasing()
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