changed is_valid in RTOS(added bool verbose when verbose is true, messages are printed

This commit is contained in:
Claudia Werner 2012-10-22 14:12:07 +00:00
parent 167984a073
commit 8ed86471dc
2 changed files with 199 additions and 566 deletions

View File

@ -78,20 +78,10 @@ public:
using Base::orientation;
using Base ::show_all;
using Base ::show_face;
using Base ::show_vertex;
using Base ::show_vertex;
using Base::clear;
#endif
private:
typedef std::list<Face_handle> Faces_around_stack;
public:
//CONSTRUCTORS
Regular_triangulation_on_sphere_2(const Gt& gt=Gt())
@ -101,8 +91,8 @@ public:
//CHECK
bool is_valid(bool verbose = false, int level = 0) const;
bool is_valid_face(Face_handle fh) const;
bool is_valid_vertex(Vertex_handle fh) const;
bool is_valid_face(Face_handle fh,bool verbose = false, int level = 0 ) const;
bool is_valid_vertex(Vertex_handle fh, bool verbose = false, int level = 0) const;
bool test_conflict(const Point &p, Face_handle fh) const;
@ -141,10 +131,8 @@ public:
//INSERTION
Vertex_handle insert(const Point &p, Locate_type lt, Face_handle loc, int li );
Vertex_handle insert(const Point &p, Face_handle f = Face_handle() );
Vertex_handle insert_first(const Point &p);
Vertex_handle insert_second(const Point &p);
Vertex_handle insert_first(const Point &p);
Vertex_handle insert_second(const Point &p);
Vertex_handle insert_outside_affine_hull_regular(const Point& p,bool plane=false);
Vertex_handle insert_hole_approach_2(const Point &p, Locate_type lt, Face_handle loc, int li) ;
Vertex_handle insert_in_plane_triangulation(const Point &p);
@ -152,7 +140,7 @@ public:
//REMOVAL
void remove_degree_3(Vertex_handle v, Face_handle f = Face_handle());
void remove(Vertex_handle v);
void remove_1D(Vertex_handle v);
void remove_1D(Vertex_handle v);
void remove_2D(Vertex_handle v);
bool test_dim_down(Vertex_handle v);
bool test_dim_up(Point p);
@ -164,7 +152,7 @@ public:
template < class InputIterator >
int insert(InputIterator first, InputIterator last)
int insert(InputIterator first, InputIterator last)
{
int n = number_of_vertices();
@ -196,8 +184,8 @@ void delete_faces(FaceIt face_begin, FaceIt face_end)
template <class Stream>
Stream &write_vertices(Stream &out,std::vector<Vertex_handle> &t)
template <class Stream>
Stream &write_vertices(Stream &out,std::vector<Vertex_handle> &t)
{
for(typename std::vector<Vertex_handle>::iterator it= t.begin(); it!= t.end(); ++it)
{if((*it)->face()==Face_handle()){
@ -343,9 +331,9 @@ get_conflicts_and_boundary(const Point &p, OutputItFaces fit, OutputItBoundaryE
private:
template <class OutputItFaces, class OutputItBoundaryEdges>
std::pair<OutputItFaces,OutputItBoundaryEdges>
propagate_conflicts (const Point &p, Face_handle fh, int i,
template <class OutputItFaces, class OutputItBoundaryEdges>
std::pair<OutputItFaces,OutputItBoundaryEdges>
propagate_conflicts (const Point &p, Face_handle fh, int i,
std::pair<OutputItFaces,OutputItBoundaryEdges> pit) const {
Face_handle fn = fh->neighbor(i);
if (fn->get_in_conflict_flag() ==1){
@ -372,24 +360,24 @@ private:
template < class Gt, class Tds >
bool
Regular_triangulation_on_sphere_2<Gt,Tds>::
is_valid(bool verbose, int ) const //int level
is_valid(bool verbose, int level ) const //int level
{
// cannot call for is_valid() of Base Triangulation class
// because 1) number of vertices of base class does not match
// tds.is_valid calls is_valid for each vertex
// and the test is not fullfilled by hidden vertices ...
//bool result = true;
//result = result && Triangulation_on_sphere_2<Gt,Tds>::is_valid(verbose);
bool result = true;
if ( !this-> _tds.is_valid(verbose,level) ) {
if (verbose)
std::cerr << "invalid data structure" << std::endl;
CGAL_triangulation_assertion(false);
return false;
}
bool result = true;
for(Faces_iterator fit = faces_begin();
fit != faces_end(); ++fit) {
result = result && is_valid_face(fit);
}CGAL_triangulation_assertion(result);
for(Vertices_iterator vit = vertices_begin();
vit != vertices_end(); ++vit) {
for(Faces_iterator fit = faces_begin(); fit != faces_end(); ++fit) {
result = result && is_valid_face(fit, verbose, level);
CGAL_triangulation_assertion(result);
}
for(Vertices_iterator vit = vertices_begin(); vit != vertices_end(); ++vit) {
result = result && is_valid_vertex(vit);
}CGAL_triangulation_assertion(result);
@ -413,20 +401,16 @@ Regular_triangulation_on_sphere_2<Gt,Tds>::
}
break;
case 2 :
for(Faces_iterator it=faces_begin();
it!=faces_end(); it++) {
Orientation s = orientation(it->vertex(0)->point(),
it->vertex(1)->point(),
it->vertex(2)->point());
CGAL_triangulation_assertion( s == LEFT_TURN || it->is_negative());
result = result && ( s == LEFT_TURN || it->is_negative());
for(Faces_iterator it=faces_begin(); it!=faces_end(); it++) {
Orientation s = orientation(it->vertex(0)->point(), it->vertex(1)->point(), it->vertex(2)->point());
result = result && ( s == LEFT_TURN || it->is_negative());
CGAL_triangulation_assertion(result);
}
// check number of faces. This cannot be done by the Tds
// which does not know the number of components nor the genus
result = result && (number_of_faces() == 2*(number_of_vertices()) - 4 );
CGAL_triangulation_assertion( result);
break;
@ -434,7 +418,7 @@ Regular_triangulation_on_sphere_2<Gt,Tds>::
// in any dimension
if(verbose) {
std::cerr << " nombres de sommets " << number_of_vertices() << "\t"
std::cerr << " number of vertices " << number_of_vertices() << "\t"
<< std::endl;
}
CGAL_triangulation_assertion( result);
@ -445,12 +429,13 @@ Regular_triangulation_on_sphere_2<Gt,Tds>::
template < class Gt, class Tds >
bool
Regular_triangulation_on_sphere_2<Gt,Tds>::
is_valid_vertex(Vertex_handle vh) const
is_valid_vertex(Vertex_handle vh, bool verbose, int level) const
{
bool result = true;
result = result && vh->face()->has_vertex(vh);
if ( !result) {
if ( !result )
if(verbose) {
std::cerr << " from is_valid_vertex " << std::endl;
std::cerr << "normal vertex " << &(*vh) << std::endl;
std::cerr << vh->point() << " " << std::endl;
@ -465,22 +450,23 @@ is_valid_vertex(Vertex_handle vh) const
template < class Gt, class Tds >
bool
Regular_triangulation_on_sphere_2<Gt,Tds>::
is_valid_face(Face_handle fh) const
is_valid_face(Face_handle fh, bool verbose, int level) const
{
bool result = true;
CGAL_triangulation_assertion(result);
result = fh->get_in_conflict_flag()==0;
typename Vertex_list::iterator vlit = fh->vertex_list().begin(),
vldone = fh->vertex_list().end();
vldone = fh->vertex_list().end();
for (; vlit != vldone; vlit++) {
result = result && Base::power_test(fh, (*vlit)->point()) == ON_NEGATIVE_SIDE;
CGAL_triangulation_assertion(result);
result = result && ((*vlit)->face() == fh);
if (!result) show_face(fh);
if (!result)
if (verbose)
show_face(fh);
CGAL_triangulation_assertion(result);
}
return result;
@ -550,47 +536,38 @@ insert_in_plane_triangulation(const Point &p){
return v;
}
}
template <class Gt, class Tds >
typename Triangulation_on_sphere_2<Gt,Tds>::Vertex_handle
Regular_triangulation_on_sphere_2<Gt,Tds>::
insert_first(const Point& p)
{
CGAL_triangulation_precondition(number_of_vertices() == 0);
Vertex_handle v =this-> _tds.insert_first();
v->set_point(p);
return v;
}
template <class Gt, class Tds >
typename Triangulation_on_sphere_2<Gt,Tds>::Vertex_handle
Regular_triangulation_on_sphere_2<Gt,Tds>::
insert_first(const Point& p)
{
CGAL_triangulation_precondition(number_of_vertices() == 0);
Vertex_handle v =this-> _tds.insert_first();
v->set_point(p);
return v;
}
template <class Gt, class Tds >
typename Triangulation_on_sphere_2<Gt,Tds>::Vertex_handle
Regular_triangulation_on_sphere_2<Gt,Tds>::
insert_second(const Point& p)
{
CGAL_triangulation_precondition(number_of_vertices() == 1);
template <class Gt, class Tds >
typename Triangulation_on_sphere_2<Gt,Tds>::Vertex_handle
Regular_triangulation_on_sphere_2<Gt,Tds>::
insert_second(const Point& p)
{
CGAL_triangulation_precondition(number_of_vertices() == 1);
Vertex_handle v =this-> _tds.insert_second();
v->set_point(p);
return v;
}
Vertex_handle v =this-> _tds.insert_second();
v->set_point(p);
return v;
}
template < class Gt, class Tds >
typename Regular_triangulation_on_sphere_2<Gt,Tds>::Vertex_handle
Regular_triangulation_on_sphere_2<Gt,Tds>::
insert_hole_approach_2(const Point &p, Locate_type lt, Face_handle loc, int li) {
insert_hole_approach_2(const Point &p, Locate_type lt, Face_handle loc, int li) {
//!!!!!!!!!!!!TODO point is valide!!!!!!!!!!!!!!!!!!!!!!!!
Vertex_handle v;
@ -633,7 +610,7 @@ Regular_triangulation_on_sphere_2<Gt,Tds>::
get_conflicts_and_boundary(p, std::back_inserter(faces), std::back_inserter(edges), loc);
v =this-> _tds.star_hole(edges.begin(), edges.end());
v =this->_tds.star_hole(edges.begin(), edges.end());
v->set_point(p);
delete_faces(faces.begin(), faces.end());
@ -644,12 +621,97 @@ Regular_triangulation_on_sphere_2<Gt,Tds>::
return v;
}
}
template <class Gt, class Tds >
typename Triangulation_on_sphere_2<Gt,Tds>::Vertex_handle
Regular_triangulation_on_sphere_2<Gt,Tds>::
insert_outside_affine_hull_regular(const Point& p,bool plane)
{
if(dimension()==0){
Vertex_handle v=vertices_begin();
Vertex_handle u=v->face()->neighbor(0)->vertex(0);
Vertex_handle nv;
//orientation is given by the 2 first points
if( collinear_between(v->point(),u->point(),p) || orientation(u->point(),v->point(),p) == LEFT_TURN )
nv=Base::tds().insert_dim_up(v,false);
else
nv=Base::tds().insert_dim_up(v,true);
nv->set_point(p);
CGAL_triangulation_assertion( orientation(edges_begin()->first->vertex(0)->point(),
edges_begin()->first->vertex(1)->point(),
edges_begin()->first->neighbor(0)->vertex(1)->point())
!= RIGHT_TURN );
//seting negative edge if needed
bool done=false;
Edges_iterator eit=edges_begin();
do{
Face_handle f=eit->first;
Face_handle fn=f->neighbor(0);
Point q=fn->vertex(1)->point();
if(collinear_between(f->vertex(0)->point(),f->vertex(1)->point(),q)){
f->negative()=true;
}
else{
f->negative()=false;
}
++eit;
}while( eit!=edges_end() && !done );
return nv;
}
else{ //dimension=1
bool conform = false;
Face_handle f = (edges_begin())->first;
if(plane){//points coplanar with geom_traits->sphere
Orientation orient = orientation( f->vertex(0)->point(),f->vertex(1)->point(),p);
conform = ( orient == COUNTERCLOCKWISE);
}
else{//three vertices non-coplanar with geom_traits->sphere
Face_handle fn=f->neighbor(0);
const Point p0=f->vertex(0)->point();
const Point p1=f->vertex(1)->point();
const Point p2=fn->vertex(1)->point();
Orientation orient = orientation(p0, p1, p2);
Orientation orient2 = Base::power_test(p0, p1, p2, p);
CGAL_triangulation_assertion(orient);
if(orient2==POSITIVE)
conform =true;
}
Vertex_handle v = this->_tds.insert_dim_up( f->vertex(0), conform);
v->set_point(p);
this->_negative=faces_begin();
//seting negative faces if needed
Faces_iterator fit;
for(fit = faces_begin(); fit != faces_end(); fit++) {
if(orientation(fit)==NEGATIVE){
fit->negative()=true;
this->_negative=fit;
}
}
return v;
}
}
template < class Gt, class Tds >
bool
Regular_triangulation_on_sphere_2<Gt,Tds>::
@ -703,7 +765,7 @@ Regular_triangulation_on_sphere_2<Gt,Tds>::
remove_degree_3(Vertex_handle v, Face_handle f)
{
if (f == Face_handle()) f=v->face();
Base::remove_degree_3(v,f);
this->_tds.remove_degree_3(v,f);
}
@ -728,14 +790,14 @@ CGAL_triangulation_precondition( v != Vertex_handle() );
template <class Gt, class Tds >
void
Regular_triangulation_on_sphere_2<Gt, Tds>::
remove_1D(Vertex_handle v)
{
this->_tds.remove_1D(v);
update_negative_faces();
}
template <class Gt, class Tds >
void
Regular_triangulation_on_sphere_2<Gt, Tds>::
remove_1D(Vertex_handle v)
{
this->_tds.remove_1D(v);
update_negative_faces();
}
@ -787,21 +849,21 @@ test_dim_down(Vertex_handle v)
}
template <class Gt, class Tds >
bool
Regular_triangulation_on_sphere_2<Gt,Tds>::
test_dim_up(Point p){
template <class Gt, class Tds >
bool
Regular_triangulation_on_sphere_2<Gt,Tds>::
test_dim_up(Point p){
// dimension of triangulation increase from 1 to 2 iff the new vertex in not coplanar with the old vertices
//first three points of triangulation
Face_handle f=edges_begin()->first;
Vertex_handle v1=f->vertex(0);
Vertex_handle v2=f->vertex(1);
Vertex_handle v3=f->neighbor(0)->vertex(1);
Face_handle f=edges_begin()->first;
Vertex_handle v1=f->vertex(0);
Vertex_handle v2=f->vertex(1);
Vertex_handle v3=f->neighbor(0)->vertex(1);
return (Base::power_test(v1->point(), v2->point(), v3->point(),p)!=ON_ORIENTED_BOUNDARY);
return (Base::power_test(v1->point(), v2->point(), v3->point(),p)!=ON_ORIENTED_BOUNDARY);
}
}
template < class Gt, class Tds >
@ -954,173 +1016,6 @@ fill_hole_regular(std::list<Edge> & first_hole)
}
//-----------------------------------------------------------------------------REGULAR------------------------------------------------------//
template <class Gt, class Tds >
typename Triangulation_on_sphere_2<Gt,Tds>::Vertex_handle
Regular_triangulation_on_sphere_2<Gt,Tds>::
insert_outside_affine_hull_regular(const Point& p,bool plane)
{
if(dimension()==0){
Vertex_handle v=vertices_begin();
Vertex_handle u=v->face()->neighbor(0)->vertex(0);
Vertex_handle nv;
//orientation is given by the 2 first points
if( collinear_between(v->point(),u->point(),p) || orientation(u->point(),v->point(),p) == LEFT_TURN )
nv=Base::tds().insert_dim_up(v,false);
else
nv=Base::tds().insert_dim_up(v,true);
nv->set_point(p);
CGAL_triangulation_assertion( orientation(edges_begin()->first->vertex(0)->point(),
edges_begin()->first->vertex(1)->point(),
edges_begin()->first->neighbor(0)->vertex(1)->point())
!= RIGHT_TURN );
//seting negative edge if needed
bool done=false;
Edges_iterator eit=edges_begin();
do{
Face_handle f=eit->first;
Face_handle fn=f->neighbor(0);
Point q=fn->vertex(1)->point();
if(collinear_between(f->vertex(0)->point(),f->vertex(1)->point(),q)){
f->negative()=true;
}
else{
f->negative()=false;
}
++eit;
}while( eit!=edges_end() && !done );
return nv;
}
else{ //dimension=1
bool conform = false;
Face_handle f = (edges_begin())->first;
if(plane){//points coplanar with geom_traits->sphere
Orientation orient = orientation( f->vertex(0)->point(),f->vertex(1)->point(),p);
conform = ( orient == COUNTERCLOCKWISE);
}
else{//three vertices non-coplanar with geom_traits->sphere
Face_handle fn=f->neighbor(0);
const Point p0=f->vertex(0)->point();
const Point p1=f->vertex(1)->point();
const Point p2=fn->vertex(1)->point();
//
Orientation orient = orientation(p0, p1, p2);
Orientation orient2 = Base::power_test(p0, p1, p2, p);
CGAL_triangulation_assertion(orient);
if(orient2==POSITIVE)
conform =true;
}
Vertex_handle v = this->_tds.insert_dim_up( f->vertex(0), conform);
v->set_point(p);
this->_negative=faces_begin();
//seting negative faces if needed
bool neg_found=false;
Faces_iterator fit;
int numb =0;
for(fit = faces_begin(); fit != faces_end(); fit++) {
if(orientation(fit)==NEGATIVE){
fit->negative()=true;
neg_found=true;
numb ++;
this->_negative=fit;
}
}
return v;
}
}
//-------------------------------------------------------------------CLASS DEFINITION--------------------------------------------------------//
/*
//Tag to distinguish Delaunay from Regular triangulations
typedef Tag_true Weighted_tag;
public:
Regular_triangulation_on_sphere_2(const Regular_triangulation_on_sphere_2 &rt);
Regular_triangulation_on_sphere_2 & operator=(const Regular_triangulation_on_sphere_2 &tr);
//
// DUAL
Bare_point dual (Face_handle f) const;
Object dual(const Edge &e) const ;
Object dual(const Edge_circulator& ec) const;
Object dual(const Finite_edges_iterator& ei) const;
Bare_point weighted_circumcenter(Face_handle f) const;
Bare_point weighted_circumcenter(const Point& p0,
const Point& p1,
const Point& p2) const;
// Insertion, Deletion and Flip
Vertex_handle push_back(const Point &p);
// Vertex_handle file_input(std::istream& is);
// void file_output(std::ostream& os) const;
public:
void copy_triangulation(const Self& tr);
private:
void copy_triangulation_();
public:
template < class Stream>
Stream& draw_dual(Stream & ps) const
{
Finite_edges_iterator eit = finite_edges_begin();
for (; eit != finite_edges_end(); ++eit) {
Object o = dual(eit);
typename Geom_traits::Line_2 l;
typename Geom_traits::Ray_2 r;
typename Geom_traits::Segment_2 s;
if (CGAL::assign(s,o)) ps << s;
if (CGAL::assign(r,o)) ps << r;
if (CGAL::assign(l,o)) ps << l;
}
return ps;
}
// nearest power vertex query
Vertex_handle nearest_power_vertex(const Bare_point& p) const;
*/
//---------------------------------------------------------------------END CLASS DEFINITION--------------------------------------------------//

View File

@ -20,7 +20,6 @@
#include <CGAL/Triangulation_data_structure_2.h>
#include <CGAL/Triangulation_vertex_base_2.h>
#include <CGAL/Triangulation_face_base_on_sphere_2.h>
//#include <CGAL/Triangulation_line_face_circulator_sphere_2.h>
#include <CGAL/Random.h>
#include <CGAL/spatial_sort.h>
@ -64,8 +63,6 @@ public:
typedef Tds Triangulation_data_structure;
typedef Gt Geom_traits;
typedef typename Geom_traits::Point_2 Point;
// typedef typename Geom_traits::Segment_2 Segment;
// typedef typename Geom_traits::Triangle_2 Triangle;
typedef typename Geom_traits::Orientation_2 Orientation_2;
typedef typename Geom_traits::Coradial_sphere_2 Coradial_sphere_2;
typedef typename Geom_traits::Compare_x_2 Compare_x;
@ -125,7 +122,6 @@ public:
// TESTS
bool is_edge(Vertex_handle va, Vertex_handle vb) const;
bool is_edge(Vertex_handle va, Vertex_handle vb, Face_handle& fr,int & i) const;
//bool includes_edge(Vertex_handle va, Vertex_handle vb,Vertex_handle& vbb, Face_handle& fr, int & i) const;
bool is_face(Vertex_handle v1, Vertex_handle v2, Vertex_handle v3) const;
bool is_face(Vertex_handle v1, Vertex_handle v2, Vertex_handle v3, Face_handle &fr) const;
@ -138,7 +134,8 @@ public:
int dimension() const { return _tds.dimension();}
size_type number_of_vertices() const {return _tds.number_of_vertices();}
size_type number_of_faces() const;
size_type number_of_faces() const{return _tds.number_of_faces();}
int number_of_negative_faces();
@ -156,35 +153,16 @@ public:
Orientation orientation(const Face_handle f) const;
//Comparison_result compare_x(const Point& p, const Point& q) const;
//Comparison_result compare_y(const Point& p, const Point& q) const;
Oriented_side power_test(const Point &p,
const Point &q,
const Point &r,
const Point &s) const;
Oriented_side power_test(const Point &p,
const Point &q,
const Point &r) const;
Oriented_side power_test(const Point &p,
const Point &r) const;
Oriented_side power_test(const Face_handle &f,
const Point &p) const;
Oriented_side power_test(const Face_handle& f, int i,
const Point &p) const;
Oriented_side
oriented_side(const Point &p0, const Point &p1,
const Point &p2, const Point &p) const;
Bounded_side
bounded_side(const Point &p0, const Point &p1,
const Point &p2, const Point &p) const;
Oriented_side
oriented_side(Face_handle f, const Point &p) const;
Oriented_side power_test(const Point &p,const Point &q, const Point &r, const Point &s) const;
Oriented_side power_test(const Point &p,const Point &q, const Point &r) const;
Oriented_side power_test(const Point &p, const Point &r) const;
Oriented_side power_test(const Face_handle &f,const Point &p) const;
Oriented_side power_test(const Face_handle& f, int i, const Point &p) const;
Oriented_side oriented_side(const Point &p0, const Point &p1,const Point &p2, const Point &p) const;
Bounded_side bounded_side(const Point &p0, const Point &p1,const Point &p2, const Point &p) const;
Oriented_side oriented_side(Face_handle f, const Point &p) const;
bool xy_equal(const Point& p, const Point& q) const;
bool
collinear_between(const Point& p, const Point& q, const Point& r) const;
bool collinear_between(const Point& p, const Point& q, const Point& r) const;
//------------------------------------------------------------------DEBUG---------------------------------------------------
@ -200,8 +178,7 @@ public:
Face_handle f2, int i2,
Face_handle f3, int i3);
Face_handle create_face(Face_handle f1, int i1,
Face_handle f2, int i2);
Face_handle create_face(Face_handle f1, int i1, Face_handle f2, int i2);
Face_handle create_face();
Face_handle create_face(Face_handle f, int i, Vertex_handle v);
@ -273,32 +250,9 @@ public:
/*---------------------------------------------------------------------TEMPLATE MEMBERS--------------------------------------*/
public:
template<class EdgeIt, class FaceIt>
Vertex_handle star_hole( const Point& p,
EdgeIt edge_begin,
EdgeIt edge_end,
FaceIt face_begin,
FaceIt face_end) {
typedef typename Triangulation_data_structure::Edge Tds_Edge;
typedef typename Triangulation_data_structure::Face Tds_Face;
Vertex_handle v = _tds.star_hole( edge_begin, edge_end,
face_begin, face_end);
v->set_point(p);
return v;
}
template<class FaceIt>
void delete_faces(FaceIt face_begin, FaceIt face_end)
{
@ -380,14 +334,6 @@ swap(Triangulation_on_sphere_2 &tr)
//ACCESS FUNCTION
//CHECKING
template <class Gt, class Tds >
typename Triangulation_on_sphere_2<Gt, Tds>::size_type
Triangulation_on_sphere_2<Gt, Tds>::
number_of_faces() const
{
size_type count = _tds.number_of_faces();
return count;
}
template <class Gt, class Tds >
bool
@ -439,54 +385,6 @@ is_edge(Vertex_handle va, Vertex_handle vb, Face_handle& fr, int & i) const
{
return _tds.is_edge(va, vb, fr, i);
}
/*
template <class Gt, class Tds >
bool
Triangulation_on_sphere_2<Gt, Tds>::
includes_edge(Vertex_handle va, Vertex_handle vb,
Vertex_handle & vbb,
Face_handle& fr, int & i) const
// returns true if the line segment ab contains an edge e of t
// incident to a, false otherwise
// if true, vbb becomes the vertex of e distinct from a
// fr is the face incident to e and e=(fr,i)
// fr is on the right side of a->b
{
Vertex_handle v;
Orientation orient;
int indv;
Edge_circulator ec = incident_edges(va), done(ec);
if (ec != 0) {
do {
//find the index of the other vertex of *ec
indv = 3 - ((*ec).first)->index(va) - (*ec).second ;
v = ((*ec).first)->vertex(indv);
if (!is_infinite(v)) {
if (v==vb) {
vbb = vb;
fr=(*ec).first;
i= (*ec).second;
return true;
}
else {
orient = orientation(va->point(),
vb->point(),
v->point());
if((orient==COLLINEAR) &&
(collinear_between (va->point(),
v->point(),
vb->point()))) {
vbb = v;
fr=(*ec).first;
i= (*ec).second;
return true;
}
}
}
} while (++ec != done);
}
return false;
}*/
template <class Gt, class Tds >
inline bool
@ -929,35 +827,7 @@ locate(const Point &p,
//--------------------------------------------------------------------ITERATORS AND CIRCULATORS--------------------------------------
/*template <class Gt, class Tds >
inline
typename Triangulation_on_sphere_2<Gt, Tds>::size_type
Triangulation_on_sphere_2<Gt, Tds>::
degree(Vertex_handle v) const
{
return _tds.degree(v);
}
template <class Gt, class Tds >
inline
int
Triangulation_on_sphere_2<Gt, Tds>::
mirror_index(Face_handle f, int i) const
{
return _tds.mirror_index(f,i);
}
template <class Gt, class Tds >
inline
typename Triangulation_on_sphere_2<Gt, Tds>::Vertex_handle
Triangulation_on_sphere_2<Gt, Tds>::
mirror_vertex(Face_handle f, int i) const
{
return _tds.mirror_vertex(f,i);
}
*/
//------------------------------------------------------------------------------PREDICATES-----------------------------------------------------------------
template <class Gt, class Tds >
inline
@ -984,32 +854,14 @@ orientation(const Face_handle f) const
{
return orientation(f->vertex(0)->point(),f->vertex(1)->point(),f->vertex(2)->point());
}
/*
template <class Gt, class Tds >
inline
Comparison_result
Triangulation_on_sphere_2<Gt, Tds>::
compare_x(const Point& p, const Point& q) const
{
return geom_traits().compare_x_2_object()(p,q);
}
template <class Gt, class Tds >
inline
Comparison_result
Triangulation_on_sphere_2<Gt, Tds>::
compare_y(const Point& p, const Point& q) const
{
return geom_traits().compare_y_2_object()(p,q);
}
*/
/*
template <class Gt, class Tds >
Oriented_side
Triangulation_on_sphere_2<Gt, Tds>::
oriented_side(const Point &p0, const Point &p1,
const Point &p2, const Point &p) const
{
oriented_side(const Point &p0, const Point &p1,const Point &p2, const Point &p) const {
Orientation o1 = orientation(p0, p1, p),
o2 = orientation(p1, p2, p),
o3 = orientation(p2, p0, p);
@ -1017,7 +869,7 @@ oriented_side(const Point &p0, const Point &p1,
if(orientation(p0, p1, p2)==POSITIVE){
if(o1==POSITIVE)
if(o2==POSITIVE)
if(o3==POSITIVE)
if(o3==POSITIVE)
return ON_POSITIVE_SIDE;
if (o1 == COLLINEAR){
@ -1046,7 +898,7 @@ oriented_side(const Point &p0, const Point &p1,
else
return ON_NEGATIVE_SIDE;
}
}
}*/
@ -1102,98 +954,8 @@ oriented_side(const Point &p0, const Point &p1,
return geom_traits().power_test_2_object()(p,q,r);
}
/*template <class Gt, class Tds >
Oriented_side
Triangulation_on_sphere_2<Gt, Tds>::
oriented_side(const Point &p0, const Point &p1,
const Point &p2, const Point &p) const
{
Orientation o1 = orientation(p0, p1, p),
o2 = orientation(p1, p2, p),
o3 = orientation(p2, p0, p);
//if(orientation(p0, p1, p2)==POSITIVE &&orientation(p0, p1, p2, p)==POSITIVE ){
//return ON_POSITIVE_SIDE;
// //}
if(orientation(p0, p1, p2)==NEGATIVE &&geom_traits().power_test_2_object()(p0, p1, p2, p)== NEGATIVE ){
return ON_POSITIVE_SIDE;
}
if (o1 == COLLINEAR){
if (o2 == COLLINEAR || o3 == COLLINEAR) return ON_ORIENTED_BOUNDARY;
return ON_NEGATIVE_SIDE;
}
if (o2 == COLLINEAR){
if (o1 == COLLINEAR || o3 == COLLINEAR) return ON_ORIENTED_BOUNDARY;
return ON_NEGATIVE_SIDE;
}
if (o3 == COLLINEAR){
if (o2 == COLLINEAR || o1 == COLLINEAR) return ON_ORIENTED_BOUNDARY;
return ON_NEGATIVE_SIDE;
}
return ON_NEGATIVE_SIDE;
//return ON_NEGATIVE_SIDE;
} */
/*template <class Gt, class Tds >
Oriented_side
Triangulation_on_sphere_2<Gt, Tds>::
oriented_side(const Point &p0, const Point &p1,
const Point &p2, const Point &p) const
{
Orientation o1 = orientation(p0, p1, p),
o2 = orientation(p1, p2, p),
o3 = orientation(p2, p0, p);
if(orientation(p0, p1, p2)==POSITIVE){
if(o1==POSITIVE)
if(o2==POSITIVE)
if(o3==POSITIVE)
return ON_POSITIVE_SIDE;
if (o1 == COLLINEAR){
if (o2 == COLLINEAR || o3 == COLLINEAR) return ON_ORIENTED_BOUNDARY;
return ON_NEGATIVE_SIDE;
}
if (o2 == COLLINEAR){
if (o1 == COLLINEAR || o3 == COLLINEAR) return ON_ORIENTED_BOUNDARY;
return ON_NEGATIVE_SIDE;
}
if (o3 == COLLINEAR){
if (o2 == COLLINEAR || o1 == COLLINEAR) return ON_ORIENTED_BOUNDARY;
return ON_NEGATIVE_SIDE;
}
return ON_NEGATIVE_SIDE;
}else{
if(o1==POSITIVE)
if(o2==POSITIVE)
if(o3==POSITIVE)
return ON_POSITIVE_SIDE;
if(o1==NEGATIVE && o2==NEGATIVE && o3==NEGATIVE )
return ON_NEGATIVE_SIDE;
else
return ON_POSITIVE_SIDE;
}
}
*/
/*
template <class Gt, class Tds >
Bounded_side
Triangulation_on_sphere_2<Gt, Tds>::
@ -1240,30 +1002,7 @@ oriented_side(Face_handle f, const Point &p) const
f->vertex(2)->point(),
p);
}
/*
template < class Gt, class Tds >
Oriented_side
Triangulation_on_sphere_2<Gt,Tds>::
side_of_oriented_circle(Face_handle f, const Point & p) const
{
if ( ! is_infinite(f) ) {
typename Gt::Side_of_oriented_circle_2
in_circle = geom_traits().side_of_oriented_circle_2_object();
return in_circle(f->vertex(0)->point(),
f->vertex(1)->point(),
f->vertex(2)->point(),p);
}
int i = f->index(infinite_vertex());
Orientation o = orientation(f->vertex(ccw(i))->point(),
f->vertex(cw(i))->point(),
p);
return (o == NEGATIVE) ? ON_NEGATIVE_SIDE :
(o == POSITIVE) ? ON_POSITIVE_SIDE :
ON_ORIENTED_BOUNDARY;
}
*/
template <class Gt, class Tds >
bool
Triangulation_on_sphere_2<Gt, Tds>::
@ -1363,7 +1102,7 @@ show_face(Face_handle fh) const
std::cerr <<"]" << std::endl;
break;
case 1:
std::cerr <<"point :" ; show_vertex(fh->vertex(0));
std::cerr <<"point :" ; show_vertex(fh->vertex(0));
std::cerr <<" / voisin " << &(*(fh->neighbor(0)));
std::cerr <<"[" ; show_vertex(fh->neighbor(0)->vertex(0));
std::cerr <<"/" ; show_vertex(fh->neighbor(0)->vertex(1));
@ -1447,7 +1186,6 @@ template <class Gt, class Tds >
inline
void
Triangulation_on_sphere_2<Gt, Tds>::
delete_face(Face_handle f)
{
_tds.delete_face(f);