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
449 tbb::parallel_for(tbb::blocked_range<size_t>(0, m_triangulations.size()),
450 [&]( tbb::blocked_range<size_t> r ) {
451 for (size_t i = r.begin(); i != r.end(); ++i) {
452 num_inconsistencies.local() += try_to_solve_inconsistencies_in_a_local_triangulation(
453 i, max_perturb, std::back_inserter(tls_updated_points.local()));
454 }
455 });
456
457 num_inconsistent_stars = num_inconsistencies.combine(std::plus<std::size_t>());
458 updated_points =
459 tls_updated_points.combine([](std::vector<std::size_t> const &x, std::vector<std::size_t> const &y) {
460 std::vector<std::size_t> res;
461 res.reserve(x.size() + y.size());
462 res.insert(res.end(), x.begin(), x.end());
463 res.insert(res.end(), y.begin(), y.end());
464 return res;
465 });
466 } else {
467#endif // GUDHI_USE_TBB
468 // Sequential
469 for (std::size_t i = 0; i < m_triangulations.size(); ++i) {
470 num_inconsistent_stars +=
471 try_to_solve_inconsistencies_in_a_local_triangulation(i, max_perturb, std::back_inserter(updated_points));
472 }
473#if defined(GUDHI_USE_TBB)
474 }
475#endif // GUDHI_USE_TBB
476
477#ifdef GUDHI_TC_PROFILING
478 t_fix_step.end();
479#endif
480
481#if defined(GUDHI_TC_SHOW_DETAILED_STATS_FOR_INCONSISTENCIES) || defined(DEBUG_TRACES)
482 std::cerr << "\nEncountered during fix:\n"
483 << " * Num stars containing inconsistent simplices: " << red << num_inconsistent_stars << white << " ("
484 << 100. * num_inconsistent_stars / m_points.size() << "%)\n";
485#endif
486
487#ifdef GUDHI_TC_PROFILING
488 std::cerr << yellow << "done in " << t_fix_step.num_seconds() << " seconds.\n" << white;
489#elif defined(DEBUG_TRACES)
490 std::cerr << yellow << "done.\n" << white;
491#endif
492
493 if (num_inconsistent_stars > 0) refresh_tangential_complex(updated_points);
494
495#ifdef GUDHI_TC_PERFORM_EXTRA_CHECKS
496 // Confirm that all stars were actually refreshed
497 std::size_t num_inc_1 = std::get<1>(number_of_inconsistent_simplices(false));
498 refresh_tangential_complex();
499 std::size_t num_inc_2 = std::get<1>(number_of_inconsistent_simplices(false));
500 if (num_inc_1 != num_inc_2)
501 std::cerr << red << "REFRESHMENT CHECK: FAILED. (" << num_inc_1 << " vs " << num_inc_2 << ")\n" << white;
502 else
503 std::cerr << green << "REFRESHMENT CHECK: PASSED.\n" << white;
504#endif
505
506#ifdef GUDHI_TC_SHOW_DETAILED_STATS_FOR_INCONSISTENCIES
507 std::tuple<std::size_t, std::size_t, std::size_t> stats_after = number_of_inconsistent_simplices(false);
508
509 std::cerr << "\nAfter fix:\n"
510 << " * Total number of simplices in stars (incl. duplicates): " << std::get<0>(stats_after) << "\n"
511 << " * Num inconsistent simplices in stars (incl. duplicates): " << red << std::get<1>(stats_after)
512 << white << " (" << 100. * std::get<1>(stats_after) / std::get<0>(stats_after) << "%)\n"
513 << " * Number of stars containing inconsistent simplices: " << red << std::get<2>(stats_after) << white
514 << " (" << 100. * std::get<2>(stats_after) / m_points.size() << "%)\n";
515
516 stats_before = stats_after;
517#endif
518
519 if (info.num_steps == 0) info.initial_num_inconsistent_stars = num_inconsistent_stars;
520
521 if (num_inconsistent_stars < info.best_num_inconsistent_stars)
522 info.best_num_inconsistent_stars = num_inconsistent_stars;
523
524 info.final_num_inconsistent_stars = num_inconsistent_stars;
525
526 done = (num_inconsistent_stars == 0);
527 if (!done) {
528 ++info.num_steps;
529 if (time_limit > 0. && t.num_seconds() > time_limit) {
530#ifdef DEBUG_TRACES
531 std::cerr << red << "Time limit reached.\n" << white;
532#endif
533 info.success = false;
534 return info;
535 }
536 }
537 }
538
539#ifdef DEBUG_TRACES
540 std::cerr << green << "Fixed!\n" << white;
541#endif
542 info.success = true;
543 return info;
544 }
545
549 std::size_t num_simplices = 0;
551 std::size_t num_inconsistent_simplices = 0;
553 std::size_t num_inconsistent_stars = 0;
554 };
555
558
560#ifdef DEBUG_TRACES
561 bool verbose = true
562#else
563 bool verbose = false
564#endif
565 ) const {
567
568 // For each triangulation
569 for (std::size_t idx = 0; idx < m_points.size(); ++idx) {
570 bool is_star_inconsistent = false;
571
572 // For each cell
573 Star::const_iterator it_inc_simplex = m_stars[idx].begin();
574 Star::const_iterator it_inc_simplex_end = m_stars[idx].end();
575 for (; it_inc_simplex != it_inc_simplex_end; ++it_inc_simplex) {
576 // Don't check infinite cells
577 if (is_infinite(*it_inc_simplex)) continue;
578
579 Simplex c = *it_inc_simplex;
580 c.insert(idx); // Add the missing index
581
582 if (!is_simplex_consistent(c)) {
584 is_star_inconsistent = true;
585 }
586
587 ++stats.num_simplices;
588 }
589 stats.num_inconsistent_stars += is_star_inconsistent;
590 }
591
592 if (verbose) {
593 std::cerr << "\n==========================================================\n"
594 << "Inconsistencies:\n"
595 << " * Total number of simplices in stars (incl. duplicates): " << stats.num_simplices << "\n"
596 << " * Number of inconsistent simplices in stars (incl. duplicates): "
597 << stats.num_inconsistent_simplices << " ("
598 << 100. * stats.num_inconsistent_simplices / stats.num_simplices << "%)\n"
599 << " * Number of stars containing inconsistent simplices: " << stats.num_inconsistent_stars << " ("
600 << 100. * stats.num_inconsistent_stars / m_points.size() << "%)\n"
601 << "==========================================================\n";
602 }
603
604 return stats;
605 }
606
617 template <typename Simplex_tree_>
618 int create_complex(Simplex_tree_ &tree,
619 bool export_inconsistent_simplices = true
621 ,
622 bool export_infinite_simplices = false, Simplex_set *p_inconsistent_simplices = NULL
624 ) const {
625#if defined(DEBUG_TRACES) || defined(GUDHI_TC_PROFILING)
626 std::cerr << yellow << "\nExporting the TC as a Simplex_tree... " << white;
627#endif
628#ifdef GUDHI_TC_PROFILING
629 Gudhi::Clock t;
630#endif
631
632 int max_dim = -1;
633
634 // Ordered vertices to be inserted first by the create_complex method to avoid quadratic complexity.
635 std::vector<typename Simplex_tree_::Vertex_handle> vertices(m_points.size());
636 std::iota(vertices.begin(), vertices.end(), 0);
637 tree.insert_batch_vertices(vertices);
638
639 // For each triangulation
640 for (std::size_t idx = 0; idx < m_points.size(); ++idx) {
641 // For each cell of the star
642 Star::const_iterator it_inc_simplex = m_stars[idx].begin();
643 Star::const_iterator it_inc_simplex_end = m_stars[idx].end();
644 for (; it_inc_simplex != it_inc_simplex_end; ++it_inc_simplex) {
645 Simplex c = *it_inc_simplex;
646
647 // Don't export infinite cells
648 if (!export_infinite_simplices && is_infinite(c)) continue;
649
650 if (static_cast<int>(c.size()) > max_dim) max_dim = static_cast<int>(c.size());
651 // Add the missing center vertex
652 c.insert(idx);
653
654 if (!export_inconsistent_simplices && !is_simplex_consistent(c)) continue;
655
656 // Try to insert the simplex
657 bool inserted = tree.insert_simplex_and_subfaces(c).second;
658
659 // Inconsistent?
660 if (p_inconsistent_simplices && inserted && !is_simplex_consistent(c)) {
661 p_inconsistent_simplices->insert(c);
662 }
663 }
664 }
665
666#ifdef GUDHI_TC_PROFILING
667 t.end();
668 std::cerr << yellow << "done in " << t.num_seconds() << " seconds.\n" << white;
669#elif defined(DEBUG_TRACES)
670 std::cerr << yellow << "done.\n" << white;
671#endif
672
673 return max_dim;
674 }
675
676 // First clears the complex then exports the TC into it
677 // Returns the max dimension of the simplices
678 // check_lower_and_higher_dim_simplices : 0 (false), 1 (true), 2 (auto)
679 // If the check is enabled, the function:
680 // - won't insert the simplex if it is already in a higher dim simplex
681 // - will erase any lower-dim simplices that are faces of the new simplex
682 // "auto" (= 2) will enable the check as a soon as it encounters a
683 // simplex whose dimension is different from the previous ones.
684 // N.B.: The check is quite expensive.
685
686 int create_complex(Simplicial_complex &complex, bool export_inconsistent_simplices = true,
687 bool export_infinite_simplices = false, int check_lower_and_higher_dim_simplices = 2,
688 Simplex_set *p_inconsistent_simplices = NULL) const {
689#if defined(DEBUG_TRACES) || defined(GUDHI_TC_PROFILING)
690 std::cerr << yellow << "\nExporting the TC as a Simplicial_complex... " << white;
691#endif
692#ifdef GUDHI_TC_PROFILING
693 Gudhi::Clock t;
694#endif
695
696 int max_dim = -1;
697 complex.clear();
698
699 // For each triangulation
700 for (std::size_t idx = 0; idx < m_points.size(); ++idx) {
701 // For each cell of the star
702 Star::const_iterator it_inc_simplex = m_stars[idx].begin();
703 Star::const_iterator it_inc_simplex_end = m_stars[idx].end();
704 for (; it_inc_simplex != it_inc_simplex_end; ++it_inc_simplex) {
705 Simplex c = *it_inc_simplex;
706
707 // Don't export infinite cells
708 if (!export_infinite_simplices && is_infinite(c)) continue;
709
710 if (static_cast<int>(c.size()) > max_dim) max_dim = static_cast<int>(c.size());
711 // Add the missing center vertex
712 c.insert(idx);
713
714 if (!export_inconsistent_simplices && !is_simplex_consistent(c)) continue;
715
716 // Unusual simplex dim?
717 if (check_lower_and_higher_dim_simplices == 2 && max_dim != -1 && static_cast<int>(c.size()) != max_dim) {
718 // Let's activate the check
719 std::cerr << red
720 << "Info: check_lower_and_higher_dim_simplices ACTIVATED. "
721 "Export might be take some time...\n"
722 << white;
723 check_lower_and_higher_dim_simplices = 1;
724 }
725
726 // Try to insert the simplex
727 bool added = complex.add_simplex(c, check_lower_and_higher_dim_simplices == 1);
728
729 // Inconsistent?
730 if (p_inconsistent_simplices && added && !is_simplex_consistent(c)) {
731 p_inconsistent_simplices->insert(c);
732 }
733 }
734 }
735
736#ifdef GUDHI_TC_PROFILING
737 t.end();
738 std::cerr << yellow << "done in " << t.num_seconds() << " seconds.\n" << white;
739#elif defined(DEBUG_TRACES)
740 std::cerr << yellow << "done.\n" << white;
741#endif
742
743 return max_dim;
744 }
745
746 template <typename ProjectionFunctor = CGAL::Identity<Point> >
747 std::ostream &export_to_off(const Simplicial_complex &complex, std::ostream &os,
748 Simplex_set const *p_simpl_to_color_in_red = NULL,
749 Simplex_set const *p_simpl_to_color_in_green = NULL,
750 Simplex_set const *p_simpl_to_color_in_blue = NULL,
751 ProjectionFunctor const &point_projection = ProjectionFunctor()) const {
752 return export_to_off(os, false, p_simpl_to_color_in_red, p_simpl_to_color_in_green, p_simpl_to_color_in_blue,
753 &complex, point_projection);
754 }
755
756 template <typename ProjectionFunctor = CGAL::Identity<Point> >
757 std::ostream &export_to_off(std::ostream &os, bool color_inconsistencies = false,
758 Simplex_set const *p_simpl_to_color_in_red = NULL,
759 Simplex_set const *p_simpl_to_color_in_green = NULL,
760 Simplex_set const *p_simpl_to_color_in_blue = NULL,
761 const Simplicial_complex *p_complex = NULL,
762 ProjectionFunctor const &point_projection = ProjectionFunctor()) const {
763 if (m_points.empty()) return os;
764
765 if (m_ambient_dim < 2) {
766 std::cerr << "Error: export_to_off => ambient dimension should be >= 2.\n";
767 os << "Error: export_to_off => ambient dimension should be >= 2.\n";
768 return os;
769 }
770 if (m_ambient_dim > 3) {
771 std::cerr << "Warning: export_to_off => ambient dimension should be "
772 "<= 3. Only the first 3 coordinates will be exported.\n";
773 }
774
775 if (m_intrinsic_dim < 1 || m_intrinsic_dim > 3) {
776 std::cerr << "Error: export_to_off => intrinsic dimension should be "
777 "between 1 and 3.\n";
778 os << "Error: export_to_off => intrinsic dimension should be "
779 "between 1 and 3.\n";
780 return os;
781 }
782
783 std::stringstream output;
784 std::size_t num_simplices, num_vertices;
785 export_vertices_to_off(output, num_vertices, false, point_projection);
786 if (p_complex) {
787 export_simplices_to_off(*p_complex, output, num_simplices, p_simpl_to_color_in_red, p_simpl_to_color_in_green,
788 p_simpl_to_color_in_blue);
789 } else {
790 export_simplices_to_off(output, num_simplices, color_inconsistencies, p_simpl_to_color_in_red,
791 p_simpl_to_color_in_green, p_simpl_to_color_in_blue);
792 }
793
794#ifdef GUDHI_TC_EXPORT_NORMALS
795 os << "N";
796#endif
797
798 os << "OFF \n"
799 << num_vertices << " " << num_simplices << " "
800 << "0 \n"
801 << output.str();
802
803 return os;
804 }
805
806 private:
807 void refresh_tangential_complex() {
808#if defined(DEBUG_TRACES) || defined(GUDHI_TC_PROFILING)
809 std::cerr << yellow << "\nRefreshing TC... " << white;
810#endif
811
812#ifdef GUDHI_TC_PROFILING
813 Gudhi::Clock t;
814#endif
815#ifdef GUDHI_USE_TBB
816 // Parallel
817 if (boost::is_convertible<Concurrency_tag, CGAL::Parallel_tag>::value) {
818 tbb::parallel_for(tbb::blocked_range<size_t>(0, m_points.size()), Compute_tangent_triangulation(*this));
819 } else {
820#endif // GUDHI_USE_TBB
821 // Sequential
822 for (std::size_t i = 0; i < m_points.size(); ++i) compute_tangent_triangulation(i);
823#ifdef GUDHI_USE_TBB
824 }
825#endif // GUDHI_USE_TBB
826
827#ifdef GUDHI_TC_PROFILING
828 t.end();
829 std::cerr << yellow << "done in " << t.num_seconds() << " seconds.\n" << white;
830#elif defined(DEBUG_TRACES)
831 std::cerr << yellow << "done.\n" << white;
832#endif
833 }
834
835 // If the list of perturbed points is provided, it is much faster
836 template <typename Point_indices_range>
837 void refresh_tangential_complex(Point_indices_range const &perturbed_points_indices) {
838#if defined(DEBUG_TRACES) || defined(GUDHI_TC_PROFILING)
839 std::cerr << yellow << "\nRefreshing TC... " << white;
840#endif
841
842#ifdef GUDHI_TC_PROFILING
843 Gudhi::Clock t;
844#endif
845
846 // ANN tree containing only the perturbed points
847 Points_ds updated_pts_ds(m_points, perturbed_points_indices);
848
849#ifdef GUDHI_USE_TBB
850 // Parallel
851 if (boost::is_convertible<Concurrency_tag, CGAL::Parallel_tag>::value) {
852 tbb::parallel_for(tbb::blocked_range<size_t>(0, m_points.size()),
853 Refresh_tangent_triangulation(*this, updated_pts_ds));
854 } else {
855#endif // GUDHI_USE_TBB
856 // Sequential
857 for (std::size_t i = 0; i < m_points.size(); ++i) refresh_tangent_triangulation(i, updated_pts_ds);
858#ifdef GUDHI_USE_TBB
859 }
860#endif // GUDHI_USE_TBB
861
862#ifdef GUDHI_TC_PROFILING
863 t.end();
864 std::cerr << yellow << "done in " << t.num_seconds() << " seconds.\n" << white;
865#elif defined(DEBUG_TRACES)
866 std::cerr << yellow << "done.\n" << white;
867#endif
868 }
869
870 void export_inconsistent_stars_to_OFF_files(std::string const &filename_base) const {
871 // For each triangulation
872 for (std::size_t idx = 0; idx < m_points.size(); ++idx) {
873 // We build a SC along the way in case it's inconsistent
874 Simplicial_complex sc;
875 // For each cell
876 bool is_inconsistent = false;
877 Star::const_iterator it_inc_simplex = m_stars[idx].begin();
878 Star::const_iterator it_inc_simplex_end = m_stars[idx].end();
879 for (; it_inc_simplex != it_inc_simplex_end; ++it_inc_simplex) {
880 // Skip infinite cells
881 if (is_infinite(*it_inc_simplex)) continue;
882
883 Simplex c = *it_inc_simplex;
884 c.insert(idx); // Add the missing index
885
886 sc.add_simplex(c);
887
888 // If we do not already know this star is inconsistent, test it
889 if (!is_inconsistent && !is_simplex_consistent(c)) is_inconsistent = true;
890 }
891
892 if (is_inconsistent) {
893 // Export star to OFF file
894 std::stringstream output_filename;
895 output_filename << filename_base << "_" << idx << ".off";
896 std::ofstream off_stream(output_filename.str().c_str());
897 export_to_off(sc, off_stream);
898 }
899 }
900 }
901
902 class Compare_distance_to_ref_point {
903 public:
904 Compare_distance_to_ref_point(Point const &ref, K const &k) : m_ref(ref), m_k(k) {}
905
906 bool operator()(Point const &p1, Point const &p2) {
907 typename K::Squared_distance_d sqdist = m_k.squared_distance_d_object();
908 return sqdist(p1, m_ref) < sqdist(p2, m_ref);
909 }
910
911 private:
912 Point const &m_ref;
913 K const &m_k;
914 };
915
916#ifdef GUDHI_USE_TBB
917 // Functor for compute_tangential_complex function
918 class Compute_tangent_triangulation {
919 Tangential_complex &m_tc;
920
921 public:
922 // Constructor
923 Compute_tangent_triangulation(Tangential_complex &tc) : m_tc(tc) {}
924
925 // Constructor
926 Compute_tangent_triangulation(const Compute_tangent_triangulation &ctt) : m_tc(ctt.m_tc) {}
927
928 // operator()
929 void operator()(const tbb::blocked_range<size_t> &r) const {
930 for (size_t i = r.begin(); i != r.end(); ++i) m_tc.compute_tangent_triangulation(i);
931 }
932 };
933
934 // Functor for refresh_tangential_complex function
935 class Refresh_tangent_triangulation {
936 Tangential_complex &m_tc;
937 Points_ds const &m_updated_pts_ds;
938
939 public:
940 // Constructor
941 Refresh_tangent_triangulation(Tangential_complex &tc, Points_ds const &updated_pts_ds)
942 : m_tc(tc), m_updated_pts_ds(updated_pts_ds) {}
943
944 // Constructor
945 Refresh_tangent_triangulation(const Refresh_tangent_triangulation &ctt)
946 : m_tc(ctt.m_tc), m_updated_pts_ds(ctt.m_updated_pts_ds) {}
947
948 // operator()
949 void operator()(const tbb::blocked_range<size_t> &r) const {
950 for (size_t i = r.begin(); i != r.end(); ++i) m_tc.refresh_tangent_triangulation(i, m_updated_pts_ds);
951 }
952 };
953#endif // GUDHI_USE_TBB
954
955 bool is_infinite(Simplex const &s) const { return *s.rbegin() == (std::numeric_limits<std::size_t>::max)(); }
956
957 // Output: "triangulation" is a Regular Triangulation containing at least the
958 // star of "center_pt"
959 // Returns the handle of the center vertex
960 Tr_vertex_handle compute_star(std::size_t i, const Point &center_pt, const Tangent_space_basis &tsb,
961 Triangulation &triangulation, bool verbose = false) {
962 int tangent_space_dim = tsb.dimension();
963 const Tr_traits &local_tr_traits = triangulation.geom_traits();
964
965 // Kernel functor & objects
966 typename K::Squared_distance_d k_sqdist = m_k.squared_distance_d_object();
967
968 // Triangulation's traits functor & objects
969 typename Tr_traits::Compute_weight_d point_weight = local_tr_traits.compute_weight_d_object();
970#if CGAL_VERSION_NR < 1050200000
971 typename Tr_traits::Power_center_d power_center = local_tr_traits.power_center_d_object();
972#else
973 typename Tr_traits::Construct_power_sphere_d power_center = local_tr_traits.construct_power_sphere_d_object();
974#endif
975
976 //***************************************************
977 // Build a minimal triangulation in the tangent space
978 // (we only need the star of p)
979 //***************************************************
980
981 // Insert p
982 Tr_point proj_wp;
983 if (i == tsb.origin()) {
984 // Insert {(0, 0, 0...), m_weights[i]}
985 proj_wp = local_tr_traits.construct_weighted_point_d_object()(
986 local_tr_traits.construct_point_d_object()(tangent_space_dim, CGAL::ORIGIN), m_weights[i]);
987 } else {
988 const Weighted_point &wp = compute_perturbed_weighted_point(i);
989 proj_wp = project_point_and_compute_weight(wp, tsb, local_tr_traits);
990 }
991
992 Tr_vertex_handle center_vertex = triangulation.insert(proj_wp);
993 center_vertex->data() = i;
994 if (verbose) std::cerr << "* Inserted point #" << i << "\n";
995
996#ifdef GUDHI_TC_VERY_VERBOSE
997 std::size_t num_attempts_to_insert_points = 1;
998 std::size_t num_inserted_points = 1;
999#endif
1000 // const int NUM_NEIGHBORS = 150;
1001 // KNS_range ins_range = m_points_ds.k_nearest_neighbors(center_pt, NUM_NEIGHBORS);
1002 INS_range ins_range = m_points_ds.incremental_nearest_neighbors(center_pt);
1003
1004 // While building the local triangulation, we keep the radius
1005 // of the sphere "star sphere" centered at "center_vertex"
1006 // and which contains all the
1007 // circumspheres of the star of "center_vertex"
1008 // If th the m_max_squared_edge_length is set the maximal radius of the "star sphere"
1009 // is at most square root of m_max_squared_edge_length
1010 std::optional<FT> squared_star_sphere_radius_plus_margin = m_max_squared_edge_length;
1011
1012 // Insert points until we find a point which is outside "star sphere"
1013 for (auto nn_it = ins_range.begin(); nn_it != ins_range.end(); ++nn_it) {
1014 std::size_t neighbor_point_idx = nn_it->first;
1015
1016 // ith point = p, which is already inserted
1017 if (neighbor_point_idx != i) {
1018 // No need to lock the Mutex_for_perturb here since this will not be
1019 // called while other threads are perturbing the positions
1020 Point neighbor_pt;
1021 FT neighbor_weight;
1022 compute_perturbed_weighted_point(neighbor_point_idx, neighbor_pt, neighbor_weight);
1023 GUDHI_CHECK(!m_max_squared_edge_length ||
1024 squared_star_sphere_radius_plus_margin.value() <= m_max_squared_edge_length.value(),
1025 std::invalid_argument("Tangential_complex::compute_star - set a bigger value with set_max_squared_edge_length."));
1026 if (squared_star_sphere_radius_plus_margin &&
1027 k_sqdist(center_pt, neighbor_pt) > squared_star_sphere_radius_plus_margin.value()) {
1028 GUDHI_CHECK(triangulation.current_dimension() >= tangent_space_dim,
1029 std::invalid_argument("Tangential_complex::compute_star - Dimension of the star is only " + \
1030 std::to_string(triangulation.current_dimension())));
1031 break;
1032 }
1033
1034 Tr_point proj_pt = project_point_and_compute_weight(neighbor_pt, neighbor_weight, tsb, local_tr_traits);
1035
1036#ifdef GUDHI_TC_VERY_VERBOSE
1037 ++num_attempts_to_insert_points;
1038#endif
1039
1040 Tr_vertex_handle vh = triangulation.insert_if_in_star(proj_pt, center_vertex);
1041 // Tr_vertex_handle vh = triangulation.insert(proj_pt);
1042 if (vh != Tr_vertex_handle() && vh->data() == (std::numeric_limits<std::size_t>::max)()) {
1043#ifdef GUDHI_TC_VERY_VERBOSE
1044 ++num_inserted_points;
1045#endif
1046 if (verbose) std::cerr << "* Inserted point #" << neighbor_point_idx << "\n";
1047
1048 vh->data() = neighbor_point_idx;
1049
1050 // Let's recompute squared_star_sphere_radius_plus_margin
1051 if (triangulation.current_dimension() >= tangent_space_dim) {
1052 squared_star_sphere_radius_plus_margin = std::nullopt;
1053 // Get the incident cells and look for the biggest circumsphere
1054 std::vector<Tr_full_cell_handle> incident_cells;
1055 triangulation.incident_full_cells(center_vertex, std::back_inserter(incident_cells));
1056 for (typename std::vector<Tr_full_cell_handle>::iterator cit = incident_cells.begin();
1057 cit != incident_cells.end(); ++cit) {
1058 Tr_full_cell_handle cell = *cit;
1059 if (triangulation.is_infinite(cell)) {
1060 squared_star_sphere_radius_plus_margin = std::nullopt;
1061 break;
1062 } else {
1063 // Note that this uses the perturbed point since it uses
1064 // the points of the local triangulation
1065 Tr_point c =
1066 power_center(boost::make_transform_iterator(cell->vertices_begin(),
1067 vertex_handle_to_point<Tr_point, Tr_vertex_handle>),
1068 boost::make_transform_iterator(cell->vertices_end(),
1069 vertex_handle_to_point<Tr_point, Tr_vertex_handle>));
1070
1071 FT sq_power_sphere_diam = 4 * point_weight(c);
1072
1073 if (!squared_star_sphere_radius_plus_margin ||
1074 sq_power_sphere_diam > squared_star_sphere_radius_plus_margin.value()) {
1075 squared_star_sphere_radius_plus_margin = sq_power_sphere_diam;
1076 }
1077 }
1078 }
1079
1080 // Let's add the margin, now
1081 // The value depends on whether we perturb weight or position
1082 if (squared_star_sphere_radius_plus_margin) {
1083 // "2*m_last_max_perturb" because both points can be perturbed
1084 squared_star_sphere_radius_plus_margin =
1085 CGAL::square(std::sqrt(squared_star_sphere_radius_plus_margin.value()) + 2 * m_last_max_perturb);
1086
1087 // Reduce the square radius to m_max_squared_edge_length if necessary
1088 if (m_max_squared_edge_length && squared_star_sphere_radius_plus_margin.value() > m_max_squared_edge_length.value()) {
1089 squared_star_sphere_radius_plus_margin = m_max_squared_edge_length.value();
1090 }
1091
1092 // Save it in `m_squared_star_spheres_radii_incl_margin`
1093 m_squared_star_spheres_radii_incl_margin[i] = squared_star_sphere_radius_plus_margin.value();
1094 } else {
1095 if (m_max_squared_edge_length) {
1096 squared_star_sphere_radius_plus_margin = m_max_squared_edge_length.value();
1097 m_squared_star_spheres_radii_incl_margin[i] = m_max_squared_edge_length.value();
1098 } else {
1099 m_squared_star_spheres_radii_incl_margin[i] = FT(-1);
1100 }
1101 }
1102 }
1103 }
1104 }
1105 }
1106
1107 return center_vertex;
1108 }
1109
1110 void refresh_tangent_triangulation(std::size_t i, Points_ds const &updated_pts_ds, bool verbose = false) {
1111 if (verbose) std::cerr << "** Refreshing tangent tri #" << i << " **\n";
1112
1113 if (m_squared_star_spheres_radii_incl_margin[i] == FT(-1)) return compute_tangent_triangulation(i, verbose);
1114
1115 Point center_point = compute_perturbed_point(i);
1116 // Among updated point, what is the closer from our center point?
1117 std::size_t closest_pt_index = updated_pts_ds.k_nearest_neighbors(center_point, 1, false).begin()->first;
1118
1119 typename K::Construct_weighted_point_d k_constr_wp = m_k.construct_weighted_point_d_object();
1120#if CGAL_VERSION_NR < 1050200000
1121 typename K::Power_distance_d k_power_dist = m_k.power_distance_d_object();
1122#else
1123 typename K::Compute_power_product_d k_power_dist = m_k.compute_power_product_d_object();
1124#endif
1125
1126 // Construct a weighted point equivalent to the star sphere
1127 Weighted_point star_sphere = k_constr_wp(compute_perturbed_point(i), m_squared_star_spheres_radii_incl_margin[i]);
1128 Weighted_point closest_updated_point = compute_perturbed_weighted_point(closest_pt_index);
1129
1130 // Is the "closest point" inside our star sphere?
1131 if (k_power_dist(star_sphere, closest_updated_point) <= FT(0)) compute_tangent_triangulation(i, verbose);
1132 }
1133
1134 void compute_tangent_triangulation(std::size_t i, bool verbose = false) {
1135 if (verbose) std::cerr << "** Computing tangent tri #" << i << " **\n";
1136 // std::cerr << "***********************************************\n";
1137
1138 // No need to lock the mutex here since this will not be called while
1139 // other threads are perturbing the positions
1140 const Point center_pt = compute_perturbed_point(i);
1141 Tangent_space_basis &tsb = m_tangent_spaces[i];
1142
1143 // Estimate the tangent space
1144 if (!m_are_tangent_spaces_computed[i]) {
1145#ifdef GUDHI_TC_EXPORT_NORMALS
1146 tsb = compute_tangent_space(center_pt, i, true /*normalize*/, &m_orth_spaces[i]);
1147#else
1148 tsb = compute_tangent_space(center_pt, i);
1149#endif
1150 }
1151
1152#if defined(GUDHI_TC_PROFILING) && defined(GUDHI_TC_VERY_VERBOSE)
1153 Gudhi::Clock t;
1154#endif
1155 int tangent_space_dim = tangent_basis_dim(i);
1156 Triangulation &local_tr = m_triangulations[i].construct_triangulation(tangent_space_dim);
1157
1158 m_triangulations[i].center_vertex() = compute_star(i, center_pt, tsb, local_tr, verbose);
1159
1160#if defined(GUDHI_TC_PROFILING) && defined(GUDHI_TC_VERY_VERBOSE)
1161 t.end();
1162 std::cerr << " - triangulation construction: " << t.num_seconds() << " s.\n";
1163 t.reset();
1164#endif
1165
1166#ifdef GUDHI_TC_VERY_VERBOSE
1167 std::cerr << "Inserted " << num_inserted_points << " points / " << num_attempts_to_insert_points
1168 << " attempts to compute the star\n";
1169#endif
1170
1171 update_star(i);
1172
1173#if defined(GUDHI_TC_PROFILING) && defined(GUDHI_TC_VERY_VERBOSE)
1174 t.end();
1175 std::cerr << " - update_star: " << t.num_seconds() << " s.\n";
1176#endif
1177 }
1178
1179 // Updates m_stars[i] directly from m_triangulations[i]
1180
1181 void update_star(std::size_t i) {
1182 Star &star = m_stars[i];
1183 star.clear();
1184 Triangulation &local_tr = m_triangulations[i].tr();
1185 Tr_vertex_handle center_vertex = m_triangulations[i].center_vertex();
1186 int cur_dim_plus_1 = local_tr.current_dimension() + 1;
1187
1188 std::vector<Tr_full_cell_handle> incident_cells;
1189 local_tr.incident_full_cells(center_vertex, std::back_inserter(incident_cells));
1190
1191 typename std::vector<Tr_full_cell_handle>::const_iterator it_c = incident_cells.begin();
1192 typename std::vector<Tr_full_cell_handle>::const_iterator it_c_end = incident_cells.end();
1193 // For each cell
1194 for (; it_c != it_c_end; ++it_c) {
1195 // Will contain all indices except center_vertex
1196 Incident_simplex incident_simplex;
1197 for (int j = 0; j < cur_dim_plus_1; ++j) {
1198 std::size_t index = (*it_c)->vertex(j)->data();
1199 if (index != i) incident_simplex.insert(index);
1200 }
1201 GUDHI_CHECK(incident_simplex.size() == cur_dim_plus_1 - 1,
1202 std::logic_error("update_star: wrong size of incident simplex"));
1203 star.push_back(incident_simplex);
1204 }
1205 }
1206
1207 // Estimates tangent subspaces using PCA
1208
1209 Tangent_space_basis compute_tangent_space(const Point &p, const std::size_t i, bool normalize_basis = true,
1210 Orthogonal_space_basis *p_orth_space_basis = NULL) {
1211 unsigned int num_pts_for_pca =
1212 (std::min)(static_cast<unsigned int>(std::pow(GUDHI_TC_BASE_VALUE_FOR_PCA, m_intrinsic_dim)),
1213 static_cast<unsigned int>(m_points.size()));
1214
1215 // Kernel functors
1216 typename K::Construct_vector_d constr_vec = m_k.construct_vector_d_object();
1217 typename K::Compute_coordinate_d coord = m_k.compute_coordinate_d_object();
1218
1219#ifdef GUDHI_TC_USE_ANOTHER_POINT_SET_FOR_TANGENT_SPACE_ESTIM
1220 KNS_range kns_range = m_points_ds_for_tse.k_nearest_neighbors(p, num_pts_for_pca, false);
1221 const Points &points_for_pca = m_points_for_tse;
1222#else
1223 KNS_range kns_range = m_points_ds.k_nearest_neighbors(p, num_pts_for_pca, false);
1224 const Points &points_for_pca = m_points;
1225#endif
1226
1227 // One row = one point
1228 Eigen::MatrixXd mat_points(num_pts_for_pca, m_ambient_dim);
1229 auto nn_it = kns_range.begin();
1230 for (unsigned int j = 0; j < num_pts_for_pca && nn_it != kns_range.end(); ++j, ++nn_it) {
1231 for (int i = 0; i < m_ambient_dim; ++i) {
1232 mat_points(j, i) = CGAL::to_double(coord(points_for_pca[nn_it->first], i));
1233 }
1234 }
1235 Eigen::MatrixXd centered = mat_points.rowwise() - mat_points.colwise().mean();
1236 Eigen::MatrixXd cov = centered.adjoint() * centered;
1237 Eigen::SelfAdjointEigenSolver<Eigen::MatrixXd> eig(cov);
1238
1239 Tangent_space_basis tsb(i); // p = compute_perturbed_point(i) here
1240
1241 // The eigenvectors are sorted in increasing order of their corresponding
1242 // eigenvalues
1243 for (int j = m_ambient_dim - 1; j >= m_ambient_dim - m_intrinsic_dim; --j) {
1244 if (normalize_basis) {
1245 Vector v = constr_vec(m_ambient_dim, eig.eigenvectors().col(j).data(),
1246 eig.eigenvectors().col(j).data() + m_ambient_dim);
1247 tsb.push_back(normalize_vector(v, m_k));
1248 } else {
1249 tsb.push_back(constr_vec(m_ambient_dim, eig.eigenvectors().col(j).data(),
1250 eig.eigenvectors().col(j).data() + m_ambient_dim));
1251 }
1252 }
1253
1254 if (p_orth_space_basis) {
1255 p_orth_space_basis->set_origin(i);
1256 for (int j = m_ambient_dim - m_intrinsic_dim - 1; j >= 0; --j) {
1257 if (normalize_basis) {
1258 Vector v = constr_vec(m_ambient_dim, eig.eigenvectors().col(j).data(),
1259 eig.eigenvectors().col(j).data() + m_ambient_dim);
1260 p_orth_space_basis->push_back(normalize_vector(v, m_k));
1261 } else {
1262 p_orth_space_basis->push_back(constr_vec(m_ambient_dim, eig.eigenvectors().col(j).data(),
1263 eig.eigenvectors().col(j).data() + m_ambient_dim));
1264 }
1265 }
1266 }
1267
1268 m_are_tangent_spaces_computed[i] = true;
1269
1270 return tsb;
1271 }
1272
1273 // Compute the space tangent to a simplex (p1, p2, ... pn)
1274 // TODO(CJ): Improve this?
1275 // Basically, it takes all the neighbor points to p1, p2... pn and runs PCA
1276 // on it. Note that most points are duplicated.
1277
1278 Tangent_space_basis compute_tangent_space(const Simplex &s, bool normalize_basis = true) {
1279 unsigned int num_pts_for_pca =
1280 (std::min)(static_cast<unsigned int>(std::pow(GUDHI_TC_BASE_VALUE_FOR_PCA, m_intrinsic_dim)),
1281 static_cast<unsigned int>(m_points.size()));
1282
1283 // Kernel functors
1284 typename K::Construct_vector_d constr_vec = m_k.construct_vector_d_object();
1285 typename K::Compute_coordinate_d coord = m_k.compute_coordinate_d_object();
1286 typename K::Squared_length_d sqlen = m_k.squared_length_d_object();
1287 typename K::Scaled_vector_d scaled_vec = m_k.scaled_vector_d_object();
1288 typename K::Scalar_product_d scalar_pdct = m_k.scalar_product_d_object();
1289 typename K::Difference_of_vectors_d diff_vec = m_k.difference_of_vectors_d_object();
1290
1291 // One row = one point
1292 Eigen::MatrixXd mat_points(s.size() * num_pts_for_pca, m_ambient_dim);
1293 unsigned int current_row = 0;
1294
1295 for (Simplex::const_iterator it_index = s.begin(); it_index != s.end(); ++it_index) {
1296 const Point &p = m_points[*it_index];
1297
1298#ifdef GUDHI_TC_USE_ANOTHER_POINT_SET_FOR_TANGENT_SPACE_ESTIM
1299 KNS_range kns_range = m_points_ds_for_tse.k_nearest_neighbors(p, num_pts_for_pca, false);
1300 const Points &points_for_pca = m_points_for_tse;
1301#else
1302 KNS_range kns_range = m_points_ds.k_nearest_neighbors(p, num_pts_for_pca, false);
1303 const Points &points_for_pca = m_points;
1304#endif
1305
1306 auto nn_it = kns_range.begin();
1307 for (; current_row < num_pts_for_pca && nn_it != kns_range.end(); ++current_row, ++nn_it) {
1308 for (int i = 0; i < m_ambient_dim; ++i) {
1309 mat_points(current_row, i) = CGAL::to_double(coord(points_for_pca[nn_it->first], i));
1310 }
1311 }
1312 }
1313 Eigen::MatrixXd centered = mat_points.rowwise() - mat_points.colwise().mean();
1314 Eigen::MatrixXd cov = centered.adjoint() * centered;
1315 Eigen::SelfAdjointEigenSolver<Eigen::MatrixXd> eig(cov);
1316
1317 Tangent_space_basis tsb;
1318
1319 // The eigenvectors are sorted in increasing order of their corresponding
1320 // eigenvalues
1321 for (int j = m_ambient_dim - 1; j >= m_ambient_dim - m_intrinsic_dim; --j) {
1322 if (normalize_basis) {
1323 Vector v = constr_vec(m_ambient_dim, eig.eigenvectors().col(j).data(),
1324 eig.eigenvectors().col(j).data() + m_ambient_dim);
1325 tsb.push_back(normalize_vector(v, m_k));
1326 } else {
1327 tsb.push_back(constr_vec(m_ambient_dim, eig.eigenvectors().col(j).data(),
1328 eig.eigenvectors().col(j).data() + m_ambient_dim));
1329 }
1330 }
1331
1332 return tsb;
1333 }
1334
1335 // Returns the dimension of the ith local triangulation
1336
1337 int tangent_basis_dim(std::size_t i) const { return m_tangent_spaces[i].dimension(); }
1338
1339 Point compute_perturbed_point(std::size_t pt_idx) const {
1340#ifdef GUDHI_TC_PERTURB_POSITION
1341 return m_k.translated_point_d_object()(m_points[pt_idx], m_translations[pt_idx]);
1342#else
1343 return m_points[pt_idx];
1344#endif
1345 }
1346
1347 void compute_perturbed_weighted_point(std::size_t pt_idx, Point &p, FT &w) const {
1348#ifdef GUDHI_TC_PERTURB_POSITION
1349 p = m_k.translated_point_d_object()(m_points[pt_idx], m_translations[pt_idx]);
1350#else
1351 p = m_points[pt_idx];
1352#endif
1353 w = m_weights[pt_idx];
1354 }
1355
1356 Weighted_point compute_perturbed_weighted_point(std::size_t pt_idx) const {
1357 typename K::Construct_weighted_point_d k_constr_wp = m_k.construct_weighted_point_d_object();
1358
1359 Weighted_point wp = k_constr_wp(
1360#ifdef GUDHI_TC_PERTURB_POSITION
1361 m_k.translated_point_d_object()(m_points[pt_idx], m_translations[pt_idx]),
1362#else
1363 m_points[pt_idx],
1364#endif
1365 m_weights[pt_idx]);
1366
1367 return wp;
1368 }
1369
1370 Point unproject_point(const Tr_point &p, const Tangent_space_basis &tsb, const Tr_traits &tr_traits) const {
1371 typename K::Translated_point_d k_transl = m_k.translated_point_d_object();
1372 typename K::Scaled_vector_d k_scaled_vec = m_k.scaled_vector_d_object();
1373 typename Tr_traits::Compute_coordinate_d coord = tr_traits.compute_coordinate_d_object();
1374
1375 Point global_point = compute_perturbed_point(tsb.origin());
1376 for (int i = 0; i < m_intrinsic_dim; ++i) global_point = k_transl(global_point, k_scaled_vec(tsb[i], coord(p, i)));
1377
1378 return global_point;
1379 }
1380
1381 // Project the point in the tangent space
1382 // Resulting point coords are expressed in tsb's space
1383 Tr_bare_point project_point(const Point &p, const Tangent_space_basis &tsb, const Tr_traits &tr_traits) const {
1384 typename K::Scalar_product_d scalar_pdct = m_k.scalar_product_d_object();
1385 typename K::Difference_of_points_d diff_points = m_k.difference_of_points_d_object();
1386
1387 Vector v = diff_points(p, compute_perturbed_point(tsb.origin()));
1388
1389 std::vector<FT> coords;
1390 // Ambiant-space coords of the projected point
1391 coords.reserve(tsb.dimension());
1392 for (std::size_t i = 0; i < m_intrinsic_dim; ++i) {
1393 // Local coords are given by the scalar product with the vectors of tsb
1394 FT coord = scalar_pdct(v, tsb[i]);
1395 coords.push_back(coord);
1396 }
1397
1398 return tr_traits.construct_point_d_object()(static_cast<int>(coords.size()), coords.begin(), coords.end());
1399 }
1400
1401 // Project the point in the tangent space
1402 // The weight will be the squared distance between p and the projection of p
1403 // Resulting point coords are expressed in tsb's space
1404
1405 Tr_point project_point_and_compute_weight(const Weighted_point &wp, const Tangent_space_basis &tsb,
1406 const Tr_traits &tr_traits) const {
1407 typename K::Point_drop_weight_d k_drop_w = m_k.point_drop_weight_d_object();
1408 typename K::Compute_weight_d k_point_weight = m_k.compute_weight_d_object();
1409 return project_point_and_compute_weight(k_drop_w(wp), k_point_weight(wp), tsb, tr_traits);
1410 }
1411
1412 // Same as above, with slightly different parameters
1413 Tr_point project_point_and_compute_weight(const Point &p, const FT w, const Tangent_space_basis &tsb,
1414 const Tr_traits &tr_traits) const {
1415 const int point_dim = m_k.point_dimension_d_object()(p);
1416
1417 typename K::Construct_point_d constr_pt = m_k.construct_point_d_object();
1418 typename K::Scalar_product_d scalar_pdct = m_k.scalar_product_d_object();
1419 typename K::Difference_of_points_d diff_points = m_k.difference_of_points_d_object();
1420 typename K::Compute_coordinate_d coord = m_k.compute_coordinate_d_object();
1421 typename K::Construct_cartesian_const_iterator_d ccci = m_k.construct_cartesian_const_iterator_d_object();
1422
1423 Point origin = compute_perturbed_point(tsb.origin());
1424 Vector v = diff_points(p, origin);
1425
1426 // Same dimension? Then the weight is 0
1427 bool same_dim = (point_dim == tsb.dimension());
1428
1429 std::vector<FT> coords;
1430 // Ambiant-space coords of the projected point
1431 std::vector<FT> p_proj(ccci(origin), ccci(origin, 0));
1432 coords.reserve(tsb.dimension());
1433 for (int i = 0; i < tsb.dimension(); ++i) {
1434 // Local coords are given by the scalar product with the vectors of tsb
1435 FT c = scalar_pdct(v, tsb[i]);
1436 coords.push_back(c);
1437
1438 // p_proj += c * tsb[i]
1439 if (!same_dim) {
1440 for (int j = 0; j < point_dim; ++j) p_proj[j] += c * coord(tsb[i], j);
1441 }
1442 }
1443
1444 // Same dimension? Then the weight is 0
1445 FT sq_dist_to_proj_pt = 0;
1446 if (!same_dim) {
1447 Point projected_pt = constr_pt(point_dim, p_proj.begin(), p_proj.end());
1448 sq_dist_to_proj_pt = m_k.squared_distance_d_object()(p, projected_pt);
1449 }
1450
1451 return tr_traits.construct_weighted_point_d_object()(
1452 tr_traits.construct_point_d_object()(static_cast<int>(coords.size()), coords.begin(), coords.end()),
1453 w - sq_dist_to_proj_pt);
1454 }
1455
1456 // Project all the points in the tangent space
1457
1458 template <typename Indexed_point_range>
1459 std::vector<Tr_point> project_points_and_compute_weights(const Indexed_point_range &point_indices,
1460 const Tangent_space_basis &tsb,
1461 const Tr_traits &tr_traits) const {
1462 std::vector<Tr_point> ret;
1463 for (typename Indexed_point_range::const_iterator it = point_indices.begin(), it_end = point_indices.end();
1464 it != it_end; ++it) {
1465 ret.push_back(project_point_and_compute_weight(compute_perturbed_weighted_point(*it), tsb, tr_traits));
1466 }
1467 return ret;
1468 }
1469
1470 // A simplex here is a local tri's full cell handle
1471
1472 bool is_simplex_consistent(Tr_full_cell_handle fch, int cur_dim) const {
1473 Simplex c;
1474 for (int i = 0; i < cur_dim + 1; ++i) {
1475 std::size_t data = fch->vertex(i)->data();
1476 c.insert(data);
1477 }
1478 return is_simplex_consistent(c);
1479 }
1480
1481 // A simplex here is a list of point indices
1482 // TODO(CJ): improve it like the other "is_simplex_consistent" below
1483
1484 bool is_simplex_consistent(Simplex const &simplex) const {
1485 // Check if the simplex is in the stars of all its vertices
1486 Simplex::const_iterator it_point_idx = simplex.begin();
1487 // For each point p of the simplex, we parse the incidents cells of p
1488 // and we check if "simplex" is among them
1489 for (; it_point_idx != simplex.end(); ++it_point_idx) {
1490 std::size_t point_idx = *it_point_idx;
1491 // Don't check infinite simplices
1492 if (point_idx == (std::numeric_limits<std::size_t>::max)()) continue;
1493
1494 Star const &star = m_stars[point_idx];
1495
1496 // What we're looking for is "simplex" \ point_idx
1497 Incident_simplex is_to_find = simplex;
1498 is_to_find.erase(point_idx);
1499
1500 // For each cell
1501 if (std::find(star.begin(), star.end(), is_to_find) == star.end()) return false;
1502 }
1503
1504 return true;
1505 }
1506
1507 // A simplex here is a list of point indices
1508 // "s" contains all the points of the simplex except "center_point"
1509 // This function returns the points whose star doesn't contain the simplex
1510 // N.B.: the function assumes that the simplex is contained in
1511 // star(center_point)
1512
1513 template <typename OutputIterator> // value_type = std::size_t
1514 bool is_simplex_consistent(std::size_t center_point,
1515 Incident_simplex const &s, // without "center_point"
1516 OutputIterator points_whose_star_does_not_contain_s,
1517 bool check_also_in_non_maximal_faces = false) const {
1518 Simplex full_simplex = s;
1519 full_simplex.insert(center_point);
1520
1521 // Check if the simplex is in the stars of all its vertices
1522 Incident_simplex::const_iterator it_point_idx = s.begin();
1523 // For each point p of the simplex, we parse the incidents cells of p
1524 // and we check if "simplex" is among them
1525 for (; it_point_idx != s.end(); ++it_point_idx) {
1526 std::size_t point_idx = *it_point_idx;
1527 // Don't check infinite simplices
1528 if (point_idx == (std::numeric_limits<std::size_t>::max)()) continue;
1529
1530 Star const &star = m_stars[point_idx];
1531
1532 // What we're looking for is full_simplex \ point_idx
1533 Incident_simplex is_to_find = full_simplex;
1534 is_to_find.erase(point_idx);
1535
1536 if (check_also_in_non_maximal_faces) {
1537 // For each simplex "is" of the star, check if ic_to_simplex is
1538 // included in "is"
1539 bool found = false;
1540 for (Star::const_iterator is = star.begin(), is_end = star.end(); !found && is != is_end; ++is) {
1541 if (std::includes(is->begin(), is->end(), is_to_find.begin(), is_to_find.end())) found = true;
1542 }
1543
1544 if (!found) *points_whose_star_does_not_contain_s++ = point_idx;
1545 } else {
1546 // Does the star contain is_to_find?
1547 if (std::find(star.begin(), star.end(), is_to_find) == star.end())
1548 *points_whose_star_does_not_contain_s++ = point_idx;
1549 }
1550 }
1551
1552 return true;
1553 }
1554
1555 // A simplex here is a list of point indices
1556 // It looks for s in star(p).
1557 // "s" contains all the points of the simplex except p.
1558 bool is_simplex_in_star(std::size_t p, Incident_simplex const &s, bool check_also_in_non_maximal_faces = true) const {
1559 Star const &star = m_stars[p];
1560
1561 if (check_also_in_non_maximal_faces) {
1562 // For each simplex "is" of the star, check if ic_to_simplex is
1563 // included in "is"
1564 bool found = false;
1565 for (Star::const_iterator is = star.begin(), is_end = star.end(); !found && is != is_end; ++is) {
1566 if (std::includes(is->begin(), is->end(), s.begin(), s.end())) found = true;
1567 }
1568
1569 return found;
1570 } else {
1571 return !(std::find(star.begin(), star.end(), s) == star.end());
1572 }
1573 }
1574
1575 void perturb(std::size_t point_idx, double max_perturb) {
1576 const Tr_traits &local_tr_traits = m_triangulations[point_idx].tr().geom_traits();
1577 typename Tr_traits::Compute_coordinate_d coord = local_tr_traits.compute_coordinate_d_object();
1578 typename K::Translated_point_d k_transl = m_k.translated_point_d_object();
1579 typename K::Construct_vector_d k_constr_vec = m_k.construct_vector_d_object();
1580 typename K::Scaled_vector_d k_scaled_vec = m_k.scaled_vector_d_object();
1581
1582 CGAL::Random_points_in_ball_d<Tr_bare_point> tr_point_in_ball_generator(
1583 m_intrinsic_dim, CGAL::get_default_random().get_double(0., max_perturb));
1584
1585 Tr_point local_random_transl =
1586 local_tr_traits.construct_weighted_point_d_object()(*tr_point_in_ball_generator++, 0);
1587 Translation_for_perturb global_transl = k_constr_vec(m_ambient_dim);
1588 const Tangent_space_basis &tsb = m_tangent_spaces[point_idx];
1589 for (int i = 0; i < m_intrinsic_dim; ++i) {
1590 global_transl = k_transl(global_transl, k_scaled_vec(tsb[i], coord(local_random_transl, i)));
1591 }
1592 // Parallel
1593#if defined(GUDHI_USE_TBB)
1594 m_p_perturb_mutexes[point_idx].lock();
1595 m_translations[point_idx] = global_transl;
1596 m_p_perturb_mutexes[point_idx].unlock();
1597 // Sequential
1598#else
1599 m_translations[point_idx] = global_transl;
1600#endif
1601 }
1602
1603 // Return true if inconsistencies were found
1604 template <typename OutputIt>
1605 bool try_to_solve_inconsistencies_in_a_local_triangulation(
1606 std::size_t tr_index, double max_perturb, OutputIt perturbed_pts_indices = CGAL::Emptyset_iterator()) {
1607 bool is_inconsistent = false;
1608
1609 Star const &star = m_stars[tr_index];
1610
1611 // For each incident simplex
1612 Star::const_iterator it_inc_simplex = star.begin();
1613 Star::const_iterator it_inc_simplex_end = star.end();
1614 for (; it_inc_simplex != it_inc_simplex_end; ++it_inc_simplex) {
1615 const Incident_simplex &incident_simplex = *it_inc_simplex;
1616
1617 // Don't check infinite cells
1618 if (is_infinite(incident_simplex)) continue;
1619
1620 Simplex c = incident_simplex;
1621 c.insert(tr_index); // Add the missing index
1622
1623 // Perturb the center point
1624 if (!is_simplex_consistent(c)) {
1625 is_inconsistent = true;
1626
1627 std::size_t idx = tr_index;
1628
1629 perturb(tr_index, max_perturb);
1630 *perturbed_pts_indices++ = idx;
1631
1632 // We will try the other cells next time
1633 break;
1634 }
1635 }
1636
1637 return is_inconsistent;
1638 }
1639
1640 // 1st line: number of points
1641 // Then one point per line
1642 std::ostream &export_point_set(std::ostream &os, bool use_perturbed_points = false,
1643 const char *coord_separator = " ") const {
1644 if (use_perturbed_points) {
1645 std::vector<Point> perturbed_points;
1646 perturbed_points.reserve(m_points.size());
1647 for (std::size_t i = 0; i < m_points.size(); ++i) perturbed_points.push_back(compute_perturbed_point(i));
1648
1649 return export_point_set(m_k, perturbed_points, os, coord_separator);
1650 } else {
1651 return export_point_set(m_k, m_points, os, coord_separator);
1652 }
1653 }
1654
1655 template <typename ProjectionFunctor = CGAL::Identity<Point> >
1656 std::ostream &export_vertices_to_off(std::ostream &os, std::size_t &num_vertices, bool use_perturbed_points = false,
1657 ProjectionFunctor const &point_projection = ProjectionFunctor()) const {
1658 if (m_points.empty()) {
1659 num_vertices = 0;
1660 return os;
1661 }
1662
1663 // If m_intrinsic_dim = 1, we output each point two times
1664 // to be able to export each segment as a flat triangle with 3 different
1665 // indices (otherwise, Meshlab detects degenerated simplices)
1666 const int N = (m_intrinsic_dim == 1 ? 2 : 1);
1667
1668 // Kernel functors
1669 typename K::Compute_coordinate_d coord = m_k.compute_coordinate_d_object();
1670
1671#ifdef GUDHI_TC_EXPORT_ALL_COORDS_IN_OFF
1672 int num_coords = m_ambient_dim;
1673#else
1674 int num_coords = (std::min)(m_ambient_dim, 3);
1675#endif
1676
1677#ifdef GUDHI_TC_EXPORT_NORMALS
1678 OS_container::const_iterator it_os = m_orth_spaces.begin();
1679#endif
1680 typename Points::const_iterator it_p = m_points.begin();
1681 typename Points::const_iterator it_p_end = m_points.end();
1682 // For each point p
1683 for (std::size_t i = 0; it_p != it_p_end; ++it_p, ++i) {
1684 Point p = point_projection(use_perturbed_points ? compute_perturbed_point(i) : *it_p);
1685 for (int ii = 0; ii < N; ++ii) {
1686 int j = 0;
1687 for (; j < num_coords; ++j) os << CGAL::to_double(coord(p, j)) << " ";
1688 if (j == 2) os << "0";
1689
1690#ifdef GUDHI_TC_EXPORT_NORMALS
1691 for (j = 0; j < num_coords; ++j) os << " " << CGAL::to_double(coord(*it_os->begin(), j));
1692#endif
1693 os << "\n";
1694 }
1695#ifdef GUDHI_TC_EXPORT_NORMALS
1696 ++it_os;
1697#endif
1698 }
1699
1700 num_vertices = N * m_points.size();
1701 return os;
1702 }
1703
1704 std::ostream &export_simplices_to_off(std::ostream &os, std::size_t &num_OFF_simplices,
1705 bool color_inconsistencies = false,
1706 Simplex_set const *p_simpl_to_color_in_red = NULL,
1707 Simplex_set const *p_simpl_to_color_in_green = NULL,
1708 Simplex_set const *p_simpl_to_color_in_blue = NULL) const {
1709 // If m_intrinsic_dim = 1, each point is output two times
1710 // (see export_vertices_to_off)
1711 num_OFF_simplices = 0;
1712 std::size_t num_maximal_simplices = 0;
1713 std::size_t num_inconsistent_maximal_simplices = 0;
1714 std::size_t num_inconsistent_stars = 0;
1715 typename Tr_container::const_iterator it_tr = m_triangulations.begin();
1716 typename Tr_container::const_iterator it_tr_end = m_triangulations.end();
1717 // For each triangulation
1718 for (std::size_t idx = 0; it_tr != it_tr_end; ++it_tr, ++idx) {
1719 bool is_star_inconsistent = false;
1720
1721 Triangulation const &tr = it_tr->tr();
1722
1723 if (tr.current_dimension() < m_intrinsic_dim) continue;
1724
1725 // Color for this star
1726 std::stringstream color;
1727 // color << rand()%256 << " " << 100+rand()%156 << " " << 100+rand()%156;
1728 color << 128 << " " << 128 << " " << 128;
1729
1730 // Gather the triangles here, with an int telling its color
1731 typedef std::vector<std::pair<Simplex, int> > Star_using_triangles;
1732 Star_using_triangles star_using_triangles;
1733
1734 // For each cell of the star
1735 Star::const_iterator it_inc_simplex = m_stars[idx].begin();
1736 Star::const_iterator it_inc_simplex_end = m_stars[idx].end();
1737 for (; it_inc_simplex != it_inc_simplex_end; ++it_inc_simplex) {
1738 Simplex c = *it_inc_simplex;
1739 c.insert(idx);
1740 std::size_t num_vertices = c.size();
1741 ++num_maximal_simplices;
1742
1743 int color_simplex = -1; // -1=no color, 0=yellow, 1=red, 2=green, 3=blue
1744 if (color_inconsistencies && !is_simplex_consistent(c)) {
1745 ++num_inconsistent_maximal_simplices;
1746 color_simplex = 0;
1747 is_star_inconsistent = true;
1748 } else {
1749 if (p_simpl_to_color_in_red && std::find(p_simpl_to_color_in_red->begin(), p_simpl_to_color_in_red->end(),
1750 c) != p_simpl_to_color_in_red->end()) {
1751 color_simplex = 1;
1752 } else if (p_simpl_to_color_in_green &&
1753 std::find(p_simpl_to_color_in_green->begin(), p_simpl_to_color_in_green->end(), c) !=
1754 p_simpl_to_color_in_green->end()) {
1755 color_simplex = 2;
1756 } else if (p_simpl_to_color_in_blue &&
1757 std::find(p_simpl_to_color_in_blue->begin(), p_simpl_to_color_in_blue->end(), c) !=
1758 p_simpl_to_color_in_blue->end()) {
1759 color_simplex = 3;
1760 }
1761 }
1762
1763 // If m_intrinsic_dim = 1, each point is output two times,
1764 // so we need to multiply each index by 2
1765 // And if only 2 vertices, add a third one (each vertex is duplicated in
1766 // the file when m_intrinsic dim = 2)
1767 if (m_intrinsic_dim == 1) {
1768 Simplex tmp_c;
1769 Simplex::iterator it = c.begin();
1770 for (; it != c.end(); ++it) tmp_c.insert(*it * 2);
1771 if (num_vertices == 2) tmp_c.insert(*tmp_c.rbegin() + 1);
1772
1773 c = tmp_c;
1774 }
1775
1776 if (num_vertices <= 3) {
1777 star_using_triangles.push_back(std::make_pair(c, color_simplex));
1778 } else {
1779 // num_vertices >= 4: decompose the simplex in triangles
1780 std::vector<bool> booleans(num_vertices, false);
1781 std::fill(booleans.begin() + num_vertices - 3, booleans.end(), true);
1782 do {
1783 Simplex triangle;
1784 Simplex::iterator it = c.begin();
1785 for (int i = 0; it != c.end(); ++i, ++it) {
1786 if (booleans[i]) triangle.insert(*it);
1787 }
1788 star_using_triangles.push_back(std::make_pair(triangle, color_simplex));
1789 } while (std::next_permutation(booleans.begin(), booleans.end()));
1790 }
1791 }
1792
1793 // For each cell
1794 Star_using_triangles::const_iterator it_simplex = star_using_triangles.begin();
1795 Star_using_triangles::const_iterator it_simplex_end = star_using_triangles.end();
1796 for (; it_simplex != it_simplex_end; ++it_simplex) {
1797 const Simplex &c = it_simplex->first;
1798
1799 // Don't export infinite cells
1800 if (is_infinite(c)) continue;
1801
1802 int color_simplex = it_simplex->second;
1803
1804 std::stringstream sstr_c;
1805
1806 Simplex::const_iterator it_point_idx = c.begin();
1807 for (; it_point_idx != c.end(); ++it_point_idx) {
1808 sstr_c << *it_point_idx << " ";
1809 }
1810
1811 os << 3 << " " << sstr_c.str();
1812 if (color_inconsistencies || p_simpl_to_color_in_red || p_simpl_to_color_in_green || p_simpl_to_color_in_blue) {
1813 switch (color_simplex) {
1814 case 0:
1815 os << " 255 255 0";
1816 break;
1817 case 1:
1818 os << " 255 0 0";
1819 break;
1820 case 2:
1821 os << " 0 255 0";
1822 break;
1823 case 3:
1824 os << " 0 0 255";
1825 break;
1826 default:
1827 os << " " << color.str();
1828 break;
1829 }
1830 }
1831 ++num_OFF_simplices;
1832 os << "\n";
1833 }
1834 if (is_star_inconsistent) ++num_inconsistent_stars;
1835 }
1836
1837#ifdef DEBUG_TRACES
1838 std::cerr << "\n==========================================================\n"
1839 << "Export from list of stars to OFF:\n"
1840 << " * Number of vertices: " << m_points.size() << "\n"
1841 << " * Total number of maximal simplices: " << num_maximal_simplices << "\n";
1842 if (color_inconsistencies) {
1843 std::cerr << " * Number of inconsistent stars: " << num_inconsistent_stars << " ("
1844 << (m_points.size() > 0 ? 100. * num_inconsistent_stars / m_points.size() : 0.) << "%)\n"
1845 << " * Number of inconsistent maximal simplices: " << num_inconsistent_maximal_simplices << " ("
1846 << (num_maximal_simplices > 0 ? 100. * num_inconsistent_maximal_simplices / num_maximal_simplices : 0.)
1847 << "%)\n";
1848 }
1849 std::cerr << "==========================================================\n";
1850#endif
1851
1852 return os;
1853 }
1854
1855 public:
1856 std::ostream &export_simplices_to_off(const Simplicial_complex &complex, std::ostream &os,
1857 std::size_t &num_OFF_simplices,
1858 Simplex_set const *p_simpl_to_color_in_red = NULL,
1859 Simplex_set const *p_simpl_to_color_in_green = NULL,
1860 Simplex_set const *p_simpl_to_color_in_blue = NULL) const {
1861 typedef Simplicial_complex::Simplex Simplex;
1862 typedef Simplicial_complex::Simplex_set Simplex_set;
1863
1864 // If m_intrinsic_dim = 1, each point is output two times
1865 // (see export_vertices_to_off)
1866 num_OFF_simplices = 0;
1867 std::size_t num_maximal_simplices = 0;
1868
1869 typename Simplex_set::const_iterator it_s = complex.simplex_range().begin();
1870 typename Simplex_set::const_iterator it_s_end = complex.simplex_range().end();
1871 // For each simplex
1872 for (; it_s != it_s_end; ++it_s) {
1873 Simplex c = *it_s;
1874 ++num_maximal_simplices;
1875
1876 int color_simplex = -1; // -1=no color, 0=yellow, 1=red, 2=green, 3=blue
1877 if (p_simpl_to_color_in_red && std::find(p_simpl_to_color_in_red->begin(), p_simpl_to_color_in_red->end(), c) !=
1878 p_simpl_to_color_in_red->end()) {
1879 color_simplex = 1;
1880 } else if (p_simpl_to_color_in_green &&
1881 std::find(p_simpl_to_color_in_green->begin(), p_simpl_to_color_in_green->end(), c) !=
1882 p_simpl_to_color_in_green->end()) {
1883 color_simplex = 2;
1884 } else if (p_simpl_to_color_in_blue &&
1885 std::find(p_simpl_to_color_in_blue->begin(), p_simpl_to_color_in_blue->end(), c) !=
1886 p_simpl_to_color_in_blue->end()) {
1887 color_simplex = 3;
1888 }
1889
1890 // Gather the triangles here
1891 typedef std::vector<Simplex> Triangles;
1892 Triangles triangles;
1893
1894 int num_vertices = static_cast<int>(c.size());
1895 // Do not export smaller dimension simplices
1896 if (num_vertices < m_intrinsic_dim + 1) continue;
1897
1898 // If m_intrinsic_dim = 1, each point is output two times,
1899 // so we need to multiply each index by 2
1900 // And if only 2 vertices, add a third one (each vertex is duplicated in
1901 // the file when m_intrinsic dim = 2)
1902 if (m_intrinsic_dim == 1) {
1903 Simplex tmp_c;
1904 Simplex::iterator it = c.begin();
1905 for (; it != c.end(); ++it) tmp_c.insert(*it * 2);
1906 if (num_vertices == 2) tmp_c.insert(*tmp_c.rbegin() + 1);
1907
1908 c = tmp_c;
1909 }
1910
1911 if (num_vertices <= 3) {
1912 triangles.push_back(c);
1913 } else {
1914 // num_vertices >= 4: decompose the simplex in triangles
1915 std::vector<bool> booleans(num_vertices, false);
1916 std::fill(booleans.begin() + num_vertices - 3, booleans.end(), true);
1917 do {
1918 Simplex triangle;
1919 Simplex::iterator it = c.begin();
1920 for (int i = 0; it != c.end(); ++i, ++it) {
1921 if (booleans[i]) triangle.insert(*it);
1922 }
1923 triangles.push_back(triangle);
1924 } while (std::next_permutation(booleans.begin(), booleans.end()));
1925 }
1926
1927 // For each cell
1928 Triangles::const_iterator it_tri = triangles.begin();
1929 Triangles::const_iterator it_tri_end = triangles.end();
1930 for (; it_tri != it_tri_end; ++it_tri) {
1931 // Don't export infinite cells
1932 if (is_infinite(*it_tri)) continue;
1933
1934 os << 3 << " ";
1935 Simplex::const_iterator it_point_idx = it_tri->begin();
1936 for (; it_point_idx != it_tri->end(); ++it_point_idx) {
1937 os << *it_point_idx << " ";
1938 }
1939
1940 if (p_simpl_to_color_in_red || p_simpl_to_color_in_green || p_simpl_to_color_in_blue) {
1941 switch (color_simplex) {
1942 case 0:
1943 os << " 255 255 0";
1944 break;
1945 case 1:
1946 os << " 255 0 0";
1947 break;
1948 case 2:
1949 os << " 0 255 0";
1950 break;
1951 case 3:
1952 os << " 0 0 255";
1953 break;
1954 default:
1955 os << " 128 128 128";
1956 break;
1957 }
1958 }
1959
1960 ++num_OFF_simplices;
1961 os << "\n";
1962 }
1963 }
1964
1965#ifdef DEBUG_TRACES
1966 std::cerr << "\n==========================================================\n"
1967 << "Export from complex to OFF:\n"
1968 << " * Number of vertices: " << m_points.size() << "\n"
1969 << " * Total number of maximal simplices: " << num_maximal_simplices << "\n"
1970 << "==========================================================\n";
1971#endif
1972
1973 return os;
1974 }
1975
1983 void set_max_squared_edge_length(FT max_squared_edge_length) { m_max_squared_edge_length = max_squared_edge_length; }
1984
1985 private:
1986 const K m_k;
1987 const int m_intrinsic_dim;
1988 const int m_ambient_dim;
1989
1990 Points m_points;
1991 Weights m_weights;
1992#ifdef GUDHI_TC_PERTURB_POSITION
1993 Translations_for_perturb m_translations;
1994#if defined(GUDHI_USE_TBB)
1995 Mutex_for_perturb *m_p_perturb_mutexes;
1996#endif
1997#endif
1998
1999 Points_ds m_points_ds;
2000 double m_last_max_perturb;
2001 std::vector<bool> m_are_tangent_spaces_computed;
2002 TS_container m_tangent_spaces;
2003#ifdef GUDHI_TC_EXPORT_NORMALS
2004 OS_container m_orth_spaces;
2005#endif
2006 Tr_container m_triangulations; // Contains the triangulations
2007 // and their center vertex
2008 Stars_container m_stars;
2009 std::vector<FT> m_squared_star_spheres_radii_incl_margin;
2010 std::optional<FT> m_max_squared_edge_length;
2011
2012#ifdef GUDHI_TC_USE_ANOTHER_POINT_SET_FOR_TANGENT_SPACE_ESTIM
2013 Points m_points_for_tse;
2014 Points_ds m_points_ds_for_tse;
2015#endif
2016}; // /class Tangential_complex
2017
2018} // end namespace tangential_complex
2019} // end namespace Gudhi
2020
2021#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
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:1983
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:559
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:618
~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:547
std::size_t num_inconsistent_simplices
Number of inconsistent simplices.
Definition: Tangential_complex.h:551
std::size_t num_inconsistent_stars
Number of stars containing at least one inconsistent simplex.
Definition: Tangential_complex.h:553
std::size_t num_simplices
Total number of simplices in stars (including duplicates that appear in several stars)
Definition: Tangential_complex.h:549