\( \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.5.1 - 3D Triangulations
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Groups Pages
Triangulation_3/info_insert_with_pair_iterator.cpp
#include <CGAL/Exact_predicates_inexact_constructions_kernel.h>
#include <CGAL/Delaunay_triangulation_3.h>
#include <CGAL/Triangulation_vertex_base_with_info_3.h>
#include <vector>
//Use the Fast_location tag. Default or Compact_location works too.
typedef Delaunay::Point Point;
int main()
{
std::vector< std::pair<Point,unsigned> > points;
points.push_back( std::make_pair(Point(0,0,0),0) );
points.push_back( std::make_pair(Point(1,0,0),1) );
points.push_back( std::make_pair(Point(0,1,0),2) );
points.push_back( std::make_pair(Point(0,0,1),3) );
points.push_back( std::make_pair(Point(2,2,2),4) );
points.push_back( std::make_pair(Point(-1,0,1),5) );
Delaunay T( points.begin(),points.end() );
CGAL_assertion( T.number_of_vertices() == 6 );
// check that the info was correctly set.
Delaunay::Finite_vertices_iterator vit;
for (vit = T.finite_vertices_begin(); vit != T.finite_vertices_end(); ++vit)
if( points[ vit->info() ].first != vit->point() ){
std::cerr << "Error different info" << std::endl;
exit(EXIT_FAILURE);
}
std::cout << "OK" << std::endl;
return 0;
}