\( \newcommand{\E}{\mathrm{E}} \) \( \newcommand{\A}{\mathrm{A}} \) \( \newcommand{\R}{\mathrm{R}} \) \( \newcommand{\N}{\mathrm{N}} \) \( \newcommand{\Q}{\mathrm{Q}} \) \( \newcommand{\Z}{\mathrm{Z}} \) \( \def\ccSum #1#2#3{ \sum_{#1}^{#2}{#3} } \def\ccProd #1#2#3{ \sum_{#1}^{#2}{#3} }\)
CGAL 4.13 - 3D Alpha Shapes
Alpha_shapes_3/ex_alpha_shapes_3.cpp
#include <CGAL/Exact_predicates_inexact_constructions_kernel.h>
#include <CGAL/Alpha_shape_3.h>
#include <CGAL/Alpha_shape_cell_base_3.h>
#include <CGAL/Alpha_shape_vertex_base_3.h>
#include <CGAL/Delaunay_triangulation_3.h>
#include <cassert>
#include <fstream>
#include <list>
typedef CGAL::Triangulation_data_structure_3<Vb,Fb> Tds;
typedef Gt::Point_3 Point;
typedef Alpha_shape_3::Alpha_iterator Alpha_iterator;
int main()
{
std::list<Point> lp;
//read input
std::ifstream is("./data/bunny_1000");
int n;
is >> n;
std::cout << "Reading " << n << " points " << std::endl;
Point p;
for( ; n>0 ; n--)
{
is >> p;
lp.push_back(p);
}
// compute alpha shape
Alpha_shape_3 as(lp.begin(),lp.end());
std::cout << "Alpha shape computed in REGULARIZED mode by default" << std::endl;
// find optimal alpha value
Alpha_iterator opt = as.find_optimal_alpha(1);
std::cout << "Optimal alpha value to get one connected component is " << *opt << std::endl;
as.set_alpha(*opt);
assert(as.number_of_solid_components() == 1);
return 0;
}