void compute_normals (Point_set& points, unsigned int neighbors)
  {
    CGAL::jet_estimate_normals<Concurrency_tag>(points.begin_or_selection_begin(), points.end(),
                                                points.point_map(),
                                                points.normal_map(),
                                                2 * neighbors);

    points.set_first_selected (CGAL::mst_orient_normals (points.begin(), points.end(),
                                                         points.point_map(),
                                                         points.normal_map(),
                                                         2 * neighbors));
    points.delete_selection();
  }
int main (int, char**)
{
  Point_set pts;

  pts.add_normal_map();
  
  bool map_added = false;
  Size_t_map echo_map;
  Color_map color_map;

  boost::tie (echo_map, map_added) = pts.add_property_map<std::size_t> ("echo");
  CGAL_assertion (map_added);
  boost::tie (color_map, map_added) = pts.add_property_map<Classification::RGB_Color> ("color");
  CGAL_assertion (map_added);

  for (std::size_t i = 0; i < 1000; ++ i)
  {
    Point_set::iterator it
      = pts.insert (Point (CGAL::get_default_random().get_double(),
                           CGAL::get_default_random().get_double(),
                           CGAL::get_default_random().get_double()),
                    Vector (CGAL::get_default_random().get_double(),
                            CGAL::get_default_random().get_double(),
                            CGAL::get_default_random().get_double()));
    echo_map[*it] = std::size_t(CGAL::get_default_random().get_int(0, 4));
    color_map[*it] = CGAL::make_array ((unsigned char)(CGAL::get_default_random().get_int(0, 255)),
                                       (unsigned char)(CGAL::get_default_random().get_int(0, 255)),
                                       (unsigned char)(CGAL::get_default_random().get_int(0, 255)));
  }

  Feature_set features;
  Feature_generator generator (features, pts, pts.point_map(),
                               5,  // using 5 scales
                               pts.normal_map(),
                               color_map, echo_map);

  CGAL_assertion (generator.number_of_scales() == 5);
  CGAL_assertion (features.size() == 80);

  Label_set labels;

  std::vector<int> training_set (pts.size(), -1);
  for (std::size_t i = 0; i < 20; ++ i)
  {
    std::ostringstream oss;
    oss << "label_" << i;
    Label_handle lh = labels.add(oss.str().c_str());

    for (std::size_t j = 0; j < 10; ++ j)
      training_set[std::size_t(CGAL::get_default_random().get_int(0, int(training_set.size())))] = int(i);
  }
  CGAL_assertion (labels.size() == 20);
  
  Classifier classifier (labels, features);
  
  classifier.train<CGAL::Sequential_tag> (training_set, 800);
#ifdef CGAL_LINKED_WITH_TBB
  classifier.train<CGAL::Parallel_tag> (training_set, 800);
#endif

  std::vector<int> label_indices(pts.size(), -1);
  
  Classification::classify<CGAL::Sequential_tag>
    (pts, labels, classifier, label_indices);

  Classification::classify_with_local_smoothing<CGAL::Sequential_tag>
    (pts, pts.point_map(), labels, classifier,
     generator.neighborhood().sphere_neighbor_query(0.01f),
     label_indices);

  Classification::classify_with_graphcut<CGAL::Sequential_tag>
    (pts, pts.point_map(), labels, classifier,
     generator.neighborhood().k_neighbor_query(12),
     0.2f, 10, label_indices);

#ifdef CGAL_LINKED_WITH_TBB
  Classification::classify<CGAL::Sequential_tag>
    (pts, labels, classifier, label_indices);

  Classification::classify_with_local_smoothing<CGAL::Sequential_tag>
    (pts, pts.point_map(), labels, classifier,
     generator.neighborhood().sphere_neighbor_query(0.01f),
     label_indices);

  Classification::classify_with_graphcut<CGAL::Sequential_tag>
    (pts, pts.point_map(), labels, classifier,
     generator.neighborhood().k_neighbor_query(12),
     0.2f, 10, label_indices);
#endif

  Classification::Evaluation evaluation (labels, training_set, label_indices);
  
  return EXIT_SUCCESS;
}
bool poisson_reconstruct(FaceGraph* graph,
                         Point_set& points,
                         typename Traits::FT sm_angle, // Min triangle angle (degrees).
                         typename Traits::FT sm_radius, // Max triangle size w.r.t. point set average spacing.
                         typename Traits::FT sm_distance, // Approximation error w.r.t. point set average spacing.
                         const QString& solver_name, // solver name
                         bool use_two_passes,
                         bool do_not_fill_holes)
{
  // Poisson implicit function
  typedef CGAL::Poisson_reconstruction_function<Traits> Poisson_reconstruction_function;

  // Surface mesher
  typedef CGAL::Surface_mesh_default_triangulation_3 STr;
  typedef CGAL::Surface_mesh_complex_2_in_triangulation_3<STr> C2t3;
  typedef CGAL::Implicit_surface_3<Traits, Poisson_reconstruction_function> Surface_3;

  // AABB tree
  typedef CGAL::AABB_face_graph_triangle_primitive<FaceGraph> Primitive;
  typedef CGAL::AABB_traits<Traits, Primitive> AABB_traits;
  typedef CGAL::AABB_tree<AABB_traits> AABB_tree;
    CGAL::Timer task_timer; task_timer.start();

    //***************************************
    // Checks requirements
    //***************************************

    if (points.size() == 0)
    {
      std::cerr << "Error: empty point set" << std::endl;
     return false;
    }

    bool points_have_normals = points.has_normal_map();
    if ( ! points_have_normals )
    {
      std::cerr << "Input point set not supported: this reconstruction method requires oriented normals" << std::endl;
     return false;
    }

    CGAL::Timer reconstruction_timer; reconstruction_timer.start();

    //***************************************
    // Computes implicit function
    //***************************************

 
    std::cerr << "Computes Poisson implicit function "
              << "using " << solver_name.toLatin1().data() << " solver...\n";
              
    
    // Creates implicit function from the point set.
    // Note: this method requires an iterator over points
    // + property maps to access each point's position and normal.
    Poisson_reconstruction_function function(points.begin_or_selection_begin(), points.end(),
                                             points.point_map(), points.normal_map());

    bool ok = false;    
    #ifdef CGAL_EIGEN3_ENABLED
    if(solver_name=="Eigen - built-in simplicial LDLt")
    {
      CGAL::Eigen_solver_traits<Eigen::SimplicialCholesky<CGAL::Eigen_sparse_matrix<double>::EigenType> > solver;
      ok = function.compute_implicit_function(solver, use_two_passes);
    }
    if(solver_name=="Eigen - built-in CG")
    {
      CGAL::Eigen_solver_traits<Eigen::ConjugateGradient<CGAL::Eigen_sparse_matrix<double>::EigenType> > solver;
      solver.solver().setTolerance(1e-6);
      solver.solver().setMaxIterations(1000);
      ok = function.compute_implicit_function(solver, use_two_passes);
    }
    #endif

    // Computes the Poisson indicator function f()
    // at each vertex of the triangulation.
    if ( ! ok )
    {
      std::cerr << "Error: cannot compute implicit function" << std::endl;
     return false;
    }

    // Prints status
    std::cerr << "Total implicit function (triangulation+refinement+solver): " << task_timer.time() << " seconds\n";
    task_timer.reset();

    //***************************************
    // Surface mesh generation
    //***************************************

    std::cerr << "Surface meshing...\n";

    // Computes average spacing
    Kernel::FT average_spacing = CGAL::compute_average_spacing<Concurrency_tag>(points.all_or_selection_if_not_empty(),
                                                                                6 /* knn = 1 ring */,
                                                                                points.parameters());

    // Gets one point inside the implicit surface
    Kernel::Point_3 inner_point = function.get_inner_point();
    Kernel::FT inner_point_value = function(inner_point);
    if(inner_point_value >= 0.0)
    {
      std::cerr << "Error: unable to seed (" << inner_point_value << " at inner_point)" << std::endl;
     return false;
    }

    // Gets implicit function's radius
    Kernel::Sphere_3 bsphere = function.bounding_sphere();
    Kernel::FT radius = std::sqrt(bsphere.squared_radius());

    // Defines the implicit surface: requires defining a
  	// conservative bounding sphere centered at inner point.
    Kernel::FT sm_sphere_radius = 5.0 * radius;
    Kernel::FT sm_dichotomy_error = sm_distance*average_spacing/1000.0; // Dichotomy error must be << sm_distance
    Surface_3 surface(function,
                      Kernel::Sphere_3(inner_point,sm_sphere_radius*sm_sphere_radius),
                      sm_dichotomy_error/sm_sphere_radius);

    // Defines surface mesh generation criteria
    CGAL::Surface_mesh_default_criteria_3<STr> criteria(sm_angle,  // Min triangle angle (degrees)
                                                        sm_radius*average_spacing,  // Max triangle size
                                                        sm_distance*average_spacing); // Approximation error

    CGAL_TRACE_STREAM << "  make_surface_mesh(sphere center=("<<inner_point << "),\n"
                      << "                    sphere radius="<<sm_sphere_radius<<",\n"
                      << "                    angle="<<sm_angle << " degrees,\n"
                      << "                    triangle size="<<sm_radius<<" * average spacing="<<sm_radius*average_spacing<<",\n"
                      << "                    distance="<<sm_distance<<" * average spacing="<<sm_distance*average_spacing<<",\n"
                      << "                    dichotomy error=distance/"<<sm_distance*average_spacing/sm_dichotomy_error<<",\n"
                      << "                    Manifold_with_boundary_tag)\n";

    // Generates surface mesh with manifold option
    STr tr; // 3D Delaunay triangulation for surface mesh generation
    C2t3 c2t3(tr); // 2D complex in 3D Delaunay triangulation
    CGAL::make_surface_mesh(c2t3,                                 // reconstructed mesh
                            surface,                              // implicit surface
                            criteria,                             // meshing criteria
                            CGAL::Manifold_with_boundary_tag());  // require manifold mesh

    // Prints status
    std::cerr << "Surface meshing: " << task_timer.time() << " seconds, "
                                     << tr.number_of_vertices() << " output vertices"
                                     << std::endl;
    task_timer.reset();

    if(tr.number_of_vertices() == 0)
     return false;

    // Converts to polyhedron
    CGAL::facets_in_complex_2_to_triangle_mesh(c2t3, *graph);

    // Prints total reconstruction duration
    std::cerr << "Total reconstruction (implicit function + meshing): " << reconstruction_timer.time() << " seconds\n";

    //***************************************
    // Computes reconstruction error
    //***************************************

    // Constructs AABB tree and computes internal KD-tree
    // data structure to accelerate distance queries
    AABB_tree tree(faces(*graph).first, faces(*graph).second, *graph);
    tree.accelerate_distance_queries();

    // Computes distance from each input point to reconstructed mesh
    double max_distance = DBL_MIN;
    double avg_distance = 0;

    std::set<typename boost::graph_traits<FaceGraph>::face_descriptor> faces_to_keep;
    
    for (Point_set::const_iterator p=points.begin_or_selection_begin(); p!=points.end(); p++)
    {
      typename AABB_traits::Point_and_primitive_id pap = tree.closest_point_and_primitive (points.point (*p));
      double distance = std::sqrt(CGAL::squared_distance (pap.first, points.point(*p)));
      
      max_distance = (std::max)(max_distance, distance);
      avg_distance += distance;

      typename boost::graph_traits<FaceGraph>::face_descriptor f = pap.second;
      faces_to_keep.insert (f);  
    }
    avg_distance /= double(points.size());

    std::cerr << "Reconstruction error:\n"
              << "  max = " << max_distance << " = " << max_distance/average_spacing << " * average spacing\n"
              << "  avg = " << avg_distance << " = " << avg_distance/average_spacing << " * average spacing\n";

    if (do_not_fill_holes)
      {
        typename boost::graph_traits<FaceGraph>::face_iterator it = faces(*graph).begin ();
        while (it != faces(*graph).end ())
	  {
            typename boost::graph_traits<FaceGraph>::face_iterator current = it ++;

            if (faces_to_keep.find (*current) == faces_to_keep.end ())
            {
              CGAL::Euler::remove_face(halfedge (*current, *graph), *graph);
            }

	  }

      }
    return true;
}