improving hierarchy, adding tests

This commit is contained in:
Manuel Caroli 2009-05-11 09:29:03 +00:00
parent 9c54b9f5da
commit 02eae34ffe
3 changed files with 121 additions and 33 deletions

1
.gitattributes vendored
View File

@ -2779,6 +2779,7 @@ Periodic_3_triangulation_3/test/Periodic_3_triangulation_3/include/CGAL/_test_cl
Periodic_3_triangulation_3/test/Periodic_3_triangulation_3/include/CGAL/_test_cls_periodic_3_offset_3.h -text
Periodic_3_triangulation_3/test/Periodic_3_triangulation_3/include/CGAL/_test_cls_periodic_3_triangulation_traits_3.h -text
Periodic_3_triangulation_3/test/Periodic_3_triangulation_3/include/CGAL/_test_periodic_3_static_filters.h -text
Periodic_3_triangulation_3/test/Periodic_3_triangulation_3/test_periodic_3_delaunay_hierarchy_3.cpp -text
Periodic_3_triangulation_3/test/Periodic_3_triangulation_3/test_periodic_3_offset_3.cpp -text
Periodic_3_triangulation_3/test/Periodic_3_triangulation_3/test_periodic_3_triangulation_traits_3.cpp -text
Point_set_2/doc_tex/Point_set_2/point_set.png -text

View File

@ -54,12 +54,13 @@ public:
using PTr_Base::number_of_vertices;
using PTr_Base::geom_traits;
using PTr_Baes::is_virtual;
using PTr_Base::is_virtual;
private:
// here is the stack of triangulations which form the hierarchy
PTr_Base* hierarchy[maxlevel];
Random random; // random generator
int level_mult_cover;
public:
Periodic_3_triangulation_hierarchy_3(
@ -73,7 +74,7 @@ public:
Periodic_3_triangulation_hierarchy_3(InputIterator first, InputIterator last,
const Iso_cuboid& domain = Iso_cuboid(0,0,0,1,1,1),
const Geom_traits& traits = Geom_traits())
: PTr_Base(domain,traits), random((long)0)
: PTr_Base(domain,traits), random((long)0), level_mult_cover(0)
{
hierarchy[0] = this;
for(int i=1; i<maxlevel; ++i)
@ -105,37 +106,36 @@ public:
template < class InputIterator >
int insert(InputIterator first, InputIterator last, bool = false)
{
int n = number_of_vertices();
int n = number_of_vertices();
std::vector<Point> points (first, last);
std::random_shuffle (points.begin(), points.end());
spatial_sort (points.begin(), points.end(), geom_traits());
std::vector<Point> points (first, last);
std::random_shuffle (points.begin(), points.end());
spatial_sort (points.begin(), points.end(), geom_traits());
// hints[i] is the cell of the previously inserted point in level i.
// Thanks to spatial sort, they are better hints than what the hierarchy
// would give us.
Cell_handle hints[maxlevel];
for (typename std::vector<Point>::const_iterator p = points.begin(), end = points.end();
p != end; ++p)
{
int vertex_level = random_level();
Vertex_handle v = hierarchy[0]->insert (*p, hints[0]);
hints[0] = v->cell();
Vertex_handle prev = v;
for (int level = 1; level <= vertex_level; ++level) {
v = hierarchy[level]->insert (*p, hints[level]);
hints[level] = v->cell();
v->set_down (prev);
prev->set_up (v);
prev = v;
}
// hints[i] is the cell of the previously inserted point in level i.
// Thanks to spatial sort, they are better hints than what the hierarchy
// would give us.
Cell_handle hints[maxlevel];
for (typename std::vector<Point>::const_iterator p = points.begin(),
end = points.end(); p != end; ++p) {
int vertex_level = random_level();
Vertex_handle v = hierarchy[0]->insert (*p, hints[0]);
hints[0] = v->cell();
Vertex_handle prev = v;
for (int level = 1; level <= vertex_level; ++level) {
v = hierarchy[level]->insert (*p, hints[level]);
hints[level] = v->cell();
v->set_down (prev);
prev->set_up (v);
prev = v;
}
return number_of_vertices() - n;
}
return number_of_vertices() - n;
}
// bool only for backward compatibility, we document void.
@ -185,7 +185,7 @@ template <class PTr >
Periodic_3_triangulation_hierarchy_3<PTr>::
Periodic_3_triangulation_hierarchy_3(
const Iso_cuboid& domain, const Geom_traits& traits)
: PTr_Base(domain, traits), random((long)0)
: PTr_Base(domain, traits), random((long)0), level_mult_cover(0)
{
hierarchy[0] = this;
for(int i=1;i<maxlevel;++i)
@ -197,7 +197,7 @@ template <class PTr>
Periodic_3_triangulation_hierarchy_3<PTr>::
Periodic_3_triangulation_hierarchy_3(
const Periodic_3_triangulation_hierarchy_3<PTr> &tr)
: PTr_Base(tr), random((long)0)
: PTr_Base(tr), random((long)0), level_mult_cover(tr.level_mult_cover)
{
hierarchy[0] = this;
for(int i=1; i<maxlevel; ++i)
@ -491,8 +491,12 @@ int
Periodic_3_triangulation_hierarchy_3<PTr>::
random_level()
{
if ( level_mult_cover < maxlevel
&& hierarchy[level_mult_cover]->number_of_sheets() == make_array(1,1,1) )
++level_mult_cover;
int l = 0;
while ( ! random(ratio) && l < maxlevel-1 )
while ( ! random(ratio) && l < level_mult_cover-1)
++l;
return l;

View File

@ -0,0 +1,83 @@
// Copyright (c) 1998 INRIA Sophia-Antipolis (France).
// All rights reserved.
//
// This file is part of CGAL (www.cgal.org); you may redistribute it under
// the terms of the Q Public License version 1.0.
// See the file LICENSE.QPL distributed with CGAL.
//
// Licensees holding a valid commercial license may use this file in
// accordance with the commercial license agreement provided with the software.
//
// This file is provided AS IS with NO WARRANTY OF ANY KIND, INCLUDING THE
// WARRANTY OF DESIGN, MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE.
//
// $URL: svn+ssh://mcaroli@scm.gforge.inria.fr/svn/cgal/trunk/Periodic_3_triangulation_3/test/Periodic_3_triangulation_3/test_periodic_3_delaunay_3.cpp $
// $Id: test_periodic_3_delaunay_3.cpp 48874 2009-04-23 13:54:38Z mcaroli $
//
//
// Author(s) : Nico Kruithof
// Manuel Caroli
#include <iostream>
#include <fstream>
#include <CGAL/Periodic_3_Delaunay_triangulation_3.h>
#include <CGAL/Periodic_3_triangulation_hierarchy_3.h>
#include <CGAL/Periodic_3_triangulation_traits_3.h>
#include <CGAL/Exact_predicates_inexact_constructions_kernel.h>
typedef CGAL::Exact_predicates_inexact_constructions_kernel K1;
typedef CGAL::Periodic_3_triangulation_traits_3<K1> PTT1;
typedef CGAL::Periodic_3_triangulation_ds_vertex_base_3<> DSVB1;
typedef CGAL::Periodic_3_triangulation_ds_cell_base_3<> DSCB1;
typedef CGAL::Triangulation_vertex_base_3<PTT1,DSVB1> VBB1;
typedef CGAL::Periodic_3_triangulation_hierarchy_vertex_base_3<VBB1> VB1;
typedef CGAL::Triangulation_cell_base_3<PTT1,DSCB1> CB1;
typedef CGAL::Triangulation_data_structure_3<VB1,CB1> TDS1;
typedef CGAL::Periodic_3_Delaunay_triangulation_3<PTT1,TDS1> PDT1;
// Explicit instantiation of the whole class :
template class CGAL::Periodic_3_triangulation_hierarchy_3<PDT1>;
#include <CGAL/Exact_predicates_exact_constructions_kernel.h>
typedef CGAL::Exact_predicates_exact_constructions_kernel K2;
typedef CGAL::Periodic_3_triangulation_traits_3<K2> PTT2;
typedef CGAL::Periodic_3_triangulation_ds_vertex_base_3<> DSVB2;
typedef CGAL::Periodic_3_triangulation_ds_cell_base_3<> DSCB2;
typedef CGAL::Triangulation_vertex_base_3<PTT2,DSVB2> VBB2;
typedef CGAL::Periodic_3_triangulation_hierarchy_vertex_base_3<VBB2> VB2;
typedef CGAL::Triangulation_cell_base_3<PTT2,DSCB2> CB2;
typedef CGAL::Triangulation_data_structure_3<VB2,CB2> TDS2;
typedef CGAL::Periodic_3_Delaunay_triangulation_3<PTT2,TDS2> PDT2;
// Explicit instantiation of the whole class :
template class CGAL::Periodic_3_triangulation_hierarchy_3<PDT2>;
#include <CGAL/MP_Float.h>
#include <CGAL/Simple_homogeneous.h>
typedef CGAL::Simple_homogeneous<CGAL::MP_Float> K3;
typedef CGAL::Periodic_3_triangulation_traits_3<K3> PTT3;
typedef CGAL::Periodic_3_triangulation_ds_vertex_base_3<> DSVB3;
typedef CGAL::Periodic_3_triangulation_ds_cell_base_3<> DSCB3;
typedef CGAL::Triangulation_vertex_base_3<PTT3,DSVB3> VBB3;
typedef CGAL::Periodic_3_triangulation_hierarchy_vertex_base_3<VBB3> VB3;
typedef CGAL::Triangulation_cell_base_3<PTT3,DSCB3> CB3;
typedef CGAL::Triangulation_data_structure_3<VB3,CB3> TDS3;
typedef CGAL::Periodic_3_Delaunay_triangulation_3<PTT3,TDS3> PDT3;
// Explicit instantiation of the whole class :
template class CGAL::Periodic_3_triangulation_hierarchy_3<PDT3>;
#include <CGAL/_test_cls_periodic_3_delaunay_3.h>
int main()
{
typedef CGAL::Periodic_3_triangulation_hierarchy_3< PDT1 > P3T3_1;
_test_cls_periodic_3_delaunay_3( P3T3_1() );
typedef CGAL::Periodic_3_triangulation_hierarchy_3< PDT2 > P3T3_2;
// this takes too much time for the test suite.
//_test_cls_periodic_3_delaunay_3( P3T3_2() );
typedef CGAL::Periodic_3_triangulation_hierarchy_3< PDT3 > P3T3_3;
// this takes too much time for the test suite.
//_test_cls_periodic_3_delaunay_3( P3T3_3(), true );
return 0;
}