mirror of https://github.com/CGAL/cgal
Added support for drawing arrangements induced by circular arcs and (linear) segment (Arr_circle_segment_traits_)
This commit is contained in:
parent
627c8d5953
commit
72235f65d5
|
|
@ -21,4 +21,7 @@ if(CGAL_Qt5_FOUND)
|
||||||
target_link_libraries(parabolas PUBLIC CGAL::CGAL_Basic_viewer)
|
target_link_libraries(parabolas PUBLIC CGAL::CGAL_Basic_viewer)
|
||||||
target_link_libraries(ellipses PUBLIC CGAL::CGAL_Basic_viewer)
|
target_link_libraries(ellipses PUBLIC CGAL::CGAL_Basic_viewer)
|
||||||
target_link_libraries(hyperbolas PUBLIC CGAL::CGAL_Basic_viewer)
|
target_link_libraries(hyperbolas PUBLIC CGAL::CGAL_Basic_viewer)
|
||||||
|
target_link_libraries(polylines PUBLIC CGAL::CGAL_Basic_viewer)
|
||||||
|
target_link_libraries(circles PUBLIC CGAL::CGAL_Basic_viewer)
|
||||||
|
target_link_libraries(circular_arcs PUBLIC CGAL::CGAL_Basic_viewer)
|
||||||
endif()
|
endif()
|
||||||
|
|
|
||||||
|
|
@ -1,6 +1,8 @@
|
||||||
//! \file examples/Arrangement_on_surface_2/circles.cpp
|
//! \file examples/Arrangement_on_surface_2/circles.cpp
|
||||||
// Constructing an arrangement of circles using the circle-segment traits.
|
// Constructing an arrangement of circles using the circle-segment traits.
|
||||||
|
|
||||||
|
#include <CGAL/draw_arrangement_2.h>
|
||||||
|
|
||||||
#include "arr_circular.h"
|
#include "arr_circular.h"
|
||||||
|
|
||||||
int main() {
|
int main() {
|
||||||
|
|
@ -27,5 +29,6 @@ int main() {
|
||||||
std::cout << "The vertex with maximal degree in the arrangement is: "
|
std::cout << "The vertex with maximal degree in the arrangement is: "
|
||||||
<< "v_max = (" << v_max->point() << ") "
|
<< "v_max = (" << v_max->point() << ") "
|
||||||
<< "with degree " << v_max->degree() << "." << std::endl;
|
<< "with degree " << v_max->degree() << "." << std::endl;
|
||||||
|
CGAL::draw(arr);
|
||||||
return 0;
|
return 0;
|
||||||
}
|
}
|
||||||
|
|
|
||||||
|
|
@ -1,6 +1,8 @@
|
||||||
//! \file examples/Arrangement_on_surface_2/circular_arc.cpp
|
//! \file examples/Arrangement_on_surface_2/circular_arc.cpp
|
||||||
// Constructing an arrangement of various circular arcs and line segments.
|
// Constructing an arrangement of various circular arcs and line segments.
|
||||||
|
|
||||||
|
#include <CGAL/draw_arrangement_2.h>
|
||||||
|
|
||||||
#include "arr_circular.h"
|
#include "arr_circular.h"
|
||||||
#include "arr_print.h"
|
#include "arr_print.h"
|
||||||
|
|
||||||
|
|
@ -56,5 +58,6 @@ int main() {
|
||||||
Arrangement arr;
|
Arrangement arr;
|
||||||
insert(arr, curves.begin(), curves.end());
|
insert(arr, curves.begin(), curves.end());
|
||||||
print_arrangement_size(arr);
|
print_arrangement_size(arr);
|
||||||
|
CGAL::draw(arr);
|
||||||
return 0;
|
return 0;
|
||||||
}
|
}
|
||||||
|
|
|
||||||
|
|
@ -4,6 +4,8 @@
|
||||||
#include <vector>
|
#include <vector>
|
||||||
#include <list>
|
#include <list>
|
||||||
|
|
||||||
|
#include <CGAL/draw_arrangement_2.h>
|
||||||
|
|
||||||
#include "arr_polylines.h"
|
#include "arr_polylines.h"
|
||||||
#include "arr_print.h"
|
#include "arr_print.h"
|
||||||
|
|
||||||
|
|
@ -44,5 +46,6 @@ int main() {
|
||||||
insert(arr, pi2);
|
insert(arr, pi2);
|
||||||
insert(arr, pi3);
|
insert(arr, pi3);
|
||||||
print_arrangement_size(arr); // print the arrangement size
|
print_arrangement_size(arr); // print the arrangement size
|
||||||
|
CGAL::draw(arr);
|
||||||
return 0;
|
return 0;
|
||||||
}
|
}
|
||||||
|
|
|
||||||
|
|
@ -25,6 +25,7 @@
|
||||||
|
|
||||||
#include <CGAL/tags.h>
|
#include <CGAL/tags.h>
|
||||||
#include <CGAL/Arr_tags.h>
|
#include <CGAL/Arr_tags.h>
|
||||||
|
#include <CGAL/Cartesian.h>
|
||||||
#include <CGAL/Arr_geometry_traits/Circle_segment_2.h>
|
#include <CGAL/Arr_geometry_traits/Circle_segment_2.h>
|
||||||
|
|
||||||
#include <fstream>
|
#include <fstream>
|
||||||
|
|
@ -378,6 +379,183 @@ public:
|
||||||
}
|
}
|
||||||
//@}
|
//@}
|
||||||
|
|
||||||
|
/// \name Functor definitions for approximations. Used by the landmarks
|
||||||
|
// point-location strategy and the drawing procedure.
|
||||||
|
//@{
|
||||||
|
typedef double Approximate_number_type;
|
||||||
|
typedef CGAL::Cartesian<Approximate_number_type> Approximate_kernel;
|
||||||
|
typedef Approximate_kernel::Point_2 Approximate_point_2;
|
||||||
|
|
||||||
|
class Approximate_2 {
|
||||||
|
protected:
|
||||||
|
using Traits = Arr_circle_segment_traits_2<Kernel, Filter>;
|
||||||
|
|
||||||
|
/*! The traits (in case it has state) */
|
||||||
|
const Traits& m_traits;
|
||||||
|
|
||||||
|
/*! Constructor
|
||||||
|
* \param traits the traits.
|
||||||
|
*/
|
||||||
|
Approximate_2(const Traits& traits) : m_traits(traits) {}
|
||||||
|
|
||||||
|
friend class Arr_circle_segment_traits_2<Kernel, Filter>;
|
||||||
|
|
||||||
|
public:
|
||||||
|
/*! Obtain an approximation of a point coordinate.
|
||||||
|
* \param p the exact point.
|
||||||
|
* \param i the coordinate index (either 0 or 1).
|
||||||
|
* \pre i is either 0 or 1.
|
||||||
|
* \return An approximation of p's x-coordinate (if i == 0), or an
|
||||||
|
* approximation of p's y-coordinate (if i == 1).
|
||||||
|
*/
|
||||||
|
Approximate_number_type operator()(const Point_2& p, int i) const {
|
||||||
|
CGAL_precondition((i == 0) || (i == 1));
|
||||||
|
return (i == 0) ? (CGAL::to_double(p.x())) : (CGAL::to_double(p.y()));
|
||||||
|
}
|
||||||
|
|
||||||
|
/*! Obtain an approximation of a point.
|
||||||
|
*/
|
||||||
|
Approximate_point_2 operator()(const Point_2& p) const
|
||||||
|
{ return Approximate_point_2(operator()(p, 0), operator()(p, 1)); }
|
||||||
|
|
||||||
|
/*! Obtain an approximation of an \f$x\f$-monotone curve.
|
||||||
|
*/
|
||||||
|
template <typename OutputIterator>
|
||||||
|
OutputIterator operator()(const X_monotone_curve_2& xcv, double error,
|
||||||
|
OutputIterator oi, bool l2r = true) const {
|
||||||
|
if (xcv.is_linear()) return approximate_segment(xcv, oi, l2r);
|
||||||
|
return approximate_arc(xcv, error, oi, l2r);;
|
||||||
|
}
|
||||||
|
|
||||||
|
private:
|
||||||
|
/*! Handle segments.
|
||||||
|
*/
|
||||||
|
template <typename OutputIterator>
|
||||||
|
OutputIterator approximate_segment(const X_monotone_curve_2& xcv,
|
||||||
|
OutputIterator oi,
|
||||||
|
bool l2r = true) const {
|
||||||
|
// std::cout << "SEGMENT\n";
|
||||||
|
auto min_vertex = m_traits.construct_min_vertex_2_object();
|
||||||
|
auto max_vertex = m_traits.construct_max_vertex_2_object();
|
||||||
|
const auto& src = (l2r) ? min_vertex(xcv) : max_vertex(xcv);
|
||||||
|
const auto& trg = (l2r) ? max_vertex(xcv) : min_vertex(xcv);
|
||||||
|
auto xs = CGAL::to_double(src.x());
|
||||||
|
auto ys = CGAL::to_double(src.y());
|
||||||
|
auto xt = CGAL::to_double(trg.x());
|
||||||
|
auto yt = CGAL::to_double(trg.y());
|
||||||
|
*oi++ = Approximate_point_2(xs, ys);
|
||||||
|
*oi++ = Approximate_point_2(xt, yt);
|
||||||
|
return oi;
|
||||||
|
}
|
||||||
|
|
||||||
|
template <typename OutputIterator, typename Op, typename Transform>
|
||||||
|
OutputIterator add_points(double x1, double y1, double t1,
|
||||||
|
double x2, double y2, double t2,
|
||||||
|
double error, OutputIterator oi,
|
||||||
|
Op op, Transform transform) const {
|
||||||
|
auto tm = (t1 + t2)*0.5;
|
||||||
|
|
||||||
|
// Compute the canocal point where the error is maximal.
|
||||||
|
double xm, ym;
|
||||||
|
op(tm, xm, ym);
|
||||||
|
|
||||||
|
auto dx = x2 - x1;
|
||||||
|
auto dy = y2 - y1;
|
||||||
|
|
||||||
|
// Compute the error; abort if it is below the threshold
|
||||||
|
auto l = std::sqrt(dx*dx + dy*dy);
|
||||||
|
auto e = std::abs((xm*dy - ym*dx + x2*y1 - x1*y2) / l);
|
||||||
|
if (e < error) return oi;
|
||||||
|
|
||||||
|
double x, y;
|
||||||
|
transform(xm, ym, x, y);
|
||||||
|
add_points(x1, y1, t1, xm, ym, tm, error, oi, op, transform);
|
||||||
|
*oi++ = Approximate_point_2(x, y);
|
||||||
|
add_points(xm, ym, tm, x2, y2, t2, error, oi, op, transform);
|
||||||
|
return oi;
|
||||||
|
}
|
||||||
|
|
||||||
|
/*! Compute the circular point given the parameter t and the transform
|
||||||
|
* data, that is, the center (translation) and the sin and cos of the
|
||||||
|
* rotation angle.
|
||||||
|
*/
|
||||||
|
void circular_point(double r, double t, double& x, double& y) const {
|
||||||
|
x = r * std::cos(t);
|
||||||
|
y = r * std::sin(t);
|
||||||
|
}
|
||||||
|
|
||||||
|
/*! Transform a point. In particular, rotate the canonical point
|
||||||
|
* (`xc`,`yc`) by an angle, the sine and cosine of which are `sint` and
|
||||||
|
* `cost`, respectively, and translate by (`cx`,`cy`).
|
||||||
|
*/
|
||||||
|
void transform_point(double xc, double yc, double cx, double cy,
|
||||||
|
double& x, double& y) const {
|
||||||
|
x = xc + cx;
|
||||||
|
y = yc + cy;
|
||||||
|
}
|
||||||
|
|
||||||
|
/*! Handle circular arcs.
|
||||||
|
*/
|
||||||
|
template <typename OutputIterator>
|
||||||
|
OutputIterator approximate_arc(const X_monotone_curve_2& xcv,
|
||||||
|
double error, OutputIterator oi,
|
||||||
|
bool l2r = true) const {
|
||||||
|
auto min_vertex = m_traits.construct_min_vertex_2_object();
|
||||||
|
auto max_vertex = m_traits.construct_max_vertex_2_object();
|
||||||
|
const auto& src = (l2r) ? min_vertex(xcv) : max_vertex(xcv);
|
||||||
|
const auto& trg = (l2r) ? max_vertex(xcv) : min_vertex(xcv);
|
||||||
|
auto xs = CGAL::to_double(src.x());
|
||||||
|
auto ys = CGAL::to_double(src.y());
|
||||||
|
auto xt = CGAL::to_double(trg.x());
|
||||||
|
auto yt = CGAL::to_double(trg.y());
|
||||||
|
|
||||||
|
const typename Kernel::Circle_2& circ = xcv.supporting_circle();
|
||||||
|
auto r_sqr = circ.squared_radius();
|
||||||
|
auto r = std::sqrt(CGAL::to_double(r_sqr));
|
||||||
|
|
||||||
|
// Obtain the center:
|
||||||
|
auto cx = CGAL::to_double(circ.center().x());
|
||||||
|
auto cy = CGAL::to_double(circ.center().y());
|
||||||
|
|
||||||
|
// Inverse transform the source and target
|
||||||
|
auto xs_t = xs - cx;
|
||||||
|
auto ys_t = ys - cy;
|
||||||
|
auto xt_t = xt - cx;
|
||||||
|
auto yt_t = yt - cy;
|
||||||
|
|
||||||
|
// Compute the parameters ts and tt such that
|
||||||
|
// source == (x(ts),y(ts)), and
|
||||||
|
// target == (x(tt),y(tt))
|
||||||
|
auto ts = std::atan2(r*ys_t, r*xs_t);
|
||||||
|
if (ts < 0) ts += 2*M_PI;
|
||||||
|
auto tt = std::atan2(r*yt_t, r*xt_t);
|
||||||
|
if (tt < 0) tt += 2*M_PI;
|
||||||
|
auto orient(xcv.orientation());
|
||||||
|
if (xcv.source() != src) orient = CGAL::opposite(orient);
|
||||||
|
if (orient == COUNTERCLOCKWISE) {
|
||||||
|
if (tt < ts) tt += 2*M_PI;
|
||||||
|
}
|
||||||
|
else {
|
||||||
|
if (ts < tt) ts += 2*M_PI;
|
||||||
|
}
|
||||||
|
|
||||||
|
*oi++ = Approximate_point_2(xs, ys);
|
||||||
|
add_points(xs_t, ys_t, ts, xt_t, yt_t, tt, error, oi,
|
||||||
|
[&](double tm, double& xm, double& ym) {
|
||||||
|
circular_point(r, tm, xm, ym);
|
||||||
|
},
|
||||||
|
[&](double xc, double& yc, double& x, double& y) {
|
||||||
|
transform_point(xc, yc, cx, cy, x, y);
|
||||||
|
});
|
||||||
|
*oi++ = Approximate_point_2(xt, yt);
|
||||||
|
return oi;
|
||||||
|
}
|
||||||
|
};
|
||||||
|
|
||||||
|
/*! Obtain an Approximate_2 functor object. */
|
||||||
|
Approximate_2 approximate_2_object() const { return Approximate_2(*this); }
|
||||||
|
//@}
|
||||||
|
|
||||||
/// \name Intersections, subdivisions, and mergings
|
/// \name Intersections, subdivisions, and mergings
|
||||||
//@{
|
//@{
|
||||||
|
|
||||||
|
|
|
||||||
|
|
@ -1674,9 +1674,9 @@ public:
|
||||||
double a;
|
double a;
|
||||||
double ts, tt;
|
double ts, tt;
|
||||||
double cx, cy;
|
double cx, cy;
|
||||||
m_traits.approximate_parabola(r_m, t_m, s_m, u_m, v_m, w_m, cost, sint,
|
m_traits.approximate_parabola(xcv,
|
||||||
xs_t, ys_t, xt_t, yt_t, a, ts, tt, cx, cy,
|
r_m, t_m, s_m, u_m, v_m, w_m, cost, sint,
|
||||||
xcv);
|
xs_t, ys_t, xt_t, yt_t, a, ts, tt, cx, cy);
|
||||||
|
|
||||||
auto ds = parabolic_arc_length(xs_t, 2.0*std::abs(ys_t));
|
auto ds = parabolic_arc_length(xs_t, 2.0*std::abs(ys_t));
|
||||||
auto dt = parabolic_arc_length(xt_t, 2.0*std::abs(yt_t));
|
auto dt = parabolic_arc_length(xt_t, 2.0*std::abs(yt_t));
|
||||||
|
|
@ -1694,9 +1694,9 @@ public:
|
||||||
double a, b;
|
double a, b;
|
||||||
double cx, cy;
|
double cx, cy;
|
||||||
double ts, tt;
|
double ts, tt;
|
||||||
m_traits.approximate_ellipse(r_m, t_m, s_m, u_m, v_m, w_m, cost, sint,
|
m_traits.approximate_ellipse(xcv, r_m, t_m, s_m, u_m, v_m, w_m, cost, sint,
|
||||||
xs_t, ys_t, ts, xt_t, yt_t, tt,
|
xs_t, ys_t, ts, xt_t, yt_t, tt,
|
||||||
a, b, cx, cy, xcv);
|
a, b, cx, cy);
|
||||||
|
|
||||||
namespace bm = boost::math;
|
namespace bm = boost::math;
|
||||||
auto ratio = b/a;
|
auto ratio = b/a;
|
||||||
|
|
@ -1746,31 +1746,29 @@ public:
|
||||||
/*! Obtain an approximation of a point.
|
/*! Obtain an approximation of a point.
|
||||||
*/
|
*/
|
||||||
Approximate_point_2 operator()(const Point_2& p) const
|
Approximate_point_2 operator()(const Point_2& p) const
|
||||||
{ return std::make_pair(operator()(p, 0), operator()(p, 1)); }
|
{ return Approximate_point_2(operator()(p, 0), operator()(p, 1)); }
|
||||||
|
|
||||||
/*! Obtain an approximation of an \f$x\f$-monotone curve.
|
/*! Obtain an approximation of an \f$x\f$-monotone curve.
|
||||||
*/
|
*/
|
||||||
template <typename OutputIterator>
|
template <typename OutputIterator>
|
||||||
OutputIterator operator()(OutputIterator oi, double error,
|
OutputIterator operator()(const X_monotone_curve_2& xcv, double error,
|
||||||
const X_monotone_curve_2& xcv,
|
OutputIterator oi, bool l2r = true) const {
|
||||||
bool l2r = true) const {
|
|
||||||
if (xcv.orientation() == COLLINEAR)
|
if (xcv.orientation() == COLLINEAR)
|
||||||
return approximate_segment(oi, xcv, l2r);
|
return approximate_segment(xcv, oi, l2r);
|
||||||
CGAL::Sign sign_conic = CGAL::sign(4*xcv.r()*xcv.s() - xcv.t()*xcv.t());
|
CGAL::Sign sign_conic = CGAL::sign(4*xcv.r()*xcv.s() - xcv.t()*xcv.t());
|
||||||
if (sign_conic == POSITIVE)
|
if (sign_conic == POSITIVE)
|
||||||
return approximate_ellipse(oi, error, xcv, l2r);
|
return approximate_ellipse(xcv, error, oi, l2r);
|
||||||
if (sign_conic == NEGATIVE)
|
if (sign_conic == NEGATIVE)
|
||||||
return approximate_hyperbola(oi, error, xcv, l2r);
|
return approximate_hyperbola(xcv, error, oi, l2r);
|
||||||
return approximate_parabola(oi, error, xcv, l2r);
|
return approximate_parabola(xcv, error, oi, l2r);
|
||||||
}
|
}
|
||||||
|
|
||||||
private:
|
private:
|
||||||
/*! Handle segments.
|
/*! Handle segments.
|
||||||
*/
|
*/
|
||||||
template <typename OutputIterator>
|
template <typename OutputIterator>
|
||||||
OutputIterator approximate_segment(OutputIterator oi,
|
OutputIterator approximate_segment(const X_monotone_curve_2& xcv,
|
||||||
const X_monotone_curve_2& xcv,
|
OutputIterator oi, bool l2r) const {
|
||||||
bool l2r) const {
|
|
||||||
// std::cout << "SEGMENT\n";
|
// std::cout << "SEGMENT\n";
|
||||||
auto min_vertex = m_traits.construct_min_vertex_2_object();
|
auto min_vertex = m_traits.construct_min_vertex_2_object();
|
||||||
auto max_vertex = m_traits.construct_max_vertex_2_object();
|
auto max_vertex = m_traits.construct_max_vertex_2_object();
|
||||||
|
|
@ -1865,9 +1863,9 @@ public:
|
||||||
* which approximates the arc.
|
* which approximates the arc.
|
||||||
*/
|
*/
|
||||||
template <typename OutputIterator>
|
template <typename OutputIterator>
|
||||||
OutputIterator approximate_ellipse(OutputIterator oi, double error,
|
OutputIterator approximate_ellipse(const X_monotone_curve_2& xcv,
|
||||||
const X_monotone_curve_2& xcv,
|
double error, OutputIterator oi,
|
||||||
bool l2r) const {
|
bool l2r = true) const {
|
||||||
// std::cout << "ELLIPSE\n";
|
// std::cout << "ELLIPSE\n";
|
||||||
auto min_vertex = m_traits.construct_min_vertex_2_object();
|
auto min_vertex = m_traits.construct_min_vertex_2_object();
|
||||||
auto max_vertex = m_traits.construct_max_vertex_2_object();
|
auto max_vertex = m_traits.construct_max_vertex_2_object();
|
||||||
|
|
@ -1887,9 +1885,9 @@ public:
|
||||||
double a, b;
|
double a, b;
|
||||||
double cx, cy;
|
double cx, cy;
|
||||||
double ts, tt;
|
double ts, tt;
|
||||||
m_traits.approximate_ellipse(r_m, t_m, s_m, u_m, v_m, w_m, cost, sint,
|
m_traits.approximate_ellipse(xcv, r_m, t_m, s_m, u_m, v_m, w_m, cost, sint,
|
||||||
xs_t, ys_t, ts, xt_t, yt_t, tt,
|
xs_t, ys_t, ts, xt_t, yt_t, tt,
|
||||||
a, b, cx, cy, xcv, l2r);
|
a, b, cx, cy, l2r);
|
||||||
// std::cout << "a, b: " << a << "," << b << std::endl;
|
// std::cout << "a, b: " << a << "," << b << std::endl;
|
||||||
|
|
||||||
*oi++ = Approximate_point_2(xs, ys);
|
*oi++ = Approximate_point_2(xs, ys);
|
||||||
|
|
@ -1954,7 +1952,7 @@ public:
|
||||||
return oi;
|
return oi;
|
||||||
}
|
}
|
||||||
|
|
||||||
/*! Compute the hyperbolic point given the parameter t and the transform
|
/*! Compute the elliptic point given the parameter t and the transform
|
||||||
* data, that is, the center (translation) and the sin and cos of the
|
* data, that is, the center (translation) and the sin and cos of the
|
||||||
* rotation angle.
|
* rotation angle.
|
||||||
*/
|
*/
|
||||||
|
|
@ -1969,8 +1967,9 @@ public:
|
||||||
* https://www.vcalc.com/wiki/vCalc/Parabola+-+arc+length
|
* https://www.vcalc.com/wiki/vCalc/Parabola+-+arc+length
|
||||||
*/
|
*/
|
||||||
template <typename OutputIterator>
|
template <typename OutputIterator>
|
||||||
OutputIterator approximate_parabola(OutputIterator oi, double error,
|
OutputIterator approximate_parabola(const X_monotone_curve_2& xcv,
|
||||||
const X_monotone_curve_2& xcv, bool l2r)
|
double error, OutputIterator oi,
|
||||||
|
bool l2r = true)
|
||||||
const {
|
const {
|
||||||
// std::cout << "PARABOLA\n";
|
// std::cout << "PARABOLA\n";
|
||||||
auto min_vertex = m_traits.construct_min_vertex_2_object();
|
auto min_vertex = m_traits.construct_min_vertex_2_object();
|
||||||
|
|
@ -1991,9 +1990,10 @@ public:
|
||||||
double a;
|
double a;
|
||||||
double ts, tt;
|
double ts, tt;
|
||||||
double cx, cy;
|
double cx, cy;
|
||||||
m_traits.approximate_parabola(r_m, t_m, s_m, u_m, v_m, w_m, cost, sint,
|
m_traits.approximate_parabola(xcv,
|
||||||
|
r_m, t_m, s_m, u_m, v_m, w_m, cost, sint,
|
||||||
xs_t, ys_t, xt_t, yt_t, a, ts, tt, cx, cy,
|
xs_t, ys_t, xt_t, yt_t, a, ts, tt, cx, cy,
|
||||||
xcv, l2r);
|
l2r);
|
||||||
// std::cout << "sint, cost: " << sint << "," << cost << std::endl;
|
// std::cout << "sint, cost: " << sint << "," << cost << std::endl;
|
||||||
// std::cout << "a: " << a << std::endl;
|
// std::cout << "a: " << a << std::endl;
|
||||||
// std::cout << "xs' = " << xs_t << "," << ys_t << std::endl;
|
// std::cout << "xs' = " << xs_t << "," << ys_t << std::endl;
|
||||||
|
|
@ -2074,9 +2074,9 @@ public:
|
||||||
/*! Handle hyperbolas.
|
/*! Handle hyperbolas.
|
||||||
*/
|
*/
|
||||||
template <typename OutputIterator>
|
template <typename OutputIterator>
|
||||||
OutputIterator approximate_hyperbola(OutputIterator oi, double error,
|
OutputIterator approximate_hyperbola(const X_monotone_curve_2& xcv,
|
||||||
const X_monotone_curve_2& xcv,
|
double error, OutputIterator oi,
|
||||||
bool l2r) const {
|
bool l2r = true) const {
|
||||||
// std::cout << "HYPERBOLA\n";
|
// std::cout << "HYPERBOLA\n";
|
||||||
auto min_vertex = m_traits.construct_min_vertex_2_object();
|
auto min_vertex = m_traits.construct_min_vertex_2_object();
|
||||||
auto max_vertex = m_traits.construct_max_vertex_2_object();
|
auto max_vertex = m_traits.construct_max_vertex_2_object();
|
||||||
|
|
@ -2093,9 +2093,10 @@ public:
|
||||||
double a, b;
|
double a, b;
|
||||||
double cx, cy;
|
double cx, cy;
|
||||||
double ts, tt;
|
double ts, tt;
|
||||||
m_traits.approximate_hyperbola(r_m, t_m, s_m, u_m, v_m, w_m, cost, sint,
|
m_traits.approximate_hyperbola(xcv, r_m, t_m, s_m, u_m, v_m, w_m,
|
||||||
|
cost, sint,
|
||||||
xs_t, ys_t, ts, xt_t, yt_t, tt,
|
xs_t, ys_t, ts, xt_t, yt_t, tt,
|
||||||
a, b, cx, cy, xcv, l2r);
|
a, b, cx, cy, l2r);
|
||||||
// std::cout << "a, b: " << a << "," << b << std::endl;
|
// std::cout << "a, b: " << a << "," << b << std::endl;
|
||||||
// std::cout << "ts, tt: " << ts << "," << tt << std::endl;
|
// std::cout << "ts, tt: " << ts << "," << tt << std::endl;
|
||||||
|
|
||||||
|
|
@ -4027,14 +4028,15 @@ public:
|
||||||
* The arc-length closed form can be found here:
|
* The arc-length closed form can be found here:
|
||||||
* https://www.vcalc.com/wiki/vCalc/Parabola+-+arc+length
|
* https://www.vcalc.com/wiki/vCalc/Parabola+-+arc+length
|
||||||
*/
|
*/
|
||||||
void approximate_parabola(double& r_m, double& t_m, double& s_m,
|
void approximate_parabola(const X_monotone_curve_2& xcv,
|
||||||
|
double& r_m, double& t_m, double& s_m,
|
||||||
double& u_m, double& v_m, double& w_m,
|
double& u_m, double& v_m, double& w_m,
|
||||||
double& cost, double& sint,
|
double& cost, double& sint,
|
||||||
double& xs_t, double& ys_t,
|
double& xs_t, double& ys_t,
|
||||||
double& xt_t, double& yt_t,
|
double& xt_t, double& yt_t,
|
||||||
double& a, double& ts, double& tt,
|
double& a, double& ts, double& tt,
|
||||||
double& cx, double& cy,
|
double& cx, double& cy,
|
||||||
const X_monotone_curve_2& xcv, bool l2r = true)
|
bool l2r = true)
|
||||||
const {
|
const {
|
||||||
auto min_vertex = construct_min_vertex_2_object();
|
auto min_vertex = construct_min_vertex_2_object();
|
||||||
auto max_vertex = construct_max_vertex_2_object();
|
auto max_vertex = construct_max_vertex_2_object();
|
||||||
|
|
@ -4095,13 +4097,14 @@ public:
|
||||||
|
|
||||||
/*! Handle ellipses.
|
/*! Handle ellipses.
|
||||||
*/
|
*/
|
||||||
void approximate_ellipse(double& r_m, double& t_m, double& s_m,
|
void approximate_ellipse(const X_monotone_curve_2& xcv,
|
||||||
|
double& r_m, double& t_m, double& s_m,
|
||||||
double& u_m, double& v_m, double& w_m,
|
double& u_m, double& v_m, double& w_m,
|
||||||
double& cost, double& sint,
|
double& cost, double& sint,
|
||||||
double& xs_t, double& ys_t, double& ts,
|
double& xs_t, double& ys_t, double& ts,
|
||||||
double& xt_t, double& yt_t, double& tt,
|
double& xt_t, double& yt_t, double& tt,
|
||||||
double& a, double& b, double& cx, double& cy,
|
double& a, double& b, double& cx, double& cy,
|
||||||
const X_monotone_curve_2& xcv, bool l2r = true)
|
bool l2r = true)
|
||||||
const {
|
const {
|
||||||
auto min_vertex = construct_min_vertex_2_object();
|
auto min_vertex = construct_min_vertex_2_object();
|
||||||
auto max_vertex = construct_max_vertex_2_object();
|
auto max_vertex = construct_max_vertex_2_object();
|
||||||
|
|
@ -4172,13 +4175,14 @@ public:
|
||||||
|
|
||||||
/*! Handle hyperbolas.
|
/*! Handle hyperbolas.
|
||||||
*/
|
*/
|
||||||
void approximate_hyperbola(double& r_m, double& t_m, double& s_m,
|
void approximate_hyperbola(const X_monotone_curve_2& xcv,
|
||||||
|
double& r_m, double& t_m, double& s_m,
|
||||||
double& u_m, double& v_m, double& w_m,
|
double& u_m, double& v_m, double& w_m,
|
||||||
double& cost, double& sint,
|
double& cost, double& sint,
|
||||||
double& xs_t, double& ys_t, double& ts,
|
double& xs_t, double& ys_t, double& ts,
|
||||||
double& xt_t, double& yt_t, double& tt,
|
double& xt_t, double& yt_t, double& tt,
|
||||||
double& a, double& b, double& cx, double& cy,
|
double& a, double& b, double& cx, double& cy,
|
||||||
const X_monotone_curve_2& xcv, bool l2r = true)
|
bool l2r = true)
|
||||||
const {
|
const {
|
||||||
auto min_vertex = construct_min_vertex_2_object();
|
auto min_vertex = construct_min_vertex_2_object();
|
||||||
auto max_vertex = construct_max_vertex_2_object();
|
auto max_vertex = construct_max_vertex_2_object();
|
||||||
|
|
|
||||||
|
|
@ -21,17 +21,19 @@
|
||||||
#include <cstdlib>
|
#include <cstdlib>
|
||||||
#include <random>
|
#include <random>
|
||||||
|
|
||||||
|
#include <experimental/type_traits>
|
||||||
|
#include <boost/hana.hpp>
|
||||||
|
|
||||||
#include <CGAL/Qt/Basic_viewer_qt.h>
|
#include <CGAL/Qt/Basic_viewer_qt.h>
|
||||||
|
|
||||||
#ifdef CGAL_USE_BASIC_VIEWER
|
#ifdef CGAL_USE_BASIC_VIEWER
|
||||||
|
|
||||||
#include <CGAL/Qt/init_ogl_context.h>
|
#include <CGAL/Qt/init_ogl_context.h>
|
||||||
#include <CGAL/Arrangement_2.h>
|
#include <CGAL/Arrangement_2.h>
|
||||||
#include <CGAL/Arr_conic_traits_2.h>
|
|
||||||
|
|
||||||
namespace CGAL {
|
namespace CGAL {
|
||||||
|
|
||||||
// Viewer class for Polygon_2
|
// Viewer class for`< Polygon_2
|
||||||
template <typename GeometryTraits_2, typename Dcel>
|
template <typename GeometryTraits_2, typename Dcel>
|
||||||
class Arr_2_basic_viewer_qt : public Basic_viewer_qt {
|
class Arr_2_basic_viewer_qt : public Basic_viewer_qt {
|
||||||
typedef GeometryTraits_2 Gt;
|
typedef GeometryTraits_2 Gt;
|
||||||
|
|
@ -45,9 +47,9 @@ class Arr_2_basic_viewer_qt : public Basic_viewer_qt {
|
||||||
typedef typename Arr::Ccb_halfedge_const_circulator
|
typedef typename Arr::Ccb_halfedge_const_circulator
|
||||||
Ccb_halfedge_const_circulator;
|
Ccb_halfedge_const_circulator;
|
||||||
|
|
||||||
// typedef double Approximate_number_type;
|
template <typename T>
|
||||||
// typedef CGAL::Cartesian<Approximate_number_type> Approximate_kernel;
|
using approximate_2_object_t =
|
||||||
// typedef Approximate_kernel::Point_2 Approximate_point_2;
|
decltype(std::declval<T&>().approximate_2_object());
|
||||||
|
|
||||||
public:
|
public:
|
||||||
/// Construct the viewer.
|
/// Construct the viewer.
|
||||||
|
|
@ -86,12 +88,33 @@ public:
|
||||||
|
|
||||||
//! Compute the bounding box
|
//! Compute the bounding box
|
||||||
CGAL::Bbox_2 bounding_box() {
|
CGAL::Bbox_2 bounding_box() {
|
||||||
|
namespace bh = boost::hana;
|
||||||
|
auto has_approximate_2_object =
|
||||||
|
bh::is_valid([](auto&& x) -> decltype(x.approximate_2_object()){});
|
||||||
|
|
||||||
// At this point we assume that the arrangement is not open, and thus the
|
// At this point we assume that the arrangement is not open, and thus the
|
||||||
// bounding box is defined by the vertices.
|
// bounding box is defined by the vertices.
|
||||||
//! The bounding box
|
//! The bounding box
|
||||||
CGAL::Bbox_2 bbox;
|
CGAL::Bbox_2 bbox;
|
||||||
for (auto it = m_arr.vertices_begin(); it != m_arr.vertices_end(); ++it)
|
for (auto it = m_arr.vertices_begin(); it != m_arr.vertices_end(); ++it) {
|
||||||
bbox += it->point().bbox();
|
bh::if_(has_approximate_2_object(Gt{}),
|
||||||
|
[&](auto& t) {
|
||||||
|
const auto* traits = this->m_arr.geometry_traits();
|
||||||
|
auto approx = traits->approximate_2_object();
|
||||||
|
auto has_operator =
|
||||||
|
bh::is_valid([](auto&& x) -> decltype(x.operator()(Point{}, int{})){});
|
||||||
|
bh::if_(has_operator(approx),
|
||||||
|
[&](auto& a) {
|
||||||
|
auto x = approx(it->point(), 0);
|
||||||
|
auto y = approx(it->point(), 1);
|
||||||
|
bbox += CGAL::Bbox_2(x, y, x, y);
|
||||||
|
},
|
||||||
|
[&](auto& a) { bbox += it->point().bbox(); }
|
||||||
|
)(approx);
|
||||||
|
},
|
||||||
|
[&](auto& t) { bbox += it->point().bbox(); }
|
||||||
|
)(*this);
|
||||||
|
}
|
||||||
return bbox;
|
return bbox;
|
||||||
}
|
}
|
||||||
|
|
||||||
|
|
@ -167,36 +190,175 @@ protected:
|
||||||
return ext;
|
return ext;
|
||||||
}
|
}
|
||||||
|
|
||||||
//!
|
//! Draw a region using aproximate coordinates.
|
||||||
virtual void draw_region(Ccb_halfedge_const_circulator circ) {
|
// Call this member function only of the geometry traits is equipped to
|
||||||
CGAL::IO::Color color(m_uni(m_rng), m_uni(m_rng), m_uni(m_rng));
|
// provide aproximate coordinates.
|
||||||
this->face_begin(color);
|
template <typename Approximate>
|
||||||
|
void draw_approximate_region(Halfedge_const_handle curr,
|
||||||
auto ext = find_smallest(circ);
|
const Approximate& approx) {
|
||||||
auto curr = ext;
|
std::vector<typename Gt::Approximate_point_2> polyline;
|
||||||
do {
|
double error(this->pixel_ratio());
|
||||||
// Skip halfedges that are "antenas":
|
bool l2r = curr->direction() == ARR_LEFT_TO_RIGHT;
|
||||||
while (curr->face() == curr->twin()->face())
|
approx(curr->curve(), error, std::back_inserter(polyline), l2r);
|
||||||
curr = curr->twin()->next();
|
if (polyline.empty()) return;
|
||||||
|
auto it = polyline.begin();
|
||||||
this->add_point_in_face(curr->source()->point());
|
auto prev = it++;
|
||||||
draw_curve(curr->curve());
|
for (; it != polyline.end(); prev = it++) {
|
||||||
curr = curr->next();
|
this->add_segment(*prev, *it);
|
||||||
} while (curr != ext);
|
this->add_point_in_face(*prev);
|
||||||
|
}
|
||||||
this->face_end();
|
|
||||||
}
|
}
|
||||||
|
|
||||||
//!
|
//! Draw an exact curve.
|
||||||
virtual void draw_curve(const X_monotone_curve& curve) {
|
template <typename XMonotoneCurve>
|
||||||
|
void draw_exact_curve(const XMonotoneCurve& curve) {
|
||||||
const auto* traits = this->m_arr.geometry_traits();
|
const auto* traits = this->m_arr.geometry_traits();
|
||||||
auto ctr_min = traits->construct_min_vertex_2_object();
|
auto ctr_min = traits->construct_min_vertex_2_object();
|
||||||
auto ctr_max = traits->construct_max_vertex_2_object();
|
auto ctr_max = traits->construct_max_vertex_2_object();
|
||||||
this->add_segment(ctr_min(curve), ctr_max(curve));
|
this->add_segment(ctr_min(curve), ctr_max(curve));
|
||||||
}
|
}
|
||||||
|
|
||||||
|
//! Draw an exact region.
|
||||||
|
void draw_exact_region(Halfedge_const_handle curr) {
|
||||||
|
this->add_point_in_face(curr->source()->point());
|
||||||
|
draw_exact_curve(curr->curve());
|
||||||
|
}
|
||||||
|
|
||||||
|
//! Utility struct
|
||||||
|
template <typename T>
|
||||||
|
struct can_call_operator_curve {
|
||||||
|
template <typename F>
|
||||||
|
constexpr auto operator()(F&& f) ->
|
||||||
|
decltype(f.template operator()<T>(X_monotone_curve{}, double{}, T{}, bool{})) {}
|
||||||
|
};
|
||||||
|
|
||||||
|
//! Draw a region.
|
||||||
|
void draw_region(Ccb_halfedge_const_circulator circ) {
|
||||||
|
/* Check whether the traits has a member function called
|
||||||
|
* approximate_2_object() and if so check whether the return type, namely
|
||||||
|
* `Approximate_2` has an appropriate operator.
|
||||||
|
*
|
||||||
|
* C++20 supports concepts and `requires` expression; see, e.g.,
|
||||||
|
* https://en.cppreference.com/w/cpp/language/constraints; thus, the first
|
||||||
|
* condition above can be elegantly verified as follows:
|
||||||
|
* constexpr bool has_approximate_2_object =
|
||||||
|
* requires(const Gt& traits) { traits.approximate_2_object(); };
|
||||||
|
*
|
||||||
|
* C++17 has experimental constructs called is_detected and
|
||||||
|
* is_detected_v that can be used to achieve the same goal.
|
||||||
|
*
|
||||||
|
* For now we use C++14 features and boost.hana instead.
|
||||||
|
*/
|
||||||
|
namespace bh = boost::hana;
|
||||||
|
auto has_approximate_2_object =
|
||||||
|
bh::is_valid([](auto&& x) -> decltype(x.approximate_2_object()){});
|
||||||
|
|
||||||
|
CGAL::IO::Color color(m_uni(m_rng), m_uni(m_rng), m_uni(m_rng));
|
||||||
|
this->face_begin(color);
|
||||||
|
|
||||||
|
const auto* traits = this->m_arr.geometry_traits();
|
||||||
|
auto ext = find_smallest(circ);
|
||||||
|
auto curr = ext;
|
||||||
|
do {
|
||||||
|
// Skip halfedges that are "antenas":
|
||||||
|
while (curr->face() == curr->twin()->face()) curr = curr->twin()->next();
|
||||||
|
|
||||||
|
bh::if_(has_approximate_2_object(Gt{}),
|
||||||
|
[&](auto& x) {
|
||||||
|
auto approx = traits->approximate_2_object();
|
||||||
|
auto has_operator = bh::is_valid(can_call_operator_curve<int>{});
|
||||||
|
bh::if_(has_operator(approx),
|
||||||
|
[&](auto& x) { x.draw_approximate_region(curr, approx); },
|
||||||
|
[&](auto& x) { x.draw_exact_region(curr); }
|
||||||
|
)(*this);
|
||||||
|
},
|
||||||
|
[&](auto& x) { x.draw_exact_region(curr); }
|
||||||
|
)(*this);
|
||||||
|
curr = curr->next();
|
||||||
|
} while (curr != ext);
|
||||||
|
|
||||||
|
this->face_end();
|
||||||
|
}
|
||||||
|
|
||||||
|
//! Draw a curve using aproximate coordinates.
|
||||||
|
// Call this member function only of the geometry traits is equipped to
|
||||||
|
// provide aproximate coordinates.
|
||||||
|
template <typename XMonotoneCurve, typename Approximate>
|
||||||
|
void draw_approximate_curve(const XMonotoneCurve& curve,
|
||||||
|
const Approximate& approx) {
|
||||||
|
std::vector<typename Gt::Approximate_point_2> polyline;
|
||||||
|
double error(this->pixel_ratio());
|
||||||
|
approx(curve, error, std::back_inserter(polyline));
|
||||||
|
if (polyline.empty()) return;
|
||||||
|
auto it = polyline.begin();
|
||||||
|
auto prev = it++;
|
||||||
|
for (; it != polyline.end(); prev = it++) this->add_segment(*prev, *it);
|
||||||
|
}
|
||||||
|
|
||||||
|
//! Draw a curve.
|
||||||
|
template <typename XMonotoneCurve>
|
||||||
|
void draw_curve(const XMonotoneCurve& curve) {
|
||||||
|
/* Check whether the traits has a member function called
|
||||||
|
* approximate_2_object() and if so check whether the return type, namely
|
||||||
|
* `Approximate_2` has an appropriate operator.
|
||||||
|
*
|
||||||
|
* C++20 supports concepts and `requires` expression; see, e.g.,
|
||||||
|
* https://en.cppreference.com/w/cpp/language/constraints; thus, the first
|
||||||
|
* condition above can be elegantly verified as follows:
|
||||||
|
* constexpr bool has_approximate_2_object =
|
||||||
|
* requires(const Gt& traits) { traits.approximate_2_object(); };
|
||||||
|
*
|
||||||
|
* C++17 has experimental constructs called is_detected and
|
||||||
|
* is_detected_v that can be used to achieve the same goal.
|
||||||
|
*
|
||||||
|
* For now we use C++14 features and boost.hana instead.
|
||||||
|
*/
|
||||||
|
#if 0
|
||||||
|
if constexpr (std::experimental::is_detected_v<approximate_2_object_t, Gt>) {
|
||||||
|
const auto* traits = this->m_arr.geometry_traits();
|
||||||
|
auto approx = traits->approximate_2_object();
|
||||||
|
draw_approximate_curve(curve, approx);
|
||||||
|
return;
|
||||||
|
}
|
||||||
|
draw_exact_curve(curve);
|
||||||
|
#else
|
||||||
|
namespace bh = boost::hana;
|
||||||
|
auto has_approximate_2_object =
|
||||||
|
bh::is_valid([](auto&& x) -> decltype(x.approximate_2_object()){});
|
||||||
|
const auto* traits = this->m_arr.geometry_traits();
|
||||||
|
bh::if_(has_approximate_2_object(Gt{}),
|
||||||
|
[&](auto& x) {
|
||||||
|
auto approx = traits->approximate_2_object();
|
||||||
|
auto has_operator = bh::is_valid(can_call_operator_curve<int>{});
|
||||||
|
bh::if_(has_operator(approx),
|
||||||
|
[&](auto& x) { x.draw_approximate_curve(curve, approx); },
|
||||||
|
[&](auto& x) { x.draw_exact_curve(curve); }
|
||||||
|
)(*this);
|
||||||
|
},
|
||||||
|
[&](auto& x) { x.draw_exact_curve(curve); }
|
||||||
|
)(*this);
|
||||||
|
#endif
|
||||||
|
}
|
||||||
|
|
||||||
//!
|
//!
|
||||||
virtual void draw_point(const Point& p) { this->add_point(p); }
|
void draw_point(const Point& p) {
|
||||||
|
namespace bh = boost::hana;
|
||||||
|
auto has_approximate_2_object =
|
||||||
|
bh::is_valid([](auto&& x) -> decltype(x.approximate_2_object()){});
|
||||||
|
bh::if_(has_approximate_2_object(Gt{}),
|
||||||
|
[&](auto& x) {
|
||||||
|
const auto* traits = x.m_arr.geometry_traits();
|
||||||
|
auto approx = traits->approximate_2_object();
|
||||||
|
auto has_operator =
|
||||||
|
bh::is_valid([](auto&& x) -> decltype(x.operator()(Point{})){});
|
||||||
|
bh::if_(has_operator(approx),
|
||||||
|
[&](auto& x) { /* x.add_point(approx(p)) */; },
|
||||||
|
[&](auto& x) { x.add_point(p); }
|
||||||
|
)(*this);
|
||||||
|
},
|
||||||
|
[&](auto& x) { x.add_point(p); }
|
||||||
|
)(*this);
|
||||||
|
}
|
||||||
|
|
||||||
//!
|
//!
|
||||||
void add_ccb(Ccb_halfedge_const_circulator circ) {
|
void add_ccb(Ccb_halfedge_const_circulator circ) {
|
||||||
|
|
@ -280,81 +442,6 @@ public:
|
||||||
{}
|
{}
|
||||||
};
|
};
|
||||||
|
|
||||||
//!
|
|
||||||
template <typename RatKernel, typename AlgKernel, typename NtTraits,
|
|
||||||
typename Dcel>
|
|
||||||
class Arr_2_viewer_qt<Arr_conic_traits_2<RatKernel, AlgKernel, NtTraits>,
|
|
||||||
Dcel> :
|
|
||||||
public Arr_2_basic_viewer_qt<Arr_conic_traits_2<RatKernel, AlgKernel,
|
|
||||||
NtTraits>, Dcel> {
|
|
||||||
public:
|
|
||||||
typedef Arr_conic_traits_2<RatKernel, AlgKernel, NtTraits> Gt;
|
|
||||||
typedef CGAL::Arrangement_2<Gt, Dcel> Arr;
|
|
||||||
typedef Arr_2_basic_viewer_qt<Gt, Dcel> Base;
|
|
||||||
typedef typename Arr::Point_2 Point;
|
|
||||||
typedef typename Arr::X_monotone_curve_2 X_monotone_curve;
|
|
||||||
typedef typename Arr::Halfedge_const_handle Halfedge_const_handle;
|
|
||||||
typedef typename Arr::Face_const_handle Face_const_handle;
|
|
||||||
typedef typename Arr::Ccb_halfedge_const_circulator
|
|
||||||
Ccb_halfedge_const_circulator;
|
|
||||||
|
|
||||||
/// Construct the viewer.
|
|
||||||
/// @param arr the arrangement to view
|
|
||||||
/// @param title the title of the window
|
|
||||||
Arr_2_viewer_qt(QWidget* parent, const Arr& arr,
|
|
||||||
const char* title = "2D Arrangement Basic Viewer") :
|
|
||||||
Base(parent, arr, title)
|
|
||||||
{}
|
|
||||||
|
|
||||||
//!
|
|
||||||
virtual void draw_region(Ccb_halfedge_const_circulator circ) {
|
|
||||||
auto& uni = this->m_uni;
|
|
||||||
auto& rng = this->m_rng;
|
|
||||||
CGAL::IO::Color color(uni(rng), uni(rng), uni(rng));
|
|
||||||
this->face_begin(color);
|
|
||||||
|
|
||||||
const auto* traits = this->m_arr.geometry_traits();
|
|
||||||
auto approx = traits->approximate_2_object();
|
|
||||||
|
|
||||||
// Find the lexicographically smallest halfedge:
|
|
||||||
auto ext = this->find_smallest(circ);
|
|
||||||
|
|
||||||
// Iterate, starting from the lexicographically smallest vertex:
|
|
||||||
auto curr = ext;
|
|
||||||
do {
|
|
||||||
// Skip halfedges that are "antenas":
|
|
||||||
while (curr->face() == curr->twin()->face())
|
|
||||||
curr = curr->twin()->next();
|
|
||||||
|
|
||||||
std::vector<typename Gt::Approximate_point_2> polyline;
|
|
||||||
double error(this->pixel_ratio());
|
|
||||||
bool l2r = curr->direction() == ARR_LEFT_TO_RIGHT;
|
|
||||||
approx(std::back_inserter(polyline), error, curr->curve(), l2r);
|
|
||||||
auto it = polyline.begin();
|
|
||||||
auto prev = it++;
|
|
||||||
for (; it != polyline.end(); prev = it++) {
|
|
||||||
this->add_segment(*prev, *it);
|
|
||||||
this->add_point_in_face(*prev);
|
|
||||||
}
|
|
||||||
curr = curr->next();
|
|
||||||
} while (curr != ext);
|
|
||||||
|
|
||||||
this->face_end();
|
|
||||||
}
|
|
||||||
|
|
||||||
//!
|
|
||||||
virtual void draw_curve(const X_monotone_curve& curve) {
|
|
||||||
const auto* traits = this->m_arr.geometry_traits();
|
|
||||||
auto approx = traits->approximate_2_object();
|
|
||||||
std::vector<typename Gt::Approximate_point_2> polyline;
|
|
||||||
double error(this->pixel_ratio());
|
|
||||||
approx(std::back_inserter(polyline), error, curve);
|
|
||||||
auto it = polyline.begin();
|
|
||||||
auto prev = it++;
|
|
||||||
for (; it != polyline.end(); prev = it++) this->add_segment(*prev, *it);
|
|
||||||
}
|
|
||||||
};
|
|
||||||
|
|
||||||
//!
|
//!
|
||||||
template <typename GeometryTraits_2, typename Dcel>
|
template <typename GeometryTraits_2, typename Dcel>
|
||||||
void draw(const Arrangement_2<GeometryTraits_2, Dcel>& arr,
|
void draw(const Arrangement_2<GeometryTraits_2, Dcel>& arr,
|
||||||
|
|
|
||||||
Loading…
Reference in New Issue