Fix broken is_hyperbolic(face/edge) framework

The previous code tries to handle finite/infinite and hyperbolic/non-hyperbolic
with a single flag, which caused errors to get stuff like finite
but non-hyperbolic edges (see e.g. https://github.com/CGAL/cgal/issues/6869).

In addition, it used tds_data(), which is something that now also
exists in the triangulation_ds_face_base_2 class.

Hence, completely re-implement the hyperbolic stack / query.

 #	../../../../Installation/CHANGES.md.orig
This commit is contained in:
Mael Rouxel-Labbé 2022-09-26 16:08:43 +02:00
parent 145a817cc1
commit cf12f90cbf
2 changed files with 108 additions and 178 deletions

View File

@ -193,10 +193,10 @@ public:
do do
{ {
_ri = cw(_iv); _ri = cw(_iv);
if (_tri.is_finite_non_hyperbolic(pos, ccw(_iv))) if (!_tri.is_Delaunay_hyperbolic(pos, ccw(_iv)))
{ {
_ri = ccw(_iv); _ri = ccw(_iv);
if (_tri.is_finite_non_hyperbolic(pos, cw(_iv))) if (!_tri.is_Delaunay_hyperbolic(pos, cw(_iv)))
{ {
pos = pos->neighbor(cw(_iv)); pos = pos->neighbor(cw(_iv));
_iv = pos->index(_v); _iv = pos->index(_v);
@ -222,10 +222,10 @@ public:
do do
{ {
_ri = cw(_iv); _ri = cw(_iv);
if (_tri.is_finite_non_hyperbolic(pos, ccw(_iv))) if (!_tri.is_Delaunay_hyperbolic(pos, ccw(_iv)))
{ {
_ri = ccw(_iv); _ri = ccw(_iv);
if (_tri.is_finite_non_hyperbolic(pos, cw(_iv))) if (!_tri.is_Delaunay_hyperbolic(pos, cw(_iv)))
{ {
pos = pos->neighbor(cw(_iv)); pos = pos->neighbor(cw(_iv));
_iv = pos->index(_v); _iv = pos->index(_v);
@ -253,10 +253,10 @@ public:
do do
{ {
_ri = ccw(_iv); _ri = ccw(_iv);
if (_tri.is_finite_non_hyperbolic(pos, cw(_iv))) if (!_tri.is_Delaunay_hyperbolic(pos, cw(_iv)))
{ {
_ri = cw(_iv); _ri = cw(_iv);
if (_tri.is_finite_non_hyperbolic(pos, ccw(_iv))) if (!_tri.is_Delaunay_hyperbolic(pos, ccw(_iv)))
{ {
pos = pos->neighbor(ccw(_iv)); pos = pos->neighbor(ccw(_iv));
_iv = pos->index(_v); _iv = pos->index(_v);
@ -396,12 +396,6 @@ public:
void clear() { Base::clear(); } void clear() { Base::clear(); }
void mark_star(Vertex_handle v) const
{
if(!is_star_bounded(v))
mark_star_faces(v);
}
template<class OutputItFaces> template<class OutputItFaces>
OutputItFaces find_conflicts(const Point& p, OutputItFaces fit, Face_handle start = Face_handle()) const OutputItFaces find_conflicts(const Point& p, OutputItFaces fit, Face_handle start = Face_handle()) const
{ {
@ -412,7 +406,7 @@ public:
Face_handle start = Face_handle()) Face_handle start = Face_handle())
{ {
Vertex_handle v = Base::insert(p, start); Vertex_handle v = Base::insert(p, start);
mark_star(v); mark_star_faces(v);
ensure_hyperbolic_face_handle(v); ensure_hyperbolic_face_handle(v);
return v; return v;
@ -427,7 +421,7 @@ public:
Face_handle loc, int li) Face_handle loc, int li)
{ {
Vertex_handle v = Base::insert(p, lt, loc, li); Vertex_handle v = Base::insert(p, lt, loc, li);
mark_star(v); mark_star_faces(v);
ensure_hyperbolic_face_handle(v); ensure_hyperbolic_face_handle(v);
return v; return v;
@ -491,15 +485,35 @@ public:
template <typename T> template <typename T>
bool is_infinite(T v) const { return Base::is_infinite(v); } bool is_infinite(T v) const { return Base::is_infinite(v); }
bool is_infinite(Face_handle f, int i) const { return Base::is_infinite(f, i); }
bool is_Delaunay_hyperbolic(Face_handle f) const bool is_Delaunay_hyperbolic(Face_handle f) const
{ {
return !Base::is_infinite(f) && !is_finite_non_hyperbolic(f); if(dimension() <= 1)
return false;
return f->hyperbolic_data().is_Delaunay_hyperbolic();
} }
bool is_Delaunay_hyperbolic(Face_handle f, int i) const bool is_Delaunay_hyperbolic(Face_handle f, int i) const
{ {
return !Base::is_infinite(f, i) && !is_finite_non_hyperbolic(f, i); if(dimension() <= 1)
return false;
if(is_infinite(f, i))
return false;
if(f->hyperbolic_data().is_Delaunay_non_hyperbolic(i))
return false;
// another incident face and corresponding index
Face_handle f2 = f->neighbor(i);
int i2 = f2->index(f);
if(f2->hyperbolic_data().is_Delaunay_non_hyperbolic(i2))
return false;
return true;
} }
bool is_Delaunay_hyperbolic(const Edge& e) const bool is_Delaunay_hyperbolic(const Edge& e) const
@ -518,38 +532,6 @@ public:
} }
private: private:
class Face_data
{
private:
// a finite face is non_hyperbolic if its circumscribing circle intersects the circle at infinity
bool _is_Delaunay_hyperbolic;
// defined only if the face is finite and non_hyperbolic
unsigned int _non_hyperbolic_edge;
public:
Face_data() : _is_Delaunay_hyperbolic(true), _non_hyperbolic_edge(UCHAR_MAX) {}
unsigned int get_non_hyperbolic_edge() const
{
CGAL_triangulation_precondition(!_is_Delaunay_hyperbolic);
CGAL_triangulation_precondition(_non_hyperbolic_edge <= 2);
return _non_hyperbolic_edge;
}
void set_non_hyperbolic_edge(unsigned int uschar)
{
CGAL_triangulation_precondition(!_is_Delaunay_hyperbolic);
CGAL_triangulation_precondition(uschar <= 2);
_non_hyperbolic_edge = uschar;
}
bool get_is_Delaunay_hyperbolic() const { return _is_Delaunay_hyperbolic; }
void set_is_Delaunay_hyperbolic(bool flag) { _is_Delaunay_hyperbolic = flag; }
};
/* /*
During the insertion of a new point in the triangulation, the added vertex points to a face. During the insertion of a new point in the triangulation, the added vertex points to a face.
This function ensures that the face to which the vertex points is hyperbolic (if there exists one). This function ensures that the face to which the vertex points is hyperbolic (if there exists one).
@ -634,173 +616,87 @@ private:
// Cannot be on the boundary here. // Cannot be on the boundary here.
lt = FACE; lt = FACE;
li = 4;
if(cs1 != cp1 || cs2 != cp2 || cs3 != cp3) if(cs1 != cp1 || cs2 != cp2 || cs3 != cp3)
return ON_NEGATIVE_SIDE; return ON_NEGATIVE_SIDE;
else else
return ON_POSITIVE_SIDE; return ON_POSITIVE_SIDE;
} }
int get_finite_non_hyperbolic_edge(Face_handle f) const
{
CGAL_triangulation_precondition(is_finite_non_hyperbolic(f));
Face_data fd = object_cast<Face_data>(f->tds_data());
return fd.get_non_hyperbolic_edge();
}
bool is_finite_non_hyperbolic(Face_handle f) const
{
if(const Face_data* td = object_cast<Face_data>(&f->tds_data()))
{
return !td->get_is_Delaunay_hyperbolic();
}
else
{
return false;
}
}
bool is_finite_non_hyperbolic(Face_handle f, int i) const
{
if(dimension() <= 1)
return false;
if(is_finite_non_hyperbolic(f) && get_finite_non_hyperbolic_edge(f) == i)
return true;
// another incident face and corresponding index
Face_handle f2 = f->neighbor(i);
int i2 = f2->index(f);
if(is_finite_non_hyperbolic(f2) && get_finite_non_hyperbolic_edge(f2) == i2)
return true;
return false;
}
bool is_finite_non_hyperbolic(const Edge& e) const
{
return is_finite_non_hyperbolic(e.first, e.second);
}
// Depth-first search (dfs) and marking the finite non_hyperbolic faces.
void mark_finite_non_hyperbolic_faces() const void mark_finite_non_hyperbolic_faces() const
{ {
if(dimension() <= 1) if(dimension() <= 1)
return; return;
std::set<Face_handle> visited_faces; for(auto fit = Base::all_faces_begin(); fit != Base::all_faces_end(); ++fit)
fit->hyperbolic_data().set_Delaunay_hyperbolic(); // finite & hyperbolic
// maintain a stack to be able to backtrack Face_handle ifh = Base::infinite_face();
// to the most recent faces which neighbors are not visited ifh->hyperbolic_data().set_infinite();
std::stack<Face_handle> backtrack;
// start from a face with infinite vertex std::stack<Face_handle> to_visit;
Face_handle current = Base::infinite_face(); to_visit.push(ifh);
// mark it as visited std::set<Face_handle> visited_faces; // @todo squat tds_data()
visited_faces.insert(current);
// put the element whose neighbors we are going to explore. while(!to_visit.empty())
backtrack.push(current);
Face_handle next;
while(!backtrack.empty())
{ {
// take a face Face_handle fh = to_visit.top();
current = backtrack.top(); to_visit.pop();
// start visiting the neighbors if(!visited_faces.insert(fh).second) // already visited previously
int i = 0; continue;
for(; i<3; ++i)
for(int i = 0; i<3; ++i)
{ {
next = current->neighbor(i); Face_handle nfh = fh->neighbor(i);
mark_face(nfh);
// if a neighbor is already visited, then stop going deeper if(is_Delaunay_hyperbolic(nfh))
if(visited_faces.find(next) != visited_faces.end())
continue; continue;
visited_faces.insert(next); to_visit.push(nfh);
mark_face(next);
// go deeper if the neighbor is non_hyperbolic
if(!is_Delaunay_hyperbolic(next))
{
backtrack.push(next);
break;
}
} }
// if all the neighbors are already visited, then remove "current" face.
if(i == 3)
backtrack.pop();
} }
} }
// check if a star is bounded by finite faces
bool is_star_bounded(Vertex_handle v) const
{
if(dimension() <= 1)
return true;
Face_handle f = v->face();
Face_handle next;
int i;
Face_handle start(f);
Face_handle opposite_face;
do
{
i = f->index(v);
next = f->neighbor(ccw(i)); // turn ccw around v
opposite_face = f->neighbor(i);
if(!is_Delaunay_hyperbolic(opposite_face))
return false;
f = next;
}
while(next != start);
return true;
}
void mark_star_faces(Vertex_handle v) const void mark_star_faces(Vertex_handle v) const
{ {
if(dimension() <= 1) if(dimension() <= 1)
return; return;
Face_handle f = v->face(); Face_handle f = v->face();
Face_handle start(f), next; Face_handle start(f);
int i;
do do
{ {
i = f->index(v);
next = f->neighbor(ccw(i)); // turn ccw around v
mark_face(f); mark_face(f);
f = next; int i = f->index(v);
} while(next != start); f = f->neighbor(ccw(i));
}
while(f != start);
} }
void mark_face(const Face_handle f) const void mark_face(const Face_handle f) const
{ {
Is_Delaunay_hyperbolic del; if(is_infinite(f))
int idx; {
bool flag = del(point(f,0), f->hyperbolic_data().set_infinite();
point(f,1), }
point(f,2), else
idx); {
int idx;
bool flag = geom_traits().is_Delaunay_hyperbolic_2_object()(point(f,0),
point(f,1),
point(f,2),
idx);
Face_data fd; if(flag)
fd.set_is_Delaunay_hyperbolic(flag); f->hyperbolic_data().set_Delaunay_hyperbolic(); // finite & hyperbolic
else
if(!flag) f->hyperbolic_data().set_Delaunay_non_hyperbolic(idx); // finite but not hyperbolic
fd.set_non_hyperbolic_edge(idx); }
f->tds_data() = make_object(fd);
} }
public: public:

View File

@ -22,6 +22,40 @@
namespace CGAL { namespace CGAL {
class Hyperbolic_data
{
typedef boost::int8_t Id;
private:
// - 2 for infinite face
// - 1 for finite, hyperbolic face
// 0 for finite, non hyperbolic with non-hyperbolic edge at index 0
// 2 for finite, non hyperbolic with non-hyperbolic edge at index 1
// 1 for finite, non hyperbolic with non-hyperbolic edge at index 2
Id _hyperbolic_tag;
public:
Hyperbolic_data(Id id = -2) : _hyperbolic_tag(id) { }
void set_infinite() { _hyperbolic_tag = -2; }
// a finite face is non_hyperbolic if its circumscribing circle intersects the circle at infinity
void set_Delaunay_hyperbolic() { _hyperbolic_tag = -1; }
bool is_Delaunay_hyperbolic() const
{
return (_hyperbolic_tag == -1);
}
// set and get the non-hyperbolic property of the edge #i
void set_Delaunay_non_hyperbolic(int i) { _hyperbolic_tag = i; }
bool is_Delaunay_non_hyperbolic(int i) const
{
return (_hyperbolic_tag == i);
}
};
template<typename Gt, template<typename Gt,
typename Fb = Triangulation_ds_face_base_2<> > typename Fb = Triangulation_ds_face_base_2<> >
class Hyperbolic_triangulation_face_base_2 class Hyperbolic_triangulation_face_base_2
@ -60,11 +94,11 @@ public:
static int ccw(int i) {return Triangulation_cw_ccw_2::ccw(i);} static int ccw(int i) {return Triangulation_cw_ccw_2::ccw(i);}
static int cw(int i) {return Triangulation_cw_ccw_2::cw(i);} static int cw(int i) {return Triangulation_cw_ccw_2::cw(i);}
CGAL::Object& tds_data() { return this->_tds_data; } Hyperbolic_data& hyperbolic_data() { return this->_hyperbolic_data; }
const CGAL::Object& tds_data() const { return this->_tds_data; } const Hyperbolic_data& hyperbolic_data() const { return this->_hyperbolic_data; }
private: private:
CGAL::Object _tds_data; Hyperbolic_data _hyperbolic_data;
}; };
} // namespace CGAL } // namespace CGAL