Loading...
Searching...
No Matches
Tangential_complex.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): Clement Jamin
4 *
5 * Copyright (C) 2016 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 TANGENTIAL_COMPLEX_H_
13#define TANGENTIAL_COMPLEX_H_
14
15#include <gudhi/Tangential_complex/config.h>
16#include <gudhi/Tangential_complex/Simplicial_complex.h>
17#include <gudhi/Tangential_complex/utilities.h>
18#include <gudhi/Kd_tree_search.h>
19#include <gudhi/console_color.h>
20#include <gudhi/Clock.h>
21#include <gudhi/Simplex_tree.h>
22#include <gudhi/Debug_utils.h>
23
24#include <CGAL/Default.h>
25#include <CGAL/Dimension.h>
26#include <CGAL/function_objects.h> // for CGAL::Identity
27#include <CGAL/Epick_d.h>
28#include <CGAL/Regular_triangulation_traits_adapter.h>
29#include <CGAL/Regular_triangulation.h>
30#include <CGAL/Delaunay_triangulation.h>
31#include <CGAL/Combination_enumerator.h>
32#include <CGAL/point_generators_d.h>
33#include <CGAL/version.h> // for CGAL_VERSION_NR
34
35#include <Eigen/Core>
36#include <Eigen/Eigen>
37#include <Eigen/src/Core/util/Macros.h> // for EIGEN_VERSION_AT_LEAST
38
39#include <boost/iterator/transform_iterator.hpp>
40#include <boost/range/adaptor/transformed.hpp>
41#include <boost/range/counting_range.hpp>
42#include <boost/math/special_functions/factorials.hpp>
43#include <boost/container/flat_set.hpp>
44
45#include <tuple>
46#include <vector>
47#include <set>
48#include <utility>
49#include <sstream>
50#include <iostream>
51#include <limits>
52#include <algorithm>
53#include <functional>
54#include <iterator>
55#include <cmath> // for std::sqrt
56#include <string>
57#include <cstddef> // for std::size_t
58#include <optional>
59#include <numeric> // for std::iota
60
61#ifdef GUDHI_USE_TBB
62#include <tbb/parallel_for.h>
63#include <tbb/combinable.h>
64#include <mutex>
65#endif
66
67// #define GUDHI_TC_EXPORT_NORMALS // Only for 3D surfaces (k=2, d=3)
68
69// Make compilation fail - required for external projects - https://github.com/GUDHI/gudhi-devel/issues/10
70#if CGAL_VERSION_NR < 1041101000
71# error Tangential_complex is only available for CGAL >= 4.11
72#endif
73
74#if !EIGEN_VERSION_AT_LEAST(3,1,0)
75# error Tangential_complex is only available for Eigen3 >= 3.1.0 installed with CGAL
76#endif
77
78namespace sps = Gudhi::spatial_searching;
79
80namespace Gudhi {
81
82namespace tangential_complex {
83
84using namespace internal;
85
86class Vertex_data {
87 public:
88 Vertex_data(std::size_t data = (std::numeric_limits<std::size_t>::max)()) : m_data(data) {}
89
90 operator std::size_t() { return m_data; }
91
92 operator std::size_t() const { return m_data; }
93
94 private:
95 std::size_t m_data;
96};
97
123template <typename Kernel_, // ambiant kernel
124 typename DimensionTag, // intrinsic dimension
125 typename Concurrency_tag = CGAL::Parallel_tag, typename Triangulation_ = CGAL::Default>
127 typedef Kernel_ K;
128 typedef typename K::FT FT;
129 typedef typename K::Point_d Point;
130 typedef typename K::Weighted_point_d Weighted_point;
131 typedef typename K::Vector_d Vector;
132
133 typedef typename CGAL::Default::Get<
134 Triangulation_,
135 CGAL::Regular_triangulation<
136 CGAL::Epick_d<DimensionTag>,
137 CGAL::Triangulation_data_structure<
138 typename CGAL::Epick_d<DimensionTag>::Dimension,
139 CGAL::Triangulation_vertex<CGAL::Regular_triangulation_traits_adapter<CGAL::Epick_d<DimensionTag> >,
140 Vertex_data>,
141 CGAL::Triangulation_full_cell<
142 CGAL::Regular_triangulation_traits_adapter<CGAL::Epick_d<DimensionTag> > > > > >::type Triangulation;
143 typedef typename Triangulation::Geom_traits Tr_traits;
144 typedef typename Triangulation::Weighted_point Tr_point;
145 typedef typename Tr_traits::Base::Point_d Tr_bare_point;
146 typedef typename Triangulation::Vertex_handle Tr_vertex_handle;
147 typedef typename Triangulation::Full_cell_handle Tr_full_cell_handle;
148 typedef typename Tr_traits::Vector_d Tr_vector;
149
150#if defined(GUDHI_USE_TBB)
151 typedef std::mutex Mutex_for_perturb;
152 typedef Vector Translation_for_perturb;
153 typedef std::vector<Atomic_wrapper<FT> > Weights;
154#else
155 typedef Vector Translation_for_perturb;
156 typedef std::vector<FT> Weights;
157#endif
158 typedef std::vector<Translation_for_perturb> Translations_for_perturb;
159
160 // Store a local triangulation and a handle to its center vertex
161
162 struct Tr_and_VH {
163 public:
164 Tr_and_VH() : m_tr(NULL) {}
165
166 Tr_and_VH(int dim) : m_tr(new Triangulation(dim)) {}
167
168 ~Tr_and_VH() { destroy_triangulation(); }
169
170 Triangulation &construct_triangulation(int dim) {
171 delete m_tr;
172 m_tr = new Triangulation(dim);
173 return tr();
174 }
175
176 void destroy_triangulation() {
177 delete m_tr;
178 m_tr = NULL;
179 }
180
181 Triangulation &tr() { return *m_tr; }
182
183 Triangulation const &tr() const { return *m_tr; }
184
185 Tr_vertex_handle const &center_vertex() const { return m_center_vertex; }
186
187 Tr_vertex_handle &center_vertex() { return m_center_vertex; }
188
189 private:
190 Triangulation *m_tr;
191 Tr_vertex_handle m_center_vertex;
192 };
193
194 public:
195 typedef Basis<K> Tangent_space_basis;
196 typedef Basis<K> Orthogonal_space_basis;
197 typedef std::vector<Tangent_space_basis> TS_container;
198 typedef std::vector<Orthogonal_space_basis> OS_container;
199
200 typedef std::vector<Point> Points;
201
202 typedef boost::container::flat_set<std::size_t> Simplex;
203 typedef std::set<Simplex> Simplex_set;
204
205 private:
207 typedef typename Points_ds::KNS_range KNS_range;
208 typedef typename Points_ds::INS_range INS_range;
209
210 typedef std::vector<Tr_and_VH> Tr_container;
211 typedef std::vector<Vector> Vectors;
212
213 // An Incident_simplex is the list of the vertex indices
214 // except the center vertex
215 typedef boost::container::flat_set<std::size_t> Incident_simplex;
216 typedef std::vector<Incident_simplex> Star;
217 typedef std::vector<Star> Stars_container;
218
219 // For transform_iterator
220
221 static const Tr_point &vertex_handle_to_point(Tr_vertex_handle vh) { return vh->point(); }
222
223 template <typename P, typename VH>
224 static const P &vertex_handle_to_point(VH vh) {
225 return vh->point();
226 }
227
228 public:
229 typedef internal::Simplicial_complex Simplicial_complex;
230
240 template <typename Point_range>
242#ifdef GUDHI_TC_USE_ANOTHER_POINT_SET_FOR_TANGENT_SPACE_ESTIM
243 InputIterator first_for_tse, InputIterator last_for_tse,
244#endif
245 const K &k = K())
246 : m_k(k),
247 m_intrinsic_dim(intrinsic_dimension),
248 m_ambient_dim(points.empty() ? 0 : k.point_dimension_d_object()(*points.begin())),
249 m_points(points.begin(), points.end()),
250 m_weights(m_points.size(), FT(0))
251#if defined(GUDHI_USE_TBB) && defined(GUDHI_TC_PERTURB_POSITION)
252 ,
253 m_p_perturb_mutexes(NULL)
254#endif
255 ,
256 m_points_ds(m_points),
257 m_last_max_perturb(0.),
258 m_are_tangent_spaces_computed(m_points.size(), false),
259 m_tangent_spaces(m_points.size(), Tangent_space_basis())
260#ifdef GUDHI_TC_EXPORT_NORMALS
261 ,
262 m_orth_spaces(m_points.size(), Orthogonal_space_basis())
263#endif
264#ifdef GUDHI_TC_USE_ANOTHER_POINT_SET_FOR_TANGENT_SPACE_ESTIM
265 ,
266 m_points_for_tse(first_for_tse, last_for_tse),
267 m_points_ds_for_tse(m_points_for_tse)
268#endif
269 {
270 }
271
274#if defined(GUDHI_USE_TBB) && defined(GUDHI_TC_PERTURB_POSITION)
275 delete[] m_p_perturb_mutexes;
276#endif
277 }
278
280 int intrinsic_dimension() const { return m_intrinsic_dim; }
281
283 int ambient_dimension() const { return m_ambient_dim; }
284
285 Points const &points() const { return m_points; }
286
292 Point get_point(std::size_t vertex) const { return m_points[vertex]; }
293
299 Point get_perturbed_point(std::size_t vertex) const { return compute_perturbed_point(vertex); }
300
302
303 std::size_t number_of_vertices() const { return m_points.size(); }
304
305 void set_weights(const Weights &weights) { m_weights = weights; }
306
307 void set_tangent_planes(const TS_container &tangent_spaces
308#ifdef GUDHI_TC_EXPORT_NORMALS
309 ,
310 const OS_container &orthogonal_spaces
311#endif
312 ) {
313#ifdef GUDHI_TC_EXPORT_NORMALS
314 GUDHI_CHECK(m_points.size() == tangent_spaces.size() && m_points.size() == orthogonal_spaces.size(),
315 std::logic_error("Wrong sizes"));
316#else
317 GUDHI_CHECK(m_points.size() == tangent_spaces.size(), std::logic_error("Wrong sizes"));
318#endif
319 m_tangent_spaces = tangent_spaces;
320#ifdef GUDHI_TC_EXPORT_NORMALS
321 m_orth_spaces = orthogonal_spaces;
322#endif
323 for (std::size_t i = 0; i < m_points.size(); ++i) m_are_tangent_spaces_computed[i] = true;
324 }
325
332#ifdef GUDHI_TC_PERFORM_EXTRA_CHECKS
333 std::cerr << red << "WARNING: GUDHI_TC_PERFORM_EXTRA_CHECKS is defined. "
334 << "Computation might be slower than usual.\n"
335 << white;
336#endif
337
338#if defined(GUDHI_TC_PROFILING) && defined(GUDHI_USE_TBB)
339 Gudhi::Clock t;
340#endif
341
342 // We need to do that because we don't want the container to copy the
343 // already-computed triangulations (while resizing) since it would
344 // invalidate the vertex handles stored beside the triangulations
345 m_triangulations.resize(m_points.size());
346 m_stars.resize(m_points.size());
347 m_squared_star_spheres_radii_incl_margin.resize(m_points.size(), FT(-1));
348#ifdef GUDHI_TC_PERTURB_POSITION
349 if (m_points.empty()) {
350 m_translations.clear();
351 } else {
352 m_translations.resize(m_points.size(), m_k.construct_vector_d_object()(m_ambient_dim));
353 }
354#if defined(GUDHI_USE_TBB)
355 delete[] m_p_perturb_mutexes;
356 m_p_perturb_mutexes = new Mutex_for_perturb[m_points.size()];
357#endif
358#endif
359
360#ifdef GUDHI_USE_TBB
361 // Parallel
362 if (boost::is_convertible<Concurrency_tag, CGAL::Parallel_tag>::value) {
363 tbb::parallel_for(tbb::blocked_range<size_t>(0, m_points.size()), Compute_tangent_triangulation(*this));
364 } else {
365#endif // GUDHI_USE_TBB
366 // Sequential
367 for (std::size_t i = 0; i < m_points.size(); ++i) compute_tangent_triangulation(i);
368#ifdef GUDHI_USE_TBB
369 }
370#endif // GUDHI_USE_TBB
371
372#if defined(GUDHI_TC_PROFILING) && defined(GUDHI_USE_TBB)
373 t.end();
374 std::cerr << "Tangential complex computed in " << t.num_seconds() << " seconds.\n";
375#endif
376 }
377
381 bool success = false;
383 unsigned int num_steps = 0;
390 };
391
397 Fix_inconsistencies_info fix_inconsistencies_using_perturbation(double max_perturb, double time_limit = -1.) {
399
400 if (time_limit == 0.) return info;
401
402 Gudhi::Clock t;
403
404#ifdef GUDHI_TC_SHOW_DETAILED_STATS_FOR_INCONSISTENCIES
405 std::tuple<std::size_t, std::size_t, std::size_t> stats_before = number_of_inconsistent_simplices(false);
406
407 if (std::get<1>(stats_before) == 0) {
408#ifdef DEBUG_TRACES
409 std::cerr << "Nothing to fix.\n";
410#endif
411 info.success = false;
412 return info;
413 }
414#endif // GUDHI_TC_SHOW_DETAILED_STATS_FOR_INCONSISTENCIES
415
416 m_last_max_perturb = max_perturb;
417
418 bool done = false;
419 info.best_num_inconsistent_stars = m_triangulations.size();
420 info.num_steps = 0;
421 while (!done) {
422#ifdef GUDHI_TC_SHOW_DETAILED_STATS_FOR_INCONSISTENCIES
423 std::cerr << "\nBefore fix step:\n"
424 << " * Total number of simplices in stars (incl. duplicates): " << std::get<0>(stats_before) << "\n"
425 << " * Num inconsistent simplices in stars (incl. duplicates): " << red << std::get<1>(stats_before)
426 << white << " (" << 100. * std::get<1>(stats_before) / std::get<0>(stats_before) << "%)\n"
427 << " * Number of stars containing inconsistent simplices: " << red << std::get<2>(stats_before)
428 << white << " (" << 100. * std::get<2>(stats_before) / m_points.size() << "%)\n";
429#endif
430
431#if defined(DEBUG_TRACES) || defined(GUDHI_TC_PROFILING)
432 std::cerr << yellow << "\nAttempt to fix inconsistencies using perturbations - step #" << info.num_steps + 1
433 << "... " << white;
434#endif
435
436 std::size_t num_inconsistent_stars = 0;
437 std::vector<std::size_t> updated_points;
438
439#ifdef GUDHI_TC_PROFILING
440 Gudhi::Clock t_fix_step;
441#endif
442
443 // Parallel
444#if defined(GUDHI_USE_TBB)
445 if (boost::is_convertible<Concurrency_tag, CGAL::Parallel_tag>::value) {
446 tbb::combinable<std::size_t> num_inconsistencies;
447 tbb::combinable<std::vector<std::size_t> > tls_updated_points;
448 tbb::parallel_for(tbb::blocked_range<size_t>(0, m_triangulations.size()),
449 Try_to_solve_inconsistencies_in_a_local_triangulation(*this, max_perturb, num_inconsistencies,
450 tls_updated_points));
451 num_inconsistent_stars = num_inconsistencies.combine(std::plus<std::size_t>());
452 updated_points =
453 tls_updated_points.combine([](std::vector<std::size_t> const &x, std::vector<std::size_t> const &y) {
454 std::vector<std::size_t> res;
455 res.reserve(x.size() + y.size());
456 res.insert(res.end(), x.begin(), x.end());
457 res.insert(res.end(), y.begin(), y.end());
458 return res;
459 });
460 } else {
461#endif // GUDHI_USE_TBB
462 // Sequential
463 for (std::size_t i = 0; i < m_triangulations.size(); ++i) {
464 num_inconsistent_stars +=
465 try_to_solve_inconsistencies_in_a_local_triangulation(i, max_perturb, std::back_inserter(updated_points));
466 }
467#if defined(GUDHI_USE_TBB)
468 }
469#endif // GUDHI_USE_TBB
470
471#ifdef GUDHI_TC_PROFILING
472 t_fix_step.end();
473#endif
474
475#if defined(GUDHI_TC_SHOW_DETAILED_STATS_FOR_INCONSISTENCIES) || defined(DEBUG_TRACES)
476 std::cerr << "\nEncountered during fix:\n"
477 << " * Num stars containing inconsistent simplices: " << red << num_inconsistent_stars << white << " ("
478 << 100. * num_inconsistent_stars / m_points.size() << "%)\n";
479#endif
480
481#ifdef GUDHI_TC_PROFILING
482 std::cerr << yellow << "done in " << t_fix_step.num_seconds() << " seconds.\n" << white;
483#elif defined(DEBUG_TRACES)
484 std::cerr << yellow << "done.\n" << white;
485#endif
486
487 if (num_inconsistent_stars > 0) refresh_tangential_complex(updated_points);
488
489#ifdef GUDHI_TC_PERFORM_EXTRA_CHECKS
490 // Confirm that all stars were actually refreshed
491 std::size_t num_inc_1 = std::get<1>(number_of_inconsistent_simplices(false));
492 refresh_tangential_complex();
493 std::size_t num_inc_2 = std::get<1>(number_of_inconsistent_simplices(false));
494 if (num_inc_1 != num_inc_2)
495 std::cerr << red << "REFRESHMENT CHECK: FAILED. (" << num_inc_1 << " vs " << num_inc_2 << ")\n" << white;
496 else
497 std::cerr << green << "REFRESHMENT CHECK: PASSED.\n" << white;
498#endif
499
500#ifdef GUDHI_TC_SHOW_DETAILED_STATS_FOR_INCONSISTENCIES
501 std::tuple<std::size_t, std::size_t, std::size_t> stats_after = number_of_inconsistent_simplices(false);
502
503 std::cerr << "\nAfter fix:\n"
504 << " * Total number of simplices in stars (incl. duplicates): " << std::get<0>(stats_after) << "\n"
505 << " * Num inconsistent simplices in stars (incl. duplicates): " << red << std::get<1>(stats_after)
506 << white << " (" << 100. * std::get<1>(stats_after) / std::get<0>(stats_after) << "%)\n"
507 << " * Number of stars containing inconsistent simplices: " << red << std::get<2>(stats_after) << white
508 << " (" << 100. * std::get<2>(stats_after) / m_points.size() << "%)\n";
509
510 stats_before = stats_after;
511#endif
512
513 if (info.num_steps == 0) info.initial_num_inconsistent_stars = num_inconsistent_stars;
514
515 if (num_inconsistent_stars < info.best_num_inconsistent_stars)
516 info.best_num_inconsistent_stars = num_inconsistent_stars;
517
518 info.final_num_inconsistent_stars = num_inconsistent_stars;
519
520 done = (num_inconsistent_stars == 0);
521 if (!done) {
522 ++info.num_steps;
523 if (time_limit > 0. && t.num_seconds() > time_limit) {
524#ifdef DEBUG_TRACES
525 std::cerr << red << "Time limit reached.\n" << white;
526#endif
527 info.success = false;
528 return info;
529 }
530 }
531 }
532
533#ifdef DEBUG_TRACES
534 std::cerr << green << "Fixed!\n" << white;
535#endif
536 info.success = true;
537 return info;
538 }
539
543 std::size_t num_simplices = 0;
547 std::size_t num_inconsistent_stars = 0;
548 };
549
552
554#ifdef DEBUG_TRACES
555 bool verbose = true
556#else
557 bool verbose = false
558#endif
559 ) const {
561
562 // For each triangulation
563 for (std::size_t idx = 0; idx < m_points.size(); ++idx) {
564 bool is_star_inconsistent = false;
565
566 // For each cell
567 Star::const_iterator it_inc_simplex = m_stars[idx].begin();
568 Star::const_iterator it_inc_simplex_end = m_stars[idx].end();
569 for (; it_inc_simplex != it_inc_simplex_end; ++it_inc_simplex) {
570 // Don't check infinite cells
571 if (is_infinite(*it_inc_simplex)) continue;
572
573 Simplex c = *it_inc_simplex;
574 c.insert(idx); // Add the missing index
575
576 if (!is_simplex_consistent(c)) {
578 is_star_inconsistent = true;
579 }
580
581 ++stats.num_simplices;
582 }
583 stats.num_inconsistent_stars += is_star_inconsistent;
584 }
585
586 if (verbose) {
587 std::cerr << "\n==========================================================\n"
588 << "Inconsistencies:\n"
589 << " * Total number of simplices in stars (incl. duplicates): " << stats.num_simplices << "\n"
590 << " * Number of inconsistent simplices in stars (incl. duplicates): "
591 << stats.num_inconsistent_simplices << " ("
592 << 100. * stats.num_inconsistent_simplices / stats.num_simplices << "%)\n"
593 << " * Number of stars containing inconsistent simplices: " << stats.num_inconsistent_stars << " ("
594 << 100. * stats.num_inconsistent_stars / m_points.size() << "%)\n"
595 << "==========================================================\n";
596 }
597
598 return stats;
599 }
600
611 template <typename Simplex_tree_>
612 int create_complex(Simplex_tree_ &tree,
613 bool export_inconsistent_simplices = true
615 ,
616 bool export_infinite_simplices = false, Simplex_set *p_inconsistent_simplices = NULL
618 ) const {
619#if defined(DEBUG_TRACES) || defined(GUDHI_TC_PROFILING)
620 std::cerr << yellow << "\nExporting the TC as a Simplex_tree... " << white;
621#endif
622#ifdef GUDHI_TC_PROFILING
623 Gudhi::Clock t;
624#endif
625
626 int max_dim = -1;
627
628 // Ordered vertices to be inserted first by the create_complex method to avoid quadratic complexity.
629 std::vector<typename Simplex_tree_::Vertex_handle> vertices(m_points.size());
630 std::iota(vertices.begin(), vertices.end(), 0);
631 tree.insert_batch_vertices(vertices);
632
633 // For each triangulation
634 for (std::size_t idx = 0; idx < m_points.size(); ++idx) {
635 // For each cell of the star
636 Star::const_iterator it_inc_simplex = m_stars[idx].begin();
637 Star::const_iterator it_inc_simplex_end = m_stars[idx].end();
638 for (; it_inc_simplex != it_inc_simplex_end; ++it_inc_simplex) {
639 Simplex c = *it_inc_simplex;
640
641 // Don't export infinite cells
642 if (!export_infinite_simplices && is_infinite(c)) continue;
643
644 if (static_cast<int>(c.size()) > max_dim) max_dim = static_cast<int>(c.size());
645 // Add the missing center vertex
646 c.insert(idx);
647
648 if (!export_inconsistent_simplices && !is_simplex_consistent(c)) continue;
649
650 // Try to insert the simplex
651 bool inserted = tree.insert_simplex_and_subfaces(c).second;
652
653 // Inconsistent?
654 if (p_inconsistent_simplices && inserted && !is_simplex_consistent(c)) {
655 p_inconsistent_simplices->insert(c);
656 }
657 }
658 }
659
660#ifdef GUDHI_TC_PROFILING
661 t.end();
662 std::cerr << yellow << "done in " << t.num_seconds() << " seconds.\n" << white;
663#elif defined(DEBUG_TRACES)
664 std::cerr << yellow << "done.\n" << white;
665#endif
666
667 return max_dim;
668 }
669
670 // First clears the complex then exports the TC into it
671 // Returns the max dimension of the simplices
672 // check_lower_and_higher_dim_simplices : 0 (false), 1 (true), 2 (auto)
673 // If the check is enabled, the function:
674 // - won't insert the simplex if it is already in a higher dim simplex
675 // - will erase any lower-dim simplices that are faces of the new simplex
676 // "auto" (= 2) will enable the check as a soon as it encounters a
677 // simplex whose dimension is different from the previous ones.
678 // N.B.: The check is quite expensive.
679
680 int create_complex(Simplicial_complex &complex, bool export_inconsistent_simplices = true,
681 bool export_infinite_simplices = false, int check_lower_and_higher_dim_simplices = 2,
682 Simplex_set *p_inconsistent_simplices = NULL) const {
683#if defined(DEBUG_TRACES) || defined(GUDHI_TC_PROFILING)
684 std::cerr << yellow << "\nExporting the TC as a Simplicial_complex... " << white;
685#endif
686#ifdef GUDHI_TC_PROFILING
687 Gudhi::Clock t;
688#endif
689
690 int max_dim = -1;
691 complex.clear();
692
693 // For each triangulation
694 for (std::size_t idx = 0; idx < m_points.size(); ++idx) {
695 // For each cell of the star
696 Star::const_iterator it_inc_simplex = m_stars[idx].begin();
697 Star::const_iterator it_inc_simplex_end = m_stars[idx].end();
698 for (; it_inc_simplex != it_inc_simplex_end; ++it_inc_simplex) {
699 Simplex c = *it_inc_simplex;
700
701 // Don't export infinite cells
702 if (!export_infinite_simplices && is_infinite(c)) continue;
703
704 if (static_cast<int>(c.size()) > max_dim) max_dim = static_cast<int>(c.size());
705 // Add the missing center vertex
706 c.insert(idx);
707
708 if (!export_inconsistent_simplices && !is_simplex_consistent(c)) continue;
709
710 // Unusual simplex dim?
711 if (check_lower_and_higher_dim_simplices == 2 && max_dim != -1 && static_cast<int>(c.size()) != max_dim) {
712 // Let's activate the check
713 std::cerr << red
714 << "Info: check_lower_and_higher_dim_simplices ACTIVATED. "
715 "Export might be take some time...\n"
716 << white;
717 check_lower_and_higher_dim_simplices = 1;
718 }
719
720 // Try to insert the simplex
721 bool added = complex.add_simplex(c, check_lower_and_higher_dim_simplices == 1);
722
723 // Inconsistent?
724 if (p_inconsistent_simplices && added && !is_simplex_consistent(c)) {
725 p_inconsistent_simplices->insert(c);
726 }
727 }
728 }
729
730#ifdef GUDHI_TC_PROFILING
731 t.end();
732 std::cerr << yellow << "done in " << t.num_seconds() << " seconds.\n" << white;
733#elif defined(DEBUG_TRACES)
734 std::cerr << yellow << "done.\n" << white;
735#endif
736
737 return max_dim;
738 }
739
740 template <typename ProjectionFunctor = CGAL::Identity<Point> >
741 std::ostream &export_to_off(const Simplicial_complex &complex, std::ostream &os,
742 Simplex_set const *p_simpl_to_color_in_red = NULL,
743 Simplex_set const *p_simpl_to_color_in_green = NULL,
744 Simplex_set const *p_simpl_to_color_in_blue = NULL,
745 ProjectionFunctor const &point_projection = ProjectionFunctor()) const {
746 return export_to_off(os, false, p_simpl_to_color_in_red, p_simpl_to_color_in_green, p_simpl_to_color_in_blue,
747 &complex, point_projection);
748 }
749
750 template <typename ProjectionFunctor = CGAL::Identity<Point> >
751 std::ostream &export_to_off(std::ostream &os, bool color_inconsistencies = false,
752 Simplex_set const *p_simpl_to_color_in_red = NULL,
753 Simplex_set const *p_simpl_to_color_in_green = NULL,
754 Simplex_set const *p_simpl_to_color_in_blue = NULL,
755 const Simplicial_complex *p_complex = NULL,
756 ProjectionFunctor const &point_projection = ProjectionFunctor()) const {
757 if (m_points.empty()) return os;
758
759 if (m_ambient_dim < 2) {
760 std::cerr << "Error: export_to_off => ambient dimension should be >= 2.\n";
761 os << "Error: export_to_off => ambient dimension should be >= 2.\n";
762 return os;
763 }
764 if (m_ambient_dim > 3) {
765 std::cerr << "Warning: export_to_off => ambient dimension should be "
766 "<= 3. Only the first 3 coordinates will be exported.\n";
767 }
768
769 if (m_intrinsic_dim < 1 || m_intrinsic_dim > 3) {
770 std::cerr << "Error: export_to_off => intrinsic dimension should be "
771 "between 1 and 3.\n";
772 os << "Error: export_to_off => intrinsic dimension should be "
773 "between 1 and 3.\n";
774 return os;
775 }
776
777 std::stringstream output;
778 std::size_t num_simplices, num_vertices;
779 export_vertices_to_off(output, num_vertices, false, point_projection);
780 if (p_complex) {
781 export_simplices_to_off(*p_complex, output, num_simplices, p_simpl_to_color_in_red, p_simpl_to_color_in_green,
782 p_simpl_to_color_in_blue);
783 } else {
784 export_simplices_to_off(output, num_simplices, color_inconsistencies, p_simpl_to_color_in_red,
785 p_simpl_to_color_in_green, p_simpl_to_color_in_blue);
786 }
787
788#ifdef GUDHI_TC_EXPORT_NORMALS
789 os << "N";
790#endif
791
792 os << "OFF \n"
793 << num_vertices << " " << num_simplices << " "
794 << "0 \n"
795 << output.str();
796
797 return os;
798 }
799
800 private:
801 void refresh_tangential_complex() {
802#if defined(DEBUG_TRACES) || defined(GUDHI_TC_PROFILING)
803 std::cerr << yellow << "\nRefreshing TC... " << white;
804#endif
805
806#ifdef GUDHI_TC_PROFILING
807 Gudhi::Clock t;
808#endif
809#ifdef GUDHI_USE_TBB
810 // Parallel
811 if (boost::is_convertible<Concurrency_tag, CGAL::Parallel_tag>::value) {
812 tbb::parallel_for(tbb::blocked_range<size_t>(0, m_points.size()), Compute_tangent_triangulation(*this));
813 } else {
814#endif // GUDHI_USE_TBB
815 // Sequential
816 for (std::size_t i = 0; i < m_points.size(); ++i) compute_tangent_triangulation(i);
817#ifdef GUDHI_USE_TBB
818 }
819#endif // GUDHI_USE_TBB
820
821#ifdef GUDHI_TC_PROFILING
822 t.end();
823 std::cerr << yellow << "done in " << t.num_seconds() << " seconds.\n" << white;
824#elif defined(DEBUG_TRACES)
825 std::cerr << yellow << "done.\n" << white;
826#endif
827 }
828
829 // If the list of perturbed points is provided, it is much faster
830 template <typename Point_indices_range>
831 void refresh_tangential_complex(Point_indices_range const &perturbed_points_indices) {
832#if defined(DEBUG_TRACES) || defined(GUDHI_TC_PROFILING)
833 std::cerr << yellow << "\nRefreshing TC... " << white;
834#endif
835
836#ifdef GUDHI_TC_PROFILING
837 Gudhi::Clock t;
838#endif
839
840 // ANN tree containing only the perturbed points
841 Points_ds updated_pts_ds(m_points, perturbed_points_indices);
842
843#ifdef GUDHI_USE_TBB
844 // Parallel
845 if (boost::is_convertible<Concurrency_tag, CGAL::Parallel_tag>::value) {
846 tbb::parallel_for(tbb::blocked_range<size_t>(0, m_points.size()),
847 Refresh_tangent_triangulation(*this, updated_pts_ds));
848 } else {
849#endif // GUDHI_USE_TBB
850 // Sequential
851 for (std::size_t i = 0; i < m_points.size(); ++i) refresh_tangent_triangulation(i, updated_pts_ds);
852#ifdef GUDHI_USE_TBB
853 }
854#endif // GUDHI_USE_TBB
855
856#ifdef GUDHI_TC_PROFILING
857 t.end();
858 std::cerr << yellow << "done in " << t.num_seconds() << " seconds.\n" << white;
859#elif defined(DEBUG_TRACES)
860 std::cerr << yellow << "done.\n" << white;
861#endif
862 }
863
864 void export_inconsistent_stars_to_OFF_files(std::string const &filename_base) const {
865 // For each triangulation
866 for (std::size_t idx = 0; idx < m_points.size(); ++idx) {
867 // We build a SC along the way in case it's inconsistent
868 Simplicial_complex sc;
869 // For each cell
870 bool is_inconsistent = false;
871 Star::const_iterator it_inc_simplex = m_stars[idx].begin();
872 Star::const_iterator it_inc_simplex_end = m_stars[idx].end();
873 for (; it_inc_simplex != it_inc_simplex_end; ++it_inc_simplex) {
874 // Skip infinite cells
875 if (is_infinite(*it_inc_simplex)) continue;
876
877 Simplex c = *it_inc_simplex;
878 c.insert(idx); // Add the missing index
879
880 sc.add_simplex(c);
881
882 // If we do not already know this star is inconsistent, test it
883 if (!is_inconsistent && !is_simplex_consistent(c)) is_inconsistent = true;
884 }
885
886 if (is_inconsistent) {
887 // Export star to OFF file
888 std::stringstream output_filename;
889 output_filename << filename_base << "_" << idx << ".off";
890 std::ofstream off_stream(output_filename.str().c_str());
891 export_to_off(sc, off_stream);
892 }
893 }
894 }
895
896 class Compare_distance_to_ref_point {
897 public:
898 Compare_distance_to_ref_point(Point const &ref, K const &k) : m_ref(ref), m_k(k) {}
899
900 bool operator()(Point const &p1, Point const &p2) {
901 typename K::Squared_distance_d sqdist = m_k.squared_distance_d_object();
902 return sqdist(p1, m_ref) < sqdist(p2, m_ref);
903 }
904
905 private:
906 Point const &m_ref;
907 K const &m_k;
908 };
909
910#ifdef GUDHI_USE_TBB
911 // Functor for compute_tangential_complex function
912 class Compute_tangent_triangulation {
913 Tangential_complex &m_tc;
914
915 public:
916 // Constructor
917 Compute_tangent_triangulation(Tangential_complex &tc) : m_tc(tc) {}
918
919 // Constructor
920 Compute_tangent_triangulation(const Compute_tangent_triangulation &ctt) : m_tc(ctt.m_tc) {}
921
922 // operator()
923 void operator()(const tbb::blocked_range<size_t> &r) const {
924 for (size_t i = r.begin(); i != r.end(); ++i) m_tc.compute_tangent_triangulation(i);
925 }
926 };
927
928 // Functor for refresh_tangential_complex function
929 class Refresh_tangent_triangulation {
930 Tangential_complex &m_tc;
931 Points_ds const &m_updated_pts_ds;
932
933 public:
934 // Constructor
935 Refresh_tangent_triangulation(Tangential_complex &tc, Points_ds const &updated_pts_ds)
936 : m_tc(tc), m_updated_pts_ds(updated_pts_ds) {}
937
938 // Constructor
939 Refresh_tangent_triangulation(const Refresh_tangent_triangulation &ctt)
940 : m_tc(ctt.m_tc), m_updated_pts_ds(ctt.m_updated_pts_ds) {}
941
942 // operator()
943 void operator()(const tbb::blocked_range<size_t> &r) const {
944 for (size_t i = r.begin(); i != r.end(); ++i) m_tc.refresh_tangent_triangulation(i, m_updated_pts_ds);
945 }
946 };
947#endif // GUDHI_USE_TBB
948
949 bool is_infinite(Simplex const &s) const { return *s.rbegin() == (std::numeric_limits<std::size_t>::max)(); }
950
951 // Output: "triangulation" is a Regular Triangulation containing at least the
952 // star of "center_pt"
953 // Returns the handle of the center vertex
954 Tr_vertex_handle compute_star(std::size_t i, const Point &center_pt, const Tangent_space_basis &tsb,
955 Triangulation &triangulation, bool verbose = false) {
956 int tangent_space_dim = tsb.dimension();
957 const Tr_traits &local_tr_traits = triangulation.geom_traits();
958
959 // Kernel functor & objects
960 typename K::Squared_distance_d k_sqdist = m_k.squared_distance_d_object();
961
962 // Triangulation's traits functor & objects
963 typename Tr_traits::Compute_weight_d point_weight = local_tr_traits.compute_weight_d_object();
964#if CGAL_VERSION_NR < 1050200000
965 typename Tr_traits::Power_center_d power_center = local_tr_traits.power_center_d_object();
966#else
967 typename Tr_traits::Construct_power_sphere_d power_center = local_tr_traits.construct_power_sphere_d_object();
968#endif
969
970 //***************************************************
971 // Build a minimal triangulation in the tangent space
972 // (we only need the star of p)
973 //***************************************************
974
975 // Insert p
976 Tr_point proj_wp;
977 if (i == tsb.origin()) {
978 // Insert {(0, 0, 0...), m_weights[i]}
979 proj_wp = local_tr_traits.construct_weighted_point_d_object()(
980 local_tr_traits.construct_point_d_object()(tangent_space_dim, CGAL::ORIGIN), m_weights[i]);
981 } else {
982 const Weighted_point &wp = compute_perturbed_weighted_point(i);
983 proj_wp = project_point_and_compute_weight(wp, tsb, local_tr_traits);
984 }
985
986 Tr_vertex_handle center_vertex = triangulation.insert(proj_wp);
987 center_vertex->data() = i;
988 if (verbose) std::cerr << "* Inserted point #" << i << "\n";
989
990#ifdef GUDHI_TC_VERY_VERBOSE
991 std::size_t num_attempts_to_insert_points = 1;
992 std::size_t num_inserted_points = 1;
993#endif
994 // const int NUM_NEIGHBORS = 150;
995 // KNS_range ins_range = m_points_ds.k_nearest_neighbors(center_pt, NUM_NEIGHBORS);
996 INS_range ins_range = m_points_ds.incremental_nearest_neighbors(center_pt);
997
998 // While building the local triangulation, we keep the radius
999 // of the sphere "star sphere" centered at "center_vertex"
1000 // and which contains all the
1001 // circumspheres of the star of "center_vertex"
1002 // If th the m_max_squared_edge_length is set the maximal radius of the "star sphere"
1003 // is at most square root of m_max_squared_edge_length
1004 std::optional<FT> squared_star_sphere_radius_plus_margin = m_max_squared_edge_length;
1005
1006 // Insert points until we find a point which is outside "star sphere"
1007 for (auto nn_it = ins_range.begin(); nn_it != ins_range.end(); ++nn_it) {
1008 std::size_t neighbor_point_idx = nn_it->first;
1009
1010 // ith point = p, which is already inserted
1011 if (neighbor_point_idx != i) {
1012 // No need to lock the Mutex_for_perturb here since this will not be
1013 // called while other threads are perturbing the positions
1014 Point neighbor_pt;
1015 FT neighbor_weight;
1016 compute_perturbed_weighted_point(neighbor_point_idx, neighbor_pt, neighbor_weight);
1017 GUDHI_CHECK(!m_max_squared_edge_length ||
1018 squared_star_sphere_radius_plus_margin.value() <= m_max_squared_edge_length.value(),
1019 std::invalid_argument("Tangential_complex::compute_star - set a bigger value with set_max_squared_edge_length."));
1020 if (squared_star_sphere_radius_plus_margin &&
1021 k_sqdist(center_pt, neighbor_pt) > squared_star_sphere_radius_plus_margin.value()) {
1022 GUDHI_CHECK(triangulation.current_dimension() >= tangent_space_dim,
1023 std::invalid_argument("Tangential_complex::compute_star - Dimension of the star is only " + \
1024 std::to_string(triangulation.current_dimension())));
1025 break;
1026 }
1027
1028 Tr_point proj_pt = project_point_and_compute_weight(neighbor_pt, neighbor_weight, tsb, local_tr_traits);
1029
1030#ifdef GUDHI_TC_VERY_VERBOSE
1031 ++num_attempts_to_insert_points;
1032#endif
1033
1034 Tr_vertex_handle vh = triangulation.insert_if_in_star(proj_pt, center_vertex);
1035 // Tr_vertex_handle vh = triangulation.insert(proj_pt);
1036 if (vh != Tr_vertex_handle() && vh->data() == (std::numeric_limits<std::size_t>::max)()) {
1037#ifdef GUDHI_TC_VERY_VERBOSE
1038 ++num_inserted_points;
1039#endif
1040 if (verbose) std::cerr << "* Inserted point #" << neighbor_point_idx << "\n";
1041
1042 vh->data() = neighbor_point_idx;
1043
1044 // Let's recompute squared_star_sphere_radius_plus_margin
1045 if (triangulation.current_dimension() >= tangent_space_dim) {
1046 squared_star_sphere_radius_plus_margin = std::nullopt;
1047 // Get the incident cells and look for the biggest circumsphere
1048 std::vector<Tr_full_cell_handle> incident_cells;
1049 triangulation.incident_full_cells(center_vertex, std::back_inserter(incident_cells));
1050 for (typename std::vector<Tr_full_cell_handle>::iterator cit = incident_cells.begin();
1051 cit != incident_cells.end(); ++cit) {
1052 Tr_full_cell_handle cell = *cit;
1053 if (triangulation.is_infinite(cell)) {
1054 squared_star_sphere_radius_plus_margin = std::nullopt;
1055 break;
1056 } else {
1057 // Note that this uses the perturbed point since it uses
1058 // the points of the local triangulation
1059 Tr_point c =
1060 power_center(boost::make_transform_iterator(cell->vertices_begin(),
1061 vertex_handle_to_point<Tr_point, Tr_vertex_handle>),
1062 boost::make_transform_iterator(cell->vertices_end(),
1063 vertex_handle_to_point<Tr_point, Tr_vertex_handle>));
1064
1065 FT sq_power_sphere_diam = 4 * point_weight(c);
1066
1067 if (!squared_star_sphere_radius_plus_margin ||
1068 sq_power_sphere_diam > squared_star_sphere_radius_plus_margin.value()) {
1069 squared_star_sphere_radius_plus_margin = sq_power_sphere_diam;
1070 }
1071 }
1072 }
1073
1074 // Let's add the margin, now
1075 // The value depends on whether we perturb weight or position
1076 if (squared_star_sphere_radius_plus_margin) {
1077 // "2*m_last_max_perturb" because both points can be perturbed
1078 squared_star_sphere_radius_plus_margin =
1079 CGAL::square(std::sqrt(squared_star_sphere_radius_plus_margin.value()) + 2 * m_last_max_perturb);
1080
1081 // Reduce the square radius to m_max_squared_edge_length if necessary
1082 if (m_max_squared_edge_length && squared_star_sphere_radius_plus_margin.value() > m_max_squared_edge_length.value()) {
1083 squared_star_sphere_radius_plus_margin = m_max_squared_edge_length.value();
1084 }
1085
1086 // Save it in `m_squared_star_spheres_radii_incl_margin`
1087 m_squared_star_spheres_radii_incl_margin[i] = squared_star_sphere_radius_plus_margin.value();
1088 } else {
1089 if (m_max_squared_edge_length) {
1090 squared_star_sphere_radius_plus_margin = m_max_squared_edge_length.value();
1091 m_squared_star_spheres_radii_incl_margin[i] = m_max_squared_edge_length.value();
1092 } else {
1093 m_squared_star_spheres_radii_incl_margin[i] = FT(-1);
1094 }
1095 }
1096 }
1097 }
1098 }
1099 }
1100
1101 return center_vertex;
1102 }
1103
1104 void refresh_tangent_triangulation(std::size_t i, Points_ds const &updated_pts_ds, bool verbose = false) {
1105 if (verbose) std::cerr << "** Refreshing tangent tri #" << i << " **\n";
1106
1107 if (m_squared_star_spheres_radii_incl_margin[i] == FT(-1)) return compute_tangent_triangulation(i, verbose);
1108
1109 Point center_point = compute_perturbed_point(i);
1110 // Among updated point, what is the closer from our center point?
1111 std::size_t closest_pt_index = updated_pts_ds.k_nearest_neighbors(center_point, 1, false).begin()->first;
1112
1113 typename K::Construct_weighted_point_d k_constr_wp = m_k.construct_weighted_point_d_object();
1114#if CGAL_VERSION_NR < 1050200000
1115 typename K::Power_distance_d k_power_dist = m_k.power_distance_d_object();
1116#else
1117 typename K::Compute_power_product_d k_power_dist = m_k.compute_power_product_d_object();
1118#endif
1119
1120 // Construct a weighted point equivalent to the star sphere
1121 Weighted_point star_sphere = k_constr_wp(compute_perturbed_point(i), m_squared_star_spheres_radii_incl_margin[i]);
1122 Weighted_point closest_updated_point = compute_perturbed_weighted_point(closest_pt_index);
1123
1124 // Is the "closest point" inside our star sphere?
1125 if (k_power_dist(star_sphere, closest_updated_point) <= FT(0)) compute_tangent_triangulation(i, verbose);
1126 }
1127
1128 void compute_tangent_triangulation(std::size_t i, bool verbose = false) {
1129 if (verbose) std::cerr << "** Computing tangent tri #" << i << " **\n";
1130 // std::cerr << "***********************************************\n";
1131
1132 // No need to lock the mutex here since this will not be called while
1133 // other threads are perturbing the positions
1134 const Point center_pt = compute_perturbed_point(i);
1135 Tangent_space_basis &tsb = m_tangent_spaces[i];
1136
1137 // Estimate the tangent space
1138 if (!m_are_tangent_spaces_computed[i]) {
1139#ifdef GUDHI_TC_EXPORT_NORMALS
1140 tsb = compute_tangent_space(center_pt, i, true /*normalize*/, &m_orth_spaces[i]);
1141#else
1142 tsb = compute_tangent_space(center_pt, i);
1143#endif
1144 }
1145
1146#if defined(GUDHI_TC_PROFILING) && defined(GUDHI_TC_VERY_VERBOSE)
1147 Gudhi::Clock t;
1148#endif
1149 int tangent_space_dim = tangent_basis_dim(i);
1150 Triangulation &local_tr = m_triangulations[i].construct_triangulation(tangent_space_dim);
1151
1152 m_triangulations[i].center_vertex() = compute_star(i, center_pt, tsb, local_tr, verbose);
1153
1154#if defined(GUDHI_TC_PROFILING) && defined(GUDHI_TC_VERY_VERBOSE)
1155 t.end();
1156 std::cerr << " - triangulation construction: " << t.num_seconds() << " s.\n";
1157 t.reset();
1158#endif
1159
1160#ifdef GUDHI_TC_VERY_VERBOSE
1161 std::cerr << "Inserted " << num_inserted_points << " points / " << num_attempts_to_insert_points
1162 << " attempts to compute the star\n";
1163#endif
1164
1165 update_star(i);
1166
1167#if defined(GUDHI_TC_PROFILING) && defined(GUDHI_TC_VERY_VERBOSE)
1168 t.end();
1169 std::cerr << " - update_star: " << t.num_seconds() << " s.\n";
1170#endif
1171 }
1172
1173 // Updates m_stars[i] directly from m_triangulations[i]
1174
1175 void update_star(std::size_t i) {
1176 Star &star = m_stars[i];
1177 star.clear();
1178 Triangulation &local_tr = m_triangulations[i].tr();
1179 Tr_vertex_handle center_vertex = m_triangulations[i].center_vertex();
1180 int cur_dim_plus_1 = local_tr.current_dimension() + 1;
1181
1182 std::vector<Tr_full_cell_handle> incident_cells;
1183 local_tr.incident_full_cells(center_vertex, std::back_inserter(incident_cells));
1184
1185 typename std::vector<Tr_full_cell_handle>::const_iterator it_c = incident_cells.begin();
1186 typename std::vector<Tr_full_cell_handle>::const_iterator it_c_end = incident_cells.end();
1187 // For each cell
1188 for (; it_c != it_c_end; ++it_c) {
1189 // Will contain all indices except center_vertex
1190 Incident_simplex incident_simplex;
1191 for (int j = 0; j < cur_dim_plus_1; ++j) {
1192 std::size_t index = (*it_c)->vertex(j)->data();
1193 if (index != i) incident_simplex.insert(index);
1194 }
1195 GUDHI_CHECK(incident_simplex.size() == cur_dim_plus_1 - 1,
1196 std::logic_error("update_star: wrong size of incident simplex"));
1197 star.push_back(incident_simplex);
1198 }
1199 }
1200
1201 // Estimates tangent subspaces using PCA
1202
1203 Tangent_space_basis compute_tangent_space(const Point &p, const std::size_t i, bool normalize_basis = true,
1204 Orthogonal_space_basis *p_orth_space_basis = NULL) {
1205 unsigned int num_pts_for_pca =
1206 (std::min)(static_cast<unsigned int>(std::pow(GUDHI_TC_BASE_VALUE_FOR_PCA, m_intrinsic_dim)),
1207 static_cast<unsigned int>(m_points.size()));
1208
1209 // Kernel functors
1210 typename K::Construct_vector_d constr_vec = m_k.construct_vector_d_object();
1211 typename K::Compute_coordinate_d coord = m_k.compute_coordinate_d_object();
1212
1213#ifdef GUDHI_TC_USE_ANOTHER_POINT_SET_FOR_TANGENT_SPACE_ESTIM
1214 KNS_range kns_range = m_points_ds_for_tse.k_nearest_neighbors(p, num_pts_for_pca, false);
1215 const Points &points_for_pca = m_points_for_tse;
1216#else
1217 KNS_range kns_range = m_points_ds.k_nearest_neighbors(p, num_pts_for_pca, false);
1218 const Points &points_for_pca = m_points;
1219#endif
1220
1221 // One row = one point
1222 Eigen::MatrixXd mat_points(num_pts_for_pca, m_ambient_dim);
1223 auto nn_it = kns_range.begin();
1224 for (unsigned int j = 0; j < num_pts_for_pca && nn_it != kns_range.end(); ++j, ++nn_it) {
1225 for (int i = 0; i < m_ambient_dim; ++i) {
1226 mat_points(j, i) = CGAL::to_double(coord(points_for_pca[nn_it->first], i));
1227 }
1228 }
1229 Eigen::MatrixXd centered = mat_points.rowwise() - mat_points.colwise().mean();
1230 Eigen::MatrixXd cov = centered.adjoint() * centered;
1231 Eigen::SelfAdjointEigenSolver<Eigen::MatrixXd> eig(cov);
1232
1233 Tangent_space_basis tsb(i); // p = compute_perturbed_point(i) here
1234
1235 // The eigenvectors are sorted in increasing order of their corresponding
1236 // eigenvalues
1237 for (int j = m_ambient_dim - 1; j >= m_ambient_dim - m_intrinsic_dim; --j) {
1238 if (normalize_basis) {
1239 Vector v = constr_vec(m_ambient_dim, eig.eigenvectors().col(j).data(),
1240 eig.eigenvectors().col(j).data() + m_ambient_dim);
1241 tsb.push_back(normalize_vector(v, m_k));
1242 } else {
1243 tsb.push_back(constr_vec(m_ambient_dim, eig.eigenvectors().col(j).data(),
1244 eig.eigenvectors().col(j).data() + m_ambient_dim));
1245 }
1246 }
1247
1248 if (p_orth_space_basis) {
1249 p_orth_space_basis->set_origin(i);
1250 for (int j = m_ambient_dim - m_intrinsic_dim - 1; j >= 0; --j) {
1251 if (normalize_basis) {
1252 Vector v = constr_vec(m_ambient_dim, eig.eigenvectors().col(j).data(),
1253 eig.eigenvectors().col(j).data() + m_ambient_dim);
1254 p_orth_space_basis->push_back(normalize_vector(v, m_k));
1255 } else {
1256 p_orth_space_basis->push_back(constr_vec(m_ambient_dim, eig.eigenvectors().col(j).data(),
1257 eig.eigenvectors().col(j).data() + m_ambient_dim));
1258 }
1259 }
1260 }
1261
1262 m_are_tangent_spaces_computed[i] = true;
1263
1264 return tsb;
1265 }
1266
1267 // Compute the space tangent to a simplex (p1, p2, ... pn)
1268 // TODO(CJ): Improve this?
1269 // Basically, it takes all the neighbor points to p1, p2... pn and runs PCA
1270 // on it. Note that most points are duplicated.
1271
1272 Tangent_space_basis compute_tangent_space(const Simplex &s, bool normalize_basis = true) {
1273 unsigned int num_pts_for_pca =
1274 (std::min)(static_cast<unsigned int>(std::pow(GUDHI_TC_BASE_VALUE_FOR_PCA, m_intrinsic_dim)),
1275 static_cast<unsigned int>(m_points.size()));
1276
1277 // Kernel functors
1278 typename K::Construct_vector_d constr_vec = m_k.construct_vector_d_object();
1279 typename K::Compute_coordinate_d coord = m_k.compute_coordinate_d_object();
1280 typename K::Squared_length_d sqlen = m_k.squared_length_d_object();
1281 typename K::Scaled_vector_d scaled_vec = m_k.scaled_vector_d_object();
1282 typename K::Scalar_product_d scalar_pdct = m_k.scalar_product_d_object();
1283 typename K::Difference_of_vectors_d diff_vec = m_k.difference_of_vectors_d_object();
1284
1285 // One row = one point
1286 Eigen::MatrixXd mat_points(s.size() * num_pts_for_pca, m_ambient_dim);
1287 unsigned int current_row = 0;
1288
1289 for (Simplex::const_iterator it_index = s.begin(); it_index != s.end(); ++it_index) {
1290 const Point &p = m_points[*it_index];
1291
1292#ifdef GUDHI_TC_USE_ANOTHER_POINT_SET_FOR_TANGENT_SPACE_ESTIM
1293 KNS_range kns_range = m_points_ds_for_tse.k_nearest_neighbors(p, num_pts_for_pca, false);
1294 const Points &points_for_pca = m_points_for_tse;
1295#else
1296 KNS_range kns_range = m_points_ds.k_nearest_neighbors(p, num_pts_for_pca, false);
1297 const Points &points_for_pca = m_points;
1298#endif
1299
1300 auto nn_it = kns_range.begin();
1301 for (; current_row < num_pts_for_pca && nn_it != kns_range.end(); ++current_row, ++nn_it) {
1302 for (int i = 0; i < m_ambient_dim; ++i) {
1303 mat_points(current_row, i) = CGAL::to_double(coord(points_for_pca[nn_it->first], i));
1304 }
1305 }
1306 }
1307 Eigen::MatrixXd centered = mat_points.rowwise() - mat_points.colwise().mean();
1308 Eigen::MatrixXd cov = centered.adjoint() * centered;
1309 Eigen::SelfAdjointEigenSolver<Eigen::MatrixXd> eig(cov);
1310
1311 Tangent_space_basis tsb;
1312
1313 // The eigenvectors are sorted in increasing order of their corresponding
1314 // eigenvalues
1315 for (int j = m_ambient_dim - 1; j >= m_ambient_dim - m_intrinsic_dim; --j) {
1316 if (normalize_basis) {
1317 Vector v = constr_vec(m_ambient_dim, eig.eigenvectors().col(j).data(),
1318 eig.eigenvectors().col(j).data() + m_ambient_dim);
1319 tsb.push_back(normalize_vector(v, m_k));
1320 } else {
1321 tsb.push_back(constr_vec(m_ambient_dim, eig.eigenvectors().col(j).data(),
1322 eig.eigenvectors().col(j).data() + m_ambient_dim));
1323 }
1324 }
1325
1326 return tsb;
1327 }
1328
1329 // Returns the dimension of the ith local triangulation
1330
1331 int tangent_basis_dim(std::size_t i) const { return m_tangent_spaces[i].dimension(); }
1332
1333 Point compute_perturbed_point(std::size_t pt_idx) const {
1334#ifdef GUDHI_TC_PERTURB_POSITION
1335 return m_k.translated_point_d_object()(m_points[pt_idx], m_translations[pt_idx]);
1336#else
1337 return m_points[pt_idx];
1338#endif
1339 }
1340
1341 void compute_perturbed_weighted_point(std::size_t pt_idx, Point &p, FT &w) const {
1342#ifdef GUDHI_TC_PERTURB_POSITION
1343 p = m_k.translated_point_d_object()(m_points[pt_idx], m_translations[pt_idx]);
1344#else
1345 p = m_points[pt_idx];
1346#endif
1347 w = m_weights[pt_idx];
1348 }
1349
1350 Weighted_point compute_perturbed_weighted_point(std::size_t pt_idx) const {
1351 typename K::Construct_weighted_point_d k_constr_wp = m_k.construct_weighted_point_d_object();
1352
1353 Weighted_point wp = k_constr_wp(
1354#ifdef GUDHI_TC_PERTURB_POSITION
1355 m_k.translated_point_d_object()(m_points[pt_idx], m_translations[pt_idx]),
1356#else
1357 m_points[pt_idx],
1358#endif
1359 m_weights[pt_idx]);
1360
1361 return wp;
1362 }
1363
1364 Point unproject_point(const Tr_point &p, const Tangent_space_basis &tsb, const Tr_traits &tr_traits) const {
1365 typename K::Translated_point_d k_transl = m_k.translated_point_d_object();
1366 typename K::Scaled_vector_d k_scaled_vec = m_k.scaled_vector_d_object();
1367 typename Tr_traits::Compute_coordinate_d coord = tr_traits.compute_coordinate_d_object();
1368
1369 Point global_point = compute_perturbed_point(tsb.origin());
1370 for (int i = 0; i < m_intrinsic_dim; ++i) global_point = k_transl(global_point, k_scaled_vec(tsb[i], coord(p, i)));
1371
1372 return global_point;
1373 }
1374
1375 // Project the point in the tangent space
1376 // Resulting point coords are expressed in tsb's space
1377 Tr_bare_point project_point(const Point &p, const Tangent_space_basis &tsb, const Tr_traits &tr_traits) const {
1378 typename K::Scalar_product_d scalar_pdct = m_k.scalar_product_d_object();
1379 typename K::Difference_of_points_d diff_points = m_k.difference_of_points_d_object();
1380
1381 Vector v = diff_points(p, compute_perturbed_point(tsb.origin()));
1382
1383 std::vector<FT> coords;
1384 // Ambiant-space coords of the projected point
1385 coords.reserve(tsb.dimension());
1386 for (std::size_t i = 0; i < m_intrinsic_dim; ++i) {
1387 // Local coords are given by the scalar product with the vectors of tsb
1388 FT coord = scalar_pdct(v, tsb[i]);
1389 coords.push_back(coord);
1390 }
1391
1392 return tr_traits.construct_point_d_object()(static_cast<int>(coords.size()), coords.begin(), coords.end());
1393 }
1394
1395 // Project the point in the tangent space
1396 // The weight will be the squared distance between p and the projection of p
1397 // Resulting point coords are expressed in tsb's space
1398
1399 Tr_point project_point_and_compute_weight(const Weighted_point &wp, const Tangent_space_basis &tsb,
1400 const Tr_traits &tr_traits) const {
1401 typename K::Point_drop_weight_d k_drop_w = m_k.point_drop_weight_d_object();
1402 typename K::Compute_weight_d k_point_weight = m_k.compute_weight_d_object();
1403 return project_point_and_compute_weight(k_drop_w(wp), k_point_weight(wp), tsb, tr_traits);
1404 }
1405
1406 // Same as above, with slightly different parameters
1407 Tr_point project_point_and_compute_weight(const Point &p, const FT w, const Tangent_space_basis &tsb,
1408 const Tr_traits &tr_traits) const {
1409 const int point_dim = m_k.point_dimension_d_object()(p);
1410
1411 typename K::Construct_point_d constr_pt = m_k.construct_point_d_object();
1412 typename K::Scalar_product_d scalar_pdct = m_k.scalar_product_d_object();
1413 typename K::Difference_of_points_d diff_points = m_k.difference_of_points_d_object();
1414 typename K::Compute_coordinate_d coord = m_k.compute_coordinate_d_object();
1415 typename K::Construct_cartesian_const_iterator_d ccci = m_k.construct_cartesian_const_iterator_d_object();
1416
1417 Point origin = compute_perturbed_point(tsb.origin());
1418 Vector v = diff_points(p, origin);
1419
1420 // Same dimension? Then the weight is 0
1421 bool same_dim = (point_dim == tsb.dimension());
1422
1423 std::vector<FT> coords;
1424 // Ambiant-space coords of the projected point
1425 std::vector<FT> p_proj(ccci(origin), ccci(origin, 0));
1426 coords.reserve(tsb.dimension());
1427 for (int i = 0; i < tsb.dimension(); ++i) {
1428 // Local coords are given by the scalar product with the vectors of tsb
1429 FT c = scalar_pdct(v, tsb[i]);
1430 coords.push_back(c);
1431
1432 // p_proj += c * tsb[i]
1433 if (!same_dim) {
1434 for (int j = 0; j < point_dim; ++j) p_proj[j] += c * coord(tsb[i], j);
1435 }
1436 }
1437
1438 // Same dimension? Then the weight is 0
1439 FT sq_dist_to_proj_pt = 0;
1440 if (!same_dim) {
1441 Point projected_pt = constr_pt(point_dim, p_proj.begin(), p_proj.end());
1442 sq_dist_to_proj_pt = m_k.squared_distance_d_object()(p, projected_pt);
1443 }
1444
1445 return tr_traits.construct_weighted_point_d_object()(
1446 tr_traits.construct_point_d_object()(static_cast<int>(coords.size()), coords.begin(), coords.end()),
1447 w - sq_dist_to_proj_pt);
1448 }
1449
1450 // Project all the points in the tangent space
1451
1452 template <typename Indexed_point_range>
1453 std::vector<Tr_point> project_points_and_compute_weights(const Indexed_point_range &point_indices,
1454 const Tangent_space_basis &tsb,
1455 const Tr_traits &tr_traits) const {
1456 std::vector<Tr_point> ret;
1457 for (typename Indexed_point_range::const_iterator it = point_indices.begin(), it_end = point_indices.end();
1458 it != it_end; ++it) {
1459 ret.push_back(project_point_and_compute_weight(compute_perturbed_weighted_point(*it), tsb, tr_traits));
1460 }
1461 return ret;
1462 }
1463
1464 // A simplex here is a local tri's full cell handle
1465
1466 bool is_simplex_consistent(Tr_full_cell_handle fch, int cur_dim) const {
1467 Simplex c;
1468 for (int i = 0; i < cur_dim + 1; ++i) {
1469 std::size_t data = fch->vertex(i)->data();
1470 c.insert(data);
1471 }
1472 return is_simplex_consistent(c);
1473 }
1474
1475 // A simplex here is a list of point indices
1476 // TODO(CJ): improve it like the other "is_simplex_consistent" below
1477
1478 bool is_simplex_consistent(Simplex const &simplex) const {
1479 // Check if the simplex is in the stars of all its vertices
1480 Simplex::const_iterator it_point_idx = simplex.begin();
1481 // For each point p of the simplex, we parse the incidents cells of p
1482 // and we check if "simplex" is among them
1483 for (; it_point_idx != simplex.end(); ++it_point_idx) {
1484 std::size_t point_idx = *it_point_idx;
1485 // Don't check infinite simplices
1486 if (point_idx == (std::numeric_limits<std::size_t>::max)()) continue;
1487
1488 Star const &star = m_stars[point_idx];
1489
1490 // What we're looking for is "simplex" \ point_idx
1491 Incident_simplex is_to_find = simplex;
1492 is_to_find.erase(point_idx);
1493
1494 // For each cell
1495 if (std::find(star.begin(), star.end(), is_to_find) == star.end()) return false;
1496 }
1497
1498 return true;
1499 }
1500
1501 // A simplex here is a list of point indices
1502 // "s" contains all the points of the simplex except "center_point"
1503 // This function returns the points whose star doesn't contain the simplex
1504 // N.B.: the function assumes that the simplex is contained in
1505 // star(center_point)
1506
1507 template <typename OutputIterator> // value_type = std::size_t
1508 bool is_simplex_consistent(std::size_t center_point,
1509 Incident_simplex const &s, // without "center_point"
1510 OutputIterator points_whose_star_does_not_contain_s,
1511 bool check_also_in_non_maximal_faces = false) const {
1512 Simplex full_simplex = s;
1513 full_simplex.insert(center_point);
1514
1515 // Check if the simplex is in the stars of all its vertices
1516 Incident_simplex::const_iterator it_point_idx = s.begin();
1517 // For each point p of the simplex, we parse the incidents cells of p
1518 // and we check if "simplex" is among them
1519 for (; it_point_idx != s.end(); ++it_point_idx) {
1520 std::size_t point_idx = *it_point_idx;
1521 // Don't check infinite simplices
1522 if (point_idx == (std::numeric_limits<std::size_t>::max)()) continue;
1523
1524 Star const &star = m_stars[point_idx];
1525
1526 // What we're looking for is full_simplex \ point_idx
1527 Incident_simplex is_to_find = full_simplex;
1528 is_to_find.erase(point_idx);
1529
1530 if (check_also_in_non_maximal_faces) {
1531 // For each simplex "is" of the star, check if ic_to_simplex is
1532 // included in "is"
1533 bool found = false;
1534 for (Star::const_iterator is = star.begin(), is_end = star.end(); !found && is != is_end; ++is) {
1535 if (std::includes(is->begin(), is->end(), is_to_find.begin(), is_to_find.end())) found = true;
1536 }
1537
1538 if (!found) *points_whose_star_does_not_contain_s++ = point_idx;
1539 } else {
1540 // Does the star contain is_to_find?
1541 if (std::find(star.begin(), star.end(), is_to_find) == star.end())
1542 *points_whose_star_does_not_contain_s++ = point_idx;
1543 }
1544 }
1545
1546 return true;
1547 }
1548
1549 // A simplex here is a list of point indices
1550 // It looks for s in star(p).
1551 // "s" contains all the points of the simplex except p.
1552 bool is_simplex_in_star(std::size_t p, Incident_simplex const &s, bool check_also_in_non_maximal_faces = true) const {
1553 Star const &star = m_stars[p];
1554
1555 if (check_also_in_non_maximal_faces) {
1556 // For each simplex "is" of the star, check if ic_to_simplex is
1557 // included in "is"
1558 bool found = false;
1559 for (Star::const_iterator is = star.begin(), is_end = star.end(); !found && is != is_end; ++is) {
1560 if (std::includes(is->begin(), is->end(), s.begin(), s.end())) found = true;
1561 }
1562
1563 return found;
1564 } else {
1565 return !(std::find(star.begin(), star.end(), s) == star.end());
1566 }
1567 }
1568
1569#ifdef GUDHI_USE_TBB
1570 // Functor for try_to_solve_inconsistencies_in_a_local_triangulation function
1571 class Try_to_solve_inconsistencies_in_a_local_triangulation {
1572 Tangential_complex &m_tc;
1573 double m_max_perturb;
1574 tbb::combinable<std::size_t> &m_num_inconsistencies;
1575 tbb::combinable<std::vector<std::size_t> > &m_updated_points;
1576
1577 public:
1578 // Constructor
1579 Try_to_solve_inconsistencies_in_a_local_triangulation(Tangential_complex &tc, double max_perturb,
1580 tbb::combinable<std::size_t> &num_inconsistencies,
1581 tbb::combinable<std::vector<std::size_t> > &updated_points)
1582 : m_tc(tc),
1583 m_max_perturb(max_perturb),
1584 m_num_inconsistencies(num_inconsistencies),
1585 m_updated_points(updated_points) {}
1586
1587 // Constructor
1588 Try_to_solve_inconsistencies_in_a_local_triangulation(
1589 const Try_to_solve_inconsistencies_in_a_local_triangulation &tsilt)
1590 : m_tc(tsilt.m_tc),
1591 m_max_perturb(tsilt.m_max_perturb),
1592 m_num_inconsistencies(tsilt.m_num_inconsistencies),
1593 m_updated_points(tsilt.m_updated_points) {}
1594
1595 // operator()
1596 void operator()(const tbb::blocked_range<size_t> &r) const {
1597 for (size_t i = r.begin(); i != r.end(); ++i) {
1598 m_num_inconsistencies.local() += m_tc.try_to_solve_inconsistencies_in_a_local_triangulation(
1599 i, m_max_perturb, std::back_inserter(m_updated_points.local()));
1600 }
1601 }
1602 };
1603#endif // GUDHI_USE_TBB
1604
1605 void perturb(std::size_t point_idx, double max_perturb) {
1606 const Tr_traits &local_tr_traits = m_triangulations[point_idx].tr().geom_traits();
1607 typename Tr_traits::Compute_coordinate_d coord = local_tr_traits.compute_coordinate_d_object();
1608 typename K::Translated_point_d k_transl = m_k.translated_point_d_object();
1609 typename K::Construct_vector_d k_constr_vec = m_k.construct_vector_d_object();
1610 typename K::Scaled_vector_d k_scaled_vec = m_k.scaled_vector_d_object();
1611
1612 CGAL::Random_points_in_ball_d<Tr_bare_point> tr_point_in_ball_generator(
1613 m_intrinsic_dim, m_random_generator.get_double(0., max_perturb));
1614
1615 Tr_point local_random_transl =
1616 local_tr_traits.construct_weighted_point_d_object()(*tr_point_in_ball_generator++, 0);
1617 Translation_for_perturb global_transl = k_constr_vec(m_ambient_dim);
1618 const Tangent_space_basis &tsb = m_tangent_spaces[point_idx];
1619 for (int i = 0; i < m_intrinsic_dim; ++i) {
1620 global_transl = k_transl(global_transl, k_scaled_vec(tsb[i], coord(local_random_transl, i)));
1621 }
1622 // Parallel
1623#if defined(GUDHI_USE_TBB)
1624 m_p_perturb_mutexes[point_idx].lock();
1625 m_translations[point_idx] = global_transl;
1626 m_p_perturb_mutexes[point_idx].unlock();
1627 // Sequential
1628#else
1629 m_translations[point_idx] = global_transl;
1630#endif
1631 }
1632
1633 // Return true if inconsistencies were found
1634 template <typename OutputIt>
1635 bool try_to_solve_inconsistencies_in_a_local_triangulation(
1636 std::size_t tr_index, double max_perturb, OutputIt perturbed_pts_indices = CGAL::Emptyset_iterator()) {
1637 bool is_inconsistent = false;
1638
1639 Star const &star = m_stars[tr_index];
1640
1641 // For each incident simplex
1642 Star::const_iterator it_inc_simplex = star.begin();
1643 Star::const_iterator it_inc_simplex_end = star.end();
1644 for (; it_inc_simplex != it_inc_simplex_end; ++it_inc_simplex) {
1645 const Incident_simplex &incident_simplex = *it_inc_simplex;
1646
1647 // Don't check infinite cells
1648 if (is_infinite(incident_simplex)) continue;
1649
1650 Simplex c = incident_simplex;
1651 c.insert(tr_index); // Add the missing index
1652
1653 // Perturb the center point
1654 if (!is_simplex_consistent(c)) {
1655 is_inconsistent = true;
1656
1657 std::size_t idx = tr_index;
1658
1659 perturb(tr_index, max_perturb);
1660 *perturbed_pts_indices++ = idx;
1661
1662 // We will try the other cells next time
1663 break;
1664 }
1665 }
1666
1667 return is_inconsistent;
1668 }
1669
1670 // 1st line: number of points
1671 // Then one point per line
1672 std::ostream &export_point_set(std::ostream &os, bool use_perturbed_points = false,
1673 const char *coord_separator = " ") const {
1674 if (use_perturbed_points) {
1675 std::vector<Point> perturbed_points;
1676 perturbed_points.reserve(m_points.size());
1677 for (std::size_t i = 0; i < m_points.size(); ++i) perturbed_points.push_back(compute_perturbed_point(i));
1678
1679 return export_point_set(m_k, perturbed_points, os, coord_separator);
1680 } else {
1681 return export_point_set(m_k, m_points, os, coord_separator);
1682 }
1683 }
1684
1685 template <typename ProjectionFunctor = CGAL::Identity<Point> >
1686 std::ostream &export_vertices_to_off(std::ostream &os, std::size_t &num_vertices, bool use_perturbed_points = false,
1687 ProjectionFunctor const &point_projection = ProjectionFunctor()) const {
1688 if (m_points.empty()) {
1689 num_vertices = 0;
1690 return os;
1691 }
1692
1693 // If m_intrinsic_dim = 1, we output each point two times
1694 // to be able to export each segment as a flat triangle with 3 different
1695 // indices (otherwise, Meshlab detects degenerated simplices)
1696 const int N = (m_intrinsic_dim == 1 ? 2 : 1);
1697
1698 // Kernel functors
1699 typename K::Compute_coordinate_d coord = m_k.compute_coordinate_d_object();
1700
1701#ifdef GUDHI_TC_EXPORT_ALL_COORDS_IN_OFF
1702 int num_coords = m_ambient_dim;
1703#else
1704 int num_coords = (std::min)(m_ambient_dim, 3);
1705#endif
1706
1707#ifdef GUDHI_TC_EXPORT_NORMALS
1708 OS_container::const_iterator it_os = m_orth_spaces.begin();
1709#endif
1710 typename Points::const_iterator it_p = m_points.begin();
1711 typename Points::const_iterator it_p_end = m_points.end();
1712 // For each point p
1713 for (std::size_t i = 0; it_p != it_p_end; ++it_p, ++i) {
1714 Point p = point_projection(use_perturbed_points ? compute_perturbed_point(i) : *it_p);
1715 for (int ii = 0; ii < N; ++ii) {
1716 int j = 0;
1717 for (; j < num_coords; ++j) os << CGAL::to_double(coord(p, j)) << " ";
1718 if (j == 2) os << "0";
1719
1720#ifdef GUDHI_TC_EXPORT_NORMALS
1721 for (j = 0; j < num_coords; ++j) os << " " << CGAL::to_double(coord(*it_os->begin(), j));
1722#endif
1723 os << "\n";
1724 }
1725#ifdef GUDHI_TC_EXPORT_NORMALS
1726 ++it_os;
1727#endif
1728 }
1729
1730 num_vertices = N * m_points.size();
1731 return os;
1732 }
1733
1734 std::ostream &export_simplices_to_off(std::ostream &os, std::size_t &num_OFF_simplices,
1735 bool color_inconsistencies = false,
1736 Simplex_set const *p_simpl_to_color_in_red = NULL,
1737 Simplex_set const *p_simpl_to_color_in_green = NULL,
1738 Simplex_set const *p_simpl_to_color_in_blue = NULL) const {
1739 // If m_intrinsic_dim = 1, each point is output two times
1740 // (see export_vertices_to_off)
1741 num_OFF_simplices = 0;
1742 std::size_t num_maximal_simplices = 0;
1743 std::size_t num_inconsistent_maximal_simplices = 0;
1744 std::size_t num_inconsistent_stars = 0;
1745 typename Tr_container::const_iterator it_tr = m_triangulations.begin();
1746 typename Tr_container::const_iterator it_tr_end = m_triangulations.end();
1747 // For each triangulation
1748 for (std::size_t idx = 0; it_tr != it_tr_end; ++it_tr, ++idx) {
1749 bool is_star_inconsistent = false;
1750
1751 Triangulation const &tr = it_tr->tr();
1752
1753 if (tr.current_dimension() < m_intrinsic_dim) continue;
1754
1755 // Color for this star
1756 std::stringstream color;
1757 // color << rand()%256 << " " << 100+rand()%156 << " " << 100+rand()%156;
1758 color << 128 << " " << 128 << " " << 128;
1759
1760 // Gather the triangles here, with an int telling its color
1761 typedef std::vector<std::pair<Simplex, int> > Star_using_triangles;
1762 Star_using_triangles star_using_triangles;
1763
1764 // For each cell of the star
1765 Star::const_iterator it_inc_simplex = m_stars[idx].begin();
1766 Star::const_iterator it_inc_simplex_end = m_stars[idx].end();
1767 for (; it_inc_simplex != it_inc_simplex_end; ++it_inc_simplex) {
1768 Simplex c = *it_inc_simplex;
1769 c.insert(idx);
1770 std::size_t num_vertices = c.size();
1771 ++num_maximal_simplices;
1772
1773 int color_simplex = -1; // -1=no color, 0=yellow, 1=red, 2=green, 3=blue
1774 if (color_inconsistencies && !is_simplex_consistent(c)) {
1775 ++num_inconsistent_maximal_simplices;
1776 color_simplex = 0;
1777 is_star_inconsistent = true;
1778 } else {
1779 if (p_simpl_to_color_in_red && std::find(p_simpl_to_color_in_red->begin(), p_simpl_to_color_in_red->end(),
1780 c) != p_simpl_to_color_in_red->end()) {
1781 color_simplex = 1;
1782 } else if (p_simpl_to_color_in_green &&
1783 std::find(p_simpl_to_color_in_green->begin(), p_simpl_to_color_in_green->end(), c) !=
1784 p_simpl_to_color_in_green->end()) {
1785 color_simplex = 2;
1786 } else if (p_simpl_to_color_in_blue &&
1787 std::find(p_simpl_to_color_in_blue->begin(), p_simpl_to_color_in_blue->end(), c) !=
1788 p_simpl_to_color_in_blue->end()) {
1789 color_simplex = 3;
1790 }
1791 }
1792
1793 // If m_intrinsic_dim = 1, each point is output two times,
1794 // so we need to multiply each index by 2
1795 // And if only 2 vertices, add a third one (each vertex is duplicated in
1796 // the file when m_intrinsic dim = 2)
1797 if (m_intrinsic_dim == 1) {
1798 Simplex tmp_c;
1799 Simplex::iterator it = c.begin();
1800 for (; it != c.end(); ++it) tmp_c.insert(*it * 2);
1801 if (num_vertices == 2) tmp_c.insert(*tmp_c.rbegin() + 1);
1802
1803 c = tmp_c;
1804 }
1805
1806 if (num_vertices <= 3) {
1807 star_using_triangles.push_back(std::make_pair(c, color_simplex));
1808 } else {
1809 // num_vertices >= 4: decompose the simplex in triangles
1810 std::vector<bool> booleans(num_vertices, false);
1811 std::fill(booleans.begin() + num_vertices - 3, booleans.end(), true);
1812 do {
1813 Simplex triangle;
1814 Simplex::iterator it = c.begin();
1815 for (int i = 0; it != c.end(); ++i, ++it) {
1816 if (booleans[i]) triangle.insert(*it);
1817 }
1818 star_using_triangles.push_back(std::make_pair(triangle, color_simplex));
1819 } while (std::next_permutation(booleans.begin(), booleans.end()));
1820 }
1821 }
1822
1823 // For each cell
1824 Star_using_triangles::const_iterator it_simplex = star_using_triangles.begin();
1825 Star_using_triangles::const_iterator it_simplex_end = star_using_triangles.end();
1826 for (; it_simplex != it_simplex_end; ++it_simplex) {
1827 const Simplex &c = it_simplex->first;
1828
1829 // Don't export infinite cells
1830 if (is_infinite(c)) continue;
1831
1832 int color_simplex = it_simplex->second;
1833
1834 std::stringstream sstr_c;
1835
1836 Simplex::const_iterator it_point_idx = c.begin();
1837 for (; it_point_idx != c.end(); ++it_point_idx) {
1838 sstr_c << *it_point_idx << " ";
1839 }
1840
1841 os << 3 << " " << sstr_c.str();
1842 if (color_inconsistencies || p_simpl_to_color_in_red || p_simpl_to_color_in_green || p_simpl_to_color_in_blue) {
1843 switch (color_simplex) {
1844 case 0:
1845 os << " 255 255 0";
1846 break;
1847 case 1:
1848 os << " 255 0 0";
1849 break;
1850 case 2:
1851 os << " 0 255 0";
1852 break;
1853 case 3:
1854 os << " 0 0 255";
1855 break;
1856 default:
1857 os << " " << color.str();
1858 break;
1859 }
1860 }
1861 ++num_OFF_simplices;
1862 os << "\n";
1863 }
1864 if (is_star_inconsistent) ++num_inconsistent_stars;
1865 }
1866
1867#ifdef DEBUG_TRACES
1868 std::cerr << "\n==========================================================\n"
1869 << "Export from list of stars to OFF:\n"
1870 << " * Number of vertices: " << m_points.size() << "\n"
1871 << " * Total number of maximal simplices: " << num_maximal_simplices << "\n";
1872 if (color_inconsistencies) {
1873 std::cerr << " * Number of inconsistent stars: " << num_inconsistent_stars << " ("
1874 << (m_points.size() > 0 ? 100. * num_inconsistent_stars / m_points.size() : 0.) << "%)\n"
1875 << " * Number of inconsistent maximal simplices: " << num_inconsistent_maximal_simplices << " ("
1876 << (num_maximal_simplices > 0 ? 100. * num_inconsistent_maximal_simplices / num_maximal_simplices : 0.)
1877 << "%)\n";
1878 }
1879 std::cerr << "==========================================================\n";
1880#endif
1881
1882 return os;
1883 }
1884
1885 public:
1886 std::ostream &export_simplices_to_off(const Simplicial_complex &complex, std::ostream &os,
1887 std::size_t &num_OFF_simplices,
1888 Simplex_set const *p_simpl_to_color_in_red = NULL,
1889 Simplex_set const *p_simpl_to_color_in_green = NULL,
1890 Simplex_set const *p_simpl_to_color_in_blue = NULL) const {
1891 typedef Simplicial_complex::Simplex Simplex;
1892 typedef Simplicial_complex::Simplex_set Simplex_set;
1893
1894 // If m_intrinsic_dim = 1, each point is output two times
1895 // (see export_vertices_to_off)
1896 num_OFF_simplices = 0;
1897 std::size_t num_maximal_simplices = 0;
1898
1899 typename Simplex_set::const_iterator it_s = complex.simplex_range().begin();
1900 typename Simplex_set::const_iterator it_s_end = complex.simplex_range().end();
1901 // For each simplex
1902 for (; it_s != it_s_end; ++it_s) {
1903 Simplex c = *it_s;
1904 ++num_maximal_simplices;
1905
1906 int color_simplex = -1; // -1=no color, 0=yellow, 1=red, 2=green, 3=blue
1907 if (p_simpl_to_color_in_red && std::find(p_simpl_to_color_in_red->begin(), p_simpl_to_color_in_red->end(), c) !=
1908 p_simpl_to_color_in_red->end()) {
1909 color_simplex = 1;
1910 } else if (p_simpl_to_color_in_green &&
1911 std::find(p_simpl_to_color_in_green->begin(), p_simpl_to_color_in_green->end(), c) !=
1912 p_simpl_to_color_in_green->end()) {
1913 color_simplex = 2;
1914 } else if (p_simpl_to_color_in_blue &&
1915 std::find(p_simpl_to_color_in_blue->begin(), p_simpl_to_color_in_blue->end(), c) !=
1916 p_simpl_to_color_in_blue->end()) {
1917 color_simplex = 3;
1918 }
1919
1920 // Gather the triangles here
1921 typedef std::vector<Simplex> Triangles;
1922 Triangles triangles;
1923
1924 int num_vertices = static_cast<int>(c.size());
1925 // Do not export smaller dimension simplices
1926 if (num_vertices < m_intrinsic_dim + 1) continue;
1927
1928 // If m_intrinsic_dim = 1, each point is output two times,
1929 // so we need to multiply each index by 2
1930 // And if only 2 vertices, add a third one (each vertex is duplicated in
1931 // the file when m_intrinsic dim = 2)
1932 if (m_intrinsic_dim == 1) {
1933 Simplex tmp_c;
1934 Simplex::iterator it = c.begin();
1935 for (; it != c.end(); ++it) tmp_c.insert(*it * 2);
1936 if (num_vertices == 2) tmp_c.insert(*tmp_c.rbegin() + 1);
1937
1938 c = tmp_c;
1939 }
1940
1941 if (num_vertices <= 3) {
1942 triangles.push_back(c);
1943 } else {
1944 // num_vertices >= 4: decompose the simplex in triangles
1945 std::vector<bool> booleans(num_vertices, false);
1946 std::fill(booleans.begin() + num_vertices - 3, booleans.end(), true);
1947 do {
1948 Simplex triangle;
1949 Simplex::iterator it = c.begin();
1950 for (int i = 0; it != c.end(); ++i, ++it) {
1951 if (booleans[i]) triangle.insert(*it);
1952 }
1953 triangles.push_back(triangle);
1954 } while (std::next_permutation(booleans.begin(), booleans.end()));
1955 }
1956
1957 // For each cell
1958 Triangles::const_iterator it_tri = triangles.begin();
1959 Triangles::const_iterator it_tri_end = triangles.end();
1960 for (; it_tri != it_tri_end; ++it_tri) {
1961 // Don't export infinite cells
1962 if (is_infinite(*it_tri)) continue;
1963
1964 os << 3 << " ";
1965 Simplex::const_iterator it_point_idx = it_tri->begin();
1966 for (; it_point_idx != it_tri->end(); ++it_point_idx) {
1967 os << *it_point_idx << " ";
1968 }
1969
1970 if (p_simpl_to_color_in_red || p_simpl_to_color_in_green || p_simpl_to_color_in_blue) {
1971 switch (color_simplex) {
1972 case 0:
1973 os << " 255 255 0";
1974 break;
1975 case 1:
1976 os << " 255 0 0";
1977 break;
1978 case 2:
1979 os << " 0 255 0";
1980 break;
1981 case 3:
1982 os << " 0 0 255";
1983 break;
1984 default:
1985 os << " 128 128 128";
1986 break;
1987 }
1988 }
1989
1990 ++num_OFF_simplices;
1991 os << "\n";
1992 }
1993 }
1994
1995#ifdef DEBUG_TRACES
1996 std::cerr << "\n==========================================================\n"
1997 << "Export from complex to OFF:\n"
1998 << " * Number of vertices: " << m_points.size() << "\n"
1999 << " * Total number of maximal simplices: " << num_maximal_simplices << "\n"
2000 << "==========================================================\n";
2001#endif
2002
2003 return os;
2004 }
2005
2013 void set_max_squared_edge_length(FT max_squared_edge_length) { m_max_squared_edge_length = max_squared_edge_length; }
2014
2015 private:
2016 const K m_k;
2017 const int m_intrinsic_dim;
2018 const int m_ambient_dim;
2019
2020 Points m_points;
2021 Weights m_weights;
2022#ifdef GUDHI_TC_PERTURB_POSITION
2023 Translations_for_perturb m_translations;
2024#if defined(GUDHI_USE_TBB)
2025 Mutex_for_perturb *m_p_perturb_mutexes;
2026#endif
2027#endif
2028
2029 Points_ds m_points_ds;
2030 double m_last_max_perturb;
2031 std::vector<bool> m_are_tangent_spaces_computed;
2032 TS_container m_tangent_spaces;
2033#ifdef GUDHI_TC_EXPORT_NORMALS
2034 OS_container m_orth_spaces;
2035#endif
2036 Tr_container m_triangulations; // Contains the triangulations
2037 // and their center vertex
2038 Stars_container m_stars;
2039 std::vector<FT> m_squared_star_spheres_radii_incl_margin;
2040 std::optional<FT> m_max_squared_edge_length;
2041
2042#ifdef GUDHI_TC_USE_ANOTHER_POINT_SET_FOR_TANGENT_SPACE_ESTIM
2043 Points m_points_for_tse;
2044 Points_ds m_points_ds_for_tse;
2045#endif
2046
2047 mutable CGAL::Random m_random_generator;
2048}; // /class Tangential_complex
2049
2050} // end namespace tangential_complex
2051} // end namespace Gudhi
2052
2053#endif // TANGENTIAL_COMPLEX_H_
Spatial tree data structure to perform (approximate) nearest and furthest neighbor search.
Definition: Kd_tree_search.h:70
K_neighbor_search KNS_range
The range returned by a k-nearest or k-furthest neighbor search. Its value type is std::pair<std::siz...
Definition: Kd_tree_search.h:104
INS_range incremental_nearest_neighbors(Point const &p, FT eps=FT(0)) const
Search incrementally for the nearest neighbors from a query point.
Definition: Kd_tree_search.h:271
KNS_range k_nearest_neighbors(Point const &p, unsigned int k, bool sorted=true, FT eps=FT(0)) const
Search for the k-nearest neighbors from a query point.
Definition: Kd_tree_search.h:245
Incremental_neighbor_search INS_range
The range returned by an incremental nearest or furthest neighbor search. Its value type is std::pair...
Definition: Kd_tree_search.h:112
Tangential complex data structure.
Definition: Tangential_complex.h:126
Tangential_complex(Point_range points, int intrinsic_dimension, const K &k=K())
Constructor from a range of points. Points are copied into the instance, and a search data structure ...
Definition: Tangential_complex.h:241
Fix_inconsistencies_info fix_inconsistencies_using_perturbation(double max_perturb, double time_limit=-1.)
Attempts to fix inconsistencies by perturbing the point positions.
Definition: Tangential_complex.h:397
void compute_tangential_complex()
Computes the tangential complex.
Definition: Tangential_complex.h:331
int intrinsic_dimension() const
Returns the intrinsic dimension of the manifold.
Definition: Tangential_complex.h:280
Point get_perturbed_point(std::size_t vertex) const
Returns the perturbed position of the point corresponding to the vertex given as parameter.
Definition: Tangential_complex.h:299
void set_max_squared_edge_length(FT max_squared_edge_length)
Sets the maximal possible squared edge length for the edges in the triangulations.
Definition: Tangential_complex.h:2013
Point get_point(std::size_t vertex) const
Returns the point corresponding to the vertex given as parameter.
Definition: Tangential_complex.h:292
Num_inconsistencies number_of_inconsistent_simplices(bool verbose=false) const
Definition: Tangential_complex.h:553
int ambient_dimension() const
Returns the ambient dimension.
Definition: Tangential_complex.h:283
std::size_t number_of_vertices() const
Returns the number of vertices.
Definition: Tangential_complex.h:303
int create_complex(Simplex_tree_ &tree, bool export_inconsistent_simplices=true) const
Exports the complex into a Simplex_tree.
Definition: Tangential_complex.h:612
~Tangential_complex()
Destructor.
Definition: Tangential_complex.h:273
Type returned by Tangential_complex::fix_inconsistencies_using_perturbation.
Definition: Tangential_complex.h:379
std::size_t final_num_inconsistent_stars
final number of inconsistent stars
Definition: Tangential_complex.h:389
bool success
true if all inconsistencies could be removed, false if the time limit has been reached before
Definition: Tangential_complex.h:381
unsigned int num_steps
number of steps performed
Definition: Tangential_complex.h:383
std::size_t best_num_inconsistent_stars
best number of inconsistent stars during the process
Definition: Tangential_complex.h:387
std::size_t initial_num_inconsistent_stars
initial number of inconsistent stars
Definition: Tangential_complex.h:385
Type returned by Tangential_complex::number_of_inconsistent_simplices.
Definition: Tangential_complex.h:541
std::size_t num_inconsistent_simplices
Number of inconsistent simplices.
Definition: Tangential_complex.h:545
std::size_t num_inconsistent_stars
Number of stars containing at least one inconsistent simplex.
Definition: Tangential_complex.h:547
std::size_t num_simplices
Total number of simplices in stars (including duplicates that appear in several stars)
Definition: Tangential_complex.h:543