Fixing the Skin surface package for kernels other than the EPIC-kernel.

This commit is contained in:
Nico Kruithof 2010-11-30 20:34:12 +00:00
parent d3a6f51f6d
commit 8364150a57
6 changed files with 48 additions and 41 deletions

1
.gitattributes vendored
View File

@ -3255,6 +3255,7 @@ Skin_surface_3/test/Skin_surface_3/data/test6.cin -text
Skin_surface_3/test/Skin_surface_3/data/test7.cin -text
Skin_surface_3/test/Skin_surface_3/data/test8.cin -text
Skin_surface_3/test/Skin_surface_3/data/test9.cin -text
Skin_surface_3/test/Skin_surface_3/degenerate_test_exact.cpp -text
Skin_surface_3/test/Skin_surface_3/union_of_balls_test.cpp -text
Snap_rounding_2/demo/Snap_rounding_2/Qt3/help/index.html svneol=native#text/html
Snap_rounding_2/doc_tex/Snap_rounding_2/isr_vs_sr.gif -text svneol=unset#image/gif

View File

@ -91,7 +91,7 @@ public:
private:
// Triangulated_mixed_complex:
typedef Simple_cartesian<Interval_nt_advanced> FK;
typedef CGAL::Exact_predicates_exact_constructions_kernel FK;
typedef Triangulation_vertex_base_with_info_3<Vertex_info, FK> Vb;
typedef Triangulation_cell_base_with_info_3<Cell_info, FK> Cb;
typedef Triangulation_data_structure_3<Vb,Cb> Tds;
@ -364,10 +364,7 @@ void
Skin_surface_base_3<MixedComplexTraits_3>::
intersect(TMC_Cell_handle ch, int i, int j,
Bare_point &p) const {
typedef typename Bare_point::R Traits;
typedef typename Traits::FT FT;
Cartesian_converter<FK,
typename Geometric_traits::Bare_point::R> converter;
Cartesian_converter<FK, Gt> converter;
Bare_point p1 = converter(ch->vertex(i)->point());
Bare_point p2 = converter(ch->vertex(j)->point());
@ -523,7 +520,7 @@ construct_bounding_box()
Bare_point mid(bbox.xmin() + dx/2,
bbox.ymin() + dy/2,
bbox.zmin() + dz/2);
double dr =
FT dr =
(dx+dy+dz+sqrt(CGAL::to_double(max_weight))+.001) / gt.get_shrink();
Weighted_point wp;

View File

@ -95,7 +95,9 @@ public:
case 0: {
vh = s;
Surface_weighted_point wp = r2s_converter(vh->point());
create_sphere(wp.point(), -wp.weight(), shrink, 1);
create_sphere(wp.point(),
-wp.weight(),
r2s_converter(shrink), 1);
break;
}
case 1: {
@ -105,14 +107,13 @@ public:
Surface_weighted_point p1 =
r2s_converter(e.first->vertex(e.third)->point());
create_hyperboloid
(typename Surface_regular_traits::
Construct_weighted_circumcenter_3()(p0,p1),
typename Surface_regular_traits::
Compute_squared_radius_smallest_orthogonal_sphere_3()(p0,p1),
p0 - p1,
shrink,
1);
create_hyperboloid(typename Surface_regular_traits::
Construct_weighted_circumcenter_3()(p0,p1),
typename Surface_regular_traits::
Compute_squared_radius_smallest_orthogonal_sphere_3()(p0,p1),
p0 - p1,
r2s_converter(shrink),
1);
break;
}
case 2: {
@ -124,35 +125,34 @@ public:
Surface_weighted_point p2 =
r2s_converter(f.first->vertex((f.second+3)&3)->point());
create_hyperboloid
(typename Surface_regular_traits::
Construct_weighted_circumcenter_3()(p0,p1,p2),
typename Surface_regular_traits::
Compute_squared_radius_smallest_orthogonal_sphere_3()(p0,p1,p2),
typename Surface_regular_traits::
Construct_orthogonal_vector_3()(p0,p1,p2),
1-shrink,
-1);
create_hyperboloid(typename Surface_regular_traits::
Construct_weighted_circumcenter_3()(p0,p1,p2),
typename Surface_regular_traits::
Compute_squared_radius_smallest_orthogonal_sphere_3()(p0,p1,p2),
typename Surface_regular_traits::
Construct_orthogonal_vector_3()(p0,p1,p2),
r2s_converter(1-shrink),
-1);
break;
}
case 3: {
ch = s;
const Surface_weighted_point pts[4] = {
r2s_converter(ch->vertex(0)->point()),
r2s_converter(ch->vertex(1)->point()),
r2s_converter(ch->vertex(2)->point()),
r2s_converter(ch->vertex(3)->point())
};
const Surface_weighted_point pts[4] =
{
r2s_converter(ch->vertex(0)->point()),
r2s_converter(ch->vertex(1)->point()),
r2s_converter(ch->vertex(2)->point()),
r2s_converter(ch->vertex(3)->point())
};
create_sphere
(typename Surface_regular_traits::
Construct_weighted_circumcenter_3()
(pts[0],pts[1],pts[2],pts[3]),
typename Surface_regular_traits::
Compute_squared_radius_smallest_orthogonal_sphere_3()
(pts[0],pts[1],pts[2],pts[3]),
1-shrink,
-1);
create_sphere(typename Surface_regular_traits::
Construct_weighted_circumcenter_3()
(pts[0],pts[1],pts[2],pts[3]),
typename Surface_regular_traits::
Compute_squared_radius_smallest_orthogonal_sphere_3()
(pts[0],pts[1],pts[2],pts[3]),
r2s_converter(1-shrink),
-1);
}
}
}

View File

@ -132,7 +132,7 @@ public:
_tmc(triangulated_mixed_complex),
observer(observer),
triangulation_incr_builder(triangulated_mixed_complex),
construct_anchor_point_3_obj(shrink),
construct_anchor_point_3_obj(r2t_converter_object(shrink)),
compute_anchor_obj(regular),
verbose(verbose) {

View File

@ -1,5 +1,6 @@
// examples/Skin_surface_3/skin_surface_simple.C
#include <CGAL/Exact_predicates_inexact_constructions_kernel.h>
#include <CGAL/Exact_predicates_exact_constructions_kernel.h>
#include <CGAL/Skin_surface_3.h>
#include <CGAL/Polyhedron_3.h>
#include <CGAL/mesh_skin_surface_3.h>
@ -9,13 +10,18 @@
#include <fstream>
#include <algorithm>
#ifdef CGAL_SKIN_SURFACE_USE_EXACT_CONSTRUCTION_KERNEL
typedef CGAL::Exact_predicates_exact_constructions_kernel K;
#else
typedef CGAL::Exact_predicates_inexact_constructions_kernel K;
#endif
typedef CGAL::Skin_surface_traits_3<K> Traits;
typedef CGAL::Skin_surface_3<Traits> Skin_surface_3;
typedef Skin_surface_3::FT FT;
typedef Skin_surface_3::Weighted_point Weighted_point;
typedef Weighted_point::Point Bare_point;
typedef CGAL::Polyhedron_3<K> Polyhedron;
typedef CGAL::Exact_predicates_inexact_constructions_kernel IK;
typedef CGAL::Polyhedron_3<IK> Polyhedron;
class Test_file {
public:

View File

@ -0,0 +1,3 @@
#define CGAL_SKIN_SURFACE_USE_EXACT_CONSTRUCTION_KERNEL 1
#include "./degenerate_test.cpp"