Working on compilation errors

This commit is contained in:
Nico Kruithof 2006-11-24 13:56:23 +00:00
parent d931482ce2
commit 8ac0951040
4 changed files with 95 additions and 82 deletions

View File

@ -49,8 +49,8 @@ public:
HDS_point intersection(Cell_iterator const ch, int i, int j) const { HDS_point intersection(Cell_iterator const ch, int i, int j) const {
// Precondition: ch is not an infinite cell: their surface is not set // Precondition: ch is not an infinite cell: their surface is not set
Skin_point p; Skin_point p;
ss_3.intersect(ch->vertex(i), ch->vertex(j), ss_3.intersect(ch, i, j, //->vertex(i), ch->vertex(j), //ch->info(),
ch->info(), p); p);
return return
Cartesian_converter<typename Skin_point::R, typename HDS_point::R>()(p); Cartesian_converter<typename Skin_point::R, typename HDS_point::R>()(p);

View File

@ -96,7 +96,9 @@ private:
typedef Triangulation_cell_base_with_info_3< typedef Triangulation_cell_base_with_info_3<
std::pair<Simplex, Quadratic_surface *>, FK> Cb; std::pair<Simplex, Quadratic_surface *>, FK> Cb;
typedef Triangulation_data_structure_3<Vb,Cb> Tds; typedef Triangulation_data_structure_3<Vb,Cb> Tds;
public:
typedef Triangulation_3<FK, Tds> TMC; typedef Triangulation_3<FK, Tds> TMC;
private:
typedef typename TMC::Vertex_iterator TMC_Vertex_iterator; typedef typename TMC::Vertex_iterator TMC_Vertex_iterator;
typedef typename TMC::Cell_iterator TMC_Cell_iterator; typedef typename TMC::Cell_iterator TMC_Cell_iterator;
typedef typename TMC::Vertex_handle TMC_Vertex_handle; typedef typename TMC::Vertex_handle TMC_Vertex_handle;
@ -104,19 +106,19 @@ private:
public: public:
template < class WP_iterator > template < class WP_iterator >
Skin_surface_3(WP_iterator begin, WP_iterator end, Skin_surface_3(WP_iterator begin, WP_iterator end,
FT shrink_factor, FT shrink,
bool grow_balls = true, bool grow_balls = true,
Gt gt_ = Gt(), Gt gt_ = Gt(),
bool _verbose = false bool _verbose = false
) )
: gt(gt_), verbose(_verbose) { : gt(gt_), verbose(_verbose) {
gt.set_shrink(shrink_factor); gt.set_shrink(shrink);
CGAL_assertion(begin != end); CGAL_assertion(begin != end);
if (grow_balls) { if (grow_balls) {
for (; begin != end; begin++) { for (; begin != end; begin++) {
regular.insert(Weighted_point(*begin, begin->weight()/gt.get_shrink())); regular.insert(Weighted_point(*begin, begin->weight()/shrink_factor()));
} }
} else { } else {
regular.insert(begin, end); regular.insert(begin, end);
@ -131,8 +133,8 @@ public:
// Construct the Triangulated_mixed_complex: // Construct the Triangulated_mixed_complex:
Triangulated_mixed_complex_observer_3<TMC, Self> observer(shrink_factor); Triangulated_mixed_complex_observer_3<TMC, Self> observer(shrink_factor());
triangulate_mixed_complex_3(regular, shrink_factor, tmc, observer, true); triangulate_mixed_complex_3(regular, shrink_factor(), tmc, observer, true);
{ // NGHK: debug code: { // NGHK: debug code:
CGAL_assertion(tmc.is_valid()); CGAL_assertion(tmc.is_valid());
@ -186,7 +188,7 @@ public:
std::string(CGAL_PRETTY_FUNCTION)); std::string(CGAL_PRETTY_FUNCTION));
Protect_FPU_rounding<false> P(CGAL_FE_TONEAREST); Protect_FPU_rounding<false> P(CGAL_FE_TONEAREST);
typedef Exact_predicates_exact_constructions_kernel EK; typedef Exact_predicates_exact_constructions_kernel EK;
Skin_surface_traits_3<EK> exact_traits(gt.get_shrink()); Skin_surface_traits_3<EK> exact_traits(shrink_factor());
typename Skin_surface_traits_3<EK>::Bare_point p_exact = typename Skin_surface_traits_3<EK>::Bare_point p_exact =
get_anchor_point(vit->info(), exact_traits); get_anchor_point(vit->info(), exact_traits);
@ -291,7 +293,7 @@ private:
} }
// CMCT_Cell locate_tet(const Bare_point &p, const Simplex &s) const { // CMCT_Cell locate_tet(const Bare_point &p, const Simplex &s) const {
// Skin_surface_traits_3<Exact_predicates_exact_constructions_kernel> // Skin_surface_traits_3<Exact_predicates_exact_constructions_kernel>
// traits(gt.get_shrink()); // traits(shrink_factor());
// std::vector<CMCT_Cell> cells; // std::vector<CMCT_Cell> cells;
@ -350,7 +352,7 @@ public:
// exact computation of the sign on a vertex of the TMC // exact computation of the sign on a vertex of the TMC
// Sign sign(const CMCT_Vertex_handle vh) const { // Sign sign(const CMCT_Vertex_handle vh) const {
// typedef Exact_predicates_exact_constructions_kernel K; // typedef Exact_predicates_exact_constructions_kernel K;
// Skin_surface_traits_3<K> traits(gt.get_shrink()); // Skin_surface_traits_3<K> traits(shrink_factor());
// typename K::Point_3 p = mc_triangulator->location(vh, traits); // typename K::Point_3 p = mc_triangulator->location(vh, traits);
@ -476,15 +478,18 @@ public:
p = midpoint(p1, p2); p = midpoint(p1, p2);
} }
void intersect(const TMC_Vertex_handle &p1, void intersect(TMC_Cell_handle ch, int i, int j,
const TMC_Vertex_handle &p2, //TMC_Vertex_handle p2,
Quadratic_surface *surf,
Bare_point &p) const { Bare_point &p) const {
typedef typename Bare_point::R Traits; typedef typename Bare_point::R Traits;
typedef typename Traits::FT FT; typedef typename Traits::FT FT;
Cartesian_converter<Traits, Cartesian_converter<FK,
typename Geometric_traits::Bare_point::R> converter; typename Geometric_traits::Bare_point::R> converter;
Quadratic_surface *surf = ch->info().second;
Bare_point p1 = converter(ch->vertex(i)->point());
Bare_point p2 = converter(ch->vertex(j)->point());
FT sq_dist = squared_distance(p1,p2); FT sq_dist = squared_distance(p1,p2);
// Use value to make the computation robust (endpoints near the surface) // Use value to make the computation robust (endpoints near the surface)
if (value(surf, p1) > value(surf, p2)) std::swap(p1, p2); if (value(surf, p1) > value(surf, p2)) std::swap(p1, p2);
@ -592,7 +597,7 @@ public:
case 0: case 0:
{ {
Vertex_handle vh = sim; Vertex_handle vh = sim;
return Quadratic_surface(conv(vh->point()), gt.get_shrink()); return Quadratic_surface(conv(vh->point()), shrink_factor());
break; break;
} }
case 1: case 1:
@ -600,7 +605,7 @@ public:
Edge e = sim; Edge e = sim;
Weighted_point p0 = conv(e.first->vertex(e.second)->point()); Weighted_point p0 = conv(e.first->vertex(e.second)->point());
Weighted_point p1 = conv(e.first->vertex(e.third)->point()); Weighted_point p1 = conv(e.first->vertex(e.third)->point());
return Quadratic_surface(p0, p1, gt.get_shrink()); return Quadratic_surface(p0, p1, shrink_factor());
break; break;
} }
case 2: case 2:
@ -609,7 +614,7 @@ public:
Weighted_point p0 = conv(f.first->vertex((f.second+1)&3)->point()); Weighted_point p0 = conv(f.first->vertex((f.second+1)&3)->point());
Weighted_point p1 = conv(f.first->vertex((f.second+2)&3)->point()); Weighted_point p1 = conv(f.first->vertex((f.second+2)&3)->point());
Weighted_point p2 = conv(f.first->vertex((f.second+3)&3)->point()); Weighted_point p2 = conv(f.first->vertex((f.second+3)&3)->point());
return Quadratic_surface(p0,p1,p2, gt.get_shrink()); return Quadratic_surface(p0,p1,p2, shrink_factor());
break; break;
} }
case 3: case 3:
@ -619,7 +624,7 @@ public:
Weighted_point p1 = conv(ch->vertex(1)->point()); Weighted_point p1 = conv(ch->vertex(1)->point());
Weighted_point p2 = conv(ch->vertex(2)->point()); Weighted_point p2 = conv(ch->vertex(2)->point());
Weighted_point p3 = conv(ch->vertex(3)->point()); Weighted_point p3 = conv(ch->vertex(3)->point());
return Quadratic_surface(p0,p1,p2,p3, gt.get_shrink()); return Quadratic_surface(p0,p1,p2,p3, shrink_factor());
break; break;
} }
} }
@ -662,7 +667,7 @@ public:
const Regular &get_regular_triangulation() const { const Regular &get_regular_triangulation() const {
return regular; return regular;
} }
FT get_shrink_factor() const { FT shrink_factor() const {
return gt.get_shrink(); return gt.get_shrink();
} }

View File

@ -107,9 +107,8 @@ private:
typedef typename Anchor_map::iterator Anchor_map_iterator; typedef typename Anchor_map::iterator Anchor_map_iterator;
public: public:
Mixed_complex_triangulator_3( Mixed_complex_triangulator_3(Regular const &regular,
Regular const &regular, Triangulated_mixed_complex const &triangulated_mixed_complex,
Triangulated_mixed_complex &triangulated_mixed_complex,
bool verbose) bool verbose)
: regular(regular), : regular(regular),
_tmc(triangulated_mixed_complex), _tmc(triangulated_mixed_complex),
@ -120,8 +119,7 @@ public:
build(); build();
} }
Mixed_complex_triangulator_3( Mixed_complex_triangulator_3(Regular const&regular,
Regular &regular,
Triangulated_mixed_complex &triangulated_mixed_complex, Triangulated_mixed_complex &triangulated_mixed_complex,
Triangulated_mixed_complex_observer &observer, Triangulated_mixed_complex_observer &observer,
bool verbose) bool verbose)
@ -777,7 +775,7 @@ template <
class TriangulatedMixedComplexObserver_3> class TriangulatedMixedComplexObserver_3>
void void
triangulate_mixed_complex_3 triangulate_mixed_complex_3
(RegularTriangulation_3 &rt, (RegularTriangulation_3 const &rt,
typename RegularTriangulation_3::Geom_traits::FT shrink, typename RegularTriangulation_3::Geom_traits::FT shrink,
TriangulatedMixedComplex_3 &tmc, TriangulatedMixedComplex_3 &tmc,
TriangulatedMixedComplexObserver_3 &observer, TriangulatedMixedComplexObserver_3 &observer,
@ -791,21 +789,21 @@ triangulate_mixed_complex_3
} }
template < // template <
class RegularTriangulation_3, // class RegularTriangulation_3,
class TriangulatedMixedComplex_3> // class TriangulatedMixedComplex_3>
void // void
triangulate_mixed_complex_3 // triangulate_mixed_complex_3
(RegularTriangulation_3 const &regular, // (RegularTriangulation_3 const &regular,
typename RegularTriangulation_3::Geom_Traits::FT shrink, // typename RegularTriangulation_3::Geom_Traits::FT shrink,
TriangulatedMixedComplex_3 &tmc, // TriangulatedMixedComplex_3 &tmc,
bool verbose) // bool verbose)
{ // {
Triangulated_mixed_complex_observer_3< // Triangulated_mixed_complex_observer_3<
TriangulatedMixedComplex_3, const RegularTriangulation_3> // TriangulatedMixedComplex_3, const RegularTriangulation_3>
observer(1); // observer(1);
triangulate_mixed_complex_3(regular, shrink, tmc, observer, verbose); // triangulate_mixed_complex_3(regular, shrink, tmc, observer, verbose);
} // }
CGAL_END_NAMESPACE CGAL_END_NAMESPACE

View File

@ -4,16 +4,17 @@
#include <CGAL/Skin_surface_3.h> #include <CGAL/Skin_surface_3.h>
#include <CGAL/Polyhedron_3.h> #include <CGAL/Polyhedron_3.h>
#include <CGAL/mesh_skin_surface_3.h> #include <CGAL/mesh_skin_surface_3.h>
#include <CGAL/triangulate_mixed_complex_3.h>
#include <list> #include <list>
#include <fstream> #include <fstream>
typedef CGAL::Exact_predicates_inexact_constructions_kernel Inexact_K; typedef CGAL::Exact_predicates_inexact_constructions_kernel Inexact_K;
typedef CGAL::Skin_surface_traits_3<Inexact_K> Traits; typedef CGAL::Skin_surface_traits_3<Inexact_K> Traits;
typedef CGAL::Skin_surface_3<Traits> Skin_surface_3; typedef CGAL::Skin_surface_3<Traits> Skin_surface;
typedef Skin_surface_3::Simplex Simplex; typedef Skin_surface::Simplex Simplex;
typedef Skin_surface_3::FT FT; typedef Skin_surface::FT FT;
typedef Skin_surface_3::Weighted_point Weighted_point; typedef Skin_surface::Weighted_point Weighted_point;
typedef Weighted_point::Point Bare_point; typedef Weighted_point::Point Bare_point;
typedef CGAL::Polyhedron_3<Traits> Polyhedron; typedef CGAL::Polyhedron_3<Traits> Polyhedron;
typedef CGAL::Exact_predicates_exact_constructions_kernel Exact_K; typedef CGAL::Exact_predicates_exact_constructions_kernel Exact_K;
@ -22,16 +23,15 @@ typedef CGAL::Skin_surface_quadratic_surface_3<Exact_K> Quadratic_surface;
CGAL::Cartesian_converter<Exact_K,Inexact_K> e2i_converter; CGAL::Cartesian_converter<Exact_K,Inexact_K> e2i_converter;
CGAL::Cartesian_converter<Inexact_K,Exact_K> i2e_converter; CGAL::Cartesian_converter<Inexact_K,Exact_K> i2e_converter;
typedef CGAL::Triangulation_vertex_base_with_info_3<Simplex, Exact_K> Vb; // Triangulated_mixed_complex:
typedef CGAL::Triangulation_cell_base_with_info_3<Quadratic_surface *, Exact_K> Cb; typedef Skin_surface::TMC TMC;
typedef CGAL::Triangulation_data_structure_3<Vb,Cb> Tds;
typedef CGAL::Triangulation_3<Exact_K, Tds> Triangulated_mixed_complex; typedef TMC::Vertex_iterator TMC_Vertex_iterator;
//typedef Skin_surface_3::Triangulated_mixed_complex Triangulated_mixed_complex; typedef TMC::Cell_iterator TMC_Cell_iterator;
typedef Triangulated_mixed_complex::Vertex_handle Tmc_Vertex_handle; typedef TMC::Vertex_handle TMC_Vertex_handle;
typedef Triangulated_mixed_complex::Finite_vertices_iterator typedef TMC::Cell_handle TMC_Cell_handle;
Tmc_Finite_vertices_iterator; typedef TMC::Finite_vertices_iterator TMC_Finite_vertices_iterator;
typedef Triangulated_mixed_complex::Finite_cells_iterator typedef TMC::Finite_cells_iterator TMC_Finite_cells_iterator;
Tmc_Finite_cells_iterator;
#include <fstream> #include <fstream>
@ -49,49 +49,59 @@ public:
while (in >> x >> y >> z >> w) while (in >> x >> y >> z >> w)
l.push_front(Weighted_point(Bare_point(x,y,z),w)); l.push_front(Weighted_point(Bare_point(x,y,z),w));
Skin_surface_3 skin_surface(l.begin(), l.end(), s); Skin_surface skin_surface(l.begin(), l.end(), s);
Triangulated_mixed_complex tmc; TMC tmc;
CGAL::Triangulated_mixed_complex_observer_3<TMC, Skin_surface>
observer(skin_surface.shrink_factor());
triangulate_mixed_complex_3(skin_surface.get_regular_triangulation(), triangulate_mixed_complex_3(skin_surface.get_regular_triangulation(),
skin_surface.get_shrink_factor(), skin_surface.shrink_factor(),
tmc, false); tmc,
observer,
false);
// CGAL::triangulate_mixed_complex_3(skin_surface.get_regular_triangulation(),
// skin_surface.shrink_factor(),
// tmc,
// false // verbose
// );
for (Tmc_Finite_vertices_iterator vit = tmc.finite_vertices_begin(); for (TMC_Finite_vertices_iterator vit = tmc.finite_vertices_begin();
vit != tmc.finite_vertices_end(); vit++) { vit != tmc.finite_vertices_end(); vit++) {
if (tmc.is_infinite(vit->cell())) { if (tmc.is_infinite(vit->cell())) {
std::cerr << "ERROR: is_infinite (main)" << std::endl; std::cerr << "ERROR: is_infinite (main)" << std::endl;
} }
Quadratic_surface::FT val = vit->cell()->info()->value(vit->point()); Quadratic_surface::FT
std::list<Triangulated_mixed_complex::Cell_handle> cells; val = vit->cell()->info().second->value(vit->point());
std::list<TMC_Cell_handle> cells;
tmc.incident_cells(vit, std::back_inserter(cells)); tmc.incident_cells(vit, std::back_inserter(cells));
for (std::list<Triangulated_mixed_complex::Cell_handle>::iterator cell = for (std::list<TMC_Cell_handle>::iterator cell =
cells.begin(); cells.begin();
cell != cells.end(); cell++) { cell != cells.end(); cell++) {
if (!tmc.is_infinite(*cell)) { if (!tmc.is_infinite(*cell)) {
Quadratic_surface::FT val2 = (*cell)->info()->value(vit->point()); Quadratic_surface::FT val2 = (*cell)->info().second->value(vit->point());
CGAL_assertion(val == val2); CGAL_assertion(val == val2);
} }
} }
} }
for (Tmc_Finite_cells_iterator cit = tmc.finite_cells_begin(); // for (TMC_Finite_cells_iterator cit = tmc.finite_cells_begin();
cit != tmc.finite_cells_end(); cit++) { // cit != tmc.finite_cells_end(); cit++) {
Bare_point baryc = // Bare_point baryc =
e2i_converter(cit->vertex(0)->point() + // e2i_converter(cit->vertex(0)->point() +
(cit->vertex(1)->point()-cit->vertex(0)->point())/4 + // (cit->vertex(1)->point()-cit->vertex(0)->point())/4 +
(cit->vertex(2)->point()-cit->vertex(0)->point())/4 + // (cit->vertex(2)->point()-cit->vertex(0)->point())/4 +
(cit->vertex(3)->point()-cit->vertex(0)->point())/4); // (cit->vertex(3)->point()-cit->vertex(0)->point())/4);
if (tmc.tetrahedron(cit).has_on_bounded_side(i2e_converter(baryc))) { // if (tmc.tetrahedron(cit).has_on_bounded_side(i2e_converter(baryc))) {
Quadratic_surface::FT val1 = cit->info()->value(i2e_converter(baryc)); // Quadratic_surface::FT val1 = cit->info()->value(i2e_converter(baryc));
Simplex s = skin_surface.locate_mixed(baryc); // Simplex s = skin_surface.locate_mixed(baryc);
Quadratic_surface::FT val2 = // Quadratic_surface::FT val2 =
skin_surface.construct_surface(s, Exact_K()).value(i2e_converter(baryc)); // skin_surface.construct_surface(s, Exact_K()).value(i2e_converter(baryc));
// std::cout << val1 << " " << val2 << " " << val2-val1 << std::endl; // // std::cout << val1 << " " << val2 << " " << val2-val1 << std::endl;
CGAL_assertion(val1==val2); // CGAL_assertion(val1==val2);
} else { // } else {
std::cout << "Barycenter on unbounded side, due to rounding errors\n"; // std::cout << "Barycenter on unbounded side, due to rounding errors\n";
} // }
} // }
} }
private: private: