No need for a geometric traits and EPEC

This commit is contained in:
Andreas Fabri 2015-06-10 12:28:22 +02:00
parent 2b5ad00c42
commit e6420cae53
10 changed files with 71 additions and 37 deletions

View File

@ -11,16 +11,16 @@ data structure that describes the surface. The vertices and facets of the 2D tr
store handles to the vertices and faces of the 3D triangulation, which enables the user to explore the
2D as well as 3D neighborhood of vertices and facets of the surface.
\tparam Gt must be a model of `Kernel`.
\tparam Dt must be a `Delaunay_triangulation_3` with
`Advancing_front_surface_reconstruction_vertex_base_3` and `Advancing_front_surface_reconstruction_cell_base_3` blended into the vertex and cell type, and the geometric traits class must be the `Exact_predicates_inexact_constructions_kernel`.
`Advancing_front_surface_reconstruction_vertex_base_3` and `Advancing_front_surface_reconstruction_cell_base_3` blended into the vertex and cell type.
The default uses the `Exact_predicates_inexact_constructions_kernel` as geometric traits class.
\tparam Filter must be a functor with `bool operator()(Gt::Point_3,Gt::Point_3,Gt::Point_3)` that allows the user to filter candidate triangles, for example based on its size.
It defaults to a functor that always returns `false`.
\tparam Filter must be a functor with `bool operator()(Point,Point,Point)` that allows the user to filter candidate triangles, for example based on its size.
The type `Point` must be the point type of the geometric traits class of the triangulation. It defaults to a functor that always returns `false`.
*/
template< typename Gt, typename Dt, typename Filter>
template< typename Dt, typename Filter>
class Advancing_front_surface_reconstruction {
public:
@ -199,7 +199,8 @@ has_on_surface(Triangulation_data_structure_2::Vertex_handle v2) const;
For a sequence of points computes a sequence of index triples
describing the faces of the reconstructed surface.
\tparam PointInputIterator must be an input iterator with 3D points from the `Exact_predicates_inexact_constructions_kernel` as value type.
\tparam PointInputIterator must be an input iterator with 3D points as value type. This point type must
be convertible to `Exact_predicates_inexact_constructions_kernel::Point_3` with the `Cartesian_converter`.
\tparam IndicesOutputIterator must be an output iterator to which
`CGAL::cpp11::tuple<std::size_t,std::size_t,std::size_t>` can be assigned.
@ -223,7 +224,8 @@ describing the faces of the reconstructed surface.
For a sequence of points computes a sequence of index triples
describing the faces of the reconstructed surface.
\tparam PointInputIterator must be an input iterator with 3D points from the `Exact_predicates_inexact_constructions_kernel` as value type.
\tparam PointInputIterator must be an input iterator with 3D points as value type. This point type must
be convertible to `Exact_predicates_inexact_constructions_kernel::Point_3` with the `Cartesian_converter`.
\tparam IndicesOutputIterator must be an output iterator to which
`CGAL::cpp11::tuple<std::size_t,std::size_t,std::size_t>` can be assigned.
\tparam Filter must be a functor with `bool operator()(Gt::Point_3,Gt::Point_3,Gt::Point_3)`.

View File

@ -5,7 +5,7 @@
#include <boost/foreach.hpp>
typedef CGAL::Exact_predicates_inexact_constructions_kernel K;
typedef CGAL::Advancing_front_surface_reconstruction<K> Reconstruction;
typedef CGAL::Advancing_front_surface_reconstruction<> Reconstruction;
typedef Reconstruction::Triangulation_3 Triangulation_3;
typedef Reconstruction::Outlier_range Outlier_range;
typedef Reconstruction::Boundary_range Boundary_range;

View File

@ -5,7 +5,7 @@
#include <CGAL/Advancing_front_surface_reconstruction.h>
typedef CGAL::Exact_predicates_inexact_constructions_kernel K;
typedef CGAL::Advancing_front_surface_reconstruction<K> Reconstruction;
typedef CGAL::Advancing_front_surface_reconstruction<> Reconstruction;
typedef Reconstruction::Triangulation_3 Triangulation_3;
typedef Reconstruction::Triangulation_data_structure_2 TDS_2;
typedef K::Point_3 Point_3;

View File

@ -1,12 +1,12 @@
#include <iostream>
#include <fstream>
#include <algorithm>
#include <CGAL/Exact_predicates_inexact_constructions_kernel.h>
#include <CGAL/Simple_cartesian.h>
#include <CGAL/Advancing_front_surface_reconstruction.h>
#include <CGAL/tuple.h>
typedef CGAL::Exact_predicates_inexact_constructions_kernel Kernel;
typedef Kernel::Point_3 Point_3;
typedef CGAL::Simple_cartesian<double> K;
typedef K::Point_3 Point_3;
typedef CGAL::cpp11::tuple<std::size_t,std::size_t,std::size_t> Facet;
@ -27,8 +27,16 @@ struct Perimeter {
: bound(bound)
{}
bool operator()(const Point_3& p, const Point_3& q, const Point_3& r) const
// The point type that will be injected here will be
// CGAL::Exact_predicates_inexact_constructions_kernel::Point_3
template <typename Point>
bool operator()(const Point& p, const Point& q, const Point& r) const
{
// bound == 0 is better than bound < infinity
// as it avoids the distance computations
if(bound == 0){
return false;
}
double d = sqrt(squared_distance(p,q));
if(d>bound) return true;
d += sqrt(squared_distance(p,r)) ;
@ -48,7 +56,7 @@ int main(int argc, char* argv[])
std::istream_iterator<Point_3>(),
std::back_inserter(points));
Perimeter perimeter(0.5);
Perimeter perimeter(0);
CGAL::advancing_front_surface_reconstruction(points.begin(),
points.end(),
std::back_inserter(facets),

View File

@ -26,6 +26,7 @@
#include <CGAL/Exact_predicates_inexact_constructions_kernel.h>
#include <CGAL/Triangulation_data_structure_3.h>
#include <CGAL/Delaunay_triangulation_3.h>
#include <CGAL/Cartesian_converter.h>
#include <cstdio>
#include <cstring>
@ -166,14 +167,15 @@ public:
};
template <class Kernel,
class Triangulation = Delaunay_triangulation_3<Kernel, Triangulation_data_structure_3<Advancing_front_surface_reconstruction_vertex_base_3<Kernel>, Advancing_front_surface_reconstruction_cell_base_3<Kernel> > >,
template <
class Triangulation = Delaunay_triangulation_3<Exact_predicates_inexact_constructions_kernel, Triangulation_data_structure_3<Advancing_front_surface_reconstruction_vertex_base_3<Exact_predicates_inexact_constructions_kernel>, Advancing_front_surface_reconstruction_cell_base_3<Exact_predicates_inexact_constructions_kernel> > >,
class Reject = Always_false>
class Advancing_front_surface_reconstruction {
public:
typedef Triangulation Triangulation_3;
typedef Advancing_front_surface_reconstruction<Kernel,Triangulation_3,Reject> Extract;
typedef typename Triangulation_3::Geom_traits Kernel;
typedef Advancing_front_surface_reconstruction<Triangulation_3,Reject> Extract;
typedef typename Triangulation_3::Geom_traits Geom_traits;
typedef typename Kernel::FT coord_type;
@ -2305,11 +2307,30 @@ namespace AFSR {
template <typename T>
struct Auto_count : public std::unary_function<const T&,std::pair<T,int> >{
mutable int i;
Auto_count() : i(0){}
Auto_count()
: i(0)
{}
std::pair<T,int> operator()(const T& p) const {
return std::make_pair(p,i++);
}
};
template <typename T, typename CC>
struct Auto_count_cc : public std::unary_function<const T&,std::pair<T,int> >{
mutable int i;
CC cc;
Auto_count_cc(CC cc)
: i(0), cc(cc)
{}
template <typename T2>
std::pair<T,int> operator()(const T2& p) const {
return std::make_pair(cc(p),i++);
}
};
}
@ -2328,7 +2349,7 @@ advancing_front_surface_reconstruction(PointIterator b,
typedef Triangulation_data_structure_3<LVb,LCb> Tds;
typedef Delaunay_triangulation_3<Kernel,Tds> Triangulation_3;
typedef Advancing_front_surface_reconstruction<Kernel,Triangulation_3> Reconstruction;
typedef Advancing_front_surface_reconstruction<Triangulation_3> Reconstruction;
typedef Kernel::Point_3 Point_3;
Triangulation_3 dt( boost::make_transform_iterator(b, AFSR::Auto_count<Point_3>()),
@ -2360,11 +2381,14 @@ advancing_front_surface_reconstruction(PointIterator b,
typedef Triangulation_data_structure_3<LVb,LCb> Tds;
typedef Delaunay_triangulation_3<Kernel,Tds> Triangulation_3;
typedef Advancing_front_surface_reconstruction<Kernel,Triangulation_3,Filter> Reconstruction;
typedef Advancing_front_surface_reconstruction<Triangulation_3,Filter> Reconstruction;
typedef std::iterator_traits<PointIterator>::value_type InputPoint;
typedef typename Kernel_traits<InputPoint>::Kernel InputKernel;
typedef Cartesian_converter<InputKernel,Kernel> CC;
typedef Kernel::Point_3 Point_3;
Triangulation_3 dt( boost::make_transform_iterator(b, AFSR::Auto_count<Point_3>()),
boost::make_transform_iterator(e, AFSR::Auto_count<Point_3>() ) );
CC cc=CC();
Triangulation_3 dt( boost::make_transform_iterator(b, AFSR::Auto_count_cc<Point_3,CC>(cc)),
boost::make_transform_iterator(e, AFSR::Auto_count_cc<Point_3,CC>(cc) ) );
AFSR_options opt;
opt.K = radius_ratio_bound;
@ -2391,7 +2415,7 @@ advancing_front_surface_reconstructionP(PointIterator b,
typedef Triangulation_data_structure_3<LVb,LCb> Tds;
typedef Delaunay_triangulation_3<Kernel,Tds> Triangulation_3;
typedef Advancing_front_surface_reconstruction<Kernel,Triangulation_3,Filter> Reconstruction;
typedef Advancing_front_surface_reconstruction<Triangulation_3,Filter> Reconstruction;
typedef typename Kernel::Point_3 Point_3;
Triangulation_3 dt( boost::make_transform_iterator(b, AFSR::Auto_count<Point_3>()),
@ -2419,7 +2443,7 @@ advancing_front_surface_reconstructionP(PointIterator b,
typedef Triangulation_data_structure_3<LVb,LCb> Tds;
typedef Delaunay_triangulation_3<Kernel,Tds> Triangulation_3;
typedef Advancing_front_surface_reconstruction<Kernel,Triangulation_3> Reconstruction;
typedef Advancing_front_surface_reconstruction<Triangulation_3> Reconstruction;
typedef typename Kernel::Point_3 Point_3;
Triangulation_3 dt( boost::make_transform_iterator(b, AFSR::Auto_count<Point_3>()),

View File

@ -31,7 +31,7 @@
namespace CGAL {
template <class A, class B, class C> class Advancing_front_surface_reconstruction;
template <class B, class C> class Advancing_front_surface_reconstruction;
template <class K, class VertexBase = Triangulation_vertex_base_3<K> >
class Advancing_front_surface_reconstruction_vertex_base_3 : public VertexBase
@ -45,7 +45,7 @@ public:
};
template <class A, class B,class C> friend class Advancing_front_surface_reconstruction;
template <class B,class C> friend class Advancing_front_surface_reconstruction;
typedef VertexBase Base;

View File

@ -25,7 +25,7 @@
namespace CGAL {
template <class Kernel, class Triangulation, class Filter>
template <class Triangulation, class Filter>
class Advancing_front_polyhedron_reconstruction;
namespace AFSR {

View File

@ -22,15 +22,15 @@
namespace CGAL {
template <class Kernel, class Triangulation, class Filter>
template <class Triangulation, class Filter>
class Advancing_front_surface_reconstruction;
namespace AFSR {
template <class Kernel, class Triangulation, class TDS, class Filter>
template <class Triangulation, class TDS, class Filter>
typename TDS::Vertex_handle
construct_surface(TDS& tds, const CGAL::Advancing_front_surface_reconstruction<Kernel,Triangulation,Filter>& surface)
construct_surface(TDS& tds, const CGAL::Advancing_front_surface_reconstruction<Triangulation,Filter>& surface)
{
typedef typename TDS::Vertex_handle Vertex_handle;

View File

@ -23,9 +23,9 @@ namespace CGAL {
namespace AFSR {
template <class Kernel, class Triangulation, class TDS, class Filter>
template <class Triangulation, class TDS, class Filter>
typename TDS::Vertex_handle
orient(TDS& tds, const Advancing_front_surface_reconstruction<Kernel,Triangulation,Filter>& surface)
orient(TDS& tds, const Advancing_front_surface_reconstruction<Triangulation,Filter>& surface)
{
typedef typename TDS::Vertex_handle Vertex_handle;

View File

@ -24,16 +24,16 @@
namespace CGAL {
template <class Kernel, class Triangulation, class Filter>
template <class Triangulation, class Filter>
class Advancing_front_surface_reconstruction;
template <class OutputIterator, class Kernel, class Triangulation, class Filter>
template <class OutputIterator, class Triangulation, class Filter>
OutputIterator
write_triple_indices(OutputIterator out, const Advancing_front_surface_reconstruction<Kernel,Triangulation,Filter>& S)
write_triple_indices(OutputIterator out, const Advancing_front_surface_reconstruction<Triangulation,Filter>& S)
{
typedef Advancing_front_surface_reconstruction<Kernel,Triangulation,Filter> Surface;
typedef Advancing_front_surface_reconstruction<Triangulation,Filter> Surface;
typedef typename Surface::TDS_2 TDS_2;
typedef typename TDS_2::Face_iterator Face_iterator;