\( \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.4 - 2D Alpha Shapes
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Groups Pages
Alpha_shapes_2/ex_alpha_shapes_2.cpp
#include <CGAL/Exact_predicates_inexact_constructions_kernel.h>
#include <CGAL/algorithm.h>
#include <CGAL/Delaunay_triangulation_2.h>
#include <CGAL/Alpha_shape_2.h>
#include <iostream>
#include <fstream>
#include <vector>
#include <list>
typedef K::FT FT;
typedef K::Point_2 Point;
typedef K::Segment_2 Segment;
typedef CGAL::Triangulation_data_structure_2<Vb,Fb> Tds;
typedef CGAL::Delaunay_triangulation_2<K,Tds> Triangulation_2;
typedef Alpha_shape_2::Alpha_shape_edges_iterator Alpha_shape_edges_iterator;
template <class OutputIterator>
void
alpha_edges( const Alpha_shape_2& A,
{
for(Alpha_shape_edges_iterator it = A.alpha_shape_edges_begin();
it != A.alpha_shape_edges_end();
++it){
*out++ = A.segment(*it);
}
}
template <class OutputIterator>
bool
file_input(OutputIterator out)
{
std::ifstream is("./data/fin", std::ios::in);
if(is.fail()){
std::cerr << "unable to open file for input" << std::endl;
return false;
}
int n;
is >> n;
std::cout << "Reading " << n << " points from file" << std::endl;
CGAL::cpp11::copy_n(std::istream_iterator<Point>(is), n, out);
return true;
}
// Reads a list of points and returns a list of segments
// corresponding to the Alpha shape.
int main()
{
std::list<Point> points;
if(! file_input(std::back_inserter(points))){
return -1;
}
Alpha_shape_2 A(points.begin(), points.end(),
FT(10000),
Alpha_shape_2::GENERAL);
std::vector<Segment> segments;
alpha_edges( A, std::back_inserter(segments));
std::cout << "Alpha Shape computed" << std::endl;
std::cout << segments.size() << " alpha shape edges" << std::endl;
std::cout << "Optimal alpha: " << *A.find_optimal_alpha(1)<<std::endl;
return 0;
}