int -> size_t / Cartesian -> ExactExact

This commit is contained in:
Andreas Fabri 2010-06-25 09:45:35 +00:00
parent 1358700d4c
commit 4e2a83ded1
7 changed files with 28 additions and 26 deletions

View File

@ -30,6 +30,7 @@ and for the triangulation data structure (Tds).
\ccTypedef{typedef Gt::Segment_2 Segment;}{the segment type}
\ccTypedef{typedef Gt::Circle_2 Circle;}{the circle type}
\ccTypedef{typedef Gt::FT Numb_type;}{the representation field number type.}
\ccNestedType{Triangulation::size_type}{the size type of the underlying triangulation.}
\ccNestedType{Triangulation::Vertex}{the vertex type of the underlying triangulation.}
\ccNestedType{Triangulation::Edge}{the edge type of the underlying triangulation.}
\ccNestedType{Triangulation::Vertex_handle }{handles to vertices.}
@ -71,7 +72,7 @@ If $v$ is the only vertex in \ccVar\ , $NULL$ is returned.
}
\ccMethod{template<class OutputIterator>
OutputIterator nearest_neighbors(Point p, int k, OutputIterator res);}
OutputIterator nearest_neighbors(Point p, size_type k, OutputIterator res);}
{ computes the $k$ nearest neighbors of $p$ in \ccVar, and places the
handles to the corresponding vertices as a sequence of objects of type
Vertex\_handle in a container of value type of $res$
@ -80,7 +81,7 @@ returns an output iterator pointing to the position beyond the end
of the sequence. }
\ccMethod{template<class OutputIterator>
OutputIterator nearest_neighbors(Vertex_handle v, int k,OutputIterator res);}
OutputIterator nearest_neighbors(Vertex_handle v, size_type k,OutputIterator res);}
{ computes the $k$ nearest neighbors of $v$, and places them as a sequence of objects of type
Vertex\_handle in a container of value type of $res$
which points to the first object in the sequence. The function

View File

@ -11,7 +11,7 @@ other taking a vertex handle.
\ccInclude{CGAL/nearest_neighbor_delaunay_2.h}
\ccFunction{template<class Dt, class OutputIterator>
OutputIterator nearest_neighbors(Dt& delau, const Dt::Point& p, int k, OutputIterator res);}
OutputIterator nearest_neighbors(Dt& delau, const Dt::Point& p, Dt::size_type k, OutputIterator res);}
{computes the $k$ nearest neighbors of $p$ in $delau$, and places the
handles to the corresponding vertices as a sequence of objects of type
Vertex\_handle in a container of value type of $res$
@ -37,7 +37,7 @@ the Delaunay triangulation data type:
\end{itemize}
\ccFunction{template<class Dt, class OutputIterator>
OutputIterator nearest_neighbors(Dt& delau, Dt::Vertex_handle v, int k, OutputIterator res);}
OutputIterator nearest_neighbors(Dt& delau, Dt::Vertex_handle v, Dt::size_type k, OutputIterator res);}
{computes the $k$ nearest neighbors of $v$ (including $v$) in $delau$, and places them as a sequence of objects of type
Vertex\_handle in a container of value type of $res$
which points to the first object in the sequence. The function

View File

@ -54,6 +54,7 @@ public:
typedef Triangulation_2<Gt,Tds> Triangulation;
typedef typename Triangulation::size_type size_type;
typedef typename Triangulation::Locate_type Locate_type;
typedef typename Triangulation::Face_handle Face_handle;
typedef typename Triangulation::Vertex_handle Vertex_handle;
@ -175,9 +176,9 @@ public:
}
template<class OutputIterator>
OutputIterator nearest_neighbors(Point p, int k, OutputIterator res)
OutputIterator nearest_neighbors(Point p, size_type k, OutputIterator res)
{
int n = number_of_vertices();
size_type n = number_of_vertices();
if ( k <= 0 ) return res;
if ( n <= k ) { // return all finite vertices ...
@ -213,9 +214,9 @@ public:
}
template<class OutputIterator>
OutputIterator nearest_neighbors(Vertex_handle v, int k,OutputIterator res)
OutputIterator nearest_neighbors(Vertex_handle v, size_type k,OutputIterator res)
{
int n = number_of_vertices();
size_type n = number_of_vertices();
if ( k <= 0 ) return res;
if ( n <= k ) { // return all (finite) vertices ...
@ -234,10 +235,10 @@ public:
void nearest_neighbors_list(Vertex_handle v,
int k,
size_type k,
std::list<Vertex_handle>& res)
{
int n = number_of_vertices();
size_type n = number_of_vertices();
if ( k <= 0 ) return;
if ( n <= k ) { vertices(std::back_inserter(res)); return; }
@ -290,9 +291,9 @@ public:
// dfs
// for marking nodes in search procedures
int cur_mark;
size_type cur_mark;
Unique_hash_map<Vertex_handle, int> mark;
Unique_hash_map<Vertex_handle, size_type> mark;
void init_vertex_marks()
{
@ -303,7 +304,7 @@ public:
void init_dfs()
{
cur_mark++;
if (cur_mark == INT_MAX) init_vertex_marks();
if (cur_mark == std::numeric_limits<size_type>::max()) init_vertex_marks();
}
void mark_vertex(Vertex_handle vh)

View File

@ -99,10 +99,11 @@ template<class Dt, class OutputIterator>
OutputIterator nearest_neighbors(Dt& delau, const typename Dt::Point& p, int k, OutputIterator res)
{
typedef typename Dt::Geom_traits Gt;
typedef typename Dt::size_type size_type;
typedef typename Dt::Vertex_handle Vertex_handle;
typedef typename Dt::Vertex_iterator Vertex_iterator;
int n = delau.number_of_vertices();
size_type n = delau.number_of_vertices();
if ( k <= 0 ) return res;
if ( n <= k ) { // return all finite vertices ...
@ -146,10 +147,11 @@ template<class Dt, class OutputIterator>
OutputIterator nearest_neighbors(const Dt& delau, typename Dt::Vertex_handle v, int k, OutputIterator res)
{
typedef typename Dt::Geom_traits Gt;
typedef typename Dt::size_type size_type;
typedef typename Dt::Vertex_handle Vertex_handle;
typedef typename Dt::Vertex_iterator Vertex_iterator;
int n = delau.number_of_vertices();
size_type n = delau.number_of_vertices();
if ( k <= 0 ) return res;
if ( n <= k ) { // return all (finite) vertices ...
@ -188,6 +190,7 @@ template<class Dt, class T2>
void nearest_neighbors_list(const Dt& delau, typename Dt::Vertex_handle v, int k, std::list<T2>& res)
{
typedef typename Dt::Geom_traits Gt;
typedef typename Dt::size_type size_type;
typedef typename Dt::Vertex_handle Vertex_handle;
typedef typename Dt::Vertex_iterator Vertex_iterator;
typedef typename Dt::Vertex_circulator Vertex_circulator;
@ -197,7 +200,7 @@ void nearest_neighbors_list(const Dt& delau, typename Dt::Vertex_handle v, int k
typedef typename Gt::Compute_squared_distance_2 Compute_squared_distance_2;
typedef Unique_hash_map<Vertex_handle, Numb_type> MAP_TYPE;
int n = delau.number_of_vertices();
size_type n = delau.number_of_vertices();
if ( k <= 0 ) return;
if ( n <= k ) {

View File

@ -1,13 +1,12 @@
#include <CGAL/basic.h>
#include <list>
#include <vector>
#include <CGAL/Cartesian.h>
#include <CGAL/Exact_predicates_inexact_constructions_kernel.h>
#include <CGAL/Point_set_2.h>
using namespace CGAL;
using namespace std;
typedef Cartesian<double> K;
typedef CGAL::Exact_predicates_inexact_constructions_kernel K;
typedef CGAL::Point_set_2<K>::Edge Edge;
typedef CGAL::Point_set_2<K>::Edge_iterator Edge_iterator;

View File

@ -1,14 +1,12 @@
#include <CGAL/basic.h>
#include <list>
#include <vector>
#include <CGAL/Cartesian.h>
#include <CGAL/Exact_predicates_inexact_constructions_kernel.h>
#include <CGAL/nearest_neighbor_delaunay_2.h>
using namespace CGAL;
using namespace std;
typedef double coord_type;
typedef CGAL::Cartesian<coord_type> K;
typedef CGAL::Exact_predicates_inexact_constructions_kernel K;
typedef CGAL::Delaunay_triangulation_2<K> Delaunay;
typedef CGAL::Delaunay_triangulation_2<K>::Edge Edge;
typedef CGAL::Delaunay_triangulation_2<K>::Edge_iterator Edge_iterator;

View File

@ -1,9 +1,9 @@
#include <CGAL/Simple_cartesian.h>
#include <CGAL/Exact_predicates_inexact_constructions_kernel.h>
#include <CGAL/range_search_delaunay_2.h>
#include <list>
typedef double coord_type;
typedef CGAL::Simple_cartesian<coord_type> Gt;
typedef CGAL::Exact_predicates_inexact_constructions_kernel Gt;
typedef CGAL::Delaunay_triangulation_2<Gt> Delaunay;
typedef CGAL::Delaunay_triangulation_2<Gt>::Edge_iterator Edge_iterator;
typedef CGAL::Delaunay_triangulation_2<Gt>::Vertex_handle Vertex_handle;