modifications OutputIterator

This commit is contained in:
Remy Thomasse 2014-04-11 14:02:29 +02:00
parent 9e6fc4f698
commit 03ef4e6789
5 changed files with 111 additions and 97 deletions

View File

@ -3,13 +3,21 @@ namespace CGAL {
/*!
\ingroup PkgGenerators
\brief Computes a random convex polygon by writing its vertices (oriented
counterclockwise) in the list `l`, as the convex hull of \f$ n \f$ random points in a disc centered in \f$0\f$ with radius \f$radius\f$.
counterclockwise) in an `OutputIterator`, as the convex hull of \f$ n \f$ random points in a disc centered in \f$0\f$ with radius \f$radius\f$.
The generated polygon will have an average number of vertices \f$ n^\frac{1}{3}(1+o(1))\f$.
\pre \f$n \geq 3 \f$
\cgalHeading{Requirements}
`Generator` has to be a Boost random generator, such as `boost::random::mt19937`.
\cgalHeading{Requires}
- `Generator` has to be a Boost random generator, such as `boost::random::mt19937`.
- `fast` is a Boolean value, `true` for a time efficiency behavior and `false` for a memory efficiency behavior.
- The `OutputIterator` must accept values of type `Point_2<Traits>`.
\cgalHeading{Implementation}
@ -27,6 +35,6 @@ The following program displays a random simple polygon made of \f$10000\f$ point
*/
template<class P,class Generator>
void random_convex_hull_in_disc_2(size_t n, typename Kernel_traits<P>::Kernel::FT radius, std::list<P> & l,Generator & g, bool fast=true );
template < class OutputIterator, class Traits, class Generator >
void random_convex_hull_in_disc_2(std::size_t n, double radius, Generator & gen, OutputIterator it, const Traits & traits, bool fast=true);
} /* namespace CGAL */

View File

@ -49,6 +49,7 @@ achieve random permutations for otherwise regular generators (
- `CGAL::random_polygon_2()`
- `CGAL::random_selection()`
- `CGAL::random_convex_hull_in_disc_2()`
## Variables ##
- `CGAL::default_random`

View File

@ -1,36 +1,27 @@
#include <CGAL/Simple_cartesian.h>
#include <CGAL/Exact_predicates_inexact_constructions_kernel.h>
#include <CGAL/random_convex_hull_in_disc_2.h>
#include <CGAL/Polygon_2_algorithms.h>
#include <boost/random.hpp>
#include <iostream>
#include <vector>
using namespace CGAL;
typedef Exact_predicates_inexact_constructions_kernel K;
typedef K::Point_2 Point;
typedef K::FT FT;
#ifdef CGAL_USE_GMP
#include <CGAL/Gmpz.h>
typedef CGAL::Gmpz FT;
#else
// NOTE: the choice of double here for a number type may cause problems
// for degenerate point sets
#include <CGAL/double.h>
typedef double FT;
#endif
#include <fstream>
#include <list>
typedef CGAL::Simple_cartesian<FT> K;
typedef K::Point_2 Point_2;
const int n=10000;
const FT radius=1.0;
const FT RADIUS=1.0;
int main( )
{
std::list<Point_2> point_set;
int N=10000;
std::vector<Point> v;
boost::random::mt19937 gen;
gen.seed(time(0));
CGAL::random_convex_hull_in_disc_2(n,radius,point_set,gen);
int size = point_set.size();
FT area=CGAL::polygon_area_2(point_set.begin(),point_set.end(),K());
std::cout<<"A random convex polygon inscribed in a disc with "<<size<<" vertices and area"<<area<<"has been generated."<<std::endl;;
random_convex_hull_in_disc_2(N,RADIUS,gen,std::back_inserter(v),K());
size_t size = v.size();
FT area=polygon_area_2(v.begin(),v.end(),K());
std::cout<<"A random convex polygon inscribed in a disc with "<<size<<" vertices and area "<<area<<" has been generated."<<std::endl;
return 0;
}

View File

@ -40,28 +40,28 @@
namespace CGAL{
namespace internal{
template<class P >
struct compare_points_angle{
bool operator()(const P& p, const P& q ){
typedef typename Kernel_traits<P>::Kernel Traits;
Traits ch_traits;
typedef typename Traits::Left_turn_2 Left_turn;
Left_turn left_turn = ch_traits.left_turn_2_object();
if ( to_double(p.y())>0 )
{
if (to_double(q.y())>0) return left_turn(ORIGIN,p,q);
else return false;
struct compare_points_angle{
bool operator()(const P& p, const P& q ){
typedef typename Kernel_traits<P>::Kernel Traits;
Traits ch_traits;
typedef typename Traits::Left_turn_2 Left_turn;
Left_turn left_turn = ch_traits.left_turn_2_object();
if ( to_double(p.y())>0 )
{
if (to_double(q.y())>0) return left_turn(ORIGIN,p,q);
else return false;
}
else {
if (to_double(q.y())>0) return true;
else return left_turn(ORIGIN,p,q);
}
}
else {
if (to_double(q.y())>0) return true;
else return left_turn(ORIGIN,p,q);
}
}
};
};
//////////////////////////////////////
template <class P, class GEN>
void generate_points_annulus(long n,double a, double b,double small_radius, double big_radius,std::list<P> & l,GEN & gen){ //generate n points between a and b
@ -105,14 +105,14 @@ namespace CGAL{
template<class P, class Traits >
void Graham_without_sort_2( std::list<P> &l, const Traits& ch_traits){
if (l.size()>3){
typedef typename Traits::Left_turn_2 Left_turn;
Left_turn left_turn = ch_traits.left_turn_2_object();
typename std::list<P>::iterator pmin=l.begin();
for (typename std::list<P>::iterator it=l.begin(); it!=l.end(); ++it) {
if ((*pmin).x()>(*it).x())
{
pmin=it;
}
typedef typename Traits::Left_turn_2 Left_turn;
Left_turn left_turn = ch_traits.left_turn_2_object();
typename std::list<P>::iterator pmin=l.begin();
for (typename std::list<P>::iterator it=l.begin(); it!=l.end(); ++it) {
if ((*pmin).x()>(*it).x())
{
pmin=it;
}
}//*pmin is the extremal point on the left
typename std::list<P>::iterator u=pmin;
typename std::list<P>::iterator u_next=u;
@ -123,7 +123,7 @@ namespace CGAL{
Cyclic_increment_iterator(u_next_next,l);
while (u_next !=pmin){
if (left_turn(*u,*u_next,*u_next_next)){
Cyclic_increment_iterator(u,l);
Cyclic_increment_iterator(u_next,l);
@ -143,32 +143,32 @@ namespace CGAL{
}
}
}
}
}
} //namespace CGAL::internal
}
//////////////////////////////////////////////////////////////////////////////
template<class P,class GEN>
void random_convex_hull_in_disc_2(size_t n, typename Kernel_traits<P>::Kernel::FT radius, std::list<P> & l,GEN & gen, bool fast=true ){
CGAL_precondition( n >= 3);
typedef typename Kernel_traits<P>::Kernel K;
typedef typename Kernel_traits<P>::Kernel::FT FT;
size_t simulated_points=0;
size_t generated_points=0;
do
void random_convex_hull_in_disc_2(std::size_t n, double radius, std::list<P> & l,GEN & gen, bool fast=true ){
CGAL_precondition( n >= 3);
typedef typename Kernel_traits<P>::Kernel K;
std::size_t simulated_points=0;
std::size_t generated_points=0;
do
{ //Initialisation
size_t init=std::min( (size_t)100,n-simulated_points );
internal::generate_points_annulus(init,-CGAL_PI, CGAL_PI,0,to_double(radius),l,gen);
//std::size_t init=std::min( (std::size_t)100,n-simulated_points );
std::size_t init=std::min( static_cast<std::size_t>(100), n-simulated_points );
generate_points_annulus(init,-CGAL_PI, CGAL_PI,0,to_double(radius),l,gen);
simulated_points+=init;
generated_points+=init;
internal::Graham_without_sort_2(l,K());
Graham_without_sort_2(l,K());
} while ((bounded_side_2(l.begin(),l.end(),P (0,0),K())!=ON_BOUNDED_SIDE)&&(simulated_points<n)); //initialisation such that 0 in P_n
size_t T=n;
if (!fast) T=(size_t)std::floor(n/std::pow(log(n),2));
std::size_t T=n;
//if (!fast) T=(size_t)std::floor(n/std::pow(log(n),2));
if (!fast) T=static_cast<std::size_t>( std::floor( n/std::pow(log(n),2) ) );
while (simulated_points<n)
{
//l is a list coming from a convex hull operation. we are moving the points s.t the angles are from -pi to pi.
@ -189,27 +189,27 @@ namespace CGAL{
}
}
FT squared_radius=radius*radius;
FT squared_small_radius=squared_radius;
double squared_radius=radius*radius;
double squared_small_radius=squared_radius;
{
P zero(0,0);
typename std::list<P>::iterator it=l.begin();
typename std::list<P>::iterator it2=++it;
for(;it!=l.end();++it,internal::Cyclic_increment_iterator(it2,l)){ //computation of annulus
for(;it!=l.end();++it,Cyclic_increment_iterator(it2,l)){ //computation of annulus
typename K::Segment_2 s(*it,*it2);
FT temp=squared_distance(s,zero);
double temp=to_double(squared_distance(s,zero));
if (squared_small_radius>temp) squared_small_radius=temp;
}
}//squared_small_radius=squared small radius of the annulus
FT p_disc=squared_small_radius/squared_radius;
size_t nb;
double p_disc=squared_small_radius/squared_radius;
std::size_t nb;
if (simulated_points< T){nb=std::min(simulated_points,n-simulated_points);}
else {nb=std::min(T,n-simulated_points); }
boost::random::binomial_distribution<long> dbin(nb,to_double(p_disc));
boost::random::binomial_distribution<long> dbin(nb,p_disc);
boost::random::variate_generator<GEN&, boost::random::binomial_distribution<long> >bin(gen,dbin);
//How many points are falling in the small disc and wont be generated:
@ -217,18 +217,32 @@ namespace CGAL{
simulated_points+=k_disc;
std::list<P> m;
internal::generate_points_annulus(nb-k_disc,-CGAL_PI, CGAL_PI,std::sqrt(to_double(squared_small_radius)),to_double(radius),m,gen);
internal::generate_points_annulus(nb-k_disc,-CGAL_PI, CGAL_PI,std::sqrt(squared_small_radius),radius,m,gen);
l.merge(m,internal::compare_points_angle<P>());
generated_points+=nb-k_disc;
simulated_points+=nb-k_disc;
m.clear();
internal::Graham_without_sort_2(l,K());
Graham_without_sort_2(l,K());
}
}
} //namespace CGAL::internal
///
template <class OutputIterator, class Traits, class Generator>
void random_convex_hull_in_disc_2(std::size_t n, double radius, Generator & gen, OutputIterator it, const Traits & traits, bool fast=true)
{
typedef Point_2<Traits> Points;
std::list<Points> l;
internal::random_convex_hull_in_disc_2(n,radius, l,gen,fast);
// for(typename std::list<Points>::iterator i=l.begin();i!=l.end();++i)
// {
// *it=*i;
// ++*it;
// }
std::copy(l.begin(),l.end(), it);
}
}//namespace CGAL
#endif

View File

@ -269,7 +269,7 @@ MainWindow::on_actionGeneratePolytopeInDisc_triggered()
typedef CGAL::Points_on_segment_2<Point_2> PG;
boost::random::mt19937 gen;
gen.seed(time(0));
std::list<Point_2> list_of_points;
std::vector<Point_2> points;
QRectF rect = CGAL::Qt::viewportsBbox(&scene);
CGAL::Qt::Converter<K> convert;
Iso_rectangle_2 isor = convert(rect);
@ -300,11 +300,11 @@ MainWindow::on_actionGeneratePolytopeInDisc_triggered()
segments.reserve(segments.size() + 100);
CGAL::random_convex_hull_in_disc_2(number_of_points,100,list_of_points,gen);
std::list<Point_2>::iterator it2=list_of_points.begin();
for(std::list<Point_2>::iterator it=list_of_points.begin();it!=list_of_points.end();it++){
CGAL::random_convex_hull_in_disc_2(number_of_points,100,gen,std::back_inserter(points),K());
std::vector<Point_2>::iterator it2=points.begin();
for(std::vector<Point_2>::iterator it=points.begin();it!=points.end();it++){
it2++;
if (it2==list_of_points.end()) it2=list_of_points.begin();
if (it2==points.end()) it2=points.begin();
Segment_2 p(*it+offset,*it2+offset);
segments.push_back(p);
}