cgal/Nef_2/noweb/Segment_overlay.lw

2407 lines
85 KiB
Plaintext

%------------------------------------------------------------------------------
%KILLSTART DISS REP
%LDEL TRACE.*?\)\;
\documentclass[a4paper]{article}
\usepackage{MyLweb}
\input{defs}
\excludeversion{ignoreindiss}
\includeversion{onlyindiss}
\excludeversion{ignore}
\begin{document}
\title{Sweeping segments in a generic framework}
\author{K. Mehlhorn \and S. Naeher \and M. Seel}
\maketitle
\tableofcontents
%KILLEND REP
\section{Introduction}
%KILLEND DISS
We describe a generic sweep algorithm of line segments along the lines
of the algorithm that is part of the LEDA library. We basically
transferred the segment sweep algorithm as described in the LEDA book
\cite{ledabook} into our generic sweep framework |generic_sweep|. We
focus on special adaptations and refer the user to the description in
\cite[chapter 10]{ledabook} for a deeper understanding.
\begin{onlyindiss}
Calculating the overlay of a set of segments is not a theoretical
problem anymore. We know its complexity and there are many existing
presentations of it. We could have described this module as a black
box and just stated the input and output properties of it. We decided
to add the implementation description because this presentation
stresses the techniques of generic programming and completes the
layered design of Nef polyhedra from the interface class down to the
sweep engine. Users who have the corresponding insights can skip this
section.
\end{onlyindiss}
\displayeps{Segment_overlay1}{The design of the segment overlay
module. |Segment_overlay_traits| implements the concept
|GenericSweepTraits|. There are actually two instances
|leda_seg_overlay_traits| and |stl_seg_overlay_traits| realizing
|Segment_overlay_traits| depending on the module configuration. The
three template parameters allow an adaptation depending on input,
output, and geometry.}{10cm}
To use our generic sweep framework we implement a traits model
|Segment_overlay_traits| that is plugged into the class
|generic_sweep<T>| (the bottom layer). For the concept of the
parameter |T| please refer to class |GenericSweepTraits| in the
appendix on page \pageref{GenericSweepTraits}. For the functionality
of class |generic_sweep| see the manual page on page
\pageref{generic_sweep}. |generic_sweep<Segment_overlay_traits<...>>|
is a generic sweep framework for the calculation of the overlay of
segments.
If you browse the original algorithm |SWEEP_SEGMENTS| with respect to
code dependencies, you find that it is hard-wired to several LEDA
modules. The wires of the original algorithm are the \emph{geometric
kernel} which is used, the \emph{input interface} which is a list of
segments, and the \emph{output interface} which is a LEDA embedded
graph. We decouple the above from concrete data types by introducing
concepts for the three units: input, output, geometry.
The \emph{input concept} is easy. We use iterators defining an
iterator range of segments (corresponding to the list of segments in
the LEDA sweep).
The \emph{output concept} is a plane map data structure like LEDA
plane maps (bidirected, embedded graphs). For an introduction refer to
\cite[chapter 8]{ledabook}. But of course there are several other
standard data structures in the literature like the \emph{Halfedge
Data Structures} (HDS), and \emph{Directed Cyclic Edge Lists}
(DCEL). See for example the textbooks \cite{mmmo:CG,prep-sham:CG} and
the CGAL HDS implementation in the manual \cite{wwwcgalhome}.
Sometimes the output has to be enriched by some additional bookkeeping
data structures (e.g. maps) to associate additional information with
the vertices and edges of the graph.
The \emph{geometric concept} contains the geometry used for the
algorithmic decisions of the sweep: geometric types like points and
segments and the primitive operations on them. Our correctness
considerations are based on affine planar geometry. However this sweep
framework works also instantiated with other geometry models.
The genericity is achieved by encapsulating the three concepts into
three template parameters of |Segment_\-overlay_\-traits<>|, which
allows a user to transport his geometry, geometric primitives and his
input and output structure into the segment sweep framework. We will
describe the concept of the geometric types, primitives and the
concept of the output graph structure below. We first give an abstract
introduction of the output produced.
The ouput of the algorithm is the result of transformations of an
output object triggered by method calls of the traits class. There are
some methods that can be used to manipulate a graph structure $G =
(V,E)$. If the implementation of those graph manipulation methods
follows the semantic description below then the output graph obtains a
certain structure. Additionally, there is also some kind of message
passing associated with the output structure. This allows a user to
refine necessary additional information from the sweep.
Assume the output contains a graph structure $G = (V,E)$ which can
represent an embedded plane map (a bidirected graph where reversal
edges are paired and where the nodes are embedded into the plane by
associating point coordinates). We state some properties of the
output production:
\begin{enumerate}
\renewcommand{\labelenumi}{O\theenumi.}
\item All end points and intersection points of segments are called
events and trigger calls to create new nodes $v \in V$ which obtain
the knowledge about their embedding via a point |p|. The creation is
done in the lexicographic order on points as specified by the user in
the geometric traits class.
\item The sweep explores the skeleton of the planar subdivision
induced by the set of input segments of the iterator range. Thus for
each segment |s| there are method calls which can be used to create
edges in $G$ such that the straight line embedding of their union
corresponds to |s|.
\item Each halfedge $e \in E$ gets to know a list of input segments
supporting its straight line embedding.
\item Each node $v \in V$ gets to know a halfedge below (vertical
ray-shooting property) where degeneracies are broken with a
left-closed perturbation scheme (all edges include their left source
node during the ray shoot). Additionally each node $v \in V$ gets to
know the input segments which start at, end at or contain its
embedding |point(v)|. Finally each node gets to know explicitly if it
is supported by a trivial input segment.
\item If the edge creation operations follow the semantic description
then each node |v| has an adjacency list such that visiting all
adjacent nodes while iterating the adjacency list corresponds to a
counterclockwise rotation around |v| (the straight line drawing of $G$
is counterclockwise \emph{order-preserving}). All adjacency lists
additionally have the property of a \emph{forward-prefix}. This means
that forward-oriented edges\footnote{an edge |e| is called
forward-oriented when $|point(source(e))| <_{lex} |point(target(e))|$.}
build a prefix in each adjacency list.
\end{enumerate}
\displayeps{Segment_overlay2}{The three concepts that allow adaptation
of the generic segment sweep. In the output concept the abbreviations
are |V| for |Vertex_handle|, |E| for |Halfedge_handle|, and |I| for
|Iterator|.}{14.5cm}
The interfaces of the three concepts are depicted in figure
\figref{Segment_overlay2}. For the actual semantics please consult the
manual pages for |SegmentOverlayOutput| on page
\pageref{SegmentOverlayOutput} and |SegmentOverlayGeometry_2| on page
\pageref{SegmentOverlayGeometry_2} in the appendix.
@ \subsection{Formalizing the sweep --- Invariants}\label{sweepinvars}
We sweep the plane from left to right by a vertical line
|SL|. Whenever we encounter an event point we have to take actions to
produce the output structure. Both endpoints of each segment and
non-degenerate intersection points of any two non-overlapping segments
define our \emph{events}. To ensure the correct actions we resort to
a list of invariants on which we can rely just before we encounter an
event and which we ensure by certain actions for the status
thereafter.
The sweep is determined by an interaction between two major data
structures, an event queue |XS| and a sweep status structure |YS|. The
event queue |XS| controls the stepping of |SL| across the plane. We
store a point as the key of each event. The order of the events
corresponds to the lexicographical order on all points as defined in
the geometry kernel. |YS| stores segments intersecting the sweep line
|SL| ordered according to their intersection points from bottom to
top. Note that for a segment |s| in |YS| $|source(s)| <_{lex}
|target(s)|$ according to our lexicographic order on points.
The above description talks about vertical sweep lines and geometry
which is left and right of the sweep line. As soon as there are events
with identical $x$-coordinates and vertical segments we have to be
more accurate. Imagine a sweep line which is slanted by an
infinitesimal angle counterclockwise thus processing the points on the
vertical line bottom-up. A more accurate intuition is created by a
sweepline consisting of an infinitesimal small step down, where the
vertical ray down is right of the current sweep position and the
vertical ray up is left of the current sweep position and the
horizontal step marks the current position while going from $-\infty$
to $+\infty$ along the $y$-axis of the coordinate system. For an
extensive treatment please refer to \cite[chapter 10]{ledabook}.
\displaylps{sweepline}{Sketching the sweepline intuition in the case
of degeneracies. The left figure shows the sweepline before the event
at $p$. The right figure shows the figure just after the event at
$p$.}
We treat the degenerate cases of several segments ending,
intersecting, and starting in one point explicitly in our code. We
also handle the possibility that several segments overlap. This
implies that just before an event the sweep line may intersect a whole
bundle of segments containing the event. Some extend through it, some
end there. And just after the event the bundle of the segments
extending though the event may be enriched by several segments
starting at the event. The segments of both bundles can be ordered
according to their points of intersection with the sweep line. In case
of overlapping segments their point of intersection with the sweep
line is identical. To break the tie we take the order on the identity
of each input segment\footnote{Internally we handle segments via
pointers. Then their identity is just the memory address.}. Note that
due to the degeneracy events can have multiple character.
\begin{invar}\label{xsinvars}
The event queue |XS| is a sorted sequence of items |<p,x>| called
events, where |p| is the embedding point of the event and |x| is an
associated information link. \emph{Starting} and \emph{ending} events
refer to the endpoints of all segments. \emph{Intersection} events are
defined by the points of intersection of two non-overlapping
segments. At any sweep position |XS| contains all starting and ending
events and all intersection events right of |SL| that are the result
of an intersection of segments that are neighbors in |YS|. The order
in |XS| is the lexicographic order on points |compare_xy()| introduced
by the geometry concept.
\end{invar}
\begin{invar}\label{ysinvars}
The sweep status structure |YS| is a sorted sequence of items |<s,y>|
where |s| is a segment intersected by the sweep line. The order is
defined by the points of intersection of the segments and the sweep
line from bottom to top.
\end{invar}
Consider any bundle of segments ending at or extending through an
event. We want to save on geometric calculations. Therefore the
information slot of the items in |YS| is used to identify such
bundles.
\begin{invar}\label{yslinks}
Let |<s,y>| and |<s',y'>| be two successing items in |YS|. The
information |y| is used as a flag how the two segments are
geometrically related. If |s| and |s'| are overlapping at the current
sweep position then |y| points to its successor item |<s',y'>| in
|YS|. If |s| and |s'| intersect right of the sweep line then |y|
points to the corresponding event in |XS|. Otherwise it is the null
handle.
Additionally for all items |it| we know an edge in the output graph
that is supported by the segment via a map |Edge_of[it]|.
\end{invar}
\begin{invar}\label{xslink}
For all intersection events and ending events |<p,x>| in |XS| the
information |x| is a link to an item |<s,y>| in |YS| such that the
segment |s| contains |p|. For ending events the invariant is
established as soon as the segment |s| enters |YS|.
\end{invar}
This construction allows us to shortcut from an event into the range
of interest within |YS| without using the standard (logarithmic)
search operations of |YS|.
\begin{lemma}
Starting from any intersection or ending event we can identify the
bundle of segments/items in |YS| ending at or extending through the
event in time proportional to the size of the bundle.
\end{lemma}
\begin{proof}
Consider an event |<p,x>| that is not only starting event. Then |x| is
a link to an item |<s,y>| in |YS| such that |s| contains |p| due to
Invariant \ref{xslink}. Iterating locally up and down starting from
|x| allows us to identify the range of items in |YS| whose segments
contain |p| by the marking links according to Invariant \ref{yslinks}.
\end{proof}
At last, at each event we ensure partial output correctness that leads
to global output correctness after the last event.
\begin{invar}\label{outputcorr}
Assume the output model's operations are defined according to our
specification. Then the overlay of all segments that are fully left
of our sweep line is correctly calculated.
\end{invar}
The following hashing idea saves again on geometric calculations and
search.
\begin{invar}\label{intersectionhash}
For each pair of segments |(s1,s2)| in |YS| and intersecting right of
|SL| which have been neighbors in |YS| once (and might have been
separated afterwards) there's a hash table short-cut to the
corresponding intersection event |it = IEvent(s1,s2)|.
\end{invar}
Note that the algorithmic correctness of this framework has two
aspects. There are global considerations and local considerations.
The global considerations concern issues like why the algorithm
terminates and why it does calculate the output we ask for. The local
considerations concern the aspects of the event handling. Namely, why
our event handling and the intialization phase of our framework
ensures the invariants stated above. Taking both issues together we
obtain the certainty that our code module actually calculates the
overlay correctly.
For the global correctness we will see that all starting and ending
events are handled in our code and inflated into the machinery in the
initialization phase. One exception that we have to incorporate into
our code is the occurance of trivial segments. The final question is
if we catch all intersection events? There's a trivial observation why
we cannot miss any such event.
\begin{lemma}
If we ensure that all our invariants hold at all events then we cannot
miss any intersection event.
\end{lemma}
\begin{proof}
Assume we miss an intersection event |ev|. If we miss several take the
lexicographically smallest one. Then the event just before |ev| was
correctly treated and the two segments that imply |ev| are neighbors
in |YS|. But Invariant \ref{xsinvars} tells us that |ev| is in |XS|
which leads to a contradiction.
\end{proof}
Note that global termination is not a big issue in sweep frameworks.
As soon as all events have gone we stop the iteration. Note finally
that if we locally keep Invariant \ref{outputcorr} just after each
event then we also know that our result is correctly calculated. We
now will link the above insights to the code.
\subsection{Two generic sweep traits models}
We use our generic plane sweep paradigm to execute the sweep. We
offer two models to plug into the |generic_sweep| framework. For the
concept see |GenericSweepTraits| in the appendix. One based on LEDA
\cite{ledabook} and including several runtime optimizations, the other
one based purely on STL data structures
\cite{MusserSaini:STLbook}. This design was necessary to allow using
the sweep module when CGAL is installed without the presence of LEDA.
\begin{ignoreindiss}
The sweep traits class wraps it all.
<<Segment_overlay_traits.h>>=
<<CGAL Header>>
#ifndef CGAL_SEGMENT_OVERLAY_TRAITS_H
#define CGAL_SEGMENT_OVERLAY_TRAITS_H
#include <cassert>
#undef _DEBUG
#define _DEBUG 23
#include <CGAL/Nef_2/debug.h>
//#define INCLUDEBOTH
#if defined(CGAL_USE_LEDA) || defined(INCLUDEBOTH)
<<leda segment overlay model>>
#endif // defined(CGAL_USE_LEDA) || defined(INCLUDEBOTH)
#if !defined(CGAL_USE_LEDA) || defined(INCLUDEBOTH)
<<stl segment overlay model>>
#endif // !defined(CGAL_USE_LEDA) || defined(INCLUDEBOTH)
namespace CGAL {
#ifdef CGAL_USE_LEDA
#define Segment_overlay_traits leda_seg_overlay_traits
static const char* sweepversion = "LEDA segment overlay sweep";
#else
#define Segment_overlay_traits stl_seg_overlay_traits
static const char* sweepversion = "STL segment overlay sweep";
#endif
} // namespace CGAL
#include <CGAL/generic_sweep.h>
#endif // CGAL_SEGMENT_OVERLAY_TRAITS_H
<<leda segment overlay model>>=
#if CGAL_LEDA_VERSION < 500
#include <LEDA/tuple.h>
#include <LEDA/slist.h>
#include <LEDA/list.h>
#include <LEDA/map.h>
#include <LEDA/map2.h>
#include <LEDA/sortseq.h>
#include <LEDA/p_queue.h>
#else
#include <LEDA/core/tuple.h>
#include <LEDA/core/slist.h>
#include <LEDA/core/list.h>
#include <LEDA/core/map.h>
#include <LEDA/core/map2.h>
#include <LEDA/core/sortseq.h>
#include <LEDA/core/p_queue.h>
#endif
#include <utility>
#include <strstream>
namespace CGAL {
<<leda debugging routines>>
<<leda segment overlay traits class>>
} // namespace CGAL
@
\end{ignoreindiss}
\subsection{The LEDA traits model}
Our class obtains three template types which have to be models
for the corresponding concepts described below.
<<leda segment overlay traits class>>=
template <typename IT, typename PMDEC, typename GEOM>
class leda_seg_overlay_traits {
public:
<<leda introducing the types from the traits>>
<<leda internal segment type>>
// types interfacing the generic sweep frame:
ITERATOR its, ite;
OUTPUT& GO;
const GEOMETRY& K;
<<leda order types for segments and points>>
<<leda sweep data structures>>
<<leda helping operations>>
<<leda operation for keeping the intersection invariant>>
<<leda initialization of the sweep>>
<<leda iteration control>>
<<leda handling the event>>
<<leda postprocessing of the sweep>>
}; // leda_seg_overlay_traits
@ The following types are introduced by the traits classes. See the
the concept descriptions |SegmentOverlayOutput| for |PMDEC| and
|SegmentOverlayGeometry_2| for |GEOM|. The iterator concept required
is that of an STL input iterator.
<<leda introducing the types from the traits>>=
typedef IT ITERATOR;
typedef std::pair<IT,IT> INPUT;
typedef PMDEC OUTPUT;
typedef typename PMDEC::Vertex_handle Vertex_handle;
typedef typename PMDEC::Halfedge_handle Halfedge_handle;
typedef GEOM GEOMETRY;
typedef typename GEOMETRY::Point_2 Point_2;
typedef typename GEOMETRY::Segment_2 Segment_2;
@ We define an internal segment type |ISegment| based on a pointer to
allow two constructions. At first we want to be able to couple the
internal segment to the input object. We achieve this by maintaining
both as a pair |two_tuple<Segment,ITERATOR>| in a list. Secondly our
internal segment type is just a pointer to such a pair. We can thus
not only compare internal segments geometrically by the pair's first
component, but also check identity by the pointer address.
<<leda internal segment type>>=
typedef leda_two_tuple<Segment_2,ITERATOR> seg_pair;
typedef seg_pair* ISegment;
typedef leda_list<seg_pair> IList;
typedef typename IList::iterator ilist_iterator;
@ The predicates below are solely based on a few geometric kernel
predicates. The clever observation is the fact, that the comparison
predicate |cmp_segs_at_sweepline| is only called with one segment
containing the sweep point. We know that assertion from the
specification of LEDA sortseqs. Compared to the LEDA sweep algorithm
we add non-geometric sentinel segments to avoid checking of boundary
cases.
<<leda order types for segments and points>>=
class cmp_segs_at_sweepline : public leda_cmp_base<ISegment>
{ const Point_2& p;
ISegment s_bottom, s_top; // sentinel segments
const GEOMETRY& K;
public:
cmp_segs_at_sweepline(const Point_2& pi,
ISegment s1, ISegment s2, const GEOMETRY& k) :
p(pi), s_bottom(s1), s_top(s2), K(k) {}
int operator()(const ISegment& is1, const ISegment& is2) const
{ // Precondition: p is identical to the left endpoint of s1 or s2.
if ( is2 == s_top || is1 == s_bottom ) return -1;
if ( is1 == s_top || is2 == s_bottom ) return 1;
if ( is1 == is2 ) return 0;
const Segment_2& s1 = is1->first();
const Segment_2& s2 = is2->first();
int s = 0;
if ( p == K.source(s1) ) s = K.orientation(s2,p);
else if ( p == K.source(s2) ) s = - K.orientation(s1,p);
else CGAL_assertion_msg(0,"compare error in sweep.");
if ( s || K.is_degenerate(s1) || K.is_degenerate(s2) )
return s;
s = K.orientation(s2,K.target(s1));
if (s==0) return ( is1 - is2 );
// overlapping segments are not equal
return s;
}
};
@ The lexicographic order on points is just transferred from the
geometric kernel.
<<leda order types for segments and points>>=
struct cmp_pnts_xy : public leda_cmp_base<Point_2>
{ const GEOMETRY& K;
public:
cmp_pnts_xy(const GEOMETRY& k) : K(k) {}
int operator()(const Point_2& p1, const Point_2& p2) const
{ return K.compare_xy(p1,p2); }
};
@ We use several LEDA data structures. The x-structure |XS| is an
ordered sequence of items based on the key type |Point_2|. For the
item concept of LEDA refer to the manual \cite{ledaman} or the LEDA
book. We associate a link into the y-structure for intersection and
ending events to save on unnecessary search operations within |YS|.
When we reach an event at |p_sweep| we can thus shortcut to find an
item in |YS| which contains a segment (as its key) containing
|p_sweep|. The y-structure |YS| is a sorted sequence of |ISegments|.
The associated |seq_item| link serves two purposes: (1) it bundles
segments in |YS| together by pointing to the (lexicographic) next
event which they contain. (2) it bundles segments which overlap. See
the LEDA book for a more elaborate description.
During the sweep we associate edges from the constructed output graph
to the items in YS. We use a hash map |map<seq_item,Halfedge_handle>|
|Edge_of| for this purpose. Finally for each event point there is a
possible sequence of input segments starting at the event. To maintain
this sequence we use a priority queue |p_queue<Point_2,ISegment> SQ|.
When two segments |s1,s2| become neighbors in |YS| we check if they
intersect right of the sweep line. If they do we calculate the
intersection point |p| and insert a corresponding event into |XS|. Now
it can happen that |s1| and |s2| get again separated by a new segment
before |p| is reached. We have several possibilities in this case. We
could remove the event at |p| again to keep the space bound implied by
the Invariant \ref{xsinvars}. However the size of the output
structure anyway comprises the space for keeping the event in |XS|. We
can even do better with respect to runtime. Instead of recalculating
the geometric intersection information when |s1| and |s2| become
neighbors again we can try to recover the previously calculated event
from the two dimensional hash $|map2<ISegment,ISegment,seq_item>
IEvent|$.
The additional members are used as a central place for data storage.
|event| provides a handle on the current event queue item. |p_sweep|
is the position of the current event which is also used from the
segment comparison object |SLcmp|.
<<leda sweep data structures>>=
typedef leda_sortseq<Point_2,seq_item> EventQueue;
typedef leda_sortseq<ISegment,seq_item> SweepStatus;
typedef leda_p_queue<Point_2,ISegment> SegQueue;
typedef leda_map<seq_item,Halfedge_handle> AssocEdgeMap;
typedef leda_slist<ITERATOR> IsoList;
typedef leda_map<seq_item, IsoList* > AssocIsoMap;
typedef leda_map2<ISegment,ISegment,seq_item> EventHash;
seq_item event;
Point_2 p_sweep;
cmp_pnts_xy cmp;
EventQueue XS;
seg_pair sl,sh;
cmp_segs_at_sweepline SLcmp;
SweepStatus YS;
SegQueue SQ;
EventHash IEvent;
IList Internal;
AssocEdgeMap Edge_of;
AssocIsoMap Isos_of;
leda_seg_overlay_traits(const INPUT& in, OUTPUT& G,
const GEOMETRY& k) :
its(in.first), ite(in.second), GO(G), K(k),
cmp(K), XS(cmp), SLcmp(p_sweep,&sl,&sh,K), YS(SLcmp), SQ(cmp),
IEvent(0), Edge_of(0), Isos_of(0) {}
@ We define some code short cuts.
\begin{ignoreindiss}
<<leda helping operations>>=
leda_string dump_structures() const
{
std::ostrstream out;
out << "SQ= ";
pq_item pqit;
forall_items(pqit,SQ) {
if (SQ.prio(pqit)==XS.key(XS.succ(XS.min_item())))
{ out << SQ.inf(pqit)->first(); }
pqit = SQ.next_item(pqit);
}
seq_item sit;
out << "\nXS=\n";
forall_items(sit,XS)
out << " " << XS.key(sit) << " " << XS.inf(sit)
<<std::endl;
out << "YS=\n";
for( sit = YS.max_item(); sit; sit=YS.pred(sit) )
out << " "<<YS.key(sit)->first()<<" "<<YS.inf(sit)<<std::endl;
out << '\0';
leda_string res(out.str()); out.freeze(0); return res;
}
@ \end{ignoreindiss}
<<leda helping operations>>=
Point_2 source(ISegment is) const
{ return K.source(is->first()); }
Point_2 target(ISegment is) const
{ return K.target(is->first()); }
ITERATOR original(ISegment s) const
{ return s->second(); }
int orientation(seq_item sit, const Point_2& p) const
{ return K.orientation(YS.key(sit)->first(),p); }
bool collinear(seq_item sit1, seq_item sit2) const
{ Point_2 ps = source(YS.key(sit2)), pt = target(YS.key(sit2));
return ( orientation(sit1,ps)==0 &&
orientation(sit1,pt)==0 );
}
@ Most events trigger changes in the segment sequence along the sweep
line. We have to reflect such changes in a test for new intersection
events right of the sweep line as soon as two segments become
neighbors. The following code ensures the Invariants \ref{ysinvars},
\ref{yslinks}, \ref{xslink} and uses the hash tuning of Invariant
\ref{intersectionhash}. |s1| is the successor of |s0| in |YS|, hence,
|s0| and |s1| intersect right or above of the event iff |target(s1)|
is not left of the line supporting |s0|, and |target(s0)| is not right
of the line supporting |s1|. In this case we intersect the underlying
lines.
<<leda operation for keeping the intersection invariant>>=
void compute_intersection(seq_item sit0)
{
seq_item sit1 = YS.succ(sit0);
if ( sit0 == YS.min_item() || sit1 == YS.max_item() ) return;
ISegment s0 = YS.key(sit0);
ISegment s1 = YS.key(sit1);
int or0 = K.orientation(s0->first(),target(s1));
int or1 = K.orientation(s1->first(),target(s0));
if ( or0 <= 0 && or1 >= 0 ) {
seq_item it = IEvent(YS.key(sit0),YS.key(sit1));
if ( it==0 ) {
Point_2 q = K.intersection(s0->first(),s1->first());
it = XS.insert(q,sit0);
}
YS.change_inf(sit0, it);
}
}
@ \subsubsection*{Event Handling}
We start with the knowledge that our invariants from Section
\ref{sweepinvars} hold. First we create a new vertex |v| in the ouput
structure. Then we work in four phases: (1) We handle the ingoing
bundle which ends at |v|. (2) We communicate all the knowledge about
the new vertex (3) We have to deal with the segments starting at
|p_sweep|. (4) We clean up to reestablish missing invariants.
<<leda handling the event>>=
void process_event()
{
TRACEN("\n\n >>> process_event: "<<p_sweep<<" "<<XS[event]<<" "<<event);
Vertex_handle v = GO.new_vertex(p_sweep);
seq_item sit = XS.inf(event);
<<leda handling ending and passing segments>>
<<leda completing additional information of the new vertex>>
<<leda inserting new segments starting at nodes>>
<<leda enforcing the invariants for YS>>
}
@ We first have to locate the bundle going through |p_sweep|. We
deviate from the implementation of the LEDA algorithm |SWEEP_SEGMENTS|
in one respect. For each segment in |YS| we store a bidirected edge
pair extending along the segment. When we reach an event point we
connect these edges to the newly created node. Note that this change
is necessary if you use halfedge data structures for the output. The
original approach used temporarily incomplete edge pairs (only forward
directed halfedges) and coupled and embedded them in a postprocessing
phase. But space minimally maintained halfedges like those of the CGAL
HDS can only exist in pairs.
If there is a non-nil item |sit = XS.inf(event)| associated with
|event|, |key(sit)| is either an ending or passing segment. We use
|sit| as an entry point to compute the bundle of segments ending at or
passing through |p_sweep|. In particular, we compute the first
(|sit_first|) and the successor (|sit_succ)|) and predecessor
(|sit_pred|) items.
<<leda handling ending and passing segments>>=
seq_item sit_succ(0), sit_pred(0), sit_pred_succ(0), sit_first(0);
if (sit == nil)
<<leda check p_sweep in YS>>
/* If no segment contains p_sweep then sit_pred and sit_succ are
correctly set after the above locate operation, if a segment
contains p_sweep sit_pred and sit_succ are set below when
determining the bundle.*/
if (sit != nil) { // key(sit) is an ending or passing segment
<<leda determine upper bundle item sit_succ>>
<<leda hash upper intersection event>>
<<leda walk ingoing bundle and trigger graph updates>>
<<leda reverse continuing bundle edges>>
} // if (sit != nil)
@ As |sit == nil| we do not know if a segment stored in |YS| does
contain |p_sweep|. We have to query |YS| with a trivial segment
|(p_sweep,p_sweep)| to find out. Two results are possible. Either a
segment referenced hereafter by |sit| constains the event point or we
determine the two segments above and below |p_sweep| in |sit_pred| and
|sit_succ|.
<<leda check p_sweep in YS>>=
{
Segment_2 s_sweep = K.construct_segment(p_sweep,p_sweep);
seg_pair sp(s_sweep,ITERATOR());
sit_succ = YS.locate( &sp );
if ( sit_succ != YS.max_item() &&
orientation(sit_succ,p_sweep) == 0 )
sit = sit_succ;
else {
sit_pred = YS.pred(sit_succ);
sit_pred_succ = sit_succ;
}
TRACEN("looked up p_sweep "<<PIS(YS.key(sit_succ)));
}
@ We first walk up as long as the event is contained in the
segment referenced via |sit|.
<<leda determine upper bundle item sit_succ>>=
TRACEN("ending/passing segs");
while ( YS.inf(sit) == event ||
YS.inf(sit) == YS.succ(sit) ) // overlapping
sit = YS.succ(sit);
sit_succ = YS.succ(sit);
seq_item sit_last = sit;
@ We hash the upper event according to Invariant \ref{intersectionhash}.
<<leda hash upper intersection event>>=
seq_item xit = YS.inf(sit_last);
if (xit) {
ISegment s1 = YS.key(sit_last);
ISegment s2 = YS.key(sit_succ);
IEvent(s1,s2) = xit;
TRACEN("hashing "<<PIS(s1)<<PIS(s2)<<xit);
}
@ We walk the ingoing bundle down again and trigger the edge closing
calls for all items in the bundle (except overlapping segments). Note
that after this code chunk we have: (i) the bundle is empty if
|succ(sit_pred) == sit_first == sit_succ|, or (ii) the bundle is not
empty if |sit_first != sit_succ|.
The actions on the bundle are easy to specify. We have to glue one
edge per segment to the newly created node except when two segments
overlap. Note that the walk top-down over the bundle implies the
order-preserving embedding of the graph. Note also how we pass the
messages about the segments supporting the event via the corresponding
methods of the output object.
<<leda walk ingoing bundle and trigger graph updates>>=
bool overlapping;
do {
ISegment s = YS.key(sit);
seq_item sit_next = YS.pred(sit);
overlapping = (YS.inf(sit_next) == sit);
Halfedge_handle e = Edge_of[sit];
if ( !overlapping ) {
TRACEN("connecting edge to node "<<PIS(s)<<" "<<sit);
GO.link_as_target_and_append(v,e);
}
GO.supporting_segment(e,original(s));
if ( target(s) == p_sweep ) { // ending segment
TRACEN("ending segment "<<PIS(s));
if ( overlapping ) YS.change_inf(sit_next,YS.inf(sit));
YS.del_item(sit);
GO.ending_segment(v,original(s));
} else { // passing segment
TRACEN("passing segment "<<PIS(s));
if ( YS.inf(sit) != YS.succ(sit) )
YS.change_inf(sit, seq_item(0));
GO.passing_segment(v,original(s));
}
sit = sit_next;
}
while ( YS.inf(sit) == event || overlapping ||
YS.inf(sit) == YS.succ(sit) );
sit_pred = sit;
sit_first = sit_pred_succ = YS.succ(sit_pred); // first item of bundle
TRACE("event bundles between\n "<<PIS(YS.key(sit_succ)));
TRACEN("\n "<<PIS(YS.key(sit_pred)));
@ We have to ensure that segments that continue through the event
point have a reversed order within |YS| when |p_sweep| has been
passed. This ensures the correct order of |YS| with respect to
Invariant \ref{ysinvars}. Some complication stems from overlapping
segments. Their order based on identity may not be changed.
<<leda reverse continuing bundle edges>>=
while ( sit != sit_succ ) {
seq_item sub_first = sit;
seq_item sub_last = sub_first;
while (YS.inf(sub_last) == YS.succ(sub_last))
sub_last = YS.succ(sub_last);
if (sub_last != sub_first)
YS.reverse_items(sub_first, sub_last);
sit = YS.succ(sub_first);
}
// reverse the entire bundle
if (sit_first != sit_succ)
YS.reverse_items(YS.succ(sit_pred),YS.pred(sit_succ));
@ For the new node |v| we pass some information to the output
structure. We post the halfedge below, we post all trivial input
segments supporting the node. We obtain that information from the hash
structure |Isos_of| filled during the initialization phase.
<<leda completing additional information of the new vertex>>=
assert(sit_pred);
GO.halfedge_below(v,Edge_of[sit_pred]);
if ( Isos_of[event] != 0 ) {
const IsoList& IL = *(Isos_of[event]);
slist_item iso_it;
for (iso_it = IL.first(); iso_it; iso_it=IL.succ(iso_it) )
GO.trivial_segment(v,IL[iso_it] );
delete (Isos_of[event]); // clean up the list
}
@ We insert all segments starting at |p_sweep| into |YS| and create
the links within |YS| to mark items with overlapping segments.
<<leda inserting new segments starting at nodes>>=
ISegment next_seg;
pq_item next_it = SQ.find_min();
while ( next_it &&
(next_seg = SQ.inf(next_it), p_sweep == source(next_seg)) ) {
seq_item s_sit = YS.locate_succ(next_seg);
seq_item p_sit = YS.pred(s_sit);
TRACEN("inserting "<<PIS(next_seg)<<" at "<<PIS(YS.key(s_sit)));
if ( YS.max_item() != s_sit &&
orientation(s_sit, source(next_seg) ) == 0 &&
orientation(s_sit, target(next_seg) ) == 0 )
sit = YS.insert_at(s_sit, next_seg, s_sit);
else
sit = YS.insert_at(s_sit, next_seg, seq_item(nil));
assert(YS.succ(sit)==s_sit);
if ( YS.min_item() != p_sit &&
orientation(p_sit, source(next_seg) ) == 0 &&
orientation(p_sit, target(next_seg) ) == 0 )
YS.change_inf(p_sit, sit);
assert(YS.succ(p_sit)==sit);
XS.insert(target(next_seg), sit);
GO.starting_segment(v,original(next_seg));
// delete minimum and assign new minimum to next_seg
SQ.del_min();
next_it = SQ.find_min();
}
@ In contrast to the original LEDA segment intersection algorithm
|SWEEP_SEGMENTS| we create ``semi-open'' edges starting at the |event|
node and supported by the input segment. The iteration again ensures
the correct order-preserving embedding at the currently handled node
|v|.
<<leda inserting new segments starting at nodes>>=
for( seq_item sitl = YS.pred(sit_succ); sitl != sit_pred;
sitl = YS.pred(sitl) ) {
if ( YS.inf(sitl) != YS.succ(sitl) ) { // non-overlapping
TRACEN("non-overlapping "<<PIS(YS.key(sitl))<<" "<<sitl);
Edge_of[sitl] = GO.new_halfedge_pair_at_source(v);
} else {
TRACEN("overlapping "<<PIS(YS.key(sitl)));
Edge_of[sitl] = Edge_of[ YS.succ(sitl) ];
}
}
sit_first = YS.succ(sit_pred);
@ Depending on the outgoing bundle we determine possible intersections
between new neighbors; if |sit_pred| is no longer adjacent to its
former successor we change its intersection event to 0. Note that the
following chunk finishes Invariant \ref{xsinvars} with the help of the
method |compute_intersection()|.
<<leda enforcing the invariants for YS>>=
assert(sit_pred); assert(sit_pred_succ);
seq_item xit = YS.inf(sit_pred);
if ( xit ) {
ISegment s1 = YS.key(sit_pred);
ISegment s2 = YS.key(sit_pred_succ);
IEvent(s1,s2) = xit;
TRACEN("hashing "<<PIS(s1)<<PIS(s2)<<xit);
YS.change_inf(sit_pred, seq_item(0));
}
compute_intersection(sit_pred);
sit = YS.pred(sit_succ);
if (sit != sit_pred)
compute_intersection(sit);
@ \subsubsection*{Initialization}
We realize the propositions of the invariants at the beginning in
our sweep initialization phase. We insert all segment endpoints into
|XS|, insert sentinels into |YS|, and exploit the fact that insert
operations into the X-structure leave previously inserted points
unchanged to achieve that any pair of endpoints |p| and |q| with |p ==
q| are identical (if the geometric point type supports this).
Degenerate segments are stored in a list associated to their events.
The knowledge about their existence is transferred to the corresponding
output object as soon as it is constructed.
<<leda initialization of the sweep>>=
void initialize_structures()
{
TRACEN("initialize_structures");
ITERATOR it_s;
for ( it_s=its; it_s != ite; ++it_s ) {
Segment_2 s = *it_s;
seq_item it1 = XS.insert( K.source(s), seq_item(nil));
seq_item it2 = XS.insert( K.target(s), seq_item(nil));
if (it1 == it2) {
if ( Isos_of[it1] == 0 ) Isos_of[it1] = new IsoList;
Isos_of[it1]->push(it_s);
continue; // ignore zero-length segments in SQ/YS
}
Point_2 p = XS.key(it1);
Point_2 q = XS.key(it2);
Segment_2 s1;
if ( K.compare_xy(p,q) < 0 )
s1 = K.construct_segment(p,q);
else
s1 = K.construct_segment(q,p);
Internal.append(seg_pair(s1,it_s));
SQ.insert(K.source(s1),&Internal[Internal.last()]);
}
// insert a lower and an upper sentinel segment
YS.insert(&sl,seq_item(nil));
YS.insert(&sh,seq_item(nil));
TRACEN("end of initialization\n"<<YS.size());
}
@ Note the invariants of the sweep loop. The |event| has to be set
before the event action.
<<leda iteration control>>=
bool event_exists()
{
if (!XS.empty()) {
// event is set at end of loop and in init
event = XS.min();
p_sweep = XS.key(event);
return true;
}
return false;
}
void procede_to_next_event()
{ XS.del_item(event); }
@ The structure is finished as soon as all events have been treated.
As we always created edges in pairs and respected the adjacency list
order we have no completion phase as in \cite[chapter 10]{ledabook}.
\begin{ignoreindiss}
<<leda postprocessing of the sweep>>=
void complete_structures() {}
void check_invariants() {TRACEN("check_invariants\n"<<dump_structures());}
void check_final() {}
<<leda debugging routines>>=
#ifdef _DEBUG
#define PIS(s) (s->first())
#endif
@ \end{ignoreindiss}
\subsubsection*{Runtime Analysis}
We shortly give a runtime analysis of the above implementation. Assume
$U(S)$ is the embedded (undirected) graph created by a proper\footnote
{The output methods have to be constant time operations of the correct
semantics.} instantiation of our sweep framework. Let $n$ be the
number of input segments of the input set $S$, $n_v$ be the number of
nodes of $U(S)$, $n_e$ be the number of (undirected) edges of $U(S)$.
In the presence of overlapping segments let $\bar{n}_e := \sum_{e
\text{\ edge of\ } U(S)} s_e$ where $s_e$ is the number of segments in
$S$ supporting edge $e$ of $U(S)$. Note that during the sweep the
whole graph is constructed and explored. Then obviously all graph
related operations like creation, and support messaging take time
$O(n_v + \bar{n}_e)$.
The sweep initialization takes time $O(n \log n)$. There are $n_v$
events and at each event the sweep status structures are manipulated a
constant number of times. Pass again through the event handling
routine. All segments are inserted into |YS| and deleted from |YS|
once (during the whole sweep) and therefore that cost adds up to $O(n
\log \Labs{\mathit{YS}})$. The removal from |SQ| adds up $O(n \log
n)$. Finally all events are inserted and deleted once and this takes
$O(\Labs{XS} \log \Labs{XS})$. When subsequences of |YS| are explored
and swapped, this cost can again be dedicated to the exploration of
the graph and is therefore subsumed in the $O(n_v + \bar{n}_e)$ from
above.
We have to bound the size of |XS| and |YS|. Natural bounds are
$\Labs{XS}=O(n_v)$ (the number of nodes) and $\Labs{YS}=O(n)$ (the
number of segments). Therefore accumulating all of the above we obtain
$O(n_v + \bar{n}_e + n \log n + n_v \log (n+n_v) ) = O(n_v + \bar{n}_e
+ (n+n_v) \log (n+n_v) )$. Note that $\bar{n}_e$ and $n_v$ can be
quadratic in $n$.
\begin{lemma}\label{generic segment sweep runtime}
Assume that $n$, $n_v$, $\bar{n}_e$ are defined as above then the
runtime of the sweep algorithm is
\[O(n_v + \bar{n}_e + (n+n_v) \log (n+n_v) ).\]
\end{lemma}
@ \subsection{The STL traits model}
\begin{onlyindiss}
We will not show the STL model. It is highly analogously, sometimes
simpler, but also less efficient due to limitations in the interface
of the STL data structures. For the complete implementation see the
report \cite{TR:nefimplementation}.
\end{onlyindiss}
\begin{ignoreindiss}
We present an alternative module purely based on STL data structures.
As the STL has no hashing support we omit all the optimization used in
the LEDA approach. Additionally we cannot use any links from events
into the y-structure. As STL maps have no order manipulating
operations like |reverse_items| of LEDA |sortseq|s we have to delete
and reinsert ranges of segments at events. Therefore no shortcuts can
be stored to minimize search operations within |YS|.
<<stl segment overlay model>>=
#include <list>
#include <map>
#include <string>
#include <strstream>
namespace CGAL {
template <typename IT, typename PMDEC, typename GEOM>
class stl_seg_overlay_traits {
public:
<<stl introducing the types from the traits>>
<<stl internal segment type>>
// types interfacing the generic sweep frame
ITERATOR its, ite;
OUTPUT& GO;
const GEOMETRY& K;
<<stl order types for segments and points>>
<<stl sweep data structures>>
<<stl helping operations>>
<<stl operation for keeping the intersection invariant>>
<<stl initialization of the sweep>>
<<stl iteration control>>
<<stl handling the event>>
<<stl postprocessing of the sweep>>
}; // stl_seg_overlay_traits
} // namespace CGAL
@ The following types are introduced by the traits classes.
<<stl introducing the types from the traits>>=
typedef IT ITERATOR;
typedef std::pair<IT,IT> INPUT;
typedef PMDEC OUTPUT;
typedef typename PMDEC::Vertex_handle Vertex_handle;
typedef typename PMDEC::Halfedge_handle Halfedge_handle;
typedef GEOM GEOMETRY;
typedef typename GEOMETRY::Point_2 Point_2;
typedef typename GEOMETRY::Segment_2 Segment_2;
@ Similar to the LEDA implementation we obtain internal segments by
pointers to |pair<Segment,ITERATOR>|.
<<stl internal segment type>>=
typedef std::pair<Segment_2,ITERATOR> seg_pair;
typedef seg_pair* ISegment;
typedef std::list<seg_pair> IList;
typedef typename IList::const_iterator ilist_iterator;
@ The STL uses strict order type objects. The predicate is implemented
in analogy to |cmp_segs_at_sweepline|.
<<stl order types for segments and points>>=
class lt_segs_at_sweepline
{ const Point_2& p;
ISegment s_bottom, s_top; // sentinel segments
const GEOMETRY& K;
public:
lt_segs_at_sweepline(const Point_2& pi,
ISegment s1, ISegment s2, const GEOMETRY& k) :
p(pi), s_bottom(s1), s_top(s2), K(k) {}
lt_segs_at_sweepline(const lt_segs_at_sweepline& lt) :
p(lt.p), s_bottom(lt.s_bottom), s_top(lt.s_top), K(lt.K) {}
bool operator()(const ISegment& is1, const ISegment& is2) const
{
if ( is2 == s_top || is1 == s_bottom ) return true;
if ( is1 == s_top || is2 == s_bottom ) return false;
if ( is1 == is2 ) return false;
// Precondition: p is contained in s1 or s2.
const Segment_2& s1 = is1->first;
const Segment_2& s2 = is2->first;
int s = 0;
if ( K.orientation(s1,p) == 0 )
s = K.orientation(s2,p);
else if ( K.orientation(s2,p) == 0 )
s = - K.orientation(s1,p);
else CGAL_assertion_msg(0,"compare error in sweep.");
if ( s || K.is_degenerate(s1) || K.is_degenerate(s2) )
return ( s < 0 );
s = K.orientation(s2,K.target(s1));
if (s==0) return ( is1 - is2 ) < 0;
// overlapping segments are not equal
return ( s < 0 );
}
};
struct lt_pnts_xy
{ const GEOMETRY& K;
public:
lt_pnts_xy(const GEOMETRY& k) : K(k) {}
lt_pnts_xy(const lt_pnts_xy& lt) : K(lt.K) {}
int operator()(const Point_2& p1, const Point_2& p2) const
{ return K.compare_xy(p1,p2) < 0; }
};
@ We need three main STL data structures. The x-structure |XS| is an ordered
set based on the key type |Point_2|. The y-structure is a sorted sequence of
|ISegments|. During the sweep we associate edges from the constructed output
graph to these segments. Finally for each event point there is a possible
sequence of input segments starting at the event. To maintain this sequence we
use a STL multimap.
<<stl sweep data structures>>=
typedef std::map<ISegment, Halfedge_handle, lt_segs_at_sweepline>
SweepStatus;
typedef typename SweepStatus::iterator ss_iterator;
typedef typename SweepStatus::value_type ss_pair;
typedef std::list<ITERATOR> IsoList;
typedef std::map<Point_2, IsoList*, lt_pnts_xy> EventQueue;
typedef typename EventQueue::iterator event_iterator;
typedef typename EventQueue::value_type event_pair;
typedef std::multimap<Point_2, ISegment, lt_pnts_xy> SegQueue;
typedef typename SegQueue::iterator seg_iterator;
typedef typename SegQueue::value_type ps_pair;
event_iterator event;
Point_2 p_sweep;
EventQueue XS;
seg_pair sl,sh;
SweepStatus YS;
SegQueue SQ;
IList Internal;
stl_seg_overlay_traits(const INPUT& in, OUTPUT& G,
const GEOMETRY& k) :
its(in.first), ite(in.second), GO(G), K(k),
XS(lt_pnts_xy(K)), YS(lt_segs_at_sweepline(p_sweep,&sl,&sh,K)),
SQ(lt_pnts_xy(K)) {}
@ We have similar helpers |source(s)|, |target(s)|, |original(s)|,
|orientation()|, |collinear()| as in the LEDA framework. We don't list them
again.
<<stl helping operations>>=
std::string dump_structures() const
{
std::ostrstream out;
out << "EventQueue:\n";
typename EventQueue::const_iterator sit1;
for(sit1 = XS.begin(); sit1 != XS.end(); ++sit1)
out << " " << sit1->first << std::endl;
out << "SegQueue:\n";
typename SegQueue::const_iterator sit2;
for(sit2 = SQ.begin(); sit2 != SQ.end(); ++sit2)
out << " " << sit2->first << " " << sit2->second
<< " " << sit2->first << std::endl;
out << "SweepStatus:\n";
typename SweepStatus::const_iterator sit3;
for( sit3 = YS.begin(); sit3 != YS.end(); ++sit3 )
out << sit3->first << " " << &*(sit3->second) << std::endl;
std::string res(out.str()); out.freeze(0); return res;
}
Point_2 source(ISegment is) const
{ return K.source(is->first); }
Point_2 target(ISegment is) const
{ return K.target(is->first); }
ITERATOR original(ISegment s) const
{ return s->second; }
int orientation(ss_iterator sit, const Point_2& p) const
{ return K.orientation(sit->first->first,p); }
bool collinear(ss_iterator sit1, ss_iterator sit2) const
{ Point_2 ps = source(sit2->first), pt = target(sit2->first);
return ( orientation(sit1,ps)==0 &&
orientation(sit1,pt)==0 );
}
@ Again there's a |compute_intersection| operations. Only the hash optimization
is left out.
<<stl operation for keeping the intersection invariant>>=
void compute_intersection(ss_iterator sit0)
{
// Given an item |sit0| in the Y-structure compute the point of
// intersection with its successor and (if existing) insert it into
// the event queue and do all necessary updates.
ss_iterator sit1 = sit0; ++sit1;
TRACEN("compute_intersection "<<sit0->first<<" "<<sit1->first);
if ( sit0 == YS.begin() || sit1 == --YS.end() ) return;
const Segment_2& s0 = sit0->first->first;
const Segment_2& s1 = sit1->first->first;
int or0 = K.orientation(s0,K.target(s1));
int or1 = K.orientation(s1,K.target(s0));
if ( or0 <= 0 && or1 >= 0 ) {
Point_2 q = K.intersection(s0,s1);
XS.insert(event_pair(q,0)); // only done if none existed!!!
}
}
@ Event handling in this version is similar to the LEDA implementation,
only that most algorithmic decisions are now based on geometric
predicate evaluation.
<<stl handling the event>>=
void process_event()
{
TRACEN("\n\n >>> process_event: "<<p_sweep);
Vertex_handle v = GO.new_vertex(p_sweep);
ss_iterator sit_succ, sit_pred, sit_first, sit;
Segment_2 s_sweep = K.construct_segment(p_sweep,p_sweep);
seg_pair sp(s_sweep,ITERATOR());
sit_succ = YS.upper_bound(&sp);
sit = sit_succ; --sit;
<<stl handling ending and passing segments>>
<<stl completing additional information of the new vertex>>
<<stl inserting new segments starting at nodes>>
<<stl enforcing the invariants for YS>>
}
@ We cannot reverse the bundle passing |p_sweep|, but have
to delete all and reinsert them.
<<stl handling ending and passing segments>>=
/* |sit| is determined by upper bounding the search for the
segment (p_sweep,p_sweep) and taking its predecessor.
if the segment associated to |sit| contains |p_sweep| then
there's a bundle of segments containing |p_sweep|.
We compute the successor (|sit_succ)|) and
predecessor (|sit_pred|) items. */
if ( sit == YS.begin() || orientation(sit,p_sweep) != 0 ) {
sit_pred = sit;
sit = YS.end();
}
/* If no segments contain p_sweep then sit_pred and sit_succ are
correctly set after the above locate operation, if a segment
contains p_sweep sit_pred and sit_succ are set below when
determining the bundle.*/
if ( sit != YS.end() ) { // sit->first is ending or passing segment
// Determine upper bundle item:
TRACEN("ending/passing segs");
/* Walk down until |sit_pred|, close edges for all segments
in the bundle, delete all segments in the bundle, but
reinsert the continuing ones */
std::list<ISegment> L_tmp;
bool overlapping;
do {
ISegment s = sit->first;
ss_iterator sit_next(sit); --sit_next;
overlapping = (sit_next != YS.begin()) && collinear(sit,sit_next);
Halfedge_handle e = sit->second;
if ( overlapping ) {
TRACEN("overlapping segment "<<s);
} else {
TRACEN("connecting edge to node "<<s);
GO.link_as_target_and_append(v,e);
/* in this case we close the output edge |e| associated to
|sit| by linking |v| as its target and by appending the
twin edge to |v|'s adjacency list. */
}
GO.supporting_segment(e,original(s));
if ( target(s) == p_sweep ) {
TRACEN("ending segment "<<s);
GO.ending_segment(v,original(s));
} else { // passing segment, take care of the node here!
TRACEN("passing segment "<<s);
L_tmp.push_back(s);
GO.passing_segment(v,original(s));
}
sit = sit_next;
}
while ( sit != YS.begin() && orientation(sit,p_sweep) == 0 );
sit_pred = sit_first = sit;
++sit_first; // first item of the bundle
TRACE("event bundles between\n "<<sit_succ->first);
TRACEN("\n "<<sit_pred->first);
/* Interfaceproposition for next chunk:
- succ(sit_pred) == sit_first == sit_succ
- bundle not empty: sit_first != sit_succ
*/
// delete and reinsert the continuing bundle
YS.erase(sit_first,sit_succ);
typename std::list<ISegment>::const_iterator lit;
for ( lit = L_tmp.begin(); lit != L_tmp.end(); ++lit ) {
YS.insert(sit_pred,ss_pair(*lit,Halfedge_handle()));
}
} // if (sit != ss_iterator() )
@ Node data association is similar again.
<<stl completing additional information of the new vertex>>=
assert( sit_pred != YS.end() );
GO.halfedge_below(v,sit_pred->second);
if ( event->second != 0 ) {
const IsoList& IL = *(event->second);
typename IsoList::const_iterator iso_it;
for (iso_it = IL.begin(); iso_it != IL.end(); ++iso_it)
GO.trivial_segment(v,*iso_it);
delete (event->second);
}
@ The insertion phase is even simpler.
<<stl inserting new segments starting at nodes>>=
ISegment next_seg;
seg_iterator next_it = SQ.begin();
while ( next_it != SQ.end() &&
( next_seg = next_it->second, p_sweep == source(next_seg)) ) {
TRACEN("inserting "<<next_seg);
YS.insert(ss_pair(next_seg,Halfedge_handle()));
GO.starting_segment(v,original(next_seg));
// delete minimum and assign new minimum to next_seg
SQ.erase(SQ.begin());
next_it = SQ.begin();
}
// we insert new edge stubs, non-linked at target
ss_iterator sit_curr = sit_succ, sit_prev = sit_succ;
for( --sit_curr; sit_curr != sit_pred;
sit_prev = sit_curr, --sit_curr ) {
TRACEN("checking outedge "<<sit_curr->first<<"\n "<<sit_prev->first);
if ( sit_curr != YS.begin() && sit_prev != --YS.end() &&
collinear(sit_curr,sit_prev) ) // overlapping
sit_curr->second = sit_prev->second;
else {
TRACEN("creating new edge");
sit_curr->second = GO.new_halfedge_pair_at_source(v);
}
}
sit_first = sit_prev;
@ For the intersection invariant we always calculate.
<<stl enforcing the invariants for YS>>=
// compute possible intersections between |sit_pred| and its
// successor and |sit_succ| and its predecessor
TRACEN("pred,succ = "<<sit_pred->first<<" "<<sit_succ->first);
compute_intersection(sit_pred);
sit = sit_succ; --sit;
if (sit != sit_pred)
compute_intersection(sit);
@ We realize the propositions of the invariants at the beginning in our
sweep initialization phase. We insert the endpoints of the segments
into |XS| and store an internal segment for each non trivial
segment. We store a geometric object and the input interator in our
|Internal| list. Then we store the internal segment associated to its
starting event in our segment queue |SQ|. A trivial segment is
associated to the corresponding event in the x-structure. To finish
the initialization we have to insert two sentinel segments into the
y-structure which avoids boundary checks in any |YS| search. Also
|SQ| obtains a trailing sentinel which is never reached during our
iteration over |XS|.
<<stl initialization of the sweep>>=
void initialize_structures()
{
/* INITIALIZATION
- insert all vertices into the x-structure
- insert sentinels into y-structure
- exploit the fact that insert operations into the x-structure
leave previously inserted points unchanged to achieve that
any pair of endpoints $p$ and $q$ with |p == q| are identical
*/
TRACEN("initialize_structures");
ITERATOR it_s;
for ( it_s=its; it_s != ite; ++it_s ) {
const Segment_2& s = *it_s;
event_iterator it1 = (XS.insert(event_pair(K.source(s),0))).first;
event_iterator it2 = (XS.insert(event_pair(K.target(s),0))).first;
// note that the STL only inserts if key is not yet in XS
if (it1 == it2) {
if ( it1->second == 0 ) it1->second = new IsoList;
it1->second->push_front(it_s);
continue; // ignore zero-length segments regarding YS
}
Point_2 p = it1->first;
Point_2 q = it2->first;
Segment_2 s1;
if ( K.compare_xy(p,q) < 0 )
s1 = K.construct_segment(p,q);
else
s1 = K.construct_segment(q,p);
Internal.push_back(seg_pair(s1,it_s));
SQ.insert(ps_pair(K.source(s1),&Internal.back()));
}
// insert a lower and an upper sentinel segment to avoid special
// cases when traversing the Y-structure
YS.insert(ss_pair(&sl,Halfedge_handle()));
YS.insert(ss_pair(&sh,Halfedge_handle()));
TRACEN("end of initialization\n");
}
@ Note the invariants of the sweep loop. The |event| has to be set
before the event action.
<<stl iteration control>>=
bool event_exists()
{
if (!XS.empty()) {
// event is set at end of loop and in init
event = XS.begin();
p_sweep = event->first;
return true;
}
return false;
}
void procede_to_next_event()
{ XS.erase(event); }
<<stl postprocessing of the sweep>>=
void complete_structures() {}
void check_invariants() {TRACEN("check_invariants\n"<<dump_structures());}
void check_final() {}
@ \section{Test Cases - Models for Geometry Kernels and Output}
The following models fulfill the requirements of the traits class
concept used in the above segment overlay framework (LEDA \emph{or}
STL). See the appendix for the concept |SegmentOverlayOutput|.
\subsection{The LEDA segment sweep traits model}
As the algorithm originally was hardwired to the LEDA types the
adaptation is not too hard. Note however that we have to integrate
some information storage into the OUTPUT structure. The OUTPUT is thus
a pair: a parameterized graph plus a node map associating edges to
nodes.
<<leda_overlay_traits.h>>=
#ifndef LEDA_OVERLAY_TRAITS_H
#define LEDA_OVERLAY_TRAITS_H
#if CGAL_LEDA_VERSION < 500
#include <LEDA/rat_kernel.h>
#include <LEDA/graph.h>
#include <LEDA/slist.h>
#else
#include <LEDA/geo/rat_kernel.h>
#include <LEDA/graph/graph.h>
#include <LEDA/core/slist.h>
#endif
template <typename T>
ostream& operator<<(ostream& os, const leda_slist<T>& s)
{ return os; }
template <typename T>
istream& operator>>(istream& is, leda_slist<T>& s)
{ return is; }
<<LEDA OUTPUT model>>
<<LEDA GEOMETRY model>>
#endif // LEDA_OVERLAY_TRAITS_H
<<LEDA OUTPUT model>>=
/*{\Moptions outfile=SegmentOverlayOutput.man}*/
/*{\Moptions section=subsection}*/
/*{\Moptions print_title=yes }*/
/*{\Msubst leda_graph_decorator#SegmentOverlayOutput}*/
/*{\Manpage {SegmentOverlayOutput}{}{Output traits for segment overlay}{G}}*/
template <typename ITERATOR>
class leda_graph_decorator {
public:
/*{\Mdefinition This is the plane map decorator concept for the |PMDEC|
template parameter of |PM_seg_overlay_traits|.}*/
/*{\Mtypes 3}*/
typedef leda_node Vertex_handle;
/*{\Mtypemember the vertex handle.}*/
typedef leda_edge Halfedge_handle;
/*{\Mtypemember the halfedge handle.}*/
typedef leda_rat_point Point_2;
/*{\Mtypemember embedding type. \precond |Point_2| equals |GEOM::Point_2|.}*/
typedef GRAPH< Point_2, leda_slist<ITERATOR> > Graph;
typedef leda_node_map<Halfedge_handle> Below_map;
Graph& G;
Below_map& M;
/*{\Mcreation}*/
leda_graph_decorator(Graph& Gi, Below_map& Mi) : G(Gi), M(Mi) {}
/*{\Mtext The output model must be copy constructible.
Let |G| be an object of type |\Mname|.}*/
/*{\Moperations}*/
Vertex_handle new_vertex(const Point_2& p)
/*{\Mop creates a new vertex in the output structure and embeds it
via the point |p|.}*/
{ return G.new_node(p); }
void link_as_target_and_append(Vertex_handle v, Halfedge_handle e)
/*{\Mop makes |v| the target of |e| and appends the twin of
|e| (its reversal edge) to |v|'s adjacency list.}*/
{ Halfedge_handle erev = G.reversal(e);
G.move_edge(e,G.cyclic_adj_pred(e,G.source(e)),v);
G.move_edge(erev,v,G.target(erev));
}
Halfedge_handle new_halfedge_pair_at_source(Vertex_handle v)
/*{\Mop returns a newly created edge inserted before the first
edge of the adjacency list of |v|. It also creates a reversal edge
whose target is |v|.}*/
{ Halfedge_handle e_res,e_rev, e_first = G.first_adj_edge(v);
if ( e_first == nil ) {
e_res = G.new_edge(v,v);
e_rev = G.new_edge(v,v);
} else {
e_res = G.new_edge(e_first,v,LEDA::before);
e_rev = G.new_edge(e_first,v,LEDA::before);
}
G.set_reversal(e_res,e_rev);
return e_res;
}
/*{\Mtext \headerline{Additional sweep information}
The iterator type |ITERATOR| has to be the same type as the first type
parameter of |Segment_overlay_traits|.}*/
void supporting_segment(Halfedge_handle e, ITERATOR it)
/*{\Mop the non-trivial segment |*it| supports the edge |e|.}*/
{ G[e].append(it); }
void trivial_segment(Vertex_handle v, ITERATOR it)
/*{\Mop the trivial segment |*it| supports vertex |v|.}*/
{ }
void halfedge_below(Vertex_handle v, Halfedge_handle e)
/*{\Mop associates the edge |e| as the edge below |v|.}*/
{ M[v] = e; }
void starting_segment(Vertex_handle v, ITERATOR it)
/*{\Mop the segment |*it| starts in |v|.}*/
{}
void passing_segment(Vertex_handle v, ITERATOR it)
/*{\Mop the segment |*it| passes |v| (contains it in its
relative interior) .}*/
{}
void ending_segment(Vertex_handle v, ITERATOR it)
/*{\Mop the segment |*it| ends in |v|.}*/
{}
}; // leda_graph_decorator
<<LEDA GEOMETRY model>>=
/*{\Moptions outfile=SegmentOverlayGeometry_2.man}*/
/*{\Moptions print_title=yes }*/
/*{\Moptions section=subsection}*/
/*{\Msubst leda_geometry#SegmentOverlayGeometry_2}*/
/*{\Manpage {SegmentOverlayGeometry_2}{}
{Geometric traits for segment overlay}{K}}*/
class leda_geometry {
public:
/*{\Mdefinition This is the geometry concept for the |GEOM|
template parameter of |PM_seg_overlay_traits|.}*/
/*{\Mtypes 3}*/
typedef leda_rat_point Point_2;
/*{\Mtypemember the point type to embed the vertices. It has to be the same
type as the }*/
typedef leda_rat_segment Segment_2;
/*{\Mtypemember the segment type (input).}*/
/*{\Mcreation 4}*/
leda_geometry() {}
/*{\Mtext The kernel must be default and copy constructible.
Let |K| be an object of type |\Mname|.}*/
/*{\Moperations 2 2}*/
Point_2 source(const Segment_2& s) const
/*{\Mop returns the source of |s|.}*/
{ return s.source(); }
Point_2 target(const Segment_2& s) const
/*{\Mop returns the target of |s|.}*/
{ return s.target(); }
Segment_2 construct_segment(const Point_2& p, const Point_2& q) const
/*{\Mop creates a segment |(p,q)|.}*/
{ return Segment_2(p,q); }
int orientation(const Segment_2& s, const Point_2& p) const
/*{\Mop returns the orientation of |p| with respect to |s|: a value $<0
(=0,>0)$ iff |p| is right of (on, left of) the line through |s| from source to
target.}*/
{ return ::orientation(s,p); }
bool is_degenerate(const Segment_2& s) const
/*{\Mop returns true iff |s| is trivial (source equals target).}*/
{ return s.is_trivial(); }
int compare_xy(const Point_2& p1, const Point_2& p2) const
/*{\Mop determines the lexicographic order of points.}*/
{ return Point_2::cmp_xy(p1,p2); }
Point_2 intersection(const Segment_2& s1, const Segment_2& s2) const
/*{\Mop determines the intersection point of the lines through |s1| and
|s2|. \precond the operation can assume that the lines do intersect.}*/
{ Point_2 p;
s1.intersection_of_lines(s2,p);
return p;
}
}; // leda_geometry
@ \subsection{The visualization model}
We visualize the sweep by means of our generic sweep observer, which
obtains messages at four events: after initialization, just before a
sweep event, just after a sweep event, and after completion. We
visualize in a LEDA window. After the intialization we draw all input
segments. Before the event we place the sweep line at the event
point. After the event we draw all created edges ending at the
event. The following class |leda_visualization| can be plugged into
the sweep observer class. For the model see |GenericSweepVisualization|
in the appendix (Section \ref{GenericSweepVisualization}).
<<leda_overlay_visualization.h>>=
#ifndef LEDA_OVERLAY_VISUALIZATION_H
#define LEDA_OVERLAY_VISUALIZATION_H
#if CGAL_LEDA_VERSION < 500
#include <LEDA/rat_window.h>
#else
#include <LEDA/graphics/rat_window.h>
#endif
template <class GENSWEEPTRAITS>
class leda_visualization {
leda_window W;
public:
typedef leda_window VDEVICE;
typedef typename GENSWEEPTRAITS::Halfedge_handle Halfedge_handle;
typedef typename GENSWEEPTRAITS::Vertex_handle Vertex_handle;
typedef typename GENSWEEPTRAITS::OUTPUT Graph;
typedef typename GENSWEEPTRAITS::Segment_2 Segment_2;
typedef typename GENSWEEPTRAITS::Point_2 Point_2;
typedef typename GENSWEEPTRAITS::ISegment ISegment;
leda_visualization() : W()
{ W.set_show_coordinates(true);
W.init(-50,50,-50,1);
W.set_node_width(3);
W.set_line_width(2);
W.display();
}
void draw(const Point_2& p, leda_color c)
{ W.draw_filled_node(p.xcoordD(),p.ycoordD(),c); }
void draw(const Segment_2& s, leda_color c)
{ W.draw_segment(s.xcoord1D(),s.ycoord1D(),
s.xcoord2D(),s.ycoord2D(),c);}
void post_init_animation(GENSWEEPTRAITS& gpst)
{ typename GENSWEEPTRAITS::ITERATOR it;
for (it = gpst.its; it != gpst.ite; ++it )
if ( (*it).is_trivial() )
draw((*it).source(),leda_blue);
else
draw(*it,leda_blue);
W.read_mouse();
}
void draw_sl(const leda_point& p)
{ leda_drawing_mode mold = W.set_mode(leda_xor_mode);
W.draw_filled_node(p,leda_red);
W.draw_vline(p.xcoord(),leda_red);
W.set_mode(mold);
}
void pre_event_animation(GENSWEEPTRAITS& gpst)
{ draw_sl(gpst.p_sweep.to_point()); }
void post_event_animation(GENSWEEPTRAITS& gpst)
{ typename GENSWEEPTRAITS::OUTPUT::Graph& G(gpst.GO.G);
Halfedge_handle e;
Vertex_handle v = G.last_node();
forall_out_edges(e, v) {
if ( source(e)!=target(e) ) {
draw(leda_rat_segment(G[source(e)],G[target(e)]),leda_black);
}
}
draw(gpst.p_sweep,leda_black);
W.read_mouse();
draw_sl(gpst.p_sweep.to_point());
}
void post_completion_animation(GENSWEEPTRAITS& gpst)
{ typename GENSWEEPTRAITS::OUTPUT::Graph& G(gpst.GO.G);
Halfedge_handle e;
forall_edges(e,G) {
draw(leda_rat_segment(G[source(e)],G[target(e)]),leda_black);
}
Vertex_handle v;
forall_nodes(v,G) {
draw(G[v],leda_black);
}
}
VDEVICE& device() { return W; }
};
#endif //LEDA_OVERLAY_VISUALIZATION_H
@ \subsection{Testing the generic LEDA sweep}
<<leda_segment_overlay-test.c>>=
#if CGAL_LEDA_VERSION < 500
#include <LEDA/rat_window.h>
#else
#include <LEDA/graphics/rat_window.h>
#endif
#include "Segment_overlay_traits.h"
#include "leda_overlay_traits.h"
#include "leda_overlay_visualization.h"
#if CGAL_LEDA_VERSION < 500
#include <LEDA/stream.h>
#include <LEDA/param_handler.h>
#else
#include <LEDA/system/stream.h>
#include <LEDA/system/param_handler.h>
#endif
using namespace CGAL;
typedef leda_list<leda_rat_segment>::iterator Seg_iterator;
typedef CGAL::Segment_overlay_traits< Seg_iterator,
leda_graph_decorator<Seg_iterator>, leda_geometry>
leda_seg_sweep_traits;
typedef leda_visualization<leda_seg_sweep_traits> leda_seg_vis;
typedef CGAL::generic_sweep<leda_seg_sweep_traits> leda_seg_sweep;
typedef CGAL::sweep_observer<leda_seg_sweep,leda_seg_vis> leda_seg_sweep_obs;
int main(int argc, char** argv)
{
// SETDTHREAD(23);
leda_param_handler P(argc,argv,".overlay");
P.add_parameter("inputfile1:-i:string:e1");
leda_param_handler::init_all();
leda_string f1;
P.get_parameter("-i",f1);
leda_seg_sweep_obs Obs;
leda_list<leda_rat_segment> L;
leda_rat_segment s;
ifstream if1(f1);
if (if1) {
if1 >> L;
forall (s,L)
if ( s.is_trivial() ) Obs.device() << s.source().to_point();
else Obs.device() << s.to_segment();
}
while (Obs.device()>>s) {
L.append(s);
}
ofstream of("input.log");
forall(s,L) of << s;
of.close();
cout << "Starting " << CGAL::sweepversion << endl;
leda_seg_sweep_traits::OUTPUT::Graph G;
leda_seg_sweep_traits::OUTPUT::Below_map E_below(G);
leda_seg_sweep_traits::OUTPUT O(G,E_below);
leda_seg_sweep SI(leda_seg_sweep_traits::INPUT(L.begin(),L.end()),O);
Obs.attach(SI);
SI.sweep();
cout << "Edges and Segments:\n";
leda_edge e;
forall_edges(e,G) {
cout <<" ("<<G[source(e)].to_point()<<G[target(e)].to_point()<< ") ";
leda_slist< leda_list<leda_rat_segment>::iterator >& SL(G[e]);
slist_item it;
for ( it = SL.first(); it; it=SL.succ(it) )
cout << *(SL[it]) << " ";
cout << endl;
}
cout << "\nNodes and Edges below\n";
leda_node v;
forall_nodes(v,G) {
cout <<" "<<G[v]<<" ";
if ( E_below[v] ) {
leda_edge e = E_below[v];
cout << "("<<G[source(e)]<<G[target(e)]<< ")\n";
} else
cout << "nil\n";
}
Obs.device().read_mouse();
return 0;
}
@ We want to compare the runtimes of the two instantiations of our
generic segment overlay and additionally the original LEDA 4.1 code
|SWEEP_SEGMENTS|. The inputs are randomly distributed segments in a
square. The runtimes are in seconds on an Ultra 2 200 MHz.
\begin{center}
\input{inputs/seg_overlay_rt.tbl}
\end{center}
<<leda_segment_overlay-rt.c>>=
#define INCLUDEBOTH // to allow comparison
#include "Segment_overlay_traits.h"
#include "leda_overlay_traits.h"
#if CGAL_LEDA_VERSION < 500
#include <LEDA/param_handler.h>
#include <LEDA/random_rat_point.h>
#include <LEDA/plane_alg.h>
#else
#include <LEDA/system/param_handler.h>
#include <LEDA/geo/random_rat_point.h>
#include <LEDA/geo/plane_alg.h>
#endif
typedef leda_list<leda_rat_segment>::iterator Seg_iterator;
typedef CGAL::leda_seg_overlay_traits< Seg_iterator,
leda_graph_decorator<Seg_iterator>,leda_geometry>
leda_seg_sweep_traits;
typedef CGAL::generic_sweep<leda_seg_sweep_traits> leda_seg_sweep;
typedef CGAL::stl_seg_overlay_traits< Seg_iterator,
leda_graph_decorator<Seg_iterator>,leda_geometry>
stl_seg_sweep_traits;
typedef CGAL::generic_sweep<stl_seg_sweep_traits> stl_seg_sweep;
template <typename T>
bool equal(const leda_list<T>& L1, const leda_list<T>& L2)
{ leda_list<T>::const_iterator it1,it2;
for (it1 = L1.begin(), it2 = L2.begin();
it1 != L1.end() && it2 != L2.end();
++it1, ++it2 )
if ( *it1 != *it2 ) return false;
return (it1 == L1.end()) && (it2 == L2.end());
}
#define OUTPUT(t) \
if ( !t ) cerr <<" "<<#t<<" = "<<(t)<< endl;
int main(int argc, char** argv)
{
//SETDTHREAD(23);
int n = 100;
if ( argc > 2 ) {
cout << "usage: " << argv[0] << " #segments\n";
return 1;
}
if ( argc == 2 ) {
n = atoi(argv[1]);
}
cout << "\\begin{tabular}[t]{llll}\n";
cout << "\\#segs & LEDA classic & LEDA generic & STL generic \\\\ \\hline\n";
for (int i = 1; i < 6; n*=2,++i ) {
leda_list<leda_rat_point> LP;
leda_list<leda_rat_segment> L;
leda_rat_segment s;
random_points_in_square(2*n,100,LP);
list_item pit1 = LP.first(),pit2 = (pit1 ? LP.succ(pit1) : 0);
while ( pit1&&pit2 ) {
L.append(leda_rat_segment(LP[pit1],LP[pit2]));
pit1 = LP.succ(pit2);
pit2 = (pit1 ? LP.succ(pit1) : 0);
}
ofstream of("input.log");
of << L;
of.close();
leda_seg_sweep_traits::OUTPUT::Graph G1,G2;
leda_seg_sweep_traits::OUTPUT::Below_map E1(G1),E2(G2);
leda_seg_sweep_traits::OUTPUT O1(G1,E1),O2(G2,E2);
GRAPH<leda_rat_point,leda_rat_segment> G0;
float t0 = used_time();
SWEEP_SEGMENTS(L,G0,true);
t0 = used_time(t0);
leda_seg_sweep SSW1(leda_seg_sweep_traits::INPUT(L.begin(),L.end()),O1);
float t1 = used_time();
SSW1.sweep();
t1 = used_time(t1);
stl_seg_sweep SSW2(leda_seg_sweep_traits::INPUT(L.begin(),L.end()),O2);
float t2 = used_time();
SSW2.sweep();
t2 = used_time(t2);
cout << n << " & " << t0 << " & " << t1 << " & " << t2 << "\\\\\n";
leda_list<leda_rat_point> EP1,EP2;
leda_node v;
forall_nodes(v,G1) EP1.append(G1[v]);
forall_nodes(v,G2) EP2.append(G2[v]);
leda_list<leda_rat_segment> ES1,ES2;
leda_edge e;
forall_edges(e,G1)
ES1.append(leda_rat_segment(G1[source(e)],G1[target(e)]));
forall_edges(e,G2)
ES2.append(leda_rat_segment(G2[source(e)],G2[target(e)]));
OUTPUT(equal(EP1,EP2));
OUTPUT(equal(ES1,ES2));
} // end of for
cout << "\\end{tabular}\n";
return 0;
}
@ \subsection{A CGAL segment sweep traits model}
We define a traits model for the segment overlay sweep which calculates
the overlay of CGAL 2-dimensional kernel segments. We have
\begin{itemize}
\item input is an iterator range of CGAL segments
\item output is a halfedge data structure decorator creating
the overlay of the segments in the decorated HDS.
\end{itemize}
<<HDS_decorator.h>>=
#ifndef CGAL_HDS_DECORATOR_H
#define CGAL_HDS_DECORATOR_H
#include <CGAL/basic.h>
#include <CGAL/PM_decorator_simple.h>
template <typename HDS, typename I>
class HDS_decorator : public CGAL::PM_decorator_simple<HDS> {
public:
typedef CGAL::PM_decorator_simple<HDS> Base;
typedef HDS Plane_map;
typedef typename Base::Point_2 Point_2;
typedef typename Base::Vertex_handle Vertex_handle;
typedef typename Base::Halfedge_handle Halfedge_handle;
typedef I ITERATOR;
HDS_decorator(HDS& H) : Base(H) {}
Vertex_handle new_vertex(const Point_2& p) const
{ Vertex_handle v = Base::new_vertex();
v->point() = p; return v; }
void link_as_target_and_append(Vertex_handle v, Halfedge_handle e) const
{ Base::link_as_target_and_append(v,e); }
Halfedge_handle new_halfedge_pair_at_source(Vertex_handle v) const
{ return Base::new_halfedge_pair_at_source(v,Base::BEFORE); }
void supporting_segment(Halfedge_handle e, ITERATOR it) const {}
void trivial_segment(Vertex_handle v, ITERATOR it) const {}
void halfedge_below(Vertex_handle v, Halfedge_handle e) const {}
void starting_segment(Vertex_handle v, ITERATOR it) const {}
void passing_segment(Vertex_handle v, ITERATOR it) const {}
void ending_segment(Vertex_handle v, ITERATOR it) const {}
}; // HDS_decorator
#endif // CGAL_HDS_DECORATOR_H
@ \subsection{Affine geometry wrapper}
<<Affine_geometry.h>>=
#ifndef CGAL_AFFINE_GEOMETRY_H
#define CGAL_AFFINE_GEOMETRY_H
#include <CGAL/basic.h>
#include <CGAL/Point_2.h>
#include <CGAL/intersections.h>
#include <CGAL/squared_distance_2.h>
namespace CGAL {
/*{\Moptions outfile=AffineGeometryTraits_2.man}*/
/*{\Msubst Affine_geometry AffineGeometryTraits_2}*/
/*{\Manpage {AffineGeometryTraits_2}{}{An affine kernel traits} {K}}*/
template <typename R>
struct Affine_geometry : public R {
/*{\Mdefinition |\Mname| is a kernel concept providing affine
geometry. The concept specifies geometric types, predicates,
and constructions.}*/
/*{\Mtypes 5}*/
/*{\Mtext Local types are |Point_2|, |Segment_2|, |Direction_2|, |Line_2|,
the ring type |RT|, and the field type |FT|. See the CGAL 2d kernel for
a description of |RT| and |FT|.}*/
typedef typename R::Point_2 Point_2;
typedef typename R::Segment_2 Segment_2;
typedef typename R::Direction_2 Direction_2;
typedef typename R::Line_2 Line_2;
typedef typename R::RT RT;
typedef typename R::FT FT;
/*{\Mcreation 5}*/
/*{\Mtext The kernel must be default and copy constructible.}*/
Affine_geometry() : R() {}
Affine_geometry(const Affine_geometry& Gi) : R(GI) {}
/*{\Moperations 2}*/
Point_2 source(const Segment_2& s) const
/*{\Mop returns the source point of |s|.}*/
{ typename R::Construct_source_point_2 _source =
construct_source_point_2_object();
return _source(s); }
Point_2 target(const Segment_2& s) const
/*{\Mop returns the target point of |s|.}*/
{ typename R::Construct_target_point_2 _target =
construct_target_point_2_object();
return _target(s); }
int orientation(const Point_2& p1, const Point_2& p2, const Point_2& p3) const
/*{\Mop returns the orientation of |p3| with respect to the line
through |p1p2|.}*/
{ typename R::Orientation_2 _orientation =
orientation_2_object();
return static_cast<int> ( _orientation(p1,p2,p3) );
}
bool left_turn(const Point_2& p1, const Point_2& p2, const Point_2& p3) const
/*{\Mop return true iff the |p3| is left of the line through |p1p2|.}*/
{ return (orientation(p1,p2,p3) > 0); }
bool is_degenerate(const Segment_2& s) const
/*{\Mop return true iff |s| is degenerate.}*/
{ typename R::Is_degenerate_2 _is_degenerate =
is_degenerate_2_object();
return _is_degenerate(s); }
int compare_xy(const Point_2& p1, const Point_2& p2) const
/*{\Mop returns the lexicographic order of |p1| and |p2|.}*/
{ typename R::Compare_xy_2 _compare_xy =
compare_xy_2_object();
return static_cast<int>( _compare_xy(p1,p2) );
}
int compare_x(const Point_2& p1, const Point_2& p2) const
/*{\Mop returns the order on the $x$-coordinates of |p1| and |p2|.}*/
{ typename R::Compare_x_2 _compare_x =
compare_x_2_object();
return static_cast<int>( _compare_x(p1,p2) );
}
int compare_y(const Point_2& p1, const Point_2& p2) const
/*{\Mop returns the order on the $y$-coordinates of |p1| and |p2|.}*/
{ typename R::Compare_y_2 _compare_y =
compare_y_2_object();
return static_cast<int>( _compare_y(p1,p2) );
}
bool first_pair_closer_than_second(const Point_2& p1,
const Point_2& p2, const Point_2& p3, const Point_2& p4) const
/*{\Mop returns |true| iff the distance |p1p2| is smaller than
the distance |p3p4|.}*/
{ return ( squared_distance(p1,p2) < squared_distance(p3,p4) ); }
bool strictly_ordered_along_line(
const Point_2& p1, const Point_2& p2, const Point_2& p3) const
/*{\Mop returns |true| iff |p2| is in the relative interior of the
segment |p1p3|.}*/
{ typename R::Are_strictly_ordered_along_line_2 _ordered =
are_strictly_ordered_along_line_2_object();
return _ordered(p1,p2,p3);
}
Segment_2 construct_segment(const Point_2& p, const Point_2& q) const
/*{\Mop constructs a segment |pq|.}*/
{ typename R::Construct_segment_2 _segment =
construct_segment_2_object();
return _segment(p,q); }
int orientation(const Segment_2& s, const Point_2& p) const
/*{\Mop returns the orientation of |p| with respect to the line
through |s|.}*/
{ typename R::Orientation_2 _orientation =
orientation_2_object();
return static_cast<int> ( _orientation(source(s),target(s),p) );
}
bool contains(const Segment_2& s, const Point_2& p) const
/*{\Mop returns true iff |s| contains |p|.}*/
{ typename R::Has_on_2 _contains = has_on_2_object();
return _contains(s,p);
}
Point_2 intersection(
const Segment_2& s1, const Segment_2& s2) const
/*{\Mop returns the point of intersection of the lines supported by |s1| and |s2|.
The algorithm asserts that this intersection point exists.}*/
{ typename R::Intersect_2 _intersect =
intersect_2_object();
typename R::Construct_line_2 _line =
construct_line_2_object();
Point_2 p;
CGAL::Object result =
_intersect(_line(s1),_line(s2));
if ( !CGAL::assign(p, result) )
CGAL_assertion_msg(false,"intersection: no intersection.");
return p;
}
// additional interface for constr triang:
Direction_2 construct_direction(
const Point_2& p1, const Point_2& p2) const
/*{\Mop returns the direction of the vector |p2| - |p1|.}*/
{ typename R::Construct_direction_of_line_2 _direction =
construct_direction_of_line_2_object();
typename R::Construct_line_2 _line =
construct_line_2_object();
return _direction(_line(p1,p2)); }
bool strictly_ordered_ccw(const Direction_2& d1,
const Direction_2& d2, const Direction_2& d3) const
/*{\Mop returns |true| iff |d2| is in the interior of the
counterclockwise angular sector between |d1| and |d3|.}*/
{ return d2.counterclockwise_in_between (d1, d3); }
// additional interface for point location:
}; // Affine_geometry<R>
} //namespace CGAL
#endif //CGAL_AFFINE_GEOMETRY_H
@ \subsection{Testing the CGAL sweep}
<<cgal_segment_overlay-test.C>>=
#include <CGAL/basic.h>
#include <fstream>
#include <list>
#include <algorithm>
#include <CGAL/Cartesian.h>
#include <CGAL/Homogeneous.h>
#include <CGAL/leda_integer.h>
#include <CGAL/point_generators_2.h>
#include <CGAL/function_objects.h>
#include <CGAL/Join_input_iterator.h>
#include <CGAL/algorithm.h>
#include <CGAL/Halfedge_data_structure_default.h>
#include "Segment_overlay_traits.h"
#include "HDS_decorator.h"
#include "Affine_geometry.h"
// GEOMETRY:
typedef CGAL::Homogeneous<leda_integer> Kernel;
//typedef CGAL::Cartesian<double> Kernel;
typedef CGAL::Affine_geometry<Kernel> Aff_kernel;
typedef Kernel::Point_2 Point_2;
typedef Kernel::Segment_2 Segment_2;
// INPUT:
typedef std::list<Segment_2>::const_iterator SIterator;
// OUTPUT:
typedef CGAL::Halfedge_data_structure_default<Point_2> HDS;
typedef HDS_decorator<HDS,SIterator> HDS_dec;
// SWEEP INSTANTIATION:
typedef CGAL::Segment_overlay_traits<SIterator,HDS_dec,Aff_kernel>
CGAL_seg_sweep_traits;
typedef CGAL::generic_sweep<CGAL_seg_sweep_traits> CGAL_seg_sweep;
typedef CGAL::Creator_uniform_2<int,Point_2> Pnt_creator;
typedef CGAL::Creator_uniform_2<Point_2, Segment_2> Seg_creator;
typedef CGAL::Random_points_in_square_2<Point_2,Pnt_creator> Pnt_iterator;
typedef CGAL::Join_input_iterator_2< Pnt_iterator,Pnt_iterator,Seg_creator>
Seg_iterator;
using std::cout;
using std::endl;
int main(int argc, char** argv)
{
SETDTHREAD(23);
CGAL::set_pretty_mode ( cout );
std::list<Segment_2> L;
std::list<Segment_2>::iterator lit;
Segment_2 s;
if (argc != 1) {
std::ifstream if1( argv[1]);
std::cout << "INPUT:\n";
if ( if1)
while (if1 >> s) {
L.push_back(s);
cout << s << ", ";
}
cout << endl;
}
Pnt_iterator p1(100),p2(100);
Seg_iterator sit(p1,p2);
for (int i = 0; i < 100; ++i,++sit)
L.push_back(*sit);
// CGAL::cpp0x::copy_n(sit,100, std::back_inserter(L));
std::ofstream of("input.log");
for (lit = L.begin(); lit!= L.end(); ++lit)
of << *lit;
of.close();
cout << "Starting " << CGAL::sweepversion << endl;
HDS G;
CGAL_seg_sweep_traits::OUTPUT P(G);
CGAL_seg_sweep SI(CGAL_seg_sweep::INPUT(L.begin(),L.end()),P);
SI.sweep();
cout << "\nVertices\n";
HDS::Halfedge_iterator eit;
HDS::Vertex_iterator vit;
for( vit=G.vertices_begin(); vit != G.vertices_end();
++vit ) {
cout <<" "<<vit->point() << endl;
}
cout << "Edges:\n";
for( eit=G.halfedges_begin(); eit != G.halfedges_end();
++(++eit) ) {
cout <<" "<< eit->opposite()->vertex()->point() << "," <<
eit->vertex()->point() << endl;
}
return 0;
}
<<random_overlay-test.c>>=
#define CGAL_USE_LEDA
#include "Segment_overlay_traits.h"
#include "leda_overlay_traits.h"
#include "leda_overlay_visualization.h"
#if CGAL_LEDA_VERSION < 500
#include <LEDA/stream.h>
#include <LEDA/rat_window.h>
#include <LEDA/param_handler.h>
#include <LEDA/random_rat_point.h>
#else
#include <LEDA/system/stream.h>
#include <LEDA/graphics/rat_window.h>
#include <LEDA/system/param_handler.h>
#include <LEDA/geo/random_rat_point.h>
#endif
typedef leda_list<leda_rat_segment>::iterator Seg_iterator;
typedef CGAL::Segment_overlay_traits< Seg_iterator,
leda_graph_decorator<Seg_iterator>, leda_geometry>
leda_seg_sweep_traits;
typedef leda_visualization<leda_seg_sweep_traits> leda_seg_vis;
typedef CGAL::generic_sweep<leda_seg_sweep_traits> leda_seg_sweep;
typedef CGAL::sweep_observer<leda_seg_sweep,leda_seg_vis> leda_seg_sweep_obs;
int main(int argc, char** argv)
{
// SETDTHREAD(23);
leda_param_handler P(argc,argv,".overlay");
P.add_parameter("inputfile:-i:string:");
P.add_parameter("segs:-n:int:100");
leda_param_handler::init_all();
leda_string f1;
int n;
P.get_parameter("-i",f1);
P.get_parameter("-n",n);
leda_window W;
W.init(-500,500,-500);
W.set_show_coordinates(true);
W.display();
leda_list<leda_rat_segment> L;
leda_list<leda_rat_point> LP;
leda_rat_segment s;
leda_rat_point p,po;
ifstream if1(f1);
if (if1) {
if1 >> L;
forall (s,L)
if ( s.is_trivial() ) W << s.source().to_point();
else W << s.to_segment();
} else {
random_points_in_square(n, n, LP);
bool first = true;
forall(p,LP) {
if (first) { po = p; first = false; }
else {
L.append(leda_rat_segment(p,po)); first=true;
W << leda_rat_segment(p,po);
}
}
}
ofstream of("input.log");
of << L;
of.close();
cout << "Starting " << CGAL::sweepversion << endl;
leda_seg_sweep_traits::OUTPUT::Graph G;
leda_seg_sweep_traits::OUTPUT::Below_map E_below(G);
leda_seg_sweep_traits::OUTPUT O(G,E_below);
leda_seg_sweep SI(leda_seg_sweep_traits::INPUT(L.begin(),L.end()),O);
//Obs.attach(SI);
SI.sweep();
cout << "\nNodes and Edges below\n";
leda_node v;
forall_nodes(v,G) {
cout <<" "<<G[v]<<" ";
if ( E_below[v] ) {
leda_edge e = E_below[v];
cout << "("<<G[source(e)]<<G[target(e)]<< ")\n";
} else
cout << "nil\n";
}
W.read_mouse();
return 0;
}
@ \begin{ignore}
<<CGAL Header>>=
// ============================================================================
//
// Copyright (c) 1997-2000 The CGAL Consortium
//
// This software and related documentation is part of an INTERNAL release
// of the Computational Geometry Algorithms Library (CGAL). It is not
// intended for general use.
//
// ----------------------------------------------------------------------------
//
// release : $CGAL_Revision$
// release_date : $CGAL_Date$
//
// file : include/CGAL/Nef_2/Segment_overlay_traits.h
// package : Nef_2
// chapter : Nef Polyhedra
//
// source : nef_2d/Segment_overlay.lw
// revision : $Id$
// revision_date : $Date$
//
// author(s) : Michael Seel <seel@mpi-sb.mpg.de>
// maintainer : Michael Seel <seel@mpi-sb.mpg.de>
// coordinator : Michael Seel <seel@mpi-sb.mpg.de>
//
// implementation: generic segment intersection sweep
// ============================================================================
@ \end{ignore}
\end{ignoreindiss}
%KILLSTART DISS REP
\bibliographystyle{alpha}
\bibliography{general,diss,comp_geo}
\newpage
\section{Appendix}\label{appendix}
\input manpages/generic_sweep.man
\input manpages/GenericSweepTraits.man
\input manpages/GenericSweepVisualization.man
\input manpages/SegmentOverlayOutput.man
\input manpages/SegmentOverlayGeometry_2.man
\end{document}
%KILLEND DISS REP