Removed Ron from maintainer list. Updated Ron's email address. Cleaned up documentation and example programs

This commit is contained in:
Efi Fogel 2014-08-28 17:05:57 +03:00
parent 7ce44a2615
commit f48ce38473
26 changed files with 270 additions and 353 deletions

View File

@ -432,7 +432,7 @@ function template computes the offset of a given polygon \f$ P\f$
by a rational radius \f$ r\f$ in an exact manner. The input
polygon $P$ must be either simple or degenerate consisting of two
vertices (representing a line segment). The `traits` argument must model
the concept `ArrangementTraits_2`; it should be capable of handling
the concept `ArrangementTraits_2` and it should be capable of handling
conic arcs in an exact manner---using an instance of the
`Arr_conic_traits_2` class template with the number
types provided by the <span class="textsc">Core</span> library is the
@ -455,69 +455,73 @@ considerably slower:
\cgalExample{Minkowski_sum_2/exact_offset.cpp}
\cgalAdvancedBegin
Both functions `approximated_offset_2()` and `offset_polygon_2()`
Both functions templates `approximated_offset_2()` and `offset_polygon_2()`
also have overloaded versions that accept a decomposition strategy
and use the polygon-decomposition approach to compute (or approximate)
the offset. These functions are less efficient than their counterparts
that employ the convolution approach, and are only included in the package
for the sake of completeness.
the offset. These functions are typically considerably slower than their
counterparts that employ the convolution approach. However, similar to the
functions that compute the general Minkowski sum, they are able to compute
the offset of polygons with holes, given a decomposition strategy that
handles polygons with holes, such as the
'Polygon_vertical_decomposition_2<Kernel>' class template.
\cgalAdvancedEnd
\subsection mink_ssecinner_offset Computing Inner Offsets
An operation closely related to the offset computation, is obtaining the
<I>inner offset</I> of a polygon, or <I>insetting</I> it by a given radius.
The inset of a polygon \f$ P\f$ with radius \f$ r\f$ is the set of points
iside \f$ P\f$ whose distance from the polygon boundary, denoted
\f$ \partial P\f$, is at least \f$ r\f$, namely:
\f$ \{ p \in P \;|\; {\rm dist}(p, \partial P) \geq r \}\f$.
Note that the resulting set may not be connected in case \f$ P\f$ is a
non-convex polygon that has some narrow components, and thus is
characterized by a (possibly empty) set of polygons whose edges are
line segments and circular arcs of radius \f$ r\f$.
An operation closely related to the (outer) offset computation, is
computing the <I>inner offset</I> of a polygon, or <I>insetting</I> it
by a given radius. The inset of a polygon \f$ P\f$ with radius
\f$ r\f$ is the set of points iside \f$ P\f$ the distance of which
from the polygon boundary, denoted \f$ \partial P\f$, is at least \f$ r\f$,
namely: \f$ \{ p \in P \;|\; {\rm dist}(p, \partial P) \geq r \}\f$.
Note that the resulting point set may be dicconnected when \f$ P\f$ is a
non-convex polygon that has some narrow components. In such a case the
resulting set is characterized by a (possibly empty) set of polygons the
edges of which are line segments and circular arcs of radius \f$ r\f$.
The offset can be computed using the convolution method if we traverse the
polygon in a <I>clockwise</I> orientation (and not in <I>counterclockwise</I>
orientation, as done in case of ofsetting a polygon).
As in case of the offset functions, the Minkowski-sum package contains two
functions for insetting a simple polygon:
The inset can be computed using the convolution method traversing the
polygon in a <I>clockwise</I> order (as oppose to the
<I>counterclockwise</I> order applied in the case of ofsetting a polygon).
As with the (outer) offset functions, the Minkowski-sum package contains two
functions for insetting a simple polygon, namely,
`approximated_inset_2(P, r, eps, oi)` and `inset_polygon_2(P, r, traits, oi)`.
The function \link approximated_inset_2() `approximated_inset_2(P, r, eps, oi)`\endlink accepts a polygon
\f$ P\f$, an inset radius \f$ r\f$ and \f$ \varepsilon > 0\f$. It constructs an
approximation for the inset of \f$ P\f$ by the radius \f$ r\f$, with the approximation
error bounded by \f$ \varepsilon\f$. The function returns its output via the
output iterator `oi`, whose value-type must be
`Gps_circle_segment_traits_2::Polygon_2` representing
the polygons that approximates the inset polygon.
The \link approximated_inset_2() `approximated_inset_2(P, r, eps, oi)`\endlink
function template accepts a polygon \f$ P\f$, an inset radius \f$ r\f$,
(a floating-point number) \f$ \epsilon > 0 \f$, and an output
iterator `oi`, whose value-type must be an instance of the class template
`Gps_circle_segment_traits_2::Polygon_2`. . It constructs an approximation
of the inset of \f$ P\f$ by the radius \f$ r\f$, where the approximation
error is bounded by \f$ \epsilon\f$. The function returns the polygons that
approximate the inset polygon through the output iterator `oi`.
\cgalExample{Minkowski_sum_2/approx_inset.cpp}
Similarly, the function
\link inset_polygon_2() `inset_polygon_2(P, r, traits, oi)` \endlink computes
the exact inset of \f$ P\f$ with radius \f$ r\f$, and returns its output via the given
output iterator `oi`. The `traits` parameter should be an instance of
an arrangement-traits class that is capable of handling conic arcs in an
exact manner, whereas `oi`'s value-type must be
`Gps_traits_2::Polygon_2`.
Similarly, the function template
\link inset_polygon_2() `inset_polygon_2(P, r, traits, oi)` \endlink
computes the exact inset of \f$ P\f$ with radius \f$ r\f$, and returns
its output through the given output iterator `oi`. The `traits` parameter
must model the concept `ArrangementTraits_2` and it should be capable of
handling conic arcs in an exact manner, whereas the value-type of `oi`
must be an instance of `Gps_traits_2::Polygon_2`.
\cgalExample{Minkowski_sum_2/exact_inset.cpp}
\cgalAdvancedBegin
Unlike the offset functions, there are no overloaded versions of the inset
functions that use convex polygon decomposition to compute insets, as this
method cannot be easily generalized for inset computations.
\cgalAdvancedEnd
Unlike the functions that compute the offsets, there are no overloaded
versions of the functions that compute the insets and use convex polygon
decomposition, as decomposition approach cannot be easily generalized for
inset computations.
In this context let us mention that there exist overloaded versions of the
functions
The package also provides overloaded versions of the functions
\link approximated_offset_2() `approximated_offset_2(P, r, eps)`\endlink
and
\link offset_polygon_2() `offset_polygon_2(P, r, traits)`\endlink
that accept a <I>polygon with holes</I>
\f$ P\f$ and computed its offset. This ofset is obtain by taking the outer offset
of \f$ P\f$'s outer boundary, and computing the inner offsets of \f$ P\f$'s holes.
The former polygon defines the output boundary of \f$ P \oplus B_r\f$, and the latter
define the holes within the result.
\f$ P\f$ and compute its offset. This ofset is obtain by computing the
outer offset of the outer boundary of \f$ P\f$'s, and computing the inner
offsets of the holes of \f$ P\f$. The former polygon defines the output
boundary of \f$ P \oplus B_r\f$, and the latter define the holes within
the result.
\section Minkowski_sum_2Acknowledgements Acknowledgements

View File

@ -1,63 +1,48 @@
//! \file examples/Minkowski_sum_2/approx_inset.cpp
// Computing the approximated inset of a polygon.
#include "ms_rational_nt.h"
#include <CGAL/Lazy_exact_nt.h>
#include <CGAL/Cartesian.h>
#include <CGAL/approximated_offset_2.h>
#include <CGAL/offset_polygon_2.h>
#include <CGAL/Timer.h>
#include <fstream>
#include <iostream>
#include <list>
#include <boost/timer.hpp>
typedef CGAL::Lazy_exact_nt<Number_type> Lazy_exact_nt;
#include <CGAL/basic.h>
#include <CGAL/approximated_offset_2.h>
struct Kernel : public CGAL::Cartesian<Lazy_exact_nt> {};
typedef CGAL::Polygon_2<Kernel> Polygon_2;
#include "bops_circular.h"
typedef CGAL::Gps_circle_segment_traits_2<Kernel> Gps_traits_2;
typedef Gps_traits_2::Polygon_2 Offset_polygon_2;
typedef Gps_traits_2::Polygon_with_holes_2 Offset_polygon_with_holes_2;
typedef CGAL::Polygon_2<Kernel> Linear_polygon;
int main ()
int main(int argc, char* argv[])
{
// Open the input file.
std::ifstream in_file ("tight.dat");
// Open the input file and read a polygon.
const char* filename = (argc > 1) ? argv[1] : "tight.dat";
std::ifstream in_file(filename);
if (! in_file.is_open())
{
if (! in_file.is_open()) {
std::cerr << "Failed to open the input file." << std::endl;
return (1);
return -1;
}
// Read the input polygon.
Polygon_2 P;
Linear_polygon P;
in_file >> P;
in_file.close();
std::cout << "Read an input polygon with "
<< P.size() << " vertices." << std::endl;
std::cout << "Read an input polygon with " << P.size() << " vertices."
<< std::endl;
// Approximate the offset polygon.
const Number_type radius = 1;
const double err_bound = 0.00001;
std::list<Offset_polygon_2> inset_polygons;
std::list<Offset_polygon_2>::iterator iit;
CGAL::Timer timer;
std::list<Polygon> inset_polygons;
boost::timer timer;
approximated_inset_2(P, 1, 0.00001, std::back_inserter(inset_polygons));
double secs = timer.elapsed();
timer.start();
approximated_inset_2 (P, radius, err_bound,
std::back_inserter (inset_polygons));
timer.stop();
std::cout << "The inset comprises "
<< inset_polygons.size() << " polygon(s)." << std::endl;
for (iit = inset_polygons.begin(); iit != inset_polygons.end(); ++iit)
{
std::cout << " Polygon with "
<< iit->size() << " vertices." << std::endl;
}
std::cout << "Inset computation took "
<< timer.time() << " seconds." << std::endl;
return (0);
std::list<Polygon>::iterator it;
std::cout << "The inset comprises " << inset_polygons.size()
<< " polygon(s)." << std::endl;
for (it = inset_polygons.begin(); it != inset_polygons.end(); ++it)
std::cout << " Polygon with " << it->size() << " vertices." << std::endl;
std::cout << "Inset computation took " << secs << " seconds." << std::endl;
return 0;
}

View File

@ -15,7 +15,7 @@ int main(int argc, char* argv[])
{
// Open the input file and read a polygon.
const char* filename = (argc > 1) ? argv[1] : "spiked.dat";
std::ifstream in_file(filename);
std::ifstream in_file(filename);
if (! in_file.is_open()) {
std::cerr << "Failed to open the input file." << std::endl;
return -1;

View File

@ -0,0 +1,25 @@
#ifndef ARR_CONICS_H
#define ARR_CONICS_H
#include <CGAL/Cartesian.h>
#include <CGAL/CORE_algebraic_number_traits.h>
#include <CGAL/Arr_conic_traits_2.h>
#include <CGAL/Arrangement_2.h>
typedef CGAL::CORE_algebraic_number_traits Nt_traits;
typedef Nt_traits::Rational Rational;
typedef CGAL::Cartesian<Rational> Rat_kernel;
typedef Rat_kernel::Point_2 Rat_point;
typedef Rat_kernel::Segment_2 Rat_segment;
typedef Rat_kernel::Circle_2 Rat_circle;
typedef Nt_traits::Algebraic Algebraic;
typedef CGAL::Cartesian<Algebraic> Alg_kernel;
typedef CGAL::Arr_conic_traits_2<Rat_kernel, Alg_kernel, Nt_traits>
Traits;
typedef Traits::Point_2 Point;
typedef Traits::Curve_2 Conic_arc;
typedef Traits::X_monotone_curve_2 X_monotone_conic_arc;
typedef CGAL::Arrangement_2<Traits> Arrangement;
#endif

View File

@ -1,80 +1,60 @@
//! \file examples/Minkowski_sum_2/exact_inset.cpp
// Computing the exact inner offset of a polygon.
#include <iostream>
#include <CGAL/basic.h>
#ifndef CGAL_USE_CORE
#include <iostream>
int main ()
int main()
{
std::cout << "Sorry, this example needs CORE ..." << std::endl;
return (0);
return 0;
}
#else
#include <CGAL/Cartesian.h>
#include <CGAL/CORE_algebraic_number_traits.h>
#include <CGAL/Arr_conic_traits_2.h>
#include <fstream>
#include <boost/timer.hpp>
#include <CGAL/Gps_traits_2.h>
#include <CGAL/offset_polygon_2.h>
#include <CGAL/Timer.h>
#include <iostream>
typedef CGAL::CORE_algebraic_number_traits Nt_traits;
typedef Nt_traits::Rational Rational;
typedef Nt_traits::Algebraic Algebraic;
#include "arr_conics.h"
struct Rat_kernel : public CGAL::Cartesian<Rational> {};
struct Alg_kernel : public CGAL::Cartesian<Algebraic> {};
struct Conic_traits_2 : public CGAL::Arr_conic_traits_2<Rat_kernel,
Alg_kernel,
Nt_traits> {};
typedef CGAL::Polygon_2<Rat_kernel> Polygon;
typedef CGAL::Gps_traits_2<Traits> Gps_traits;
typedef Gps_traits::Polygon_2 Offset_polygon;
typedef CGAL::Polygon_2<Rat_kernel> Polygon_2;
typedef CGAL::Gps_traits_2<Conic_traits_2> Gps_traits_2;
typedef Gps_traits_2::Polygon_2 Offset_polygon_2;
int main ()
int main(int argc, char* argv[])
{
// Open the input file.
std::ifstream in_file ("tight.dat");
if (! in_file.is_open())
{
// Open the input file and read the input polygon.
const char* filename = (argc > 1) ? argv[1] : "tight.dat";
std::ifstream in_file(filename);
if (! in_file.is_open()) {
std::cerr << "Failed to open the input file." << std::endl;
return (1);
return -1;
}
// Read the input polygon.
Polygon_2 P;
Polygon P;
in_file >> P;
in_file.close();
std::cout << "Read an input polygon with "
<< P.size() << " vertices." << std::endl;
std::cout << "Read an input polygon with " << P.size() << " vertices."
<< std::endl;
// Compute the inner offset of the polygon.
Conic_traits_2 traits;
const Rational radius = 1;
std::list<Offset_polygon_2> inset_polygons;
std::list<Offset_polygon_2>::iterator iit;
CGAL::Timer timer;
Traits traits;
std::list<Offset_polygon> inset_polygons;
boost::timer timer;
inset_polygon_2(P, 1, traits, std::back_inserter(inset_polygons));
double secs = timer.elapsed();
timer.start();
inset_polygon_2 (P, radius, traits,
std::back_inserter (inset_polygons));
timer.stop();
std::cout << "The inset comprises "
std::list<Offset_polygon>::iterator it;
std::cout << "The inset comprises "
<< inset_polygons.size() << " polygon(s)." << std::endl;
for (iit = inset_polygons.begin(); iit != inset_polygons.end(); ++iit)
{
std::cout << " Polygon with "
<< iit->size() << " vertices." << std::endl;
}
std::cout << "Inset computation took "
<< timer.time() << " seconds." << std::endl;
return (0);
for (it = inset_polygons.begin(); it != inset_polygons.end(); ++it)
std::cout << " Polygon with " << it->size() << " vertices."
<< std::endl;
std::cout << "Inset computation took " << secs << " seconds." << std::endl;
return 0;
}
#endif

View File

@ -1,75 +1,56 @@
//! \file examples/Minkowski_sum_2/exact_offset.cpp
// Computing the exact offset of a polygon.
#include <iostream>
#include <CGAL/basic.h>
#ifndef CGAL_USE_CORE
#include <iostream>
int main ()
int main()
{
std::cout << "Sorry, this example needs CORE ..." << std::endl;
return (0);
return 0;
}
#else
#include <CGAL/Cartesian.h>
#include <CGAL/CORE_algebraic_number_traits.h>
#include <CGAL/Arr_conic_traits_2.h>
#include <fstream>
#include <boost/timer.hpp>
#include <CGAL/Gps_traits_2.h>
#include <CGAL/offset_polygon_2.h>
#include <CGAL/Timer.h>
#include <iostream>
typedef CGAL::CORE_algebraic_number_traits Nt_traits;
typedef Nt_traits::Rational Rational;
typedef Nt_traits::Algebraic Algebraic;
#include "arr_conics.h"
struct Rat_kernel : public CGAL::Cartesian<Rational> {};
struct Alg_kernel : public CGAL::Cartesian<Algebraic> {};
struct Conic_traits_2 : public CGAL::Arr_conic_traits_2<Rat_kernel,
Alg_kernel,
Nt_traits> {};
typedef CGAL::Polygon_2<Rat_kernel> Polygon;
typedef CGAL::Gps_traits_2<Traits> Gps_traits;
typedef Gps_traits::Polygon_with_holes_2 Offset_polygon_with_holes;
typedef CGAL::Polygon_2<Rat_kernel> Polygon_2;
typedef CGAL::Gps_traits_2<Conic_traits_2> Gps_traits_2;
typedef Gps_traits_2::Polygon_2 Offset_polygon_2;
typedef Gps_traits_2::Polygon_with_holes_2 Offset_polygon_with_holes_2;
int main ()
int main(int argc, char* argv[])
{
// Open the input file.
std::ifstream in_file ("spiked.dat");
if (! in_file.is_open())
{
// Open the input file and read the input polygon.
const char* filename = (argc > 1) ? argv[1] : "spiked.dat";
std::ifstream in_file(filename);
if (! in_file.is_open()) {
std::cerr << "Failed to open the input file." << std::endl;
return (1);
return -1;
}
// Read the input polygon.
Polygon_2 P;
Polygon P;
in_file >> P;
in_file.close();
std::cout << "Read an input polygon with "
<< P.size() << " vertices." << std::endl;
std::cout << "Read an input polygon with " << P.size() << " vertices."
<< std::endl;
// Compute the offset polygon.
Conic_traits_2 traits;
const Rational radius = 5;
Offset_polygon_with_holes_2 offset;
CGAL::Timer timer;
Traits traits;
boost::timer timer;
Offset_polygon_with_holes offset = CGAL::offset_polygon_2(P, 5, traits);
double secs = timer.elapsed();
timer.start();
offset = offset_polygon_2 (P, radius, traits);
timer.stop();
std::cout << "The offset polygon has "
<< offset.outer_boundary().size() << " vertices, "
<< offset.number_of_holes() << " holes." << std::endl;
std::cout << "Offset computation took "
<< timer.time() << " seconds." << std::endl;
return (0);
std::cout << "The offset polygon has " << offset.outer_boundary().size()
<< " vertices, " << offset.number_of_holes() << " holes."
<< std::endl;
std::cout << "Offset computation took " << secs << " seconds." << std::endl;
return 0;
}
#endif

View File

@ -12,10 +12,7 @@
// This file is provided AS IS with NO WARRANTY OF ANY KIND, INCLUDING THE
// WARRANTY OF DESIGN, MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE.
//
// $URL$
// $Id$
//
// Author(s) : Ron Wein <wein@post.tau.ac.il>
// Author(s) : Ron Wein <wein_r@yahoo.com>
// Andreas Fabri <Andreas.Fabri@geometryfactory.com>
// Laurent Rineau <Laurent.Rineau@geometryfactory.com>
// Efi Fogel <efif@post.tau.ac.il>
@ -52,7 +49,7 @@ protected:
typedef CGAL::Polygon_with_holes_2<Kernel, Container_> Polygon_with_holes_2;
private:
// Kernel types:
typedef typename Kernel::Point_2 Point_2;
typedef typename Kernel::Line_2 Line_2;
@ -70,8 +67,8 @@ private:
protected:
typedef typename Traits_2::Polygon_2 Offset_polygon_2;
typedef Arr_labeled_traits_2<Traits_2> Labeled_traits_2;
typedef Arr_labeled_traits_2<Traits_2> Labeled_traits_2;
typedef typename Labeled_traits_2::X_monotone_curve_2 Labeled_curve_2;
// Data members:
@ -92,7 +89,7 @@ public:
_inv_sqrt_eps = static_cast<int> (1.0 / CGAL::sqrt (_eps));
if (_inv_sqrt_eps <= 0)
_inv_sqrt_eps = 1;
}
}
protected:
@ -119,7 +116,7 @@ protected:
Vertex_circulator first, curr, next, prev;
first = pgn.vertices_circulator();
curr = first;
curr = first;
next = first;
prev = first;
@ -312,19 +309,19 @@ protected:
// In case of overflow of denom
numer = 1;
denom = max_int;
}
}
}
}
else {// if numer < 0 (overflow)
else {// if numer < 0 (overflow)
numer = max_int;
denom = 1;
}
app_d = NT (numer) / NT (denom);
app_err = sqr_d - CGAL::square (app_d);
app_err = sqr_d - CGAL::square (app_d);
while (CGAL::compare (CGAL::abs (app_err),
while (CGAL::compare (CGAL::abs (app_err),
err_bound) == CGAL::LARGER ||
CGAL::compare (app_d, abs_delta_x) != LARGER ||
CGAL::compare (app_d, abs_delta_y) != LARGER)
@ -452,7 +449,7 @@ protected:
// Compute the line l1 tangent to the circle centered at (x1, y1)
// with radius r at the approximated point op1.
l1 = f_perp_line (f_line (*curr, op1), op1);
// Compute the line l2 tangent to the circle centered at (x2, y2)
// with radius r at the approximated point op2.
l2 = f_perp_line (f_line (*next, op2), op2);
@ -460,7 +457,7 @@ protected:
// Intersect the two lines. The intersection point serves as a common
// end point for the two line segments we are about to introduce.
obj = f_intersect (l1, l2);
assign_success = CGAL::assign (mid_p, obj);
CGAL_assertion (assign_success);
CGAL_USE(assign_success);
@ -508,7 +505,7 @@ protected:
}
if (res1 != res2) orient = CGAL::LEFT_TURN;
}
// Connect the offset target point of the previous edge to the
// offset source of the current edge.
if (orient == CGAL::LEFT_TURN) {
@ -557,12 +554,12 @@ protected:
*oi++ = Labeled_curve_2 (seg2, X_curve_label (dir_right2,
cycle_id, curve_index++));
}
// Proceed to the next polygon vertex.
prev_op = op2;
prev = curr;
curr = next;
} while (curr != first);
// Close the convolution cycle by creating the final circular arc,
@ -582,7 +579,7 @@ protected:
assign_success = CGAL::assign (xarc, *xobj_it);
CGAL_assertion (assign_success);
CGAL_USE(assign_success);
++xobj_it;
bool is_last = (xobj_it == xobjs.end());

View File

@ -12,10 +12,7 @@
// This file is provided AS IS with NO WARRANTY OF ANY KIND, INCLUDING THE
// WARRANTY OF DESIGN, MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE.
//
// $URL$
// $Id$
//
// Author(s) : Ron Wein <wein@post.tau.ac.il>
// Author(s) : Ron Wein <wein_r@yahoo.com>
#ifndef CGAL_ARR_LABELED_TRAITS_2_H
#define CGAL_ARR_LABELED_TRAITS_2_H
@ -29,7 +26,7 @@ namespace CGAL {
* A meta-traits class that adds lables to points and to x-monotone curves,
* such that the comparison of two points, as well as the computation of the
* intersections between two segments can be easily filtered.
*/
*/
template <class Traits_>
class Arr_labeled_traits_2 : public Traits_
{
@ -97,7 +94,7 @@ public:
{}
/*! Constructor from an x-monotone curve an a label. */
X_monotone_curve_2 (const Base_x_monotone_curve_2& p,
X_monotone_curve_2 (const Base_x_monotone_curve_2& p,
const X_curve_label& label) :
Base_x_monotone_curve_2 (p),
_label (label)
@ -133,7 +130,7 @@ public:
// Inherited functors:
typedef typename Base_traits_2::Is_vertical_2 Is_vertical_2;
typedef typename Base_traits_2::Compare_y_at_x_2 Compare_y_at_x_2;
typedef typename Base_traits_2::Compare_y_at_x_right_2
typedef typename Base_traits_2::Compare_y_at_x_right_2
Compare_y_at_x_right_2;
typedef typename Base_traits_2::Equal_2 Equal_2;
@ -235,7 +232,7 @@ public:
else if (cv.label().right_count() == 0 && cv.label().left_count() == 1)
{
// A curve directed from right to left:
Point_label label (cv.label().component(),
Point_label label (cv.label().component(),
cv.label().is_last() ? 0 : cv.label().index()+1);
return (Point_2 (pt, label));
@ -370,7 +367,7 @@ public:
std::list<CGAL::Object> base_objs;
base->intersect_2_object() (cv1, cv2, std::back_inserter (base_objs));
if (base_objs.empty())
return (oi);
@ -387,7 +384,7 @@ public:
if (base_pt != NULL)
{
// Attach an invalid label to an itersection point.
*oi = CGAL::make_object
*oi = CGAL::make_object
(std::make_pair (Point_2 (base_pt->first), base_pt->second));
++oi;
}
@ -397,8 +394,8 @@ public:
CGAL_assertion (base_xcv != NULL);
// Attach a merged label to the overlapping curve.
*oi = CGAL::make_object
(X_monotone_curve_2 (*base_xcv,
*oi = CGAL::make_object
(X_monotone_curve_2 (*base_xcv,
X_curve_label (cv1.label(), cv2.label())));
++oi;
}

View File

@ -12,10 +12,7 @@
// This file is provided AS IS with NO WARRANTY OF ANY KIND, INCLUDING THE
// WARRANTY OF DESIGN, MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE.
//
// $URL$
// $Id$
//
// Author(s) : Ron Wein <wein@post.tau.ac.il>
// Author(s) : Ron Wein <wein_r@yahoo.com>
// Efi Fogel <efifogel@gmail.com>
#ifndef CGAL_POLYGON_DECOMPOSITION_STRATEGY_ADAPTER_H

View File

@ -12,10 +12,7 @@
// This file is provided AS IS with NO WARRANTY OF ANY KIND, INCLUDING THE
// WARRANTY OF DESIGN, MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE.
//
// $URL$
// $Id$
//
// Author(s) : Ron Wein <wein@post.tau.ac.il>
// Author(s) : Ron Wein <wein_r@yahoo.com>
#ifndef CGAL_EXACT_OFFSET_BASE_H
#define CGAL_EXACT_OFFSET_BASE_H
@ -37,9 +34,9 @@ template <class Traits_, class Container_>
class Exact_offset_base_2
{
private:
typedef Traits_ Traits_2;
// Rational kernel types:
typedef typename Traits_2::Rat_kernel Rat_kernel;
typedef typename Rat_kernel::FT Rational;
@ -66,14 +63,14 @@ private:
typedef typename Traits_2::X_monotone_curve_2 X_monotone_curve_2;
typedef CGAL::Gps_traits_2<Traits_2> Gps_traits_2;
protected:
typedef CGAL::Polygon_2<Rat_kernel, Container_> Polygon_2;
typedef CGAL::Polygon_with_holes_2<Rat_kernel,
Container_> Polygon_with_holes_2;
typedef typename Gps_traits_2::Polygon_2 Offset_polygon_2;
private:
// Polygon-related types:
@ -81,7 +78,7 @@ private:
protected:
typedef Arr_labeled_traits_2<Traits_2> Labeled_traits_2;
typedef Arr_labeled_traits_2<Traits_2> Labeled_traits_2;
typedef typename Labeled_traits_2::X_monotone_curve_2 Labeled_curve_2;
@ -89,7 +86,7 @@ public:
/*! Default constructor. */
Exact_offset_base_2 ()
{}
{}
protected:
@ -114,9 +111,9 @@ protected:
// Prepare circulators over the polygon vertices.
const bool forward = (pgn.orientation() == orient);
Vertex_circulator first, curr, next;
first = pgn.vertices_circulator();
curr = first;
curr = first;
next = first;
// Traverse the polygon vertices and edges and construct the arcs that
@ -163,7 +160,7 @@ protected:
delta_x = x2 - x1;
delta_y = y2 - y1;
len = nt_traits.sqrt (nt_traits.convert (CGAL::square (delta_x) +
len = nt_traits.sqrt (nt_traits.convert (CGAL::square (delta_x) +
CGAL::square (delta_y)));
// The angle theta between the vector v and the x-axis is given by:
@ -184,7 +181,7 @@ protected:
// Construct the first offset vertex, which corresponds to the
// source vertex of the current polygon edge.
op1 = Alg_point_2 (nt_traits.convert (x1) + trans_x,
op1 = Alg_point_2 (nt_traits.convert (x1) + trans_x,
nt_traits.convert (y1) + trans_y);
if (curr == first)
@ -201,12 +198,12 @@ protected:
arc = Curve_2 (Rat_circle_2 (*curr, sqr_r),
CGAL::COUNTERCLOCKWISE,
op2, op1);
// Subdivide the arc into x-monotone subarcs and append them to the
// convolution cycle.
xobjs.clear();
f_make_x_monotone (arc, std::back_inserter(xobjs));
for (xobj_it = xobjs.begin(); xobj_it != xobjs.end(); ++xobj_it)
{
assign_success = CGAL::assign (xarc, *xobj_it);
@ -225,7 +222,7 @@ protected:
// Construct the second offset vertex, which corresponds to the
// target vertex of the current polygon edge.
op2 = Alg_point_2 (nt_traits.convert (x2) + trans_x,
op2 = Alg_point_2 (nt_traits.convert (x2) + trans_x,
nt_traits.convert (y2) + trans_y);
// The equation of the line connecting op1 and op2 is given by:
@ -248,7 +245,7 @@ protected:
// Proceed to the next polygon vertex.
curr = next;
} while (curr != first);
if (! f_equal (op2, first_op))
@ -258,21 +255,21 @@ protected:
arc = Curve_2 (Rat_circle_2 (*first, sqr_r),
CGAL::COUNTERCLOCKWISE,
op2, first_op);
// Subdivide the arc into x-monotone subarcs and append them to the
// convolution cycle.
bool is_last;
xobjs.clear();
f_make_x_monotone (arc, std::back_inserter(xobjs));
xobj_it = xobjs.begin();
while (xobj_it != xobjs.end())
{
assign_success = CGAL::assign (xarc, *xobj_it);
CGAL_assertion (assign_success);
CGAL_USE(assign_success);
++xobj_it;
is_last = (xobj_it == xobjs.end());

View File

@ -12,10 +12,7 @@
// This file is provided AS IS with NO WARRANTY OF ANY KIND, INCLUDING THE
// WARRANTY OF DESIGN, MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE.
//
// $URL$
// $Id$
//
// Author(s) : Ron Wein <wein@post.tau.ac.il>
// Author(s) : Ron Wein <wein_r@yahoo.com>
#ifndef CGAL_MINKOWSKI_SUM_LABELS_H
#define CGAL_MINKOWSKI_SUM_LABELS_H
@ -46,7 +43,7 @@ public:
* \param component The index of the component.
* \param index Index of the point within the component.
*/
Point_label (unsigned int component,
Point_label (unsigned int component,
unsigned int index) :
_component (component),
_index (index)
@ -86,7 +83,7 @@ public:
void set_component (unsigned int component)
{
CGAL_precondition (component != 0);
_component = component;
return;
}
@ -134,7 +131,7 @@ public:
* \param is_last Is this the last curve of the component.
*/
X_curve_label (bool is_directed_right,
unsigned int component,
unsigned int component,
unsigned int index,
bool is_last = false) :
_component (component),
@ -252,13 +249,13 @@ public:
void set_component (unsigned int component)
{
CGAL_precondition (component != 0);
_component = component;
return;
}
/*!
* Set the curve index within the component.
* Set the curve index within the component.
* \param index The index in the component.
* \param is_last Is this the last curve of the component.
*/

View File

@ -12,10 +12,7 @@
// This file is provided AS IS with NO WARRANTY OF ANY KIND, INCLUDING THE
// WARRANTY OF DESIGN, MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE.
//
// $URL$
// $Id$
//
// Author(s) : Ron Wein <wein@post.tau.ac.il>
// Author(s) : Ron Wein <wein_r@yahoo.com>
// Efi Fogel <efifogel@gmail.com>
#ifndef CGAL_MINKOWSKI_SUM_CONV_H

View File

@ -12,10 +12,7 @@
// This file is provided AS IS with NO WARRANTY OF ANY KIND, INCLUDING THE
// WARRANTY OF DESIGN, MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE.
//
// $URL$
// $Id$
//
// Author(s) : Ron Wein <wein@post.tau.ac.il>
// Author(s) : Ron Wein <wein_r@yahoo.com>
#ifndef CGAL_MINKOWSKI_SUM_DECOMP_2_H
#define CGAL_MINKOWSKI_SUM_DECOMP_2_H

View File

@ -12,10 +12,7 @@
// This file is provided AS IS with NO WARRANTY OF ANY KIND, INCLUDING THE
// WARRANTY OF DESIGN, MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE.
//
// $URL$
// $Id$
//
// Author(s) : Ron Wein <wein@post.tau.ac.il>
// Author(s) : Ron Wein <wein_r@yahoo.com>
#ifndef CGAL_OFFSET_CONV_H
#define CGAL_OFFSET_CONV_H
@ -43,10 +40,10 @@ public:
typedef typename Base::Polygon_2 Polygon_2;
typedef typename Base::Polygon_with_holes_2 Polygon_with_holes_2;
typedef typename Base::Offset_polygon_2 Offset_polygon_2;
private:
typedef typename Base::Labeled_traits_2 Labeled_traits_2;
typedef typename Base::Labeled_traits_2 Labeled_traits_2;
typedef typename Base::Labeled_curve_2 Labeled_curve_2;
typedef std::list<Labeled_curve_2> Curves_list;
@ -59,11 +56,11 @@ public:
/*! Constructor. */
Offset_by_convolution_2 (const Base_& base) :
Base (base)
{}
{}
/*!
* Compute the offset of a simple polygon by a given radius.
* Note that as the input polygon may not be convex, its offset may not be
* Note that as the input polygon may not be convex, its offset may not be
* simply connected. The result is therefore represented as the outer
* boundary of the Minkowski sum (which is always a simple offset polygon)
* and a container of offset polygons, representing the holes in this "outer"
@ -126,7 +123,7 @@ public:
// that forms the outer boundary.
Curves_list cycle;
unsigned int cycle_id = 1;
_offset_polygon (pwh.outer_boundary(),
CGAL::COUNTERCLOCKWISE,
r,
@ -159,7 +156,7 @@ public:
/*!
* Compute the inset of a simple polygon by a given radius.
* Note that as the input polygon may not be convex, its offset may not be
* Note that as the input polygon may not be convex, its offset may not be
* simply connected. The result is therefore represented as a sequence of
* polygons (which may also be empty).
* \param pgn The polygon.

View File

@ -12,10 +12,7 @@
// This file is provided AS IS with NO WARRANTY OF ANY KIND, INCLUDING THE
// WARRANTY OF DESIGN, MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE.
//
// $URL$
// $Id$
//
// Author(s) : Ron Wein <wein@post.tau.ac.il>
// Author(s) : Ron Wein <wein_r@yahoo.com>
#ifndef CGAL_OFFSET_DECOMP_H
#define CGAL_OFFSET_DECOMP_H
@ -43,16 +40,16 @@ public:
typedef typename Base::Basic_kernel Kernel;
typedef typename Base::Basic_NT NT;
typedef typename Base::Polygon_2 Polygon_2;
typedef typename Base::Polygon_2 Polygon_2;
typedef typename Base::Offset_polygon_2 Offset_polygon_2;
typedef DecompStrategy_ Decomposition_strategy;
private:
typedef std::list<Polygon_2> Polygons_list;
typedef typename Polygons_list::iterator Polygons_iterator;
typedef typename Base::Labeled_traits_2 Labeled_traits_2;
typedef typename Base::Labeled_traits_2 Labeled_traits_2;
typedef typename Base::Labeled_curve_2 Labeled_curve_2;
typedef std::list<Labeled_curve_2> Curves_list;
@ -64,16 +61,16 @@ public:
/*! Constructor. */
Offset_by_decomposition_2 (const Base_& base) :
Base (base)
{}
{}
/*!
* Compute the offset of a simple polygon by a given radius.
* Note that as the input polygon may not be convex, its offset may not be
* Note that as the input polygon may not be convex, its offset may not be
* simply connected. The result is therefore represented as the outer
* boundary of the Minkowski sum (which is always a simple offset polygon)
* and a container of offset polygons, representing the holes in this "outer"
* polygon.
* \param traits Arrangement traits that can deal with line segments and
* \param traits Arrangement traits that can deal with line segments and
* circular arcs.
* \param pgn The polygon.
* \param r The offset radius.

View File

@ -12,10 +12,7 @@
// This file is provided AS IS with NO WARRANTY OF ANY KIND, INCLUDING THE
// WARRANTY OF DESIGN, MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE.
//
// $URL$
// $Id$
//
// Author(s) : Ron Wein <wein@post.tau.ac.il>
// Author(s) : Ron Wein <wein_r@yahoo.com>
#ifndef CGAL_DECOMPOSITION_STRATEGY_ADAPTER_H
#define CGAL_DECOMPOSITION_STRATEGY_ADAPTER_H
@ -48,7 +45,7 @@ template <class Kernel_, class Container_, class StrategyTag_>
class Decomposition_strategy_adapter
{
public:
typedef Kernel_ Kernel;
typedef Polygon_2<Kernel, Container_> Polygon_2;
typedef StrategyTag_ Strategy_tag;

View File

@ -12,10 +12,7 @@
// This file is provided AS IS with NO WARRANTY OF ANY KIND, INCLUDING THE
// WARRANTY OF DESIGN, MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE.
//
// $URL$
// $Id$
//
// Author(s) : Ron Wein <wein@post.tau.ac.il>
// Author(s) : Ron Wein <wein_r@yahoo.com>
#ifndef CGAL_UNION_OF_CURVE_CYCLES_2_H
#define CGAL_UNION_OF_CURVE_CYCLES_2_H
@ -151,7 +148,7 @@ public:
return (temp);
}
};
/*! Default constructor. */
@ -179,7 +176,7 @@ public:
this->_construct_arrangement (begin, end, arr);
// Produce the result. First set the outer boundary of the union, given
// as the inner boundary of the single hole in the unbounded face.
// as the inner boundary of the single hole in the unbounded face.
Face_iterator fit;
const Face_handle uf = arr.unbounded_face();
Inner_ccb_iterator iccb_it = uf->inner_ccbs_begin();
@ -193,7 +190,7 @@ public:
for (fit = arr.faces_begin(); fit != arr.faces_end(); ++fit)
{
CGAL_assertion (fit->data() != this->UNVISITED);
// If a bounded face has an inside count that equals 0, it forms a hole
// in the union.
if (! fit->is_unbounded() && fit->data() == 0)

View File

@ -12,10 +12,7 @@
// This file is provided AS IS with NO WARRANTY OF ANY KIND, INCLUDING THE
// WARRANTY OF DESIGN, MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE.
//
// $URL$
// $Id$
//
// Author(s) : Ron Wein <wein@post.tau.ac.il>
// Author(s) : Ron Wein <wein_r@yahoo.com>
#ifndef CGAL_UNION_OF_CYCLES_2_H
#define CGAL_UNION_OF_CYCLES_2_H
@ -111,7 +108,7 @@ protected:
{
he = circ;
f_next = he->twin()->face();
if (f_next->data() == UNVISITED)
{
next_count = _boundary_count (he->twin());
@ -124,9 +121,9 @@ protected:
}
++circ;
} while (circ != first);
++iccb_it;
// Make sure that there is a single hole in the unbounded face.
@ -158,11 +155,11 @@ protected:
}
else if (f_curr != f_next)
{
CGAL_assertion (f_next->data() ==
CGAL_assertion (f_next->data() ==
curr_count + _boundary_count (he->twin()));
}
++circ;
} while (circ != first);
// Go over the holes (inner CCBs) of the current face.
@ -184,7 +181,7 @@ protected:
}
else if (f_curr != f_next)
{
CGAL_assertion (f_next->data() ==
CGAL_assertion (f_next->data() ==
curr_count + _boundary_count (he->twin()));
}
++circ;
@ -192,7 +189,7 @@ protected:
} while (circ != first);
}
}
return;
}
@ -208,13 +205,13 @@ private:
if ((Arr_halfedge_direction) he->direction() == ARR_LEFT_TO_RIGHT)
{
// Halfedge is directed from left to right:
return (he->curve().label().right_count() -
return (he->curve().label().right_count() -
he->curve().label().left_count());
}
else
{
// Halfedge is directed from right to left:
return (he->curve().label().left_count() -
return (he->curve().label().left_count() -
he->curve().label().right_count());
}
}

View File

@ -12,10 +12,7 @@
// This file is provided AS IS with NO WARRANTY OF ANY KIND, INCLUDING THE
// WARRANTY OF DESIGN, MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE.
//
// $URL$
// $Id$
//
// Author(s) : Ron Wein <wein@post.tau.ac.il>
// Author(s) : Ron Wein <wein_r@yahoo.com>
#ifndef CGAL_UNION_OF_SEGMENT_CYCLES_2_H
#define CGAL_UNION_OF_SEGMENT_CYCLES_2_H
@ -82,13 +79,13 @@ public:
this->_construct_arrangement (begin, end, arr);
// Produce the result. First set the outer boundary of the union, given
// as the inner boundary of the single hole in the unbounded face.
// as the inner boundary of the single hole in the unbounded face.
Face_iterator fit;
const Face_handle uf = arr.unbounded_face();
Inner_ccb_iterator iccb_it = uf->inner_ccbs_begin();
Ccb_halfedge_circulator first, circ;
Halfedge_handle he;
out_bound.erase (out_bound.vertices_begin(), out_bound.vertices_end());
circ = first = *iccb_it;
@ -96,7 +93,7 @@ public:
{
out_bound.push_back (circ->source()->point());
--circ;
} while (circ != first);
++iccb_it;
@ -104,19 +101,19 @@ public:
for (fit = arr.faces_begin(); fit != arr.faces_end(); ++fit)
{
CGAL_assertion (fit->data() != this->UNVISITED);
// If a bounded face has an inside count that equals 0, it forms a hole
// in the union.
if (! fit->is_unbounded() && fit->data() == 0)
{
Polygon_2 pgn_hole;
circ = first = fit->outer_ccb();
do
{
pgn_hole.push_back (circ->source()->point());
--circ;
} while (circ != first);
// Insert it to the containers of holes in the Minkowski sum.

View File

@ -12,10 +12,7 @@
// This file is provided AS IS with NO WARRANTY OF ANY KIND, INCLUDING THE
// WARRANTY OF DESIGN, MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE.
//
// $URL$
// $Id$
//
// Author(s) : Ron Wein <wein@post.tau.ac.il>
// Author(s) : Ron Wein <wein_r@yahoo.com>
#ifndef CGAL_POLYGON_CONVEX_DECOMPOSITION_H
#define CGAL_POLYGON_CONVEX_DECOMPOSITION_H

View File

@ -12,9 +12,6 @@
// This file is provided AS IS with NO WARRANTY OF ANY KIND, INCLUDING THE
// WARRANTY OF DESIGN, MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE.
//
// $URL$
// $Id$
//
// Author(s) : Efi Fogel <efifogel@gmail.com>
#ifndef CGAL_POLYGON_VERTICAL_DECOMPOSITION_2_H

View File

@ -12,10 +12,7 @@
// This file is provided AS IS with NO WARRANTY OF ANY KIND, INCLUDING THE
// WARRANTY OF DESIGN, MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE.
//
// $URL$
// $Id$
//
// Author(s) : Ron Wein <wein@post.tau.ac.il>
// Author(s) : Ron Wein <wein_r@yahoo.com>
// (based on an old version by Eyal Flato)
// Efi Fogel <efifogel@gmail.com>

View File

@ -12,10 +12,7 @@
// This file is provided AS IS with NO WARRANTY OF ANY KIND, INCLUDING THE
// WARRANTY OF DESIGN, MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE.
//
// $URL$
// $Id$
//
// Author(s) : Ron Wein <wein@post.tau.ac.il>
// Author(s) : Ron Wein <wein_r@yahoo.com>
#ifndef CGAL_APPROXIMATED_OFFSET_H
#define CGAL_APPROXIMATED_OFFSET_H
@ -29,7 +26,7 @@ namespace CGAL {
/*!
* Approximate the offset of a given simple polygon by a given radius,
* using the convolution method.
* Note that as the input polygon may not be convex, its offset may not be
* Note that as the input polygon may not be convex, its offset may not be
* simply connected. The result is therefore represented as a polygon with
* holes.
* \param pgn The polygon.
@ -96,7 +93,7 @@ approximated_offset_2 (const Polygon_with_holes_2<Kernel, Container>& pwh,
* Approximate the offset of a given simple polygon by a given radius,
* by decomposing it to convex sub-polygons and computing the union of their
* offsets.
* Note that as the input polygon may not be convex, its offset may not be
* Note that as the input polygon may not be convex, its offset may not be
* simply connected. The result is therefore represented as a polygon with
* holes.
* \param pgn The polygon.
@ -132,7 +129,7 @@ approximated_offset_2 (const Polygon_2<Kernel, Container>& pgn,
/*!
* Approximate the inset of a given simple polygon by a given radius, using
* the convolution method.
* Note that as the input polygon may not be convex, its inset may not be
* Note that as the input polygon may not be convex, its inset may not be
* simply connected. The result is therefore represented as a set of polygons.
* \param pgn The polygon.
* \param r The inset radius.

View File

@ -12,10 +12,7 @@
// This file is provided AS IS with NO WARRANTY OF ANY KIND, INCLUDING THE
// WARRANTY OF DESIGN, MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE.
//
// $URL$
// $Id$
//
// Author(s) : Ron Wein <wein@post.tau.ac.il>
// Author(s) : Ron Wein <wein_r@yahoo.com>
// Efi Fogel <efifogel@gmail.com>
#ifndef CGAL_MINKOWSKI_SUM_2_H

View File

@ -12,10 +12,7 @@
// This file is provided AS IS with NO WARRANTY OF ANY KIND, INCLUDING THE
// WARRANTY OF DESIGN, MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE.
//
// $URL$
// $Id$
//
// Author(s) : Ron Wein <wein@post.tau.ac.il>
// Author(s) : Ron Wein <wein_r@yahoo.com>
#ifndef CGAL_OFFSET_POLYGON_H
#define CGAL_OFFSET_POLYGON_H
@ -29,7 +26,7 @@ namespace CGAL {
/*!
* Compute the offset of a given simple polygon by a given radius,
* using the convolution method.
* Note that as the input polygon may not be convex, its offset may not be
* Note that as the input polygon may not be convex, its offset may not be
* simply connected. The result is therefore represented as a polygon with
* holes.
* \param pgn The polygon.
@ -52,7 +49,7 @@ offset_polygon_2 (const Polygon_2<typename ConicTraits::Rat_kernel,
Offset_polygon_2 offset_bound;
std::list<Offset_polygon_2> offset_holes;
exact_offset (pgn, r,
exact_offset (pgn, r,
offset_bound, std::back_inserter(offset_holes));
return (typename Gps_traits_2<ConicTraits>::Polygon_with_holes_2
@ -85,7 +82,7 @@ offset_polygon_2 (const Polygon_with_holes_2<typename ConicTraits::Rat_kernel,
Offset_polygon_2 offset_bound;
std::list<Offset_polygon_2> offset_holes;
exact_offset (pwh, r,
exact_offset (pwh, r,
offset_bound, std::back_inserter(offset_holes));
return (typename Gps_traits_2<ConicTraits>::Polygon_with_holes_2
@ -96,7 +93,7 @@ offset_polygon_2 (const Polygon_with_holes_2<typename ConicTraits::Rat_kernel,
* Compute the offset of a given simple polygon by a given radius,
* by decomposing it to convex sub-polygons and computing the union of their
* offsets.
* Note that as the input polygon may not be convex, its offset may not be
* Note that as the input polygon may not be convex, its offset may not be
* simply connected. The result is therefore represented as a polygon with
* holes.
* \param pgn The polygon.
@ -132,7 +129,7 @@ offset_polygon_2 (const Polygon_2<typename ConicTraits::Rat_kernel,
/*!
* Compute the inset of a given simple polygon by a given radius, using the
* convolution method.
* Note that as the input polygon may not be convex, its inset may not be
* Note that as the input polygon may not be convex, its inset may not be
* simply connected. The result is therefore represented as a set of polygons.
* \param pgn The polygon.
* \param r The inset radius.

View File

@ -1,2 +1 @@
Ron Wein <wein@post.tau.ac.il>
Efi Fogel <efif@post.tau.ac.il>