Loading...
Searching...
No Matches
MEB_filtration.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): Marc Glisse
4 *
5 * Copyright (C) 2023 Inria
6 *
7 * Modification(s):
8 * - 2024/10 Vincent Rouvreau: Add Output_squared_values argument to enable/disable squared radii computation
9 * - YYYY/MM Author: Description of the modification
10 */
11
12#ifndef MEB_FILTRATION_H_
13#define MEB_FILTRATION_H_
14
15#include <CGAL/Lazy_exact_nt.h> // for CGAL::exact
16#include <CGAL/NT_converter.h>
17
18#include <vector>
19#include <utility> // for std::pair
20#include <cmath> // for std::sqrt
21
22namespace Gudhi::cech_complex {
23
46template<bool Output_squared_values = true, typename Kernel, typename SimplicialComplexForMEB, typename PointRange>
47void assign_MEB_filtration(Kernel&&k, SimplicialComplexForMEB& complex, PointRange const& points, bool exact = false) {
48 using Point_d = typename Kernel::Point_d;
49 using FT = typename Kernel::FT;
50 using Sphere = std::pair<Point_d, FT>;
51
52 using Vertex_handle = typename SimplicialComplexForMEB::Vertex_handle;
53 using Simplex_handle = typename SimplicialComplexForMEB::Simplex_handle;
54 using Filtration_value = typename SimplicialComplexForMEB::Filtration_value;
55
56 // For users to be able to define their own sqrt function on their desired Filtration_value type
57 using std::sqrt;
58
59 std::vector<Sphere> cache_;
60 std::vector<Point_d> pts;
61 CGAL::NT_converter<FT, Filtration_value> cvt;
62
63 // This block is only needed to get ambient_dim
64 if(std::begin(points) == std::end(points)) {
65 // assert(complex.is_empty());
66 return;
67 }
68 int ambient_dim = k.point_dimension_d_object()(*std::begin(points));
69
70 auto fun = [&](Simplex_handle sh, int dim){
71 using std::max;
72 if (dim == 0) complex.assign_filtration(sh, 0);
73 else if (dim == 1) {
74 // For a Simplex_tree, this would be a bit faster, but that's probably negligible
75 // Vertex_handle u = sh->first; Vertex_handle v = self_siblings(sh)->parent();
76 auto verts = complex.simplex_vertex_range(sh);
77 auto vert_it = verts.begin();
78 Vertex_handle u = *vert_it;
79 Vertex_handle v = *++vert_it;
80 auto&& pu = points[u];
81 Point_d m = k.midpoint_d_object()(pu, points[v]);
82 FT r = k.squared_distance_d_object()(m, pu);
83 if (exact) CGAL::exact(r);
84 complex.assign_key(sh, cache_.size());
85 Filtration_value filt{max(cvt(r), Filtration_value(0))};
86 if constexpr (!Output_squared_values)
87 filt = sqrt(filt);
88 complex.assign_filtration(sh, filt);
89 cache_.emplace_back(std::move(m), std::move(r));
90 } else if (dim > ambient_dim) {
91 // The sphere is always defined by at most d+1 points
92 Filtration_value maxf = 0; // max filtration of the faces
93 for (auto face : complex.boundary_simplex_range(sh)) {
94 maxf = max(maxf, complex.filtration(face));
95 }
96 complex.assign_filtration(sh, maxf);
97 } else {
98 Filtration_value maxf = 0; // max filtration of the faces
99 bool found = false;
100 for (auto face_opposite_vertex : complex.boundary_opposite_vertex_simplex_range(sh)) {
101 maxf = max(maxf, complex.filtration(face_opposite_vertex.first));
102 if (!found) {
103 auto key = complex.key(face_opposite_vertex.first);
104 Sphere const& sph = cache_[key];
105 if (k.squared_distance_d_object()(sph.first, points[face_opposite_vertex.second]) > sph.second) continue;
106 found = true;
107 complex.assign_key(sh, key);
108 // With exact computations, we could stop here
109 // complex.assign_filtration(sh, complex.filtration(face_opposite_vertex.first)); return;
110 // but because of possible rounding errors, we continue with the equivalent of make_filtration_non_decreasing
111 }
112 }
113 if (!found) {
114 // None of the faces are good enough, MEB must be the circumsphere.
115 pts.clear();
116 for (auto vertex : complex.simplex_vertex_range(sh))
117 pts.push_back(points[vertex]);
118 Point_d c = k.construct_circumcenter_d_object()(pts.begin(), pts.end());
119 FT r = k.squared_distance_d_object()(c, pts.front());
120 if (exact) CGAL::exact(r);
121 // For Epick_d, if the circumcenter computation is too unstable, we could compute
122 // int d2 = dim * dim;
123 // Filtration_value max_sanity = maxf * d2 / (d2 - 1);
124 // and use min(max_sanity, ...), which would limit how bad numerical errors can be.
125 Filtration_value filt{cvt(r)};
126 if constexpr (!Output_squared_values)
127 filt = sqrt(max(filt, Filtration_value(0)));
128 // maxf = filt except for rounding errors
129 maxf = max(maxf, filt);
130 complex.assign_key(sh, cache_.size());
131 // We could check if the simplex is maximal and avoiding adding it to the cache in that case.
132 cache_.emplace_back(std::move(c), std::move(r));
133 }
134 complex.assign_filtration(sh, maxf);
135 }
136 };
137 complex.for_each_simplex(fun);
138
139 // We could avoid computing maxf, but when !exact rounding errors may cause
140 // the filtration values to be non-monotonous, so we would need to call
141 // if (!exact) complex.make_filtration_non_decreasing();
142 // which is way more costly than computing maxf. The exact case is already so
143 // costly that it isn't worth maintaining code without maxf just for it.
144 // Cech_complex has "free" access to the max of the faces, because
145 // expansion_with_blockers computes it before the callback.
146
147 // TODO: use a map if complex does not provide key?
148}
149} // namespace Gudhi::cech_complex
150
151#endif // MEB_FILTRATION_H_
void assign_MEB_filtration(Kernel &&k, SimplicialComplexForMEB &complex, PointRange const &points, bool exact=false)
Given a simplicial complex and an embedding of its vertices, this assigns to each simplex a filtratio...
Definition MEB_filtration.h:47
Definition SimplicialComplexForMEB.h:22
unspecified Simplex_handle
Handle for a simplex.
Definition SimplicialComplexForMEB.h:24
Simplex_key key(Simplex_handle simplex)
Returns the key assigned to the 'simplex' with assign_key().
Filtration_value filtration(Simplex_handle simplex)
Returns the filtration value to the 'simplex'.
Boundary_opposite_vertex_simplex_range boundary_opposite_vertex_simplex_range(Simplex_handle simplex)
Returns a range of the pairs (simplex, opposite vertex) of the boundary of the 'simplex'.
unspecified Filtration_value
Type of filtration values.
Definition SimplicialComplexForMEB.h:29
void assign_key(Simplex_handle simplex, Simplex_key key)
Assigns this 'key' to the 'simplex'.
Simplex_vertex_range simplex_vertex_range(Simplex_handle simplex)
Returns a range over vertices (as Vertex_handle) of a given simplex.
unspecified Vertex_handle
Handle for a vertex. Must be a non-negative integer, it is also used as an index into the input list ...
Definition SimplicialComplexForMEB.h:27
int assign_filtration(Simplex_handle simplex, Filtration_value filtration)
Assigns this 'filtration' value to the 'simplex'.
void for_each_simplex(auto callback)
Calls callback(simplex, dim) for every simplex of the complex, with the guarantee that faces are visi...