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>   
   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 
   78namespace sps = Gudhi::spatial_searching;
 
   82namespace tangential_complex {
 
   84using 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; }
 
  123template <
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,
 
  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";
 
  400    if (time_limit == 0.) 
return info;
 
  404#ifdef GUDHI_TC_SHOW_DETAILED_STATS_FOR_INCONSISTENCIES 
  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;
 
  448        tbb::parallel_for(tbb::blocked_range<size_t>(0, m_triangulations.size()),
 
  449                          Try_to_solve_inconsistencies_in_a_local_triangulation(*
this, max_perturb, num_inconsistencies,
 
  450                                                                                tls_updated_points));
 
  451        num_inconsistent_stars = num_inconsistencies.combine(std::plus<std::size_t>());
 
  453            tls_updated_points.combine([](std::vector<std::size_t> 
const &x, std::vector<std::size_t> 
const &y) {
 
  454              std::vector<std::size_t> res;
 
  455              res.reserve(x.size() + y.size());
 
  456              res.insert(res.end(), x.begin(), x.end());
 
  457              res.insert(res.end(), y.begin(), y.end());
 
  463        for (std::size_t i = 0; i < m_triangulations.size(); ++i) {
 
  464          num_inconsistent_stars +=
 
  465              try_to_solve_inconsistencies_in_a_local_triangulation(i, max_perturb, std::back_inserter(updated_points));
 
  467#if defined(GUDHI_USE_TBB) 
  471#ifdef GUDHI_TC_PROFILING 
  475#if defined(GUDHI_TC_SHOW_DETAILED_STATS_FOR_INCONSISTENCIES) || defined(DEBUG_TRACES) 
  476      std::cerr << 
"\nEncountered during fix:\n" 
  477                << 
"  * Num stars containing inconsistent simplices: " << red << num_inconsistent_stars << white << 
" (" 
  478                << 100. * num_inconsistent_stars / m_points.size() << 
"%)\n";
 
  481#ifdef GUDHI_TC_PROFILING 
  482      std::cerr << yellow << 
"done in " << t_fix_step.num_seconds() << 
" seconds.\n" << white;
 
  483#elif defined(DEBUG_TRACES) 
  484      std::cerr << yellow << 
"done.\n" << white;
 
  487      if (num_inconsistent_stars > 0) refresh_tangential_complex(updated_points);
 
  489#ifdef GUDHI_TC_PERFORM_EXTRA_CHECKS 
  492      refresh_tangential_complex();
 
  494      if (num_inc_1 != num_inc_2)
 
  495        std::cerr << red << 
"REFRESHMENT CHECK: FAILED. (" << num_inc_1 << 
" vs " << num_inc_2 << 
")\n" << white;
 
  497        std::cerr << green << 
"REFRESHMENT CHECK: PASSED.\n" << white;
 
  500#ifdef GUDHI_TC_SHOW_DETAILED_STATS_FOR_INCONSISTENCIES 
  503      std::cerr << 
"\nAfter fix:\n" 
  504                << 
"  * Total number of simplices in stars (incl. duplicates): " << std::get<0>(stats_after) << 
"\n" 
  505                << 
"  * Num inconsistent simplices in stars (incl. duplicates): " << red << std::get<1>(stats_after)
 
  506                << white << 
" (" << 100. * std::get<1>(stats_after) / std::get<0>(stats_after) << 
"%)\n" 
  507                << 
"  * Number of stars containing inconsistent simplices: " << red << std::get<2>(stats_after) << white
 
  508                << 
" (" << 100. * std::get<2>(stats_after) / m_points.size() << 
"%)\n";
 
  510      stats_before = stats_after;
 
  520      done = (num_inconsistent_stars == 0);
 
  523        if (time_limit > 0. && t.num_seconds() > time_limit) {
 
  525          std::cerr << red << 
"Time limit reached.\n" << white;
 
  534    std::cerr << green << 
"Fixed!\n" << white;
 
  563    for (std::size_t idx = 0; idx < m_points.size(); ++idx) {
 
  564      bool is_star_inconsistent = 
false;
 
  567      Star::const_iterator it_inc_simplex = m_stars[idx].begin();
 
  568      Star::const_iterator it_inc_simplex_end = m_stars[idx].end();
 
  569      for (; it_inc_simplex != it_inc_simplex_end; ++it_inc_simplex) {
 
  571        if (is_infinite(*it_inc_simplex)) 
continue;
 
  573        Simplex c = *it_inc_simplex;
 
  576        if (!is_simplex_consistent(c)) {
 
  578          is_star_inconsistent = 
true;
 
  587      std::cerr << 
"\n==========================================================\n" 
  588                << 
"Inconsistencies:\n" 
  589                << 
"  * Total number of simplices in stars (incl. duplicates): " << stats.
num_simplices << 
"\n" 
  590                << 
"  * Number of inconsistent simplices in stars (incl. duplicates): " 
  595                << 
"==========================================================\n";
 
  611  template <
typename Simplex_tree_>
 
  613                     bool export_inconsistent_simplices = 
true 
  616                     bool export_infinite_simplices = 
false, Simplex_set *p_inconsistent_simplices = NULL
 
  619#if defined(DEBUG_TRACES) || defined(GUDHI_TC_PROFILING) 
  620    std::cerr << yellow << 
"\nExporting the TC as a Simplex_tree... " << white;
 
  622#ifdef GUDHI_TC_PROFILING 
  629    std::vector<typename Simplex_tree_::Vertex_handle> vertices(m_points.size());
 
  630    std::iota(vertices.begin(), vertices.end(), 0);
 
  631    tree.insert_batch_vertices(vertices);
 
  634    for (std::size_t idx = 0; idx < m_points.size(); ++idx) {
 
  636      Star::const_iterator it_inc_simplex = m_stars[idx].begin();
 
  637      Star::const_iterator it_inc_simplex_end = m_stars[idx].end();
 
  638      for (; it_inc_simplex != it_inc_simplex_end; ++it_inc_simplex) {
 
  639        Simplex c = *it_inc_simplex;
 
  642        if (!export_infinite_simplices && is_infinite(c)) 
continue;
 
  644        if (
static_cast<int>(c.size()) > max_dim) max_dim = 
static_cast<int>(c.size());
 
  648        if (!export_inconsistent_simplices && !is_simplex_consistent(c)) 
continue;
 
  651        bool inserted = tree.insert_simplex_and_subfaces(c).second;
 
  654        if (p_inconsistent_simplices && inserted && !is_simplex_consistent(c)) {
 
  655          p_inconsistent_simplices->insert(c);
 
  660#ifdef GUDHI_TC_PROFILING 
  662    std::cerr << yellow << 
"done in " << t.num_seconds() << 
" seconds.\n" << white;
 
  663#elif defined(DEBUG_TRACES) 
  664    std::cerr << yellow << 
"done.\n" << white;
 
  680  int create_complex(Simplicial_complex &complex, 
bool export_inconsistent_simplices = 
true,
 
  681                     bool export_infinite_simplices = 
false, 
int check_lower_and_higher_dim_simplices = 2,
 
  682                     Simplex_set *p_inconsistent_simplices = NULL)
 const {
 
  683#if defined(DEBUG_TRACES) || defined(GUDHI_TC_PROFILING) 
  684    std::cerr << yellow << 
"\nExporting the TC as a Simplicial_complex... " << white;
 
  686#ifdef GUDHI_TC_PROFILING 
  694    for (std::size_t idx = 0; idx < m_points.size(); ++idx) {
 
  696      Star::const_iterator it_inc_simplex = m_stars[idx].begin();
 
  697      Star::const_iterator it_inc_simplex_end = m_stars[idx].end();
 
  698      for (; it_inc_simplex != it_inc_simplex_end; ++it_inc_simplex) {
 
  699        Simplex c = *it_inc_simplex;
 
  702        if (!export_infinite_simplices && is_infinite(c)) 
continue;
 
  704        if (
static_cast<int>(c.size()) > max_dim) max_dim = 
static_cast<int>(c.size());
 
  708        if (!export_inconsistent_simplices && !is_simplex_consistent(c)) 
continue;
 
  711        if (check_lower_and_higher_dim_simplices == 2 && max_dim != -1 && 
static_cast<int>(c.size()) != max_dim) {
 
  714                    << 
"Info: check_lower_and_higher_dim_simplices ACTIVATED. " 
  715                       "Export might be take some time...\n" 
  717          check_lower_and_higher_dim_simplices = 1;
 
  721        bool added = complex.add_simplex(c, check_lower_and_higher_dim_simplices == 1);
 
  724        if (p_inconsistent_simplices && added && !is_simplex_consistent(c)) {
 
  725          p_inconsistent_simplices->insert(c);
 
  730#ifdef GUDHI_TC_PROFILING 
  732    std::cerr << yellow << 
"done in " << t.num_seconds() << 
" seconds.\n" << white;
 
  733#elif defined(DEBUG_TRACES) 
  734    std::cerr << yellow << 
"done.\n" << white;
 
  740  template <
typename ProjectionFunctor = CGAL::Identity<Po
int> >
 
  741  std::ostream &export_to_off(
const Simplicial_complex &complex, std::ostream &os,
 
  742                              Simplex_set 
const *p_simpl_to_color_in_red = NULL,
 
  743                              Simplex_set 
const *p_simpl_to_color_in_green = NULL,
 
  744                              Simplex_set 
const *p_simpl_to_color_in_blue = NULL,
 
  745                              ProjectionFunctor 
const &point_projection = ProjectionFunctor())
 const {
 
  746    return export_to_off(os, 
false, p_simpl_to_color_in_red, p_simpl_to_color_in_green, p_simpl_to_color_in_blue,
 
  747                         &complex, point_projection);
 
  750  template <
typename ProjectionFunctor = CGAL::Identity<Po
int> >
 
  751  std::ostream &export_to_off(std::ostream &os, 
bool color_inconsistencies = 
false,
 
  752                              Simplex_set 
const *p_simpl_to_color_in_red = NULL,
 
  753                              Simplex_set 
const *p_simpl_to_color_in_green = NULL,
 
  754                              Simplex_set 
const *p_simpl_to_color_in_blue = NULL,
 
  755                              const Simplicial_complex *p_complex = NULL,
 
  756                              ProjectionFunctor 
const &point_projection = ProjectionFunctor())
 const {
 
  757    if (m_points.empty()) 
return os;
 
  759    if (m_ambient_dim < 2) {
 
  760      std::cerr << 
"Error: export_to_off => ambient dimension should be >= 2.\n";
 
  761      os << 
"Error: export_to_off => ambient dimension should be >= 2.\n";
 
  764    if (m_ambient_dim > 3) {
 
  765      std::cerr << 
"Warning: export_to_off => ambient dimension should be " 
  766                   "<= 3. Only the first 3 coordinates will be exported.\n";
 
  769    if (m_intrinsic_dim < 1 || m_intrinsic_dim > 3) {
 
  770      std::cerr << 
"Error: export_to_off => intrinsic dimension should be " 
  771                   "between 1 and 3.\n";
 
  772      os << 
"Error: export_to_off => intrinsic dimension should be " 
  773            "between 1 and 3.\n";
 
  777    std::stringstream output;
 
  778    std::size_t num_simplices, num_vertices;
 
  779    export_vertices_to_off(output, num_vertices, 
false, point_projection);
 
  781      export_simplices_to_off(*p_complex, output, num_simplices, p_simpl_to_color_in_red, p_simpl_to_color_in_green,
 
  782                              p_simpl_to_color_in_blue);
 
  784      export_simplices_to_off(output, num_simplices, color_inconsistencies, p_simpl_to_color_in_red,
 
  785                              p_simpl_to_color_in_green, p_simpl_to_color_in_blue);
 
  788#ifdef GUDHI_TC_EXPORT_NORMALS 
  793       << num_vertices << 
" " << num_simplices << 
" " 
  801  void refresh_tangential_complex() {
 
  802#if defined(DEBUG_TRACES) || defined(GUDHI_TC_PROFILING) 
  803    std::cerr << yellow << 
"\nRefreshing TC... " << white;
 
  806#ifdef GUDHI_TC_PROFILING 
  811    if (boost::is_convertible<Concurrency_tag, CGAL::Parallel_tag>::value) {
 
  812      tbb::parallel_for(tbb::blocked_range<size_t>(0, m_points.size()), Compute_tangent_triangulation(*
this));
 
  816      for (std::size_t i = 0; i < m_points.size(); ++i) compute_tangent_triangulation(i);
 
  821#ifdef GUDHI_TC_PROFILING 
  823    std::cerr << yellow << 
"done in " << t.num_seconds() << 
" seconds.\n" << white;
 
  824#elif defined(DEBUG_TRACES) 
  825    std::cerr << yellow << 
"done.\n" << white;
 
  830  template <
typename Po
int_indices_range>
 
  831  void refresh_tangential_complex(Point_indices_range 
const &perturbed_points_indices) {
 
  832#if defined(DEBUG_TRACES) || defined(GUDHI_TC_PROFILING) 
  833    std::cerr << yellow << 
"\nRefreshing TC... " << white;
 
  836#ifdef GUDHI_TC_PROFILING 
  841    Points_ds updated_pts_ds(m_points, perturbed_points_indices);
 
  845    if (boost::is_convertible<Concurrency_tag, CGAL::Parallel_tag>::value) {
 
  846      tbb::parallel_for(tbb::blocked_range<size_t>(0, m_points.size()),
 
  847                        Refresh_tangent_triangulation(*
this, updated_pts_ds));
 
  851      for (std::size_t i = 0; i < m_points.size(); ++i) refresh_tangent_triangulation(i, updated_pts_ds);
 
  856#ifdef GUDHI_TC_PROFILING 
  858    std::cerr << yellow << 
"done in " << t.num_seconds() << 
" seconds.\n" << white;
 
  859#elif defined(DEBUG_TRACES) 
  860    std::cerr << yellow << 
"done.\n" << white;
 
  864  void export_inconsistent_stars_to_OFF_files(std::string 
const &filename_base)
 const {
 
  866    for (std::size_t idx = 0; idx < m_points.size(); ++idx) {
 
  868      Simplicial_complex sc;
 
  870      bool is_inconsistent = 
false;
 
  871      Star::const_iterator it_inc_simplex = m_stars[idx].begin();
 
  872      Star::const_iterator it_inc_simplex_end = m_stars[idx].end();
 
  873      for (; it_inc_simplex != it_inc_simplex_end; ++it_inc_simplex) {
 
  875        if (is_infinite(*it_inc_simplex)) 
continue;
 
  877        Simplex c = *it_inc_simplex;
 
  883        if (!is_inconsistent && !is_simplex_consistent(c)) is_inconsistent = 
true;
 
  886      if (is_inconsistent) {
 
  888        std::stringstream output_filename;
 
  889        output_filename << filename_base << 
"_" << idx << 
".off";
 
  890        std::ofstream off_stream(output_filename.str().c_str());
 
  891        export_to_off(sc, off_stream);
 
  896  class Compare_distance_to_ref_point {
 
  898    Compare_distance_to_ref_point(Point 
const &ref, K 
const &k) : m_ref(ref), m_k(k) {}
 
  900    bool operator()(Point 
const &p1, Point 
const &p2) {
 
  901      typename K::Squared_distance_d sqdist = m_k.squared_distance_d_object();
 
  902      return sqdist(p1, m_ref) < sqdist(p2, m_ref);
 
  912  class Compute_tangent_triangulation {
 
  913    Tangential_complex &m_tc;
 
  917    Compute_tangent_triangulation(Tangential_complex &tc) : m_tc(tc) {}
 
  920    Compute_tangent_triangulation(
const Compute_tangent_triangulation &ctt) : m_tc(ctt.m_tc) {}
 
  923    void operator()(
const tbb::blocked_range<size_t> &r)
 const {
 
  924      for (
size_t i = r.begin(); i != r.end(); ++i) m_tc.compute_tangent_triangulation(i);
 
  929  class Refresh_tangent_triangulation {
 
  930    Tangential_complex &m_tc;
 
  931    Points_ds 
const &m_updated_pts_ds;
 
  935    Refresh_tangent_triangulation(Tangential_complex &tc, Points_ds 
const &updated_pts_ds)
 
  936        : m_tc(tc), m_updated_pts_ds(updated_pts_ds) {}
 
  939    Refresh_tangent_triangulation(
const Refresh_tangent_triangulation &ctt)
 
  940        : m_tc(ctt.m_tc), m_updated_pts_ds(ctt.m_updated_pts_ds) {}
 
  943    void operator()(
const tbb::blocked_range<size_t> &r)
 const {
 
  944      for (
size_t i = r.begin(); i != r.end(); ++i) m_tc.refresh_tangent_triangulation(i, m_updated_pts_ds);
 
  949  bool is_infinite(Simplex 
const &s)
 const { 
return *s.rbegin() == (std::numeric_limits<std::size_t>::max)(); }
 
  954  Tr_vertex_handle compute_star(std::size_t i, 
const Point ¢er_pt, 
const Tangent_space_basis &tsb,
 
  955                                Triangulation &triangulation, 
bool verbose = 
false) {
 
  956    int tangent_space_dim = tsb.dimension();
 
  957    const Tr_traits &local_tr_traits = triangulation.geom_traits();
 
  960    typename K::Squared_distance_d k_sqdist = m_k.squared_distance_d_object();
 
  963    typename Tr_traits::Compute_weight_d point_weight = local_tr_traits.compute_weight_d_object();
 
  964#if CGAL_VERSION_NR < 1050200000 
  965    typename Tr_traits::Power_center_d power_center = local_tr_traits.power_center_d_object();
 
  967    typename Tr_traits::Construct_power_sphere_d power_center = local_tr_traits.construct_power_sphere_d_object();
 
  977    if (i == tsb.origin()) {
 
  979      proj_wp = local_tr_traits.construct_weighted_point_d_object()(
 
  980          local_tr_traits.construct_point_d_object()(tangent_space_dim, CGAL::ORIGIN), m_weights[i]);
 
  982      const Weighted_point &wp = compute_perturbed_weighted_point(i);
 
  983      proj_wp = project_point_and_compute_weight(wp, tsb, local_tr_traits);
 
  986    Tr_vertex_handle center_vertex = triangulation.insert(proj_wp);
 
  987    center_vertex->data() = i;
 
  988    if (verbose) std::cerr << 
"* Inserted point #" << i << 
"\n";
 
  990#ifdef GUDHI_TC_VERY_VERBOSE 
  991    std::size_t num_attempts_to_insert_points = 1;
 
  992    std::size_t num_inserted_points = 1;
 
 1004    std::optional<FT> squared_star_sphere_radius_plus_margin = m_max_squared_edge_length;
 
 1007    for (
auto nn_it = ins_range.begin(); nn_it != ins_range.end(); ++nn_it) {
 
 1008      std::size_t neighbor_point_idx = nn_it->first;
 
 1011      if (neighbor_point_idx != i) {
 
 1016        compute_perturbed_weighted_point(neighbor_point_idx, neighbor_pt, neighbor_weight);
 
 1017        GUDHI_CHECK(!m_max_squared_edge_length ||
 
 1018                    squared_star_sphere_radius_plus_margin.value() <= m_max_squared_edge_length.value(),
 
 1019                    std::invalid_argument(
"Tangential_complex::compute_star - set a bigger value with set_max_squared_edge_length."));
 
 1020        if (squared_star_sphere_radius_plus_margin &&
 
 1021            k_sqdist(center_pt, neighbor_pt) > squared_star_sphere_radius_plus_margin.value()) {
 
 1022          GUDHI_CHECK(triangulation.current_dimension() >= tangent_space_dim,
 
 1023                      std::invalid_argument(
"Tangential_complex::compute_star - Dimension of the star is only " + \
 
 1024                                            std::to_string(triangulation.current_dimension())));
 
 1028        Tr_point proj_pt = project_point_and_compute_weight(neighbor_pt, neighbor_weight, tsb, local_tr_traits);
 
 1030#ifdef GUDHI_TC_VERY_VERBOSE 
 1031        ++num_attempts_to_insert_points;
 
 1034        Tr_vertex_handle vh = triangulation.insert_if_in_star(proj_pt, center_vertex);
 
 1036        if (vh != Tr_vertex_handle() && vh->data() == (std::numeric_limits<std::size_t>::max)()) {
 
 1037#ifdef GUDHI_TC_VERY_VERBOSE 
 1038          ++num_inserted_points;
 
 1040          if (verbose) std::cerr << 
"* Inserted point #" << neighbor_point_idx << 
"\n";
 
 1042          vh->data() = neighbor_point_idx;
 
 1045          if (triangulation.current_dimension() >= tangent_space_dim) {
 
 1046            squared_star_sphere_radius_plus_margin = std::nullopt;
 
 1048            std::vector<Tr_full_cell_handle> incident_cells;
 
 1049            triangulation.incident_full_cells(center_vertex, std::back_inserter(incident_cells));
 
 1050            for (
typename std::vector<Tr_full_cell_handle>::iterator cit = incident_cells.begin();
 
 1051                 cit != incident_cells.end(); ++cit) {
 
 1052              Tr_full_cell_handle cell = *cit;
 
 1053              if (triangulation.is_infinite(cell)) {
 
 1054                squared_star_sphere_radius_plus_margin = std::nullopt;
 
 1060                    power_center(boost::make_transform_iterator(cell->vertices_begin(),
 
 1061                                                                vertex_handle_to_point<Tr_point, Tr_vertex_handle>),
 
 1062                                 boost::make_transform_iterator(cell->vertices_end(),
 
 1063                                                                vertex_handle_to_point<Tr_point, Tr_vertex_handle>));
 
 1065                FT sq_power_sphere_diam = 4 * point_weight(c);
 
 1067                if (!squared_star_sphere_radius_plus_margin ||
 
 1068                    sq_power_sphere_diam > squared_star_sphere_radius_plus_margin.value()) {
 
 1069                  squared_star_sphere_radius_plus_margin = sq_power_sphere_diam;
 
 1076            if (squared_star_sphere_radius_plus_margin) {
 
 1078              squared_star_sphere_radius_plus_margin =
 
 1079                  CGAL::square(std::sqrt(squared_star_sphere_radius_plus_margin.value()) + 2 * m_last_max_perturb);
 
 1082              if (m_max_squared_edge_length && squared_star_sphere_radius_plus_margin.value() > m_max_squared_edge_length.value()) {
 
 1083                squared_star_sphere_radius_plus_margin = m_max_squared_edge_length.value();
 
 1087              m_squared_star_spheres_radii_incl_margin[i] = squared_star_sphere_radius_plus_margin.value();
 
 1089              if (m_max_squared_edge_length) {
 
 1090                squared_star_sphere_radius_plus_margin = m_max_squared_edge_length.value();
 
 1091                m_squared_star_spheres_radii_incl_margin[i] = m_max_squared_edge_length.value();
 
 1093                m_squared_star_spheres_radii_incl_margin[i] = FT(-1);
 
 1101    return center_vertex;
 
 1104  void refresh_tangent_triangulation(std::size_t i, Points_ds 
const &updated_pts_ds, 
bool verbose = 
false) {
 
 1105    if (verbose) std::cerr << 
"** Refreshing tangent tri #" << i << 
" **\n";
 
 1107    if (m_squared_star_spheres_radii_incl_margin[i] == FT(-1)) 
return compute_tangent_triangulation(i, verbose);
 
 1109    Point center_point = compute_perturbed_point(i);
 
 1111    std::size_t closest_pt_index = updated_pts_ds.k_nearest_neighbors(center_point, 1, 
false).begin()->first;
 
 1113    typename K::Construct_weighted_point_d k_constr_wp = m_k.construct_weighted_point_d_object();
 
 1114#if CGAL_VERSION_NR < 1050200000 
 1115    typename K::Power_distance_d k_power_dist = m_k.power_distance_d_object();
 
 1117    typename K::Compute_power_product_d k_power_dist = m_k.compute_power_product_d_object();
 
 1121    Weighted_point star_sphere = k_constr_wp(compute_perturbed_point(i), m_squared_star_spheres_radii_incl_margin[i]);
 
 1122    Weighted_point closest_updated_point = compute_perturbed_weighted_point(closest_pt_index);
 
 1125    if (k_power_dist(star_sphere, closest_updated_point) <= FT(0)) compute_tangent_triangulation(i, verbose);
 
 1128  void compute_tangent_triangulation(std::size_t i, 
bool verbose = 
false) {
 
 1129    if (verbose) std::cerr << 
"** Computing tangent tri #" << i << 
" **\n";
 
 1134    const Point center_pt = compute_perturbed_point(i);
 
 1135    Tangent_space_basis &tsb = m_tangent_spaces[i];
 
 1138    if (!m_are_tangent_spaces_computed[i]) {
 
 1139#ifdef GUDHI_TC_EXPORT_NORMALS 
 1140      tsb = compute_tangent_space(center_pt, i, 
true , &m_orth_spaces[i]);
 
 1142      tsb = compute_tangent_space(center_pt, i);
 
 1146#if defined(GUDHI_TC_PROFILING) && defined(GUDHI_TC_VERY_VERBOSE) 
 1149    int tangent_space_dim = tangent_basis_dim(i);
 
 1150    Triangulation &local_tr = m_triangulations[i].construct_triangulation(tangent_space_dim);
 
 1152    m_triangulations[i].center_vertex() = compute_star(i, center_pt, tsb, local_tr, verbose);
 
 1154#if defined(GUDHI_TC_PROFILING) && defined(GUDHI_TC_VERY_VERBOSE) 
 1156    std::cerr << 
"  - triangulation construction: " << t.num_seconds() << 
" s.\n";
 
 1160#ifdef GUDHI_TC_VERY_VERBOSE 
 1161    std::cerr << 
"Inserted " << num_inserted_points << 
" points / " << num_attempts_to_insert_points
 
 1162              << 
" attempts to compute the star\n";
 
 1167#if defined(GUDHI_TC_PROFILING) && defined(GUDHI_TC_VERY_VERBOSE) 
 1169    std::cerr << 
"  - update_star: " << t.num_seconds() << 
" s.\n";
 
 1175  void update_star(std::size_t i) {
 
 1176    Star &star = m_stars[i];
 
 1178    Triangulation &local_tr = m_triangulations[i].tr();
 
 1179    Tr_vertex_handle center_vertex = m_triangulations[i].center_vertex();
 
 1180    int cur_dim_plus_1 = local_tr.current_dimension() + 1;
 
 1182    std::vector<Tr_full_cell_handle> incident_cells;
 
 1183    local_tr.incident_full_cells(center_vertex, std::back_inserter(incident_cells));
 
 1185    typename std::vector<Tr_full_cell_handle>::const_iterator it_c = incident_cells.begin();
 
 1186    typename std::vector<Tr_full_cell_handle>::const_iterator it_c_end = incident_cells.end();
 
 1188    for (; it_c != it_c_end; ++it_c) {
 
 1190      Incident_simplex incident_simplex;
 
 1191      for (
int j = 0; j < cur_dim_plus_1; ++j) {
 
 1192        std::size_t index = (*it_c)->vertex(j)->data();
 
 1193        if (index != i) incident_simplex.insert(index);
 
 1195      GUDHI_CHECK(incident_simplex.size() == cur_dim_plus_1 - 1,
 
 1196                  std::logic_error(
"update_star: wrong size of incident simplex"));
 
 1197      star.push_back(incident_simplex);
 
 1203  Tangent_space_basis compute_tangent_space(
const Point &p, 
const std::size_t i, 
bool normalize_basis = 
true,
 
 1204                                            Orthogonal_space_basis *p_orth_space_basis = NULL) {
 
 1205    unsigned int num_pts_for_pca =
 
 1206        (std::min)(
static_cast<unsigned int>(std::pow(GUDHI_TC_BASE_VALUE_FOR_PCA, m_intrinsic_dim)),
 
 1207                   static_cast<unsigned int>(m_points.size()));
 
 1210    typename K::Construct_vector_d constr_vec = m_k.construct_vector_d_object();
 
 1211    typename K::Compute_coordinate_d coord = m_k.compute_coordinate_d_object();
 
 1213#ifdef GUDHI_TC_USE_ANOTHER_POINT_SET_FOR_TANGENT_SPACE_ESTIM 
 1214    KNS_range kns_range = m_points_ds_for_tse.k_nearest_neighbors(p, num_pts_for_pca, 
false);
 
 1215    const Points &points_for_pca = m_points_for_tse;
 
 1218    const Points &points_for_pca = m_points;
 
 1222    Eigen::MatrixXd mat_points(num_pts_for_pca, m_ambient_dim);
 
 1223    auto nn_it = kns_range.begin();
 
 1224    for (
unsigned int j = 0; j < num_pts_for_pca && nn_it != kns_range.end(); ++j, ++nn_it) {
 
 1225      for (
int i = 0; i < m_ambient_dim; ++i) {
 
 1226        mat_points(j, i) = CGAL::to_double(coord(points_for_pca[nn_it->first], i));
 
 1229    Eigen::MatrixXd centered = mat_points.rowwise() - mat_points.colwise().mean();
 
 1230    Eigen::MatrixXd cov = centered.adjoint() * centered;
 
 1231    Eigen::SelfAdjointEigenSolver<Eigen::MatrixXd> eig(cov);
 
 1233    Tangent_space_basis tsb(i);  
 
 1237    for (
int j = m_ambient_dim - 1; j >= m_ambient_dim - m_intrinsic_dim; --j) {
 
 1238      if (normalize_basis) {
 
 1239        Vector v = constr_vec(m_ambient_dim, eig.eigenvectors().col(j).data(),
 
 1240                              eig.eigenvectors().col(j).data() + m_ambient_dim);
 
 1241        tsb.push_back(normalize_vector(v, m_k));
 
 1243        tsb.push_back(constr_vec(m_ambient_dim, eig.eigenvectors().col(j).data(),
 
 1244                                 eig.eigenvectors().col(j).data() + m_ambient_dim));
 
 1248    if (p_orth_space_basis) {
 
 1249      p_orth_space_basis->set_origin(i);
 
 1250      for (
int j = m_ambient_dim - m_intrinsic_dim - 1; j >= 0; --j) {
 
 1251        if (normalize_basis) {
 
 1252          Vector v = constr_vec(m_ambient_dim, eig.eigenvectors().col(j).data(),
 
 1253                                eig.eigenvectors().col(j).data() + m_ambient_dim);
 
 1254          p_orth_space_basis->push_back(normalize_vector(v, m_k));
 
 1256          p_orth_space_basis->push_back(constr_vec(m_ambient_dim, eig.eigenvectors().col(j).data(),
 
 1257                                                   eig.eigenvectors().col(j).data() + m_ambient_dim));
 
 1262    m_are_tangent_spaces_computed[i] = 
true;
 
 1272  Tangent_space_basis compute_tangent_space(
const Simplex &s, 
bool normalize_basis = 
true) {
 
 1273    unsigned int num_pts_for_pca =
 
 1274        (std::min)(
static_cast<unsigned int>(std::pow(GUDHI_TC_BASE_VALUE_FOR_PCA, m_intrinsic_dim)),
 
 1275                   static_cast<unsigned int>(m_points.size()));
 
 1278    typename K::Construct_vector_d constr_vec = m_k.construct_vector_d_object();
 
 1279    typename K::Compute_coordinate_d coord = m_k.compute_coordinate_d_object();
 
 1280    typename K::Squared_length_d sqlen = m_k.squared_length_d_object();
 
 1281    typename K::Scaled_vector_d scaled_vec = m_k.scaled_vector_d_object();
 
 1282    typename K::Scalar_product_d scalar_pdct = m_k.scalar_product_d_object();
 
 1283    typename K::Difference_of_vectors_d diff_vec = m_k.difference_of_vectors_d_object();
 
 1286    Eigen::MatrixXd mat_points(s.size() * num_pts_for_pca, m_ambient_dim);
 
 1287    unsigned int current_row = 0;
 
 1289    for (Simplex::const_iterator it_index = s.begin(); it_index != s.end(); ++it_index) {
 
 1290      const Point &p = m_points[*it_index];
 
 1292#ifdef GUDHI_TC_USE_ANOTHER_POINT_SET_FOR_TANGENT_SPACE_ESTIM 
 1293      KNS_range kns_range = m_points_ds_for_tse.k_nearest_neighbors(p, num_pts_for_pca, 
false);
 
 1294      const Points &points_for_pca = m_points_for_tse;
 
 1297      const Points &points_for_pca = m_points;
 
 1300      auto nn_it = kns_range.begin();
 
 1301      for (; current_row < num_pts_for_pca && nn_it != kns_range.end(); ++current_row, ++nn_it) {
 
 1302        for (
int i = 0; i < m_ambient_dim; ++i) {
 
 1303          mat_points(current_row, i) = CGAL::to_double(coord(points_for_pca[nn_it->first], i));
 
 1307    Eigen::MatrixXd centered = mat_points.rowwise() - mat_points.colwise().mean();
 
 1308    Eigen::MatrixXd cov = centered.adjoint() * centered;
 
 1309    Eigen::SelfAdjointEigenSolver<Eigen::MatrixXd> eig(cov);
 
 1311    Tangent_space_basis tsb;
 
 1315    for (
int j = m_ambient_dim - 1; j >= m_ambient_dim - m_intrinsic_dim; --j) {
 
 1316      if (normalize_basis) {
 
 1317        Vector v = constr_vec(m_ambient_dim, eig.eigenvectors().col(j).data(),
 
 1318                              eig.eigenvectors().col(j).data() + m_ambient_dim);
 
 1319        tsb.push_back(normalize_vector(v, m_k));
 
 1321        tsb.push_back(constr_vec(m_ambient_dim, eig.eigenvectors().col(j).data(),
 
 1322                                 eig.eigenvectors().col(j).data() + m_ambient_dim));
 
 1331  int tangent_basis_dim(std::size_t i)
 const { 
return m_tangent_spaces[i].dimension(); }
 
 1333  Point compute_perturbed_point(std::size_t pt_idx)
 const {
 
 1334#ifdef GUDHI_TC_PERTURB_POSITION 
 1335    return m_k.translated_point_d_object()(m_points[pt_idx], m_translations[pt_idx]);
 
 1337    return m_points[pt_idx];
 
 1341  void compute_perturbed_weighted_point(std::size_t pt_idx, Point &p, FT &w)
 const {
 
 1342#ifdef GUDHI_TC_PERTURB_POSITION 
 1343    p = m_k.translated_point_d_object()(m_points[pt_idx], m_translations[pt_idx]);
 
 1345    p = m_points[pt_idx];
 
 1347    w = m_weights[pt_idx];
 
 1350  Weighted_point compute_perturbed_weighted_point(std::size_t pt_idx)
 const {
 
 1351    typename K::Construct_weighted_point_d k_constr_wp = m_k.construct_weighted_point_d_object();
 
 1353    Weighted_point wp = k_constr_wp(
 
 1354#ifdef GUDHI_TC_PERTURB_POSITION
 
 1355        m_k.translated_point_d_object()(m_points[pt_idx], m_translations[pt_idx]),
 
 1364  Point unproject_point(
const Tr_point &p, 
const Tangent_space_basis &tsb, 
const Tr_traits &tr_traits)
 const {
 
 1365    typename K::Translated_point_d k_transl = m_k.translated_point_d_object();
 
 1366    typename K::Scaled_vector_d k_scaled_vec = m_k.scaled_vector_d_object();
 
 1367    typename Tr_traits::Compute_coordinate_d coord = tr_traits.compute_coordinate_d_object();
 
 1369    Point global_point = compute_perturbed_point(tsb.origin());
 
 1370    for (
int i = 0; i < m_intrinsic_dim; ++i) global_point = k_transl(global_point, k_scaled_vec(tsb[i], coord(p, i)));
 
 1372    return global_point;
 
 1377  Tr_bare_point project_point(
const Point &p, 
const Tangent_space_basis &tsb, 
const Tr_traits &tr_traits)
 const {
 
 1378    typename K::Scalar_product_d scalar_pdct = m_k.scalar_product_d_object();
 
 1379    typename K::Difference_of_points_d diff_points = m_k.difference_of_points_d_object();
 
 1381    Vector v = diff_points(p, compute_perturbed_point(tsb.origin()));
 
 1383    std::vector<FT> coords;
 
 1385    coords.reserve(tsb.dimension());
 
 1386    for (std::size_t i = 0; i < m_intrinsic_dim; ++i) {
 
 1388      FT coord = scalar_pdct(v, tsb[i]);
 
 1389      coords.push_back(coord);
 
 1392    return tr_traits.construct_point_d_object()(
static_cast<int>(coords.size()), coords.begin(), coords.end());
 
 1399  Tr_point project_point_and_compute_weight(
const Weighted_point &wp, 
const Tangent_space_basis &tsb,
 
 1400                                            const Tr_traits &tr_traits)
 const {
 
 1401    typename K::Point_drop_weight_d k_drop_w = m_k.point_drop_weight_d_object();
 
 1402    typename K::Compute_weight_d k_point_weight = m_k.compute_weight_d_object();
 
 1403    return project_point_and_compute_weight(k_drop_w(wp), k_point_weight(wp), tsb, tr_traits);
 
 1407  Tr_point project_point_and_compute_weight(
const Point &p, 
const FT w, 
const Tangent_space_basis &tsb,
 
 1408                                            const Tr_traits &tr_traits)
 const {
 
 1409    const int point_dim = m_k.point_dimension_d_object()(p);
 
 1411    typename K::Construct_point_d constr_pt = m_k.construct_point_d_object();
 
 1412    typename K::Scalar_product_d scalar_pdct = m_k.scalar_product_d_object();
 
 1413    typename K::Difference_of_points_d diff_points = m_k.difference_of_points_d_object();
 
 1414    typename K::Compute_coordinate_d coord = m_k.compute_coordinate_d_object();
 
 1415    typename K::Construct_cartesian_const_iterator_d ccci = m_k.construct_cartesian_const_iterator_d_object();
 
 1417    Point origin = compute_perturbed_point(tsb.origin());
 
 1418    Vector v = diff_points(p, origin);
 
 1421    bool same_dim = (point_dim == tsb.dimension());
 
 1423    std::vector<FT> coords;
 
 1425    std::vector<FT> p_proj(ccci(origin), ccci(origin, 0));
 
 1426    coords.reserve(tsb.dimension());
 
 1427    for (
int i = 0; i < tsb.dimension(); ++i) {
 
 1429      FT c = scalar_pdct(v, tsb[i]);
 
 1430      coords.push_back(c);
 
 1434        for (
int j = 0; j < point_dim; ++j) p_proj[j] += c * coord(tsb[i], j);
 
 1439    FT sq_dist_to_proj_pt = 0;
 
 1441      Point projected_pt = constr_pt(point_dim, p_proj.begin(), p_proj.end());
 
 1442      sq_dist_to_proj_pt = m_k.squared_distance_d_object()(p, projected_pt);
 
 1445    return tr_traits.construct_weighted_point_d_object()(
 
 1446        tr_traits.construct_point_d_object()(
static_cast<int>(coords.size()), coords.begin(), coords.end()),
 
 1447        w - sq_dist_to_proj_pt);
 
 1452  template <
typename Indexed_po
int_range>
 
 1453  std::vector<Tr_point> project_points_and_compute_weights(
const Indexed_point_range &point_indices,
 
 1454                                                           const Tangent_space_basis &tsb,
 
 1455                                                           const Tr_traits &tr_traits)
 const {
 
 1456    std::vector<Tr_point> ret;
 
 1457    for (
typename Indexed_point_range::const_iterator it = point_indices.begin(), it_end = point_indices.end();
 
 1458         it != it_end; ++it) {
 
 1459      ret.push_back(project_point_and_compute_weight(compute_perturbed_weighted_point(*it), tsb, tr_traits));
 
 1466  bool is_simplex_consistent(Tr_full_cell_handle fch, 
int cur_dim)
 const {
 
 1468    for (
int i = 0; i < cur_dim + 1; ++i) {
 
 1469      std::size_t data = fch->vertex(i)->data();
 
 1472    return is_simplex_consistent(c);
 
 1478  bool is_simplex_consistent(Simplex 
const &simplex)
 const {
 
 1480    Simplex::const_iterator it_point_idx = simplex.begin();
 
 1483    for (; it_point_idx != simplex.end(); ++it_point_idx) {
 
 1484      std::size_t point_idx = *it_point_idx;
 
 1486      if (point_idx == (std::numeric_limits<std::size_t>::max)()) 
continue;
 
 1488      Star 
const &star = m_stars[point_idx];
 
 1491      Incident_simplex is_to_find = simplex;
 
 1492      is_to_find.erase(point_idx);
 
 1495      if (std::find(star.begin(), star.end(), is_to_find) == star.end()) 
return false;
 
 1507  template <
typename OutputIterator>  
 
 1508  bool is_simplex_consistent(std::size_t center_point,
 
 1509                             Incident_simplex 
const &s,  
 
 1510                             OutputIterator points_whose_star_does_not_contain_s,
 
 1511                             bool check_also_in_non_maximal_faces = 
false)
 const {
 
 1512    Simplex full_simplex = s;
 
 1513    full_simplex.insert(center_point);
 
 1516    Incident_simplex::const_iterator it_point_idx = s.begin();
 
 1519    for (; it_point_idx != s.end(); ++it_point_idx) {
 
 1520      std::size_t point_idx = *it_point_idx;
 
 1522      if (point_idx == (std::numeric_limits<std::size_t>::max)()) 
continue;
 
 1524      Star 
const &star = m_stars[point_idx];
 
 1527      Incident_simplex is_to_find = full_simplex;
 
 1528      is_to_find.erase(point_idx);
 
 1530      if (check_also_in_non_maximal_faces) {
 
 1534        for (Star::const_iterator is = star.begin(), is_end = star.end(); !found && is != is_end; ++is) {
 
 1535          if (std::includes(is->begin(), is->end(), is_to_find.begin(), is_to_find.end())) found = 
true;
 
 1538        if (!found) *points_whose_star_does_not_contain_s++ = point_idx;
 
 1541        if (std::find(star.begin(), star.end(), is_to_find) == star.end())
 
 1542          *points_whose_star_does_not_contain_s++ = point_idx;
 
 1552  bool is_simplex_in_star(std::size_t p, Incident_simplex 
const &s, 
bool check_also_in_non_maximal_faces = 
true)
 const {
 
 1553    Star 
const &star = m_stars[p];
 
 1555    if (check_also_in_non_maximal_faces) {
 
 1559      for (Star::const_iterator is = star.begin(), is_end = star.end(); !found && is != is_end; ++is) {
 
 1560        if (std::includes(is->begin(), is->end(), s.begin(), s.end())) found = 
true;
 
 1565      return !(std::find(star.begin(), star.end(), s) == star.end());
 
 1571  class Try_to_solve_inconsistencies_in_a_local_triangulation {
 
 1572    Tangential_complex &m_tc;
 
 1573    double m_max_perturb;
 
 1574    tbb::combinable<std::size_t> &m_num_inconsistencies;
 
 1575    tbb::combinable<std::vector<std::size_t> > &m_updated_points;
 
 1579    Try_to_solve_inconsistencies_in_a_local_triangulation(Tangential_complex &tc, 
double max_perturb,
 
 1580                                                          tbb::combinable<std::size_t> &num_inconsistencies,
 
 1581                                                          tbb::combinable<std::vector<std::size_t> > &updated_points)
 
 1583          m_max_perturb(max_perturb),
 
 1584          m_num_inconsistencies(num_inconsistencies),
 
 1585          m_updated_points(updated_points) {}
 
 1588    Try_to_solve_inconsistencies_in_a_local_triangulation(
 
 1589        const Try_to_solve_inconsistencies_in_a_local_triangulation &tsilt)
 
 1591          m_max_perturb(tsilt.m_max_perturb),
 
 1592          m_num_inconsistencies(tsilt.m_num_inconsistencies),
 
 1593          m_updated_points(tsilt.m_updated_points) {}
 
 1596    void operator()(
const tbb::blocked_range<size_t> &r)
 const {
 
 1597      for (
size_t i = r.begin(); i != r.end(); ++i) {
 
 1598        m_num_inconsistencies.local() += m_tc.try_to_solve_inconsistencies_in_a_local_triangulation(
 
 1599            i, m_max_perturb, std::back_inserter(m_updated_points.local()));
 
 1605  void perturb(std::size_t point_idx, 
double max_perturb) {
 
 1606    const Tr_traits &local_tr_traits = m_triangulations[point_idx].tr().geom_traits();
 
 1607    typename Tr_traits::Compute_coordinate_d coord = local_tr_traits.compute_coordinate_d_object();
 
 1608    typename K::Translated_point_d k_transl = m_k.translated_point_d_object();
 
 1609    typename K::Construct_vector_d k_constr_vec = m_k.construct_vector_d_object();
 
 1610    typename K::Scaled_vector_d k_scaled_vec = m_k.scaled_vector_d_object();
 
 1612    CGAL::Random_points_in_ball_d<Tr_bare_point> tr_point_in_ball_generator(
 
 1613        m_intrinsic_dim, m_random_generator.get_double(0., max_perturb));
 
 1615    Tr_point local_random_transl =
 
 1616        local_tr_traits.construct_weighted_point_d_object()(*tr_point_in_ball_generator++, 0);
 
 1617    Translation_for_perturb global_transl = k_constr_vec(m_ambient_dim);
 
 1618    const Tangent_space_basis &tsb = m_tangent_spaces[point_idx];
 
 1619    for (
int i = 0; i < m_intrinsic_dim; ++i) {
 
 1620      global_transl = k_transl(global_transl, k_scaled_vec(tsb[i], coord(local_random_transl, i)));
 
 1623#if defined(GUDHI_USE_TBB) 
 1624    m_p_perturb_mutexes[point_idx].lock();
 
 1625    m_translations[point_idx] = global_transl;
 
 1626    m_p_perturb_mutexes[point_idx].unlock();
 
 1629    m_translations[point_idx] = global_transl;
 
 1634  template <
typename OutputIt>
 
 1635  bool try_to_solve_inconsistencies_in_a_local_triangulation(
 
 1636      std::size_t tr_index, 
double max_perturb, OutputIt perturbed_pts_indices = CGAL::Emptyset_iterator()) {
 
 1637    bool is_inconsistent = 
false;
 
 1639    Star 
const &star = m_stars[tr_index];
 
 1642    Star::const_iterator it_inc_simplex = star.begin();
 
 1643    Star::const_iterator it_inc_simplex_end = star.end();
 
 1644    for (; it_inc_simplex != it_inc_simplex_end; ++it_inc_simplex) {
 
 1645      const Incident_simplex &incident_simplex = *it_inc_simplex;
 
 1648      if (is_infinite(incident_simplex)) 
continue;
 
 1650      Simplex c = incident_simplex;
 
 1654      if (!is_simplex_consistent(c)) {
 
 1655        is_inconsistent = 
true;
 
 1657        std::size_t idx = tr_index;
 
 1659        perturb(tr_index, max_perturb);
 
 1660        *perturbed_pts_indices++ = idx;
 
 1667    return is_inconsistent;
 
 1672  std::ostream &export_point_set(std::ostream &os, 
bool use_perturbed_points = 
false,
 
 1673                                 const char *coord_separator = 
" ")
 const {
 
 1674    if (use_perturbed_points) {
 
 1675      std::vector<Point> perturbed_points;
 
 1676      perturbed_points.reserve(m_points.size());
 
 1677      for (std::size_t i = 0; i < m_points.size(); ++i) perturbed_points.push_back(compute_perturbed_point(i));
 
 1679      return export_point_set(m_k, perturbed_points, os, coord_separator);
 
 1681      return export_point_set(m_k, m_points, os, coord_separator);
 
 1685  template <
typename ProjectionFunctor = CGAL::Identity<Po
int> >
 
 1686  std::ostream &export_vertices_to_off(std::ostream &os, std::size_t &num_vertices, 
bool use_perturbed_points = 
false,
 
 1687                                       ProjectionFunctor 
const &point_projection = ProjectionFunctor())
 const {
 
 1688    if (m_points.empty()) {
 
 1696    const int N = (m_intrinsic_dim == 1 ? 2 : 1);
 
 1699    typename K::Compute_coordinate_d coord = m_k.compute_coordinate_d_object();
 
 1701#ifdef GUDHI_TC_EXPORT_ALL_COORDS_IN_OFF 
 1702    int num_coords = m_ambient_dim;
 
 1704    int num_coords = (std::min)(m_ambient_dim, 3);
 
 1707#ifdef GUDHI_TC_EXPORT_NORMALS 
 1708    OS_container::const_iterator it_os = m_orth_spaces.begin();
 
 1710    typename Points::const_iterator it_p = m_points.begin();
 
 1711    typename Points::const_iterator it_p_end = m_points.end();
 
 1713    for (std::size_t i = 0; it_p != it_p_end; ++it_p, ++i) {
 
 1714      Point p = point_projection(use_perturbed_points ? compute_perturbed_point(i) : *it_p);
 
 1715      for (
int ii = 0; ii < N; ++ii) {
 
 1717        for (; j < num_coords; ++j) os << CGAL::to_double(coord(p, j)) << 
" ";
 
 1718        if (j == 2) os << 
"0";
 
 1720#ifdef GUDHI_TC_EXPORT_NORMALS 
 1721        for (j = 0; j < num_coords; ++j) os << 
" " << CGAL::to_double(coord(*it_os->begin(), j));
 
 1725#ifdef GUDHI_TC_EXPORT_NORMALS 
 1730    num_vertices = N * m_points.size();
 
 1734  std::ostream &export_simplices_to_off(std::ostream &os, std::size_t &num_OFF_simplices,
 
 1735                                        bool color_inconsistencies = 
false,
 
 1736                                        Simplex_set 
const *p_simpl_to_color_in_red = NULL,
 
 1737                                        Simplex_set 
const *p_simpl_to_color_in_green = NULL,
 
 1738                                        Simplex_set 
const *p_simpl_to_color_in_blue = NULL)
 const {
 
 1741    num_OFF_simplices = 0;
 
 1742    std::size_t num_maximal_simplices = 0;
 
 1743    std::size_t num_inconsistent_maximal_simplices = 0;
 
 1744    std::size_t num_inconsistent_stars = 0;
 
 1745    typename Tr_container::const_iterator it_tr = m_triangulations.begin();
 
 1746    typename Tr_container::const_iterator it_tr_end = m_triangulations.end();
 
 1748    for (std::size_t idx = 0; it_tr != it_tr_end; ++it_tr, ++idx) {
 
 1749      bool is_star_inconsistent = 
false;
 
 1751      Triangulation 
const &tr = it_tr->tr();
 
 1753      if (tr.current_dimension() < m_intrinsic_dim) 
continue;
 
 1756      std::stringstream color;
 
 1758      color << 128 << 
" " << 128 << 
" " << 128;
 
 1761      typedef std::vector<std::pair<Simplex, int> > Star_using_triangles;
 
 1762      Star_using_triangles star_using_triangles;
 
 1765      Star::const_iterator it_inc_simplex = m_stars[idx].begin();
 
 1766      Star::const_iterator it_inc_simplex_end = m_stars[idx].end();
 
 1767      for (; it_inc_simplex != it_inc_simplex_end; ++it_inc_simplex) {
 
 1768        Simplex c = *it_inc_simplex;
 
 1770        std::size_t num_vertices = c.size();
 
 1771        ++num_maximal_simplices;
 
 1773        int color_simplex = -1;  
 
 1774        if (color_inconsistencies && !is_simplex_consistent(c)) {
 
 1775          ++num_inconsistent_maximal_simplices;
 
 1777          is_star_inconsistent = 
true;
 
 1779          if (p_simpl_to_color_in_red && std::find(p_simpl_to_color_in_red->begin(), p_simpl_to_color_in_red->end(),
 
 1780                                                   c) != p_simpl_to_color_in_red->end()) {
 
 1782          } 
else if (p_simpl_to_color_in_green &&
 
 1783                     std::find(p_simpl_to_color_in_green->begin(), p_simpl_to_color_in_green->end(), c) !=
 
 1784                         p_simpl_to_color_in_green->end()) {
 
 1786          } 
else if (p_simpl_to_color_in_blue &&
 
 1787                     std::find(p_simpl_to_color_in_blue->begin(), p_simpl_to_color_in_blue->end(), c) !=
 
 1788                         p_simpl_to_color_in_blue->end()) {
 
 1797        if (m_intrinsic_dim == 1) {
 
 1799          Simplex::iterator it = c.begin();
 
 1800          for (; it != c.end(); ++it) tmp_c.insert(*it * 2);
 
 1801          if (num_vertices == 2) tmp_c.insert(*tmp_c.rbegin() + 1);
 
 1806        if (num_vertices <= 3) {
 
 1807          star_using_triangles.push_back(std::make_pair(c, color_simplex));
 
 1810          std::vector<bool> booleans(num_vertices, 
false);
 
 1811          std::fill(booleans.begin() + num_vertices - 3, booleans.end(), 
true);
 
 1814            Simplex::iterator it = c.begin();
 
 1815            for (
int i = 0; it != c.end(); ++i, ++it) {
 
 1816              if (booleans[i]) triangle.insert(*it);
 
 1818            star_using_triangles.push_back(std::make_pair(triangle, color_simplex));
 
 1819          } 
while (std::next_permutation(booleans.begin(), booleans.end()));
 
 1824      Star_using_triangles::const_iterator it_simplex = star_using_triangles.begin();
 
 1825      Star_using_triangles::const_iterator it_simplex_end = star_using_triangles.end();
 
 1826      for (; it_simplex != it_simplex_end; ++it_simplex) {
 
 1827        const Simplex &c = it_simplex->first;
 
 1830        if (is_infinite(c)) 
continue;
 
 1832        int color_simplex = it_simplex->second;
 
 1834        std::stringstream sstr_c;
 
 1836        Simplex::const_iterator it_point_idx = c.begin();
 
 1837        for (; it_point_idx != c.end(); ++it_point_idx) {
 
 1838          sstr_c << *it_point_idx << 
" ";
 
 1841        os << 3 << 
" " << sstr_c.str();
 
 1842        if (color_inconsistencies || p_simpl_to_color_in_red || p_simpl_to_color_in_green || p_simpl_to_color_in_blue) {
 
 1843          switch (color_simplex) {
 
 1857              os << 
" " << color.str();
 
 1861        ++num_OFF_simplices;
 
 1864      if (is_star_inconsistent) ++num_inconsistent_stars;
 
 1868    std::cerr << 
"\n==========================================================\n" 
 1869              << 
"Export from list of stars to OFF:\n" 
 1870              << 
"  * Number of vertices: " << m_points.size() << 
"\n" 
 1871              << 
"  * Total number of maximal simplices: " << num_maximal_simplices << 
"\n";
 
 1872    if (color_inconsistencies) {
 
 1873      std::cerr << 
"  * Number of inconsistent stars: " << num_inconsistent_stars << 
" (" 
 1874                << (m_points.size() > 0 ? 100. * num_inconsistent_stars / m_points.size() : 0.) << 
"%)\n" 
 1875                << 
"  * Number of inconsistent maximal simplices: " << num_inconsistent_maximal_simplices << 
" (" 
 1876                << (num_maximal_simplices > 0 ? 100. * num_inconsistent_maximal_simplices / num_maximal_simplices : 0.)
 
 1879    std::cerr << 
"==========================================================\n";
 
 1886  std::ostream &export_simplices_to_off(
const Simplicial_complex &complex, std::ostream &os,
 
 1887                                        std::size_t &num_OFF_simplices,
 
 1888                                        Simplex_set 
const *p_simpl_to_color_in_red = NULL,
 
 1889                                        Simplex_set 
const *p_simpl_to_color_in_green = NULL,
 
 1890                                        Simplex_set 
const *p_simpl_to_color_in_blue = NULL)
 const {
 
 1891    typedef Simplicial_complex::Simplex Simplex;
 
 1892    typedef Simplicial_complex::Simplex_set Simplex_set;
 
 1896    num_OFF_simplices = 0;
 
 1897    std::size_t num_maximal_simplices = 0;
 
 1899    typename Simplex_set::const_iterator it_s = complex.simplex_range().begin();
 
 1900    typename Simplex_set::const_iterator it_s_end = complex.simplex_range().end();
 
 1902    for (; it_s != it_s_end; ++it_s) {
 
 1904      ++num_maximal_simplices;
 
 1906      int color_simplex = -1;  
 
 1907      if (p_simpl_to_color_in_red && std::find(p_simpl_to_color_in_red->begin(), p_simpl_to_color_in_red->end(), c) !=
 
 1908                                         p_simpl_to_color_in_red->end()) {
 
 1910      } 
else if (p_simpl_to_color_in_green &&
 
 1911                 std::find(p_simpl_to_color_in_green->begin(), p_simpl_to_color_in_green->end(), c) !=
 
 1912                     p_simpl_to_color_in_green->end()) {
 
 1914      } 
else if (p_simpl_to_color_in_blue &&
 
 1915                 std::find(p_simpl_to_color_in_blue->begin(), p_simpl_to_color_in_blue->end(), c) !=
 
 1916                     p_simpl_to_color_in_blue->end()) {
 
 1921      typedef std::vector<Simplex> Triangles;
 
 1922      Triangles triangles;
 
 1924      int num_vertices = 
static_cast<int>(c.size());
 
 1926      if (num_vertices < m_intrinsic_dim + 1) 
continue;
 
 1932      if (m_intrinsic_dim == 1) {
 
 1934        Simplex::iterator it = c.begin();
 
 1935        for (; it != c.end(); ++it) tmp_c.insert(*it * 2);
 
 1936        if (num_vertices == 2) tmp_c.insert(*tmp_c.rbegin() + 1);
 
 1941      if (num_vertices <= 3) {
 
 1942        triangles.push_back(c);
 
 1945        std::vector<bool> booleans(num_vertices, 
false);
 
 1946        std::fill(booleans.begin() + num_vertices - 3, booleans.end(), 
true);
 
 1949          Simplex::iterator it = c.begin();
 
 1950          for (
int i = 0; it != c.end(); ++i, ++it) {
 
 1951            if (booleans[i]) triangle.insert(*it);
 
 1953          triangles.push_back(triangle);
 
 1954        } 
while (std::next_permutation(booleans.begin(), booleans.end()));
 
 1958      Triangles::const_iterator it_tri = triangles.begin();
 
 1959      Triangles::const_iterator it_tri_end = triangles.end();
 
 1960      for (; it_tri != it_tri_end; ++it_tri) {
 
 1962        if (is_infinite(*it_tri)) 
continue;
 
 1965        Simplex::const_iterator it_point_idx = it_tri->begin();
 
 1966        for (; it_point_idx != it_tri->end(); ++it_point_idx) {
 
 1967          os << *it_point_idx << 
" ";
 
 1970        if (p_simpl_to_color_in_red || p_simpl_to_color_in_green || p_simpl_to_color_in_blue) {
 
 1971          switch (color_simplex) {
 
 1985              os << 
" 128 128 128";
 
 1990        ++num_OFF_simplices;
 
 1996    std::cerr << 
"\n==========================================================\n" 
 1997              << 
"Export from complex to OFF:\n" 
 1998              << 
"  * Number of vertices: " << m_points.size() << 
"\n" 
 1999              << 
"  * Total number of maximal simplices: " << num_maximal_simplices << 
"\n" 
 2000              << 
"==========================================================\n";
 
 2017  const int m_intrinsic_dim;
 
 2018  const int m_ambient_dim;
 
 2022#ifdef GUDHI_TC_PERTURB_POSITION 
 2023  Translations_for_perturb m_translations;
 
 2024#if defined(GUDHI_USE_TBB) 
 2025  Mutex_for_perturb *m_p_perturb_mutexes;
 
 2029  Points_ds m_points_ds;
 
 2030  double m_last_max_perturb;
 
 2031  std::vector<bool> m_are_tangent_spaces_computed;
 
 2032  TS_container m_tangent_spaces;
 
 2033#ifdef GUDHI_TC_EXPORT_NORMALS 
 2034  OS_container m_orth_spaces;
 
 2036  Tr_container m_triangulations;  
 
 2038  Stars_container m_stars;
 
 2039  std::vector<FT> m_squared_star_spheres_radii_incl_margin;
 
 2040  std::optional<FT> m_max_squared_edge_length;
 
 2042#ifdef GUDHI_TC_USE_ANOTHER_POINT_SET_FOR_TANGENT_SPACE_ESTIM 
 2043  Points m_points_for_tse;
 
 2044  Points_ds m_points_ds_for_tse;
 
 2047  mutable CGAL::Random m_random_generator;
 
Spatial tree data structure to perform (approximate) nearest and furthest neighbor search.
Definition: Kd_tree_search.h:70
K_neighbor_search KNS_range
The range returned by a k-nearest or k-furthest neighbor search. Its value type is std::pair<std::siz...
Definition: Kd_tree_search.h:104
INS_range incremental_nearest_neighbors(Point const &p, FT eps=FT(0)) const
Search incrementally for the nearest neighbors from a query point.
Definition: Kd_tree_search.h:271
KNS_range k_nearest_neighbors(Point const &p, unsigned int k, bool sorted=true, FT eps=FT(0)) const
Search for the k-nearest neighbors from a query point.
Definition: Kd_tree_search.h:245
Incremental_neighbor_search INS_range
The range returned by an incremental nearest or furthest neighbor search. Its value type is std::pair...
Definition: Kd_tree_search.h:112
Tangential complex data structure.
Definition: Tangential_complex.h:126
Tangential_complex(Point_range points, int intrinsic_dimension, const K &k=K())
Constructor from a range of points. Points are copied into the instance, and a search data structure ...
Definition: Tangential_complex.h:241
Fix_inconsistencies_info fix_inconsistencies_using_perturbation(double max_perturb, double time_limit=-1.)
Attempts to fix inconsistencies by perturbing the point positions.
Definition: Tangential_complex.h:397
void compute_tangential_complex()
Computes the tangential complex.
Definition: Tangential_complex.h:331
int intrinsic_dimension() const
Returns the intrinsic dimension of the manifold.
Definition: Tangential_complex.h:280
Point get_perturbed_point(std::size_t vertex) const
Returns the perturbed position of the point corresponding to the vertex given as parameter.
Definition: Tangential_complex.h:299
void set_max_squared_edge_length(FT max_squared_edge_length)
Sets the maximal possible squared edge length for the edges in the triangulations.
Definition: Tangential_complex.h:2013
Point get_point(std::size_t vertex) const
Returns the point corresponding to the vertex given as parameter.
Definition: Tangential_complex.h:292
Num_inconsistencies number_of_inconsistent_simplices(bool verbose=false) const
Definition: Tangential_complex.h:553
int ambient_dimension() const
Returns the ambient dimension.
Definition: Tangential_complex.h:283
std::size_t number_of_vertices() const
Returns the number of vertices.
Definition: Tangential_complex.h:303
int create_complex(Simplex_tree_ &tree, bool export_inconsistent_simplices=true) const
Exports the complex into a Simplex_tree.
Definition: Tangential_complex.h:612
~Tangential_complex()
Destructor.
Definition: Tangential_complex.h:273
Type returned by Tangential_complex::fix_inconsistencies_using_perturbation.
Definition: Tangential_complex.h:379
std::size_t final_num_inconsistent_stars
final number of inconsistent stars
Definition: Tangential_complex.h:389
bool success
true if all inconsistencies could be removed, false if the time limit has been reached before
Definition: Tangential_complex.h:381
unsigned int num_steps
number of steps performed
Definition: Tangential_complex.h:383
std::size_t best_num_inconsistent_stars
best number of inconsistent stars during the process
Definition: Tangential_complex.h:387
std::size_t initial_num_inconsistent_stars
initial number of inconsistent stars
Definition: Tangential_complex.h:385
Type returned by Tangential_complex::number_of_inconsistent_simplices.
Definition: Tangential_complex.h:541
std::size_t num_inconsistent_simplices
Number of inconsistent simplices.
Definition: Tangential_complex.h:545
std::size_t num_inconsistent_stars
Number of stars containing at least one inconsistent simplex.
Definition: Tangential_complex.h:547
std::size_t num_simplices
Total number of simplices in stars (including duplicates that appear in several stars)
Definition: Tangential_complex.h:543