mirror of https://github.com/CGAL/cgal
Adding the test where the Delaunay triangulation and a (unweighted)
regular triangulation are compared. Also added vertices_in_conflict to Regular_triangulation_3.h for that.
This commit is contained in:
parent
a42b691374
commit
6841f56078
|
|
@ -288,6 +288,17 @@ Returns the \ccc{Triple} composed of the resulting output iterators.
|
|||
with \ccc{p}.}
|
||||
}
|
||||
|
||||
\ccMethod{template <class OutputIterator>
|
||||
OutputIterator
|
||||
vertices_in_conflict(Point p, Cell_handle c,
|
||||
OutputIterator res);}
|
||||
{Similar to \ccc{find_conflicts()}, but reports the vertices which are on the
|
||||
boundary of the conflict hole of \ccc{p}, in the output iterator \ccc{res}.
|
||||
Returns the resulting output iterator.
|
||||
\ccPrecond{\ccVar.\ccc{dimension()} $\geq 2$, and \ccc{c} is a cell containing
|
||||
\ccc{p}.}
|
||||
}
|
||||
|
||||
In the weighted setting, a face (cell, facet, edge or vertex)
|
||||
is said to be a Gabriel face iff
|
||||
the smallest sphere orthogonal to the weighted
|
||||
|
|
|
|||
|
|
@ -206,6 +206,39 @@ public:
|
|||
return std::make_pair(t.first, t.second);
|
||||
}
|
||||
|
||||
// Returns the vertices on the boundary of the conflict hole.
|
||||
template <class OutputIterator>
|
||||
OutputIterator
|
||||
vertices_in_conflict(const Weighted_point&p, Cell_handle c,
|
||||
OutputIterator res) const
|
||||
{
|
||||
CGAL_triangulation_precondition(dimension() >= 2);
|
||||
|
||||
// Get the facets on the boundary of the hole.
|
||||
std::vector<Facet> facets;
|
||||
find_conflicts(p, c, std::back_inserter(facets),
|
||||
Emptyset_iterator(), Emptyset_iterator());
|
||||
|
||||
// Then extract uniquely the vertices.
|
||||
std::set<Vertex_handle> vertices;
|
||||
if (dimension() == 3) {
|
||||
for (typename std::vector<Facet>::const_iterator i = facets.begin();
|
||||
i != facets.end(); ++i) {
|
||||
vertices.insert(i->first->vertex((i->second+1)&3));
|
||||
vertices.insert(i->first->vertex((i->second+2)&3));
|
||||
vertices.insert(i->first->vertex((i->second+3)&3));
|
||||
}
|
||||
} else {
|
||||
for (typename std::vector<Facet>::const_iterator i = facets.begin();
|
||||
i != facets.end(); ++i) {
|
||||
vertices.insert(i->first->vertex(cw(i->second)));
|
||||
vertices.insert(i->first->vertex(ccw(i->second)));
|
||||
}
|
||||
}
|
||||
|
||||
return std::copy(vertices.begin(), vertices.end(), res);
|
||||
}
|
||||
|
||||
void remove (Vertex_handle v);
|
||||
|
||||
template < typename InputIterator >
|
||||
|
|
|
|||
|
|
@ -94,22 +94,8 @@ struct If <CGAL::Tag_false, Then, Else>
|
|||
typedef Else type;
|
||||
};
|
||||
|
||||
|
||||
// test_conflicts() has to be distingushed between Delaunay and Regular,
|
||||
// since Regular does not have it at the moment, and it will probably
|
||||
// have a different interface once it will be implemented anyway.
|
||||
|
||||
template < typename T, typename P >
|
||||
void test_conflicts(T&t, const P *q, CGAL::Tag_true)
|
||||
{
|
||||
bool WARNING_REGULAR_DOES_NOT_HAVE_FIND_CONFLICTS_AND_CO;
|
||||
std::cerr << "--- WARNING : find_conflicts not tested for Regular ---" << std::endl;
|
||||
for (int i=0; i<22; ++i)
|
||||
t.insert(q[i]);
|
||||
}
|
||||
|
||||
template < typename T, typename P >
|
||||
void test_conflicts(T& T3_13, const P *q, CGAL::Tag_false)
|
||||
void test_conflicts(T& T3_13, const P *q)
|
||||
{
|
||||
typedef T Cls;
|
||||
typedef typename Cls::Locate_type Locate_type;
|
||||
|
|
@ -151,12 +137,6 @@ void test_conflicts(T& T3_13, const P *q, CGAL::Tag_false)
|
|||
}
|
||||
}
|
||||
|
||||
template < typename T, typename P >
|
||||
void test_conflicts(T&t, const P *q)
|
||||
{
|
||||
test_conflicts(t, q, typename T::Weighted_tag());
|
||||
}
|
||||
|
||||
template <class Triangulation>
|
||||
void
|
||||
_test_cls_delaunay_3(const Triangulation &)
|
||||
|
|
|
|||
|
|
@ -32,10 +32,7 @@ int main()
|
|||
{
|
||||
typedef CGAL::Regular_triangulation_3<Traits> Cls;
|
||||
|
||||
if (false)
|
||||
_test_cls_delaunay_3( Cls() );
|
||||
else
|
||||
std::cout << "NOT READY TO BE TESTED" << std::endl;
|
||||
_test_cls_delaunay_3( Cls() );
|
||||
|
||||
return 0;
|
||||
}
|
||||
|
|
|
|||
Loading…
Reference in New Issue