#include <gudhi/Simplex_tree.h>
#include <gudhi/Points_3D_off_io.h>
 
#include <boost/variant.hpp>
 
#include <CGAL/Exact_predicates_inexact_constructions_kernel.h>
#include <CGAL/Delaunay_triangulation_3.h>
#include <CGAL/Alpha_shape_3.h>
#include <CGAL/Alpha_shape_vertex_base_3.h>
#include <CGAL/Alpha_shape_cell_base_3.h>
#include <CGAL/iterator.h>
 
#include <fstream>
#include <cmath>
#include <string>
#include <tuple>  
#include <map>
#include <utility>  
#include <list>
#include <vector>
 
typedef CGAL::Exact_predicates_inexact_constructions_kernel Kernel;
typedef CGAL::Alpha_shape_vertex_base_3<Kernel> Vb;
typedef CGAL::Alpha_shape_cell_base_3<Kernel> Fb;
typedef CGAL::Triangulation_data_structure_3<Vb, Fb> Tds;
typedef CGAL::Delaunay_triangulation_3<Kernel, Tds> Triangulation_3;
typedef CGAL::Alpha_shape_3<Triangulation_3> Alpha_shape_3;
 
typedef Kernel::Point_3 Point;
 
typedef Alpha_shape_3::FT Alpha_value_type;
typedef CGAL::Object Object;
typedef CGAL::Dispatch_output_iterator<
CGAL::cpp11::tuple<Object, Alpha_value_type>,
CGAL::cpp11::tuple<std::back_insert_iterator< std::vector<Object> >,
    std::back_insert_iterator< std::vector<Alpha_value_type> > > > Dispatch;
typedef Alpha_shape_3::Cell_handle Cell_handle;
typedef Alpha_shape_3::Facet Facet;
typedef Alpha_shape_3::Edge Edge;
typedef std::list<Alpha_shape_3::Vertex_handle> Vertex_list;
 
typedef std::map<Alpha_shape_3::Vertex_handle, Simplex_tree_vertex > Alpha_shape_simplex_tree_map;
typedef std::pair<Alpha_shape_3::Vertex_handle, Simplex_tree_vertex> Alpha_shape_simplex_tree_pair;
typedef std::vector< Simplex_tree_vertex > Simplex_tree_vector_vertex;
 
Vertex_list from(const Cell_handle& ch) {
  Vertex_list the_list;
  for (auto i = 0; i < 4; i++) {
#ifdef DEBUG_TRACES
    std::clog << "from cell[" << i << "]=" << ch->vertex(i)->point() << std::endl;
#endif  
    the_list.push_back(ch->vertex(i));
  }
  return the_list;
}
 
Vertex_list from(const Facet& fct) {
  Vertex_list the_list;
  for (auto i = 0; i < 4; i++) {
    if (fct.second != i) {
#ifdef DEBUG_TRACES
      std::clog << "from facet=[" << i << "]" << fct.first->vertex(i)->point() << std::endl;
#endif  
      the_list.push_back(fct.first->vertex(i));
    }
  }
  return the_list;
}
 
Vertex_list from(const Edge& edg) {
  Vertex_list the_list;
  for (auto i = 0; i < 4; i++) {
    if ((edg.second == i) || (edg.third == i)) {
#ifdef DEBUG_TRACES
      std::clog << "from edge[" << i << "]=" << edg.first->vertex(i)->point() << std::endl;
#endif  
      the_list.push_back(edg.first->vertex(i));
    }
  }
  return the_list;
}
 
Vertex_list from(const Alpha_shape_3::Vertex_handle& vh) {
  Vertex_list the_list;
#ifdef DEBUG_TRACES
  std::clog << "from vertex=" << vh->point() << std::endl;
#endif  
  the_list.push_back(vh);
  return the_list;
}
 
int main(int argc, char * const argv[]) {
  
  if (argc != 2) {
    std::cerr << "Usage: " << argv[0]
        << " path_to_off_file \n";
    return 0;
  }
 
  
  std::string offInputFile(argv[1]);
  
  
  if (!off_reader.is_valid()) {
    std::cerr << "Unable to read file " << argv[1] << std::endl;
    return 0;
  }
  
  std::vector<Point> lp = off_reader.get_point_cloud();
 
  
  Alpha_shape_3 as(lp.begin(), lp.end(), 0, Alpha_shape_3::GENERAL);
#ifdef DEBUG_TRACES
  std::clog << "Alpha shape computed in GENERAL mode" << std::endl;
#endif  
 
  
  std::vector<Object> the_objects;
  std::vector<Alpha_value_type> the_alpha_values;
 
  Dispatch disp = CGAL::dispatch_output<Object, Alpha_value_type>(std::back_inserter(the_objects),
      std::back_inserter(the_alpha_values));
 
  as.filtration_with_alpha_values(disp);
#ifdef DEBUG_TRACES
  std::clog << "filtration_with_alpha_values returns : " << the_objects.size() << " objects" << std::endl;
#endif  
 
  Alpha_shape_3::size_type count_vertices = 0;
  Alpha_shape_3::size_type count_edges = 0;
  Alpha_shape_3::size_type count_facets = 0;
  Alpha_shape_3::size_type count_cells = 0;
 
  
  Vertex_list vertex_list;
  Alpha_shape_simplex_tree_map map_cgal_simplex_tree;
  std::vector<Alpha_value_type>::iterator the_alpha_value_iterator = the_alpha_values.begin();
  for (auto object_iterator : the_objects) {
    
    if (const Cell_handle * cell = CGAL::object_cast<Cell_handle>(&object_iterator)) {
      vertex_list = from(*cell);
      count_cells++;
    } else if (const Facet * facet = CGAL::object_cast<Facet>(&object_iterator)) {
      vertex_list = from(*facet);
      count_facets++;
    } else if (const Edge * edge = CGAL::object_cast<Edge>(&object_iterator)) {
      vertex_list = from(*edge);
      count_edges++;
    } else if (const Alpha_shape_3::Vertex_handle * vertex =
              CGAL::object_cast<Alpha_shape_3::Vertex_handle>(&object_iterator)) {
      count_vertices++;
      vertex_list = from(*vertex);
    }
    
    Simplex_tree_vector_vertex the_simplex_tree;
    for (auto the_alpha_shape_vertex : vertex_list) {
      Alpha_shape_simplex_tree_map::iterator the_map_iterator = map_cgal_simplex_tree.find(the_alpha_shape_vertex);
      if (the_map_iterator == map_cgal_simplex_tree.end()) {
        
        Simplex_tree_vertex vertex = map_cgal_simplex_tree.size();
#ifdef DEBUG_TRACES
        std::clog << "vertex [" << the_alpha_shape_vertex->point() << "] not found - insert_simplex " << vertex << "\n";
#endif  
        the_simplex_tree.push_back(vertex);
        map_cgal_simplex_tree.insert(Alpha_shape_simplex_tree_pair(the_alpha_shape_vertex, vertex));
      } else {
        
        Simplex_tree_vertex vertex = the_map_iterator->second;
#ifdef DEBUG_TRACES
        std::clog << "vertex [" << the_alpha_shape_vertex->point() << "] found in " << vertex << std::endl;
#endif  
        the_simplex_tree.push_back(vertex);
      }
    }
    
#ifdef DEBUG_TRACES
    std::clog << "filtration = " << *the_alpha_value_iterator << std::endl;
#endif  
    simplex_tree.
insert_simplex(the_simplex_tree, std::sqrt(*the_alpha_value_iterator));
    if (the_alpha_value_iterator != the_alpha_values.end())
      ++the_alpha_value_iterator;
    else
      std::cerr << "This shall not happen" << std::endl;
  }
#ifdef DEBUG_TRACES
  std::clog << "vertices \t\t" << count_vertices << std::endl;
  std::clog << "edges \t\t" << count_edges << std::endl;
  std::clog << "facets \t\t" << count_facets << std::endl;
  std::clog << "cells \t\t" << count_cells << std::endl;
 
 
  std::clog << "Information of the Simplex Tree:\n";
  std::clog << 
"  Number of vertices = " << simplex_tree.
num_vertices() << 
" ";
  std::clog << 
"  Number of simplices = " << simplex_tree.
num_simplices() << std::endl << std::endl;
#endif  
 
#ifdef DEBUG_TRACES
  std::clog << "Iterator on vertices: \n";
    std::clog << vertex << " ";
  }
#endif  
 
  std::clog << simplex_tree << std::endl;
 
#ifdef DEBUG_TRACES
  std::clog << std::endl << std::endl << "Iterator on simplices:\n";
    std::clog << "   ";
      std::clog << vertex << " ";
    }
    std::clog << std::endl;
  }
#endif  
#ifdef DEBUG_TRACES
  std::clog << std::endl << std::endl << "Iterator on Simplices in the filtration, with [filtration value]:\n";
    std::clog << 
"   " << 
"[" << simplex_tree.
filtration(f_simplex) << 
"] ";
      std::clog << vertex << " ";
    }
    std::clog << std::endl;
  }
#endif  
#ifdef DEBUG_TRACES
  std::clog << std::endl << std::endl << "Iterator on Simplices in the filtration, and their boundary simplices:\n";
    std::clog << 
"   " << 
"[" << simplex_tree.
filtration(f_simplex) << 
"] ";
      std::clog << vertex << " ";
    }
    std::clog << std::endl;
 
      std::clog << 
"      " << 
"[" << simplex_tree.
filtration(b_simplex) << 
"] ";
        std::clog << vertex << " ";
      }
      std::clog << std::endl;
    }
  }
#endif  
 
  return 0;
}
Definition: Points_3D_off_io.h:140
Simplex Tree data structure for representing simplicial complexes.
Definition: Simplex_tree.h:81
Filtration_simplex_range const & filtration_simplex_range(Indexing_tag=Indexing_tag())
Returns a range over the simplices of the simplicial complex, in the order of the filtration.
Definition: Simplex_tree.h:273
std::pair< Simplex_handle, bool > insert_simplex(const InputVertexRange &simplex, Filtration_value filtration=0)
Insert a simplex, represented by a range of Vertex_handles, in the simplicial complex.
Definition: Simplex_tree.h:782
Simplex_vertex_range simplex_vertex_range(Simplex_handle sh) const
Returns a range over the vertices of a simplex.
Definition: Simplex_tree.h:284
static Filtration_value filtration(Simplex_handle sh)
Returns the filtration value of a simplex.
Definition: Simplex_tree.h:537
Options::Vertex_handle Vertex_handle
Type for the vertex handle.
Definition: Simplex_tree.h:96
size_t num_vertices() const
Returns the number of vertices in the complex.
Definition: Simplex_tree.h:574
Boundary_simplex_range boundary_simplex_range(SimplexHandle sh)
Returns a range over the simplices of the boundary of a simplex.
Definition: Simplex_tree.h:305
size_t num_simplices()
returns the number of simplices in the simplex_tree.
Definition: Simplex_tree.h:580
Complex_simplex_range complex_simplex_range()
Returns a range over the simplices of the simplicial complex.
Definition: Simplex_tree.h:239
Complex_vertex_range complex_vertex_range()
Returns a range over the vertices of the simplicial complex. The order is increasing according to < o...
Definition: Simplex_tree.h:228