12 #ifndef TANGENTIAL_COMPLEX_H_
13 #define TANGENTIAL_COMPLEX_H_
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>
24 #include <CGAL/Default.h>
25 #include <CGAL/Dimension.h>
26 #include <CGAL/function_objects.h>
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>
36 #include <Eigen/Eigen>
37 #include <Eigen/src/Core/util/Macros.h>
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>
62 #include <tbb/parallel_for.h>
63 #include <tbb/combinable.h>
70 #if CGAL_VERSION_NR < 1041101000
71 # error Tangential_complex is only available for CGAL >= 4.11
74 #if !EIGEN_VERSION_AT_LEAST(3,1,0)
75 # error Tangential_complex is only available for Eigen3 >= 3.1.0 installed with CGAL
78 namespace sps = Gudhi::spatial_searching;
82 namespace tangential_complex {
84 using namespace internal;
88 Vertex_data(std::size_t data = (std::numeric_limits<std::size_t>::max)()) : m_data(data) {}
90 operator std::size_t() {
return m_data; }
92 operator std::size_t()
const {
return m_data; }
123 template <
typename Kernel_,
124 typename DimensionTag,
125 typename Concurrency_tag = CGAL::Parallel_tag,
typename Triangulation_ = CGAL::Default>
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;
133 typedef typename CGAL::Default::Get<
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> >,
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;
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;
155 typedef Vector Translation_for_perturb;
156 typedef std::vector<FT> Weights;
158 typedef std::vector<Translation_for_perturb> Translations_for_perturb;
164 Tr_and_VH() : m_tr(NULL) {}
166 Tr_and_VH(
int dim) : m_tr(
new Triangulation(dim)) {}
168 ~Tr_and_VH() { destroy_triangulation(); }
170 Triangulation &construct_triangulation(
int dim) {
172 m_tr =
new Triangulation(dim);
176 void destroy_triangulation() {
181 Triangulation &tr() {
return *m_tr; }
183 Triangulation
const &tr()
const {
return *m_tr; }
185 Tr_vertex_handle
const ¢er_vertex()
const {
return m_center_vertex; }
187 Tr_vertex_handle ¢er_vertex() {
return m_center_vertex; }
191 Tr_vertex_handle m_center_vertex;
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;
200 typedef std::vector<Point> Points;
202 typedef boost::container::flat_set<std::size_t> Simplex;
203 typedef std::set<Simplex> Simplex_set;
210 typedef std::vector<Tr_and_VH> Tr_container;
211 typedef std::vector<Vector> Vectors;
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;
221 static const Tr_point &vertex_handle_to_point(Tr_vertex_handle vh) {
return vh->point(); }
223 template <
typename P,
typename VH>
224 static const P &vertex_handle_to_point(VH vh) {
229 typedef internal::Simplicial_complex Simplicial_complex;
240 template <
typename Po
int_range>
242 #ifdef GUDHI_TC_USE_ANOTHER_POINT_SET_FOR_TANGENT_SPACE_ESTIM
243 InputIterator first_for_tse, InputIterator last_for_tse,
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)
253 m_p_perturb_mutexes(NULL)
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
262 m_orth_spaces(m_points.size(), Orthogonal_space_basis())
264 #ifdef GUDHI_TC_USE_ANOTHER_POINT_SET_FOR_TANGENT_SPACE_ESTIM
266 m_points_for_tse(first_for_tse, last_for_tse),
267 m_points_ds_for_tse(m_points_for_tse)
274 #if defined(GUDHI_USE_TBB) && defined(GUDHI_TC_PERTURB_POSITION)
275 delete[] m_p_perturb_mutexes;
285 Points
const &points()
const {
return m_points; }
292 Point
get_point(std::size_t vertex)
const {
return m_points[vertex]; }
305 void set_weights(
const Weights &weights) { m_weights = weights; }
307 void set_tangent_planes(
const TS_container &tangent_spaces
308 #ifdef GUDHI_TC_EXPORT_NORMALS
310 const OS_container &orthogonal_spaces
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"));
317 GUDHI_CHECK(m_points.size() == tangent_spaces.size(), std::logic_error(
"Wrong sizes"));
319 m_tangent_spaces = tangent_spaces;
320 #ifdef GUDHI_TC_EXPORT_NORMALS
321 m_orth_spaces = orthogonal_spaces;
323 for (std::size_t i = 0; i < m_points.size(); ++i) m_are_tangent_spaces_computed[i] =
true;
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"
338 #if defined(GUDHI_TC_PROFILING) && defined(GUDHI_USE_TBB)
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();
352 m_translations.resize(m_points.size(), m_k.construct_vector_d_object()(m_ambient_dim));
354 #if defined(GUDHI_USE_TBB)
355 delete[] m_p_perturb_mutexes;
356 m_p_perturb_mutexes =
new Mutex_for_perturb[m_points.size()];
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));
367 for (std::size_t i = 0; i < m_points.size(); ++i) compute_tangent_triangulation(i);
372 #if defined(GUDHI_TC_PROFILING) && defined(GUDHI_USE_TBB)
374 std::cerr <<
"Tangential complex computed in " << t.num_seconds() <<
" seconds.\n";
381 bool success =
false;
383 unsigned int num_steps = 0;
385 std::size_t initial_num_inconsistent_stars = 0;
387 std::size_t best_num_inconsistent_stars = 0;
389 std::size_t final_num_inconsistent_stars = 0;
400 if (time_limit == 0.)
return info;
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);
407 if (std::get<1>(stats_before) == 0) {
409 std::cerr <<
"Nothing to fix.\n";
416 m_last_max_perturb = max_perturb;
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";
431 #if defined(DEBUG_TRACES) || defined(GUDHI_TC_PROFILING)
432 std::cerr << yellow <<
"\nAttempt to fix inconsistencies using perturbations - step #" << info.
num_steps + 1
436 std::size_t num_inconsistent_stars = 0;
437 std::vector<std::size_t> updated_points;
439 #ifdef GUDHI_TC_PROFILING
440 Gudhi::Clock t_fix_step;
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;
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()));
457 num_inconsistent_stars = num_inconsistencies.combine(std::plus<std::size_t>());
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());
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));
473 #if defined(GUDHI_USE_TBB)
477 #ifdef GUDHI_TC_PROFILING
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";
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;
493 if (num_inconsistent_stars > 0) refresh_tangential_complex(updated_points);
495 #ifdef GUDHI_TC_PERFORM_EXTRA_CHECKS
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;
503 std::cerr << green <<
"REFRESHMENT CHECK: PASSED.\n" << white;
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);
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";
516 stats_before = stats_after;
526 done = (num_inconsistent_stars == 0);
529 if (time_limit > 0. && t.num_seconds() > time_limit) {
531 std::cerr << red <<
"Time limit reached.\n" << white;
540 std::cerr << green <<
"Fixed!\n" << white;
549 std::size_t num_simplices = 0;
551 std::size_t num_inconsistent_simplices = 0;
553 std::size_t num_inconsistent_stars = 0;
569 for (std::size_t idx = 0; idx < m_points.size(); ++idx) {
570 bool is_star_inconsistent =
false;
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) {
577 if (is_infinite(*it_inc_simplex))
continue;
579 Simplex c = *it_inc_simplex;
582 if (!is_simplex_consistent(c)) {
584 is_star_inconsistent =
true;
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): "
601 <<
"==========================================================\n";
617 template <
typename Simplex_tree_>
619 bool export_inconsistent_simplices =
true
622 bool export_infinite_simplices =
false, Simplex_set *p_inconsistent_simplices = NULL
625 #if defined(DEBUG_TRACES) || defined(GUDHI_TC_PROFILING)
626 std::cerr << yellow <<
"\nExporting the TC as a Simplex_tree... " << white;
628 #ifdef GUDHI_TC_PROFILING
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);
640 for (std::size_t idx = 0; idx < m_points.size(); ++idx) {
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;
648 if (!export_infinite_simplices && is_infinite(c))
continue;
650 if (
static_cast<int>(c.size()) > max_dim) max_dim =
static_cast<int>(c.size());
654 if (!export_inconsistent_simplices && !is_simplex_consistent(c))
continue;
657 bool inserted = tree.insert_simplex_and_subfaces(c).second;
660 if (p_inconsistent_simplices && inserted && !is_simplex_consistent(c)) {
661 p_inconsistent_simplices->insert(c);
666 #ifdef GUDHI_TC_PROFILING
668 std::cerr << yellow <<
"done in " << t.num_seconds() <<
" seconds.\n" << white;
669 #elif defined(DEBUG_TRACES)
670 std::cerr << yellow <<
"done.\n" << white;
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;
692 #ifdef GUDHI_TC_PROFILING
700 for (std::size_t idx = 0; idx < m_points.size(); ++idx) {
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;
708 if (!export_infinite_simplices && is_infinite(c))
continue;
710 if (
static_cast<int>(c.size()) > max_dim) max_dim =
static_cast<int>(c.size());
714 if (!export_inconsistent_simplices && !is_simplex_consistent(c))
continue;
717 if (check_lower_and_higher_dim_simplices == 2 && max_dim != -1 &&
static_cast<int>(c.size()) != max_dim) {
720 <<
"Info: check_lower_and_higher_dim_simplices ACTIVATED. "
721 "Export might be take some time...\n"
723 check_lower_and_higher_dim_simplices = 1;
727 bool added = complex.add_simplex(c, check_lower_and_higher_dim_simplices == 1);
730 if (p_inconsistent_simplices && added && !is_simplex_consistent(c)) {
731 p_inconsistent_simplices->insert(c);
736 #ifdef GUDHI_TC_PROFILING
738 std::cerr << yellow <<
"done in " << t.num_seconds() <<
" seconds.\n" << white;
739 #elif defined(DEBUG_TRACES)
740 std::cerr << yellow <<
"done.\n" << white;
746 template <
typename ProjectionFunctor = CGAL::Identity<Po
int> >
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);
756 template <
typename ProjectionFunctor = CGAL::Identity<Po
int> >
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;
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";
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";
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";
783 std::stringstream output;
784 std::size_t num_simplices, num_vertices;
785 export_vertices_to_off(output, num_vertices,
false, point_projection);
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);
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);
794 #ifdef GUDHI_TC_EXPORT_NORMALS
799 << num_vertices <<
" " << num_simplices <<
" "
807 void refresh_tangential_complex() {
808 #if defined(DEBUG_TRACES) || defined(GUDHI_TC_PROFILING)
809 std::cerr << yellow <<
"\nRefreshing TC... " << white;
812 #ifdef GUDHI_TC_PROFILING
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));
822 for (std::size_t i = 0; i < m_points.size(); ++i) compute_tangent_triangulation(i);
827 #ifdef GUDHI_TC_PROFILING
829 std::cerr << yellow <<
"done in " << t.num_seconds() <<
" seconds.\n" << white;
830 #elif defined(DEBUG_TRACES)
831 std::cerr << yellow <<
"done.\n" << white;
836 template <
typename Po
int_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;
842 #ifdef GUDHI_TC_PROFILING
847 Points_ds updated_pts_ds(m_points, perturbed_points_indices);
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));
857 for (std::size_t i = 0; i < m_points.size(); ++i) refresh_tangent_triangulation(i, updated_pts_ds);
862 #ifdef GUDHI_TC_PROFILING
864 std::cerr << yellow <<
"done in " << t.num_seconds() <<
" seconds.\n" << white;
865 #elif defined(DEBUG_TRACES)
866 std::cerr << yellow <<
"done.\n" << white;
870 void export_inconsistent_stars_to_OFF_files(std::string
const &filename_base)
const {
872 for (std::size_t idx = 0; idx < m_points.size(); ++idx) {
874 Simplicial_complex sc;
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) {
881 if (is_infinite(*it_inc_simplex))
continue;
883 Simplex c = *it_inc_simplex;
889 if (!is_inconsistent && !is_simplex_consistent(c)) is_inconsistent =
true;
892 if (is_inconsistent) {
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);
902 class Compare_distance_to_ref_point {
904 Compare_distance_to_ref_point(Point
const &ref, K
const &k) : m_ref(ref), m_k(k) {}
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);
918 class Compute_tangent_triangulation {
919 Tangential_complex &m_tc;
923 Compute_tangent_triangulation(Tangential_complex &tc) : m_tc(tc) {}
926 Compute_tangent_triangulation(
const Compute_tangent_triangulation &ctt) : m_tc(ctt.m_tc) {}
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);
935 class Refresh_tangent_triangulation {
936 Tangential_complex &m_tc;
937 Points_ds
const &m_updated_pts_ds;
941 Refresh_tangent_triangulation(Tangential_complex &tc, Points_ds
const &updated_pts_ds)
942 : m_tc(tc), m_updated_pts_ds(updated_pts_ds) {}
945 Refresh_tangent_triangulation(
const Refresh_tangent_triangulation &ctt)
946 : m_tc(ctt.m_tc), m_updated_pts_ds(ctt.m_updated_pts_ds) {}
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);
955 bool is_infinite(Simplex
const &s)
const {
return *s.rbegin() == (std::numeric_limits<std::size_t>::max)(); }
960 Tr_vertex_handle compute_star(std::size_t i,
const Point ¢er_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();
966 typename K::Squared_distance_d k_sqdist = m_k.squared_distance_d_object();
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();
973 typename Tr_traits::Construct_power_sphere_d power_center = local_tr_traits.construct_power_sphere_d_object();
983 if (i == tsb.origin()) {
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]);
988 const Weighted_point &wp = compute_perturbed_weighted_point(i);
989 proj_wp = project_point_and_compute_weight(wp, tsb, local_tr_traits);
992 Tr_vertex_handle center_vertex = triangulation.insert(proj_wp);
993 center_vertex->data() = i;
994 if (verbose) std::cerr <<
"* Inserted point #" << i <<
"\n";
996 #ifdef GUDHI_TC_VERY_VERBOSE
997 std::size_t num_attempts_to_insert_points = 1;
998 std::size_t num_inserted_points = 1;
1002 INS_range ins_range = m_points_ds.incremental_nearest_neighbors(center_pt);
1010 std::optional<FT> squared_star_sphere_radius_plus_margin = m_max_squared_edge_length;
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;
1017 if (neighbor_point_idx != i) {
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())));
1034 Tr_point proj_pt = project_point_and_compute_weight(neighbor_pt, neighbor_weight, tsb, local_tr_traits);
1036 #ifdef GUDHI_TC_VERY_VERBOSE
1037 ++num_attempts_to_insert_points;
1040 Tr_vertex_handle vh = triangulation.insert_if_in_star(proj_pt, center_vertex);
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;
1046 if (verbose) std::cerr <<
"* Inserted point #" << neighbor_point_idx <<
"\n";
1048 vh->data() = neighbor_point_idx;
1051 if (triangulation.current_dimension() >= tangent_space_dim) {
1052 squared_star_sphere_radius_plus_margin = std::nullopt;
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;
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>));
1071 FT sq_power_sphere_diam = 4 * point_weight(c);
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;
1082 if (squared_star_sphere_radius_plus_margin) {
1084 squared_star_sphere_radius_plus_margin =
1085 CGAL::square(std::sqrt(squared_star_sphere_radius_plus_margin.value()) + 2 * m_last_max_perturb);
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();
1093 m_squared_star_spheres_radii_incl_margin[i] = squared_star_sphere_radius_plus_margin.value();
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();
1099 m_squared_star_spheres_radii_incl_margin[i] = FT(-1);
1107 return center_vertex;
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";
1113 if (m_squared_star_spheres_radii_incl_margin[i] == FT(-1))
return compute_tangent_triangulation(i, verbose);
1115 Point center_point = compute_perturbed_point(i);
1117 std::size_t closest_pt_index = updated_pts_ds.k_nearest_neighbors(center_point, 1,
false).begin()->first;
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();
1123 typename K::Compute_power_product_d k_power_dist = m_k.compute_power_product_d_object();
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);
1131 if (k_power_dist(star_sphere, closest_updated_point) <= FT(0)) compute_tangent_triangulation(i, verbose);
1134 void compute_tangent_triangulation(std::size_t i,
bool verbose =
false) {
1135 if (verbose) std::cerr <<
"** Computing tangent tri #" << i <<
" **\n";
1140 const Point center_pt = compute_perturbed_point(i);
1141 Tangent_space_basis &tsb = m_tangent_spaces[i];
1144 if (!m_are_tangent_spaces_computed[i]) {
1145 #ifdef GUDHI_TC_EXPORT_NORMALS
1146 tsb = compute_tangent_space(center_pt, i,
true , &m_orth_spaces[i]);
1148 tsb = compute_tangent_space(center_pt, i);
1152 #if defined(GUDHI_TC_PROFILING) && defined(GUDHI_TC_VERY_VERBOSE)
1155 int tangent_space_dim = tangent_basis_dim(i);
1156 Triangulation &local_tr = m_triangulations[i].construct_triangulation(tangent_space_dim);
1158 m_triangulations[i].center_vertex() = compute_star(i, center_pt, tsb, local_tr, verbose);
1160 #if defined(GUDHI_TC_PROFILING) && defined(GUDHI_TC_VERY_VERBOSE)
1162 std::cerr <<
" - triangulation construction: " << t.num_seconds() <<
" s.\n";
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";
1173 #if defined(GUDHI_TC_PROFILING) && defined(GUDHI_TC_VERY_VERBOSE)
1175 std::cerr <<
" - update_star: " << t.num_seconds() <<
" s.\n";
1181 void update_star(std::size_t i) {
1182 Star &star = m_stars[i];
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;
1188 std::vector<Tr_full_cell_handle> incident_cells;
1189 local_tr.incident_full_cells(center_vertex, std::back_inserter(incident_cells));
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();
1194 for (; it_c != it_c_end; ++it_c) {
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);
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);
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()));
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();
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;
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;
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));
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);
1239 Tangent_space_basis tsb(i);
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));
1249 tsb.push_back(constr_vec(m_ambient_dim, eig.eigenvectors().col(j).data(),
1250 eig.eigenvectors().col(j).data() + m_ambient_dim));
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));
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));
1268 m_are_tangent_spaces_computed[i] =
true;
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()));
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();
1292 Eigen::MatrixXd mat_points(s.size() * num_pts_for_pca, m_ambient_dim);
1293 unsigned int current_row = 0;
1295 for (Simplex::const_iterator it_index = s.begin(); it_index != s.end(); ++it_index) {
1296 const Point &p = m_points[*it_index];
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;
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;
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));
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);
1317 Tangent_space_basis tsb;
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));
1327 tsb.push_back(constr_vec(m_ambient_dim, eig.eigenvectors().col(j).data(),
1328 eig.eigenvectors().col(j).data() + m_ambient_dim));
1337 int tangent_basis_dim(std::size_t i)
const {
return m_tangent_spaces[i].dimension(); }
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]);
1343 return m_points[pt_idx];
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]);
1351 p = m_points[pt_idx];
1353 w = m_weights[pt_idx];
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();
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]),
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();
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)));
1378 return global_point;
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();
1387 Vector v = diff_points(p, compute_perturbed_point(tsb.origin()));
1389 std::vector<FT> coords;
1391 coords.reserve(tsb.dimension());
1392 for (std::size_t i = 0; i < m_intrinsic_dim; ++i) {
1394 FT coord = scalar_pdct(v, tsb[i]);
1395 coords.push_back(coord);
1398 return tr_traits.construct_point_d_object()(
static_cast<int>(coords.size()), coords.begin(), coords.end());
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);
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);
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();
1423 Point origin = compute_perturbed_point(tsb.origin());
1424 Vector v = diff_points(p, origin);
1427 bool same_dim = (point_dim == tsb.dimension());
1429 std::vector<FT> coords;
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) {
1435 FT c = scalar_pdct(v, tsb[i]);
1436 coords.push_back(c);
1440 for (
int j = 0; j < point_dim; ++j) p_proj[j] += c * coord(tsb[i], j);
1445 FT sq_dist_to_proj_pt = 0;
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);
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);
1458 template <
typename Indexed_po
int_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));
1472 bool is_simplex_consistent(Tr_full_cell_handle fch,
int cur_dim)
const {
1474 for (
int i = 0; i < cur_dim + 1; ++i) {
1475 std::size_t data = fch->vertex(i)->data();
1478 return is_simplex_consistent(c);
1484 bool is_simplex_consistent(Simplex
const &simplex)
const {
1486 Simplex::const_iterator it_point_idx = simplex.begin();
1489 for (; it_point_idx != simplex.end(); ++it_point_idx) {
1490 std::size_t point_idx = *it_point_idx;
1492 if (point_idx == (std::numeric_limits<std::size_t>::max)())
continue;
1494 Star
const &star = m_stars[point_idx];
1497 Incident_simplex is_to_find = simplex;
1498 is_to_find.erase(point_idx);
1501 if (std::find(star.begin(), star.end(), is_to_find) == star.end())
return false;
1513 template <
typename OutputIterator>
1514 bool is_simplex_consistent(std::size_t center_point,
1515 Incident_simplex
const &s,
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);
1522 Incident_simplex::const_iterator it_point_idx = s.begin();
1525 for (; it_point_idx != s.end(); ++it_point_idx) {
1526 std::size_t point_idx = *it_point_idx;
1528 if (point_idx == (std::numeric_limits<std::size_t>::max)())
continue;
1530 Star
const &star = m_stars[point_idx];
1533 Incident_simplex is_to_find = full_simplex;
1534 is_to_find.erase(point_idx);
1536 if (check_also_in_non_maximal_faces) {
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;
1544 if (!found) *points_whose_star_does_not_contain_s++ = point_idx;
1547 if (std::find(star.begin(), star.end(), is_to_find) == star.end())
1548 *points_whose_star_does_not_contain_s++ = point_idx;
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];
1561 if (check_also_in_non_maximal_faces) {
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;
1571 return !(std::find(star.begin(), star.end(), s) == star.end());
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();
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));
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)));
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();
1599 m_translations[point_idx] = global_transl;
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;
1609 Star
const &star = m_stars[tr_index];
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;
1618 if (is_infinite(incident_simplex))
continue;
1620 Simplex c = incident_simplex;
1624 if (!is_simplex_consistent(c)) {
1625 is_inconsistent =
true;
1627 std::size_t idx = tr_index;
1629 perturb(tr_index, max_perturb);
1630 *perturbed_pts_indices++ = idx;
1637 return is_inconsistent;
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));
1649 return export_point_set(m_k, perturbed_points, os, coord_separator);
1651 return export_point_set(m_k, m_points, os, coord_separator);
1655 template <
typename ProjectionFunctor = CGAL::Identity<Po
int> >
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()) {
1666 const int N = (m_intrinsic_dim == 1 ? 2 : 1);
1669 typename K::Compute_coordinate_d coord = m_k.compute_coordinate_d_object();
1671 #ifdef GUDHI_TC_EXPORT_ALL_COORDS_IN_OFF
1672 int num_coords = m_ambient_dim;
1674 int num_coords = (std::min)(m_ambient_dim, 3);
1677 #ifdef GUDHI_TC_EXPORT_NORMALS
1678 OS_container::const_iterator it_os = m_orth_spaces.begin();
1680 typename Points::const_iterator it_p = m_points.begin();
1681 typename Points::const_iterator it_p_end = m_points.end();
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) {
1687 for (; j < num_coords; ++j) os << CGAL::to_double(coord(p, j)) <<
" ";
1688 if (j == 2) os <<
"0";
1690 #ifdef GUDHI_TC_EXPORT_NORMALS
1691 for (j = 0; j < num_coords; ++j) os <<
" " << CGAL::to_double(coord(*it_os->begin(), j));
1695 #ifdef GUDHI_TC_EXPORT_NORMALS
1700 num_vertices = N * m_points.size();
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 {
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();
1718 for (std::size_t idx = 0; it_tr != it_tr_end; ++it_tr, ++idx) {
1719 bool is_star_inconsistent =
false;
1721 Triangulation
const &tr = it_tr->tr();
1723 if (tr.current_dimension() < m_intrinsic_dim)
continue;
1726 std::stringstream color;
1728 color << 128 <<
" " << 128 <<
" " << 128;
1731 typedef std::vector<std::pair<Simplex, int> > Star_using_triangles;
1732 Star_using_triangles star_using_triangles;
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;
1740 std::size_t num_vertices = c.size();
1741 ++num_maximal_simplices;
1743 int color_simplex = -1;
1744 if (color_inconsistencies && !is_simplex_consistent(c)) {
1745 ++num_inconsistent_maximal_simplices;
1747 is_star_inconsistent =
true;
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()) {
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()) {
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()) {
1767 if (m_intrinsic_dim == 1) {
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);
1776 if (num_vertices <= 3) {
1777 star_using_triangles.push_back(std::make_pair(c, color_simplex));
1780 std::vector<bool> booleans(num_vertices,
false);
1781 std::fill(booleans.begin() + num_vertices - 3, booleans.end(),
true);
1784 Simplex::iterator it = c.begin();
1785 for (
int i = 0; it != c.end(); ++i, ++it) {
1786 if (booleans[i]) triangle.insert(*it);
1788 star_using_triangles.push_back(std::make_pair(triangle, color_simplex));
1789 }
while (std::next_permutation(booleans.begin(), booleans.end()));
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;
1800 if (is_infinite(c))
continue;
1802 int color_simplex = it_simplex->second;
1804 std::stringstream sstr_c;
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 <<
" ";
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) {
1827 os <<
" " << color.str();
1831 ++num_OFF_simplices;
1834 if (is_star_inconsistent) ++num_inconsistent_stars;
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.)
1849 std::cerr <<
"==========================================================\n";
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;
1866 num_OFF_simplices = 0;
1867 std::size_t num_maximal_simplices = 0;
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();
1872 for (; it_s != it_s_end; ++it_s) {
1874 ++num_maximal_simplices;
1876 int color_simplex = -1;
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()) {
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()) {
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()) {
1891 typedef std::vector<Simplex> Triangles;
1892 Triangles triangles;
1894 int num_vertices =
static_cast<int>(c.size());
1896 if (num_vertices < m_intrinsic_dim + 1)
continue;
1902 if (m_intrinsic_dim == 1) {
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);
1911 if (num_vertices <= 3) {
1912 triangles.push_back(c);
1915 std::vector<bool> booleans(num_vertices,
false);
1916 std::fill(booleans.begin() + num_vertices - 3, booleans.end(),
true);
1919 Simplex::iterator it = c.begin();
1920 for (
int i = 0; it != c.end(); ++i, ++it) {
1921 if (booleans[i]) triangle.insert(*it);
1923 triangles.push_back(triangle);
1924 }
while (std::next_permutation(booleans.begin(), booleans.end()));
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) {
1932 if (is_infinite(*it_tri))
continue;
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 <<
" ";
1940 if (p_simpl_to_color_in_red || p_simpl_to_color_in_green || p_simpl_to_color_in_blue) {
1941 switch (color_simplex) {
1955 os <<
" 128 128 128";
1960 ++num_OFF_simplices;
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";
1987 const int m_intrinsic_dim;
1988 const int m_ambient_dim;
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;
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;
2006 Tr_container m_triangulations;
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;
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;
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
Gudhi namespace.
Definition: SimplicialComplexForAlpha.h:14
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