%------------------------------------------------------------------------------ %KILLSTART DISS REP %LDEL TRACE.*?\)\; %LDEL CGAL_BEGIN_NAMESPACE.* %LDEL CGAL_END_NAMESPACE.* \documentclass[a4paper]{article} \usepackage{MyLweb} \input{defs} \excludeversion{ignoreindiss} \excludeversion{ignore} \begin{document} \title{The constrained triangulation of a plane map} \author{M. Seel} \date{} \maketitle \begin{abstract} We implement a plane sweep algorithm to compute the constrained triangulation of a plane map. We use a generic sweep approach which is based on a geometric kernel concept and on a plane map decorator concept. \end{abstract} \tableofcontents %KILLEND REP \section{Introduction} We first introduce the notions or triangulation and constraining segments. Then we present the algorithmic ideas to construct a constrained triangulation by a plane sweep algorithm. In the next section the concrete implementation is shown. \begin{deff}[Triangulation \cite{prep-sham:CG}] A planar subdivision is a triangulation $T$ if all its bounded regions are triangles. A triangulation of a finite set $S$ of points is a planar graph $T(S)$ with the maximum number of edges, which is equivalent to saying that $T(S)$ is obtained by joining the points of $S$ by non-intersecting straight line segments so that every region internal to the convex hull of $S$ is a triangle. \end{deff} \begin{deff}[Constrained Triangulation] Let $S$ be a finite set of segments which have no pairwise common points except in their endpoints. A constrained triangulation $CT(S)$ is a triangulation of the endpoints of the segments in $S$ such that each non-trivial segment of $S$ corresponds to an edge in $CT(S)$. Note that we allow trivial segments. \end{deff} \newcommand{\seg}{\mathit{segment}} \newcommand{\point}{\mathit{point}} \newcommand{\rev}{\mathit{rev}} \newcommand{\source}{\mathit{source}} \newcommand{\target}{\mathit{target}} We sweep the edges and vertices of of a plane map $G=(V,E)$ to compute the constrained triangulation $CT(G) := CT(S = \seg(E)\cup \point(V))$\footnote{We slopily write $\point(V)$ for the set of all embedding points of the vertices in $V$ and $\seg(E)$ for the set of segments spanned by the end vertices of the edges in $E$.}. Actually we extend $G_{in}$ until finally $G_{out} \equiv CT(G_{in})$. What does this mean? We want to obtain a triangulation of the point set $\point(V)$ which additionally obeys the constraints given by the segments $\seg(E)$ of the graph. In the following description we identify the vertices with their positions and the edges with the segments spanned by the positions of their end points. If the graph only contains isolated vertices we calculate the standard triangulation given by the standard sweep triangulation algorithm. (see for example a description in the LEDA book \cite{ledabook}). This triangulation has the property that any vertical line $l$ subdivides the triangulation into a correctly constructed part left of it and the rest. (Fixing $l$ just remove all triangles which are intersected by it in their interior or which are right of it). This is just the result of the invariant kept during the sweep. So what is different if we calculate the constrained triangulation? Our sweepline |SL| contains all edges (segments) that currently intersect it ordered by their intersection points from bottom to top. We use a sorted sequence to store the edges. Let's conceptually identify the line with the data structure. Note that edges might touch in their endpoints as we work on a correctly embedded plane map. Also the vertices do not come to lie in the interior of a non-adjacent edge. We sweep the graph $G$ and extend it by additional edges to create the triangulation. Each input edge encountered during the sweep from $-\infty$ to $\infty$ (along the x-axis) enters |SL| that stores its position in the sweepline. At the top and the bottom of |SL| we use sentinel segments to avoid the handling of boundary cases. Note that with the above idea all edges in the sweepline are already connected via face cycles (paths) in the extended graph. At the end of the sweep we remove the sentinel edges. What is our basic invariant? We have already calculated a correct constrained triangulation of all vertices and edges that are fully left of the sweepline including a closure implied by the sweepline. Let's first assume that we don't face degeneracies like equal $x$-coordinates or vertical segments. We will comment on the removal of this restriction later. \begin{deff}[Restricted Constrained Triangulation] Let |SL| be the vertical sweepline defined by an event point $p$ in the plane. Let $S[<|SL|] \subseteq S$ be the finite set of segments already handled completely (all vertices and edges of our input graph fully left of |SL|). Let $S[SL]$ be the set of segments spanned by edges intersecting |SL|, but restricted to the closed halfplane left of |SL|. We define the restricted contrained triangulation of the edges encountered so far to be $CT[\leq|SL|] := CT(S[<|SL|] \cup S[|SL|])$. Accordingly we talk about the plane map $G[\leq|SL|]$ restricted to the closed halfplane left of |SL|, consisting of all vertices in the closed halfplane and edges connecting them. \end{deff} Note that $CT[\leq|SL|]$ is a valid constrained triangulation that has a part of its hull on |SL|. We store this triangulation in a special way. All triangulation edges which are finished are already stored in our output graph $G[\leq|SL|]$. Finished here means that they connect vertices in $G[\leq|SL|]$ totally left of |SL|. The missing edges completing $G[\leq|SL|]$ to $CT(S[\leq|SL|])$ can be seen as stored implicitly in our structures. We will explain this in a moment. \displayeps{ct1}{The triangles already finished during a sweep. Constraining segments are black. Additional triangulation edges are grey.}{6cm} \displaylps{chain}{The chain between two segments in |SL| at a point |p|.} When looking at the actions at a sweep event we'll see that we extend $G[\leq|SL|]$ in a way that it always represents a maximal subgraph of the final triangulation. It holds for example that $G[\leq|SL|]$ is connected. Each edge $e$ in |SL| is connected in $G[\leq|SL|]$ to the edge $e_s$ representing the successor edge in |SL|. Actually there's a face cycle $C = e_1, \ldots , e_k$ where $e_1 = |next|(|twin|(e_s))$ and $e_k = $|previous|$(e)$ in our bidirected representation of $G$ which proves connectivity. See figure \figref{chain}. This chain of edges has the additional property that it can be split into two x-monotone parts $C_1 = e_1,\ldots,e_i$ and $C_2 = e_{i+1},\ldots,e_k$ where the vertex $v_{vis} = source(e_{i+1})$ is a point of maximal x-coordinate visible by any point on |SL| between $e$ and its successor $e_s$ in |SL|. Moreover the slope of the lines that support the edges in $C_1$ and $C_2$ is non-increasing. We associate the edge $e_{vis} = e_{i+1}$ with source $v_{vis}$ to $e$ to allow fast access into the chain for any pair of neigbored segments in |SL|. By $e_{vis}$ we have a starting point for a visibility search along the x-monotone chains in which $v_{vis}$ is an extreme vertex. We come back to our representation of $CT(S[\leq|SL|])$. $CT(S[\leq|SL|])$ consists of the explicitly constructed part $G[\leq|SL|]$ and an implicitly possible triangulation between each pair of edges $e$, $e_s$ as depicted dashed in figure \figref{chain}. Just connect $|target|(|segment|(e)[\leq|SL|])$ to all vertices of $C_2$ and symmetrically $|target|(|segment|(e_s)[\leq|SL|])$ to all vertices of $C_1$ and both via |SL|. We will show for the event handling procedure that we keep the explicit and implicit structure consistent. Thus when reaching the end of our scenery where only our global sentinel edges are present (which frame the whole scenery) then the explicit part of $G$ contains the correct constrained triangulation of the input. \subsection*{Invariants}\label{invariants} We recapitulate our implementation invariants: \begin{enumerate} \item All edges intersecting the vertical line through the current event point are stored in |SL| ordered according to their intersection points bottom-up. We create hash links from these edges to their item in |SL| realized by a hash map |SLItem|. This map serves also as a flag of the edges intersected by the sweepline. We set this flag when entering and reset it to the default when leaving |SL|. \item For two edges $e, e_s$ that are neighbors in |SL| we maintain the chain of edges $C = e_1, \ldots, e_k$ stored as a face cycle in the output graph where $|source|(e_1) == |source|(e_s)$ and $|target|(e_k) == |source|(e)$. When closing $G[\leq|SL|]$ with our implicit construction, then we obtain a valid constrained triangulation of the objects in the closed halfspace left of |SL|. \item Each edge $e$ in |SL| knows a vertex $v_{vis}$ visible by any point on |SL| between $e$ and its successor $e_s$ in |SL|. We realize this knowledge by associating an edge $e_{vis} \in C$ to it with source $v_{vis}$. By $e_{vis}$ we have a starting point for a visibility search along the chain $C$ in which $v_{vis}$ is an extreme vertex. \end{enumerate} All kinds of degeneracies like equal x-coordinates and vertical segments can be integrated into the code if you imagine to twist the whole scenery clockwise by an infinitesimal angle. This implies that vertices are ordered lexicographically and vertical segments belong to a starting bundle of edges. Segments touching in one point are handled explicitly by an iteration over all edges starting in a vertex of the input graph. \section{Implementation} We use our generic sweep framework but open the algorithm furtheron by three template parameters of the sweep traits class |Constrained_triang_traits<>|. See the concept |GenericSweepTraits| that has to be implemented in the appendix section \ref{GenericSweepTraits}. We want to factor out the manipulation of the plane map, the geometric types and predicates used and the action which can be invoked on the additionally created triangulation edges in |CT|. Thus we add a concept for a plane map decorator |PMDEC|, a concept for a geometric kernel |GEOM| and a concept for a new-edge data accessor |NEWEDGE|. Note that the class provides the types and members according to the generic plane sweep traits concept |GenericSweepTraits|. \displayeps{Constrained_triang}{The design of the constrained triangulation module. |Constrained_triang_traits<>| implements the concept |GenericSweepTraits|. The three template parameters allow an adaptation of the input/output plane map, the geometry, and the processing of new edges.}{10cm} <>= template class Constrained_triang_traits : public PMDEC { public: <> <> <> <> <> <> <> <> }; // Constrained_triang_traits @ \begin{ignoreindiss} <>= <> // file : include/CGAL/Nef_2/Constrained_triang_traits.h <> #ifndef CGAL_PM_CONSTR_TRIANG_TRAITS_H #define CGAL_PM_CONSTR_TRIANG_TRAITS_H #include #include #include #include #include #include #include #undef _DEBUG #define _DEBUG 19 #include CGAL_BEGIN_NAMESPACE <> <> CGAL_END_NAMESPACE #endif // CGAL_PM_CONSTR_TRIANG_TRAITS_H @ \end{ignoreindiss} The geometry concept |GEOM| allows us to import the geometric types into |Constrained_triang_traits| and offers the necessary predicates as methods. See |AffineGeometryTraits_2| in the appendix for the concept description. The plane map decorator |PMDEC| exports the handle, iterator and circulator types into the class scope and provides all plane map exploration and manipulation methods. See |PMDecorator| in the appendix for the concept description. We inherit from |PMDEC| to obtain its methods into the class scope. <>= typedef Constrained_triang_traits Self; typedef PMDEC Base; // the types interfacing the sweep: typedef NEWEDGE INPUT; typedef typename PMDEC::Plane_map OUTPUT; typedef GEOM GEOMETRY; typedef typename GEOM::Point_2 Point; typedef typename GEOM::Segment_2 Segment; typedef typename GEOM::Direction_2 Direction; typedef typename Base::Halfedge_handle Halfedge_handle; typedef typename Base::Vertex_handle Vertex_handle; typedef typename Base::Face_handle Face_handle; typedef typename Base::Halfedge_iterator Halfedge_iterator; typedef typename Base::Vertex_iterator Vertex_iterator; typedef typename Base::Face_iterator Face_iterator; typedef typename Base::Halfedge_base Halfedge_base; typedef typename Base::Halfedge_around_vertex_circulator Halfedge_around_vertex_circulator; @ We have two sentinel edges which are minimum and maximum in our order by identity, we store them in a reference. Note that we need access to our plane map decorator and to our geometric kernel. The point |p| determines the sweepline position. <>= class lt_edges_in_sweepline : public PMDEC { const Point& p; const Halfedge_handle& e_bottom; const Halfedge_handle& e_top; const GEOMETRY& K; public: lt_edges_in_sweepline(const Point& pi, const Halfedge_handle& e1, const Halfedge_handle& e2, const PMDEC& D, const GEOMETRY& k) : PMDEC(D), p(pi), e_bottom(e1), e_top(e2), K(k) {} lt_edges_in_sweepline(const lt_edges_in_sweepline& lt) : PMDEC(lt), p(lt.p), e_bottom(lt.e_bottom), e_top(lt.e_top), K(lt.K) {} Segment seg(const Halfedge_handle& e) const { return K.construct_segment(point(source(e)),point(target(e))); } int orientation(Halfedge_handle e, const Point& p) const { return K.orientation(point(source(e)),point(target(e)),p); } <> }; // lt_edges_in_sweepline @ The order predicate on edges is only based on point equality and the orientation predicate on points. We have two sentinel edges which are minimum and maximum in our order by identity. For all geometric cases we have the precondition that the source vertex of one of the edges is equal to the current event point |p_sweep|. The order predicate is only used in search operation at |p_sweep| and for the insertion of new edges starting there, thus our requirement is legal. <>= bool operator()(const Halfedge_handle& e1, const Halfedge_handle& e2) const { // Precondition: // [[p]] is identical to the source of either [[e1]] or [[e2]]. if (e1 == e_bottom || e2 == e_top) return true; if (e2 == e_bottom || e1 == e_top) return false; if ( e1 == e2 ) return 0; int s = 0; if ( p == point(source(e1)) ) s = orientation(e2,p); else if ( p == point(source(e2)) ) s = - orientation(e1,p); else error_handler(1,"compare error in sweep."); if ( s || source(e1) == target(e1) || source(e2) == target(e2) ) return ( s < 0 ); s = orientation(e2,point(target(e1))); if (s==0) error_handler(1,"parallel edges not allowed."); return ( s < 0 ); } @ The order predicate on vertices maps to the lexicographic order on their embedding. <>= class lt_pnts_xy : public PMDEC { const GEOMETRY& K; public: lt_pnts_xy(const PMDEC& D, const GEOMETRY& k) : PMDEC(D), K(k) {} lt_pnts_xy(const lt_pnts_xy& lt) : PMDEC(lt), K(lt.K) {} int operator()(const Vertex_handle& v1, const Vertex_handle& v2) const { return K.compare_xy(point(v1),point(v2)) < 0; } }; // lt_pnts_xy @ We use an STL map for the |SL| and a set for the |event_Q|. Note that we use the iterators for local movement in |SL| but also as handles to the items therein which store the halfedge objects. Whenever we talk about items we misuse iterators for that concept. We store always the forward oriented halfedges. <>= typedef std::map Sweep_status_structure; typedef typename Sweep_status_structure::iterator ss_iterator; typedef typename Sweep_status_structure::value_type ss_pair; typedef std::set Event_Q; typedef typename Event_Q::const_iterator event_iterator; const GEOMETRY& K; Event_Q event_Q; event_iterator event_it; Vertex_handle event; Point p_sweep; Sweep_status_structure SL; CGAL::Unique_hash_map SLItem; const NEWEDGE& Treat_new_edge; Halfedge_handle e_low,e_high; // framing edges ! Halfedge_handle e_search; Constrained_triang_traits(const INPUT& in, OUTPUT& out, const GEOMETRY& k) : Base(out), K(k), event_Q(lt_pnts_xy(*this,K)), SL(lt_edges_in_sweepline(p_sweep,e_low,e_high,*this,K)), SLItem(SL.end()), Treat_new_edge(in) { TRACEN("Constrained Triangulation Sweep"); } @ \subsection{Event Handling} \newcommand{\pred}{\mathit{pred}} \renewcommand{\succ}{\mathit{succ}} \newcommand{\vis}{\mathit{vis}} We have to traverse a vertex $v$ of our input graph by the sweepline. Let's first assume this is an \emph{isolated vertex} appearing between two edges in |SL|. Note that between the two edges $e, e_s$ we have a local chain of edges connecting $v_1 = \source(e)$ and $v_2 = \source(e_s)$ which is similar to the global convex hull chain of the unconstrained triangulation problem. Note that this chain can be even empty if $v_1 = v_2$. In between the segments in the sweepline we know a vertex $v_{\vis}$ of maximal coordinates which is visible from any point on |SL| between $e$ and $e_s$. What happens if we encounter a new event $v$ (a vertex at |p_sweep|)? We locate the edge $e$ below $v$, and its successor $e_s$. We obtain the correct $v_{\vis}$ and an edge $e_{\vis}$ out of $v_{\vis}$ which is in the face cycle partly visible from $v$. Then we determine the visible chain of edges starting in $e_{\vis}$ as seen from $v$. Both ends are determined either by a non-visible edge or by an edge intersecting the sweepline (thus equal to $e$ or $e_s$). For all the edges in the chain we have to produce new triangles with apex $v$. \displaylps{ct_conf}{The four sweep configurations} Now the general case where the vertex $v$ can be the \emph{center of a star of edges}. Note that all edges are embedded around $v$ in counterclockwise order. This means that some subsequence of the adjacency list is a sequence of starting edges and the corresponding complement is the sequence of ending edges. For the actions we have to consider four configurations as shown in figure \figref{ct_conf}. Note that in case \textbf{A} we encounter an ending bundle. We have to triangulate up and down along the limiting edges along the arrows as long as the edges are visible or until they are crossing the sweepline. We also have to remove the ending edges from the sweepline. For the edge in |SL| below the sweep point we have to update the event vertex $v$ as its visible entry vertex into the graph. In case \textbf{B} we have only a starting bundle of edges. We first have to link $v$ (at the sweep point) to the visible vertex $v_{\vis}$ in the edge chain connecting the segments above and below the sweep point. Then we start the same actions as in case \textbf{A} along the chain of visible edges. Finally we insert all starting edges into |SL|. The cases \textbf{C} and \textbf{D} are combinations of the above. In \textbf{C} we first act as in \textbf{A} and then as in \textbf{B}. The isolated vertex in \textbf{D} is handled as in \textbf{B} but there are no edges starting. \displaylps{preevent}{The three geometric configurations of |v| at |p_sweep| with respect to $e$ and $e_s$ and the face cycle in between before the triangular extension.} What are we doing concerning our invariants? \begin{itemize} \item Between each pair of edges $e, e_s$ in |SL| we extend the chain where the target of $e$ or $e_s$ is the event vertex or $e$ is just below |p_sweep|, by triangles with apex $v$. Now $v$ becomes the visible vertex in at least one face cycle between two edges in |SL|. See figure \figref{preevent}. $v$ can be isolated as in \textbf{A}, or connected to $G$ by constraining edges as in \textbf{B} or \textbf{B}. There's a symmetric case for \textbf{B} when $|target|(e_s) == v$. The latter cases can occur combined if the ending bundle of edges consists of more than one edge. If you look into the code below, convince yourself that the chains are afterwards again x-monotone and follow our slope criterion. Note also that in case \textbf{A} we connect $v$ to $G[\leq|SL|]$ by new triangles and thus the whole graph is again connected. Now look at figure \figref{postevent} for the situation after the event. In all three cases our implicit structure $CT[\leq|SL|]$ is again consistently defined between $e$ and $e_s$. Note that combinations of the cases \textbf{B} and \textbf{C} are possible if the outgoing bundle of edges consists of more than one edge. \item We remove the edges ending at |v| from |SL| and reset the |SLItem| link to the default; we insert some edges starting at |v| into |SL| and set the |SLItem| link accordingly. Thus invariant 1 holds. \item We have to adjust the visibility information of those edges in |SL| that are below or above $v$ or that have just been inserted into |SL|. Thus invariant 3 is ensured. \end{itemize} \displaylps{postevent}{The three geometric configurations of |v| at |p_sweep| with respect to $e$ and $e_s$ and the visible vertex invariant after the insertion of new edges.} Now let's do some coding. For the visibility predicate we use the geometric orientation test provided by the geometric kernel. <>= bool edge_is_visible_from(Vertex_handle v, Halfedge_handle e) { Point p = point(v); Point p1 = point(source(e)); Point p2 = point(target(e)); return ( K.orientation(p1,p2,p)>0 ); // left_turn } @ The following operations are used to enrich the output plane map by all the edges completing the triangulation. We only create triangles looking backward from an event. We always start from an edge |e_apex| linking the event vertex to the up to now triangulated plane map. We triangulate away from |e_apex| until we can't see the examined face cycle edge or until the edge crosses the sweepline from left to right. For the visibility test we use the above predicate, for the "edge-crosses-sweepline" test we can use our map |SLItem| as a flag. <>= void triangulate_up(Halfedge_handle& e_apex) { TRACEN("triangulate_up "<>= void triangulate_between(Halfedge_handle e_upper, Halfedge_handle e_lower) { // we triangulate the interior of the whole chain between // target(e_upper) and target(e_lower) assert(source(e_upper)==source(e_lower)); TRACE("triangulate_between\n "<>= void process_event() { TRACEN("\nPROCESS_EVENT " << p_sweep); Halfedge_handle e, ep, eb_low, eb_high, e_end; if ( !is_isolated(event) ) { e = last_out_edge(event); ep = first_out_edge(event); } ss_iterator sit_pred, sit; /* PRECONDITION: only ingoing => e is lowest in ingoing bundle only outgoing => e is highest in outgoing bundle ingoing and outgoing => e is lowest in ingoing bundle */ eb_high = e_end = ep; eb_low = e; <> <> triangulate_up(eb_high); triangulate_down(eb_low); sit_pred->second = eb_low; } @ Given |e| an edge in the adjacency list of |v| |SLItem[e]| is an iterator pointing into |SL| (and not past the end |== SL.end()|) iff |e| is an edge in the bundle of edges ending at |v| (with respect to the sweep). In the latter case |--SLItem[e]| is an iterator pointing to the edge below |event|. If |SLItem[e] == SL.end()| then we have no entry point into |SL| and have to query |SL| by a call to |upper_bound|. The edge |e_search| is a loop edge with the property |source(e_search) == target(e_search)|. We use it for the geometric search in |SL|. The search is only executed when the |event| is not connected to $G[\leq|SL|]$. <>= TRACEN("determining handle in SL"); if ( e != Halfedge_handle() ) { point(target(e_search)) = p_sweep; // degenerate loop edge sit_pred = SLItem[e]; if ( sit_pred != SL.end()) sit = --sit_pred; else sit = sit_pred = --SL.upper_bound(e_search); } else { // event is isolated vertex point(target(e_search)) = p_sweep; // degenerate loop edge sit_pred = --SL.upper_bound(e_search); } @ We iterate the adjacency list \emph{clockwise} starting at |e|, thus we first encounter the ending bundle (if existing), then the starting bundle (if existing). <>= bool ending_edges(0), starting_edges(0); while ( e != Halfedge_handle() ) { // walk adjacency list clockwise if ( SLItem[e] != SL.end() ) <> else <> if (e == e_end) break; e = cyclic_adj_pred(e); } if (!ending_edges) <> @ When |e| is directed backwards, then |target(twin(e))==event|. For each wedge between two ending edges we triangulate the whole face. We also delete all ending edges from |SL| and mark the edges as such. <>= { TRACEN("ending " << seg(e)); if (ending_edges) triangulate_between(e,cyclic_adj_succ(e)); ending_edges = true; SL.erase(SLItem[e]); link_bi_edge_to(e,SL.end()); // not in SL anymore } @ For starting edges we insert them into |SL| and keep track of the last edge |eb_high| of the ending bundle. For the newly inserted edge their source is the visible vertex. <>= { TRACEN("starting "<>= { Halfedge_handle e_vis = sit_pred->second; Halfedge_handle e_vis_n = cyclic_adj_succ(e_vis); eb_low = eb_high = new_bi_edge(event,e_vis_n); TRACEN(" producing link "<>= void link_bi_edge_to(Halfedge_handle e, ss_iterator sit) { SLItem[e] = SLItem[twin(e)] = sit; } void initialize_structures() { TRACEN("initialize_structures "); for ( event=vertices_begin(); event != vertices_end(); ++event ) event_Q.insert(event); // sorted order of vertices event_it = event_Q.begin(); if ( event_Q.empty() ) return; event = *event_it; p_sweep = point(event); <> <> // we move to the second vertex: procede_to_next_event(); event_exists(); // sets p_sweep for check invariants TRACEN("EOF initialization"); } @ We insert all edges in the adjacency list of |event| into |SL|. All are marked to be in |SL| by a call to |link_bi_edge_to|. Thereby we obtain also a hashed shortcut into |SL|. In the wedge between two edges $e$ and $e_s$ which are neighbors in |SL| each edge $e$ is also the entry point for a visibility search in the wedge between $e$ and $e_s$. <>= if ( !is_isolated(event) ) { Halfedge_around_vertex_circulator e(first_out_edge(event)), eend(e); CGAL_For_all(e,eend) { TRACEN("init with "<>= Vertex_handle v_tmp = new_vertex(); point(v_tmp) = Point(); e_high = Base::new_halfedge_pair(event,v_tmp); e_low = Base::new_halfedge_pair(event,v_tmp); // this are two symbolic edges just accessed as sentinels // they carry no geometric information e_search = Base::new_halfedge_pair(v_tmp,v_tmp); // this is just a loop used for searches in SL ss_iterator sit_high = SL.insert(ss_pair(e_high,e_high)).first; ss_iterator sit_low = SL.insert(ss_pair(e_low,e_low)).first; // inserting sentinels into SL link_bi_edge_to(e_high, sit_high); link_bi_edge_to(e_low , sit_low); // we mark them being in the sweepline, which they will never leave @ Now for the iteration control. We iterate over all vertices, in the order given by the coordinates assigned to the vertices. We set |event| in the initialization and at the end of each event handling phase. We stop when |event_Q| is empty. <>= bool event_exists() { if ( event_it != event_Q.end() ) { // event is set at end of loop and in init event = *event_it; p_sweep = point(event); return true; } return false; } void procede_to_next_event() { ++event_it; } @ At the end we remove the frame from our structure. This is only the pair of edges spanning the initial wedge. <>= void complete_structures() { if (e_low != Halfedge_handle()) { delete_vertex(target(e_search)); } // removing sentinels and e_search } @ During the development we check if the adjacency lists of all vertices have the correct adjacency list embedding. The final check can do the test if the result of our sweep is a correct triangulation. As we don't delete edges or vertices from the output and don't add any vertex\footnote{apart from the temporary one which we delete at the end}, all constraining edges and all vertices are in triangulation. If we check the output structure to be a triangulation according to \cite{checking-cgta99} we have a checking module for our constrained triangulation sweep. <>= void check_ccw_local_embedding() const { PM_checker C(*this,K); C.check_order_preserving_embedding(event); } void check_invariants() { #ifdef CGAL_CHECK_EXPENSIVE if ( event_it == event_Q.end() ) return; check_ccw_local_embedding(); #endif } void check_final() { #ifdef CGAL_CHECK_EXPENSIVE PM_checker C(*this,K); C.check_is_triangulation(); #endif } @ \displayeps{NewEdge}{The data accessor concept |NewEdge| for the treatment of newly created edges.}{5cm} The following operations interface the plane map decorator. Note also that the edge data accessor allows to treat the newly created edges of the triangulation. <>= Halfedge_handle new_bi_edge(Vertex_handle v1, Vertex_handle v2) { // appended at v1 and v2 adj list Halfedge_handle e = Base::new_halfedge_pair(v1,v2); Treat_new_edge(e); return e; } Halfedge_handle new_bi_edge(Halfedge_handle e_bf, Halfedge_handle e_af) { // ccw before e_bf and after e_af Halfedge_handle e = Base::new_halfedge_pair(e_bf,e_af,Halfedge_base(), Base::BEFORE, Base::AFTER); Treat_new_edge(e); return e; } Halfedge_handle new_bi_edge(Vertex_handle v, Halfedge_handle e_bf) { // appended at v's adj list and before e_bf Halfedge_handle e = Base::new_halfedge_pair(v,e_bf,Halfedge_base(), Base::BEFORE); Treat_new_edge(e); return e; } Segment seg(Halfedge_handle e) const { return K.construct_segment(point(source(e)),point(target(e))); } Direction dir(Halfedge_handle e) const { return K.construct_direction(point(source(e)),point(target(e))); } bool is_forward(Halfedge_handle e) const { return K.compare_xy(point(source(e)),point(target(e))) < 0; } @ \subsection{Correctness and Running Time} At the end only the sentinels are in |SL|. $G[\leq|SL|]$ already consists of the constrained triangulation of the input structure. The two parts of the chain $C$ between the source of the sentinels just consist of the upper and lower convex hull chain between the lexicographic smallest and the lexicographic largest vertex. The convexity follows from the slope property. During the sweep we once encountered any vertex and any edge and integrated it into the constrained triangulation according to our invariants. Thus completeness is trivial. The size of the constrained triangulation of a set of segments is of the same order as the unconstrained triangulation of the segment end points. The sweep procedure takes time for the production of the output which is known to be linear in size. The only additional cost at each event is the insertion of the segments starting at an event point where we use a tree based dictionary to store these. This costs logarithmic time per segment to insert. We end up with the standard $O(n \log n)$ time bound where $n = \vert S \vert$. The space is dominated by the size of the produced output $O(n)$. @ \begin{ignoreindiss} <>= #ifndef LEDA_ERROR_H static void error_handler(int n, const char* s) { std::cerr << s << std::endl; exit(n); } #endif struct Do_nothing { Do_nothing() {} template void operator()(ARG&) const {} }; @ \subsection{Visualization via the generic sweep observer} <>= <> // file : include/CGAL/Nef_2/Constrained_triang_anim.h <> #ifndef CGAL_PM_CONSTR_TRIANG_ANIM_H #define CGAL_PM_CONSTR_TRIANG_ANIM_H #include CGAL_BEGIN_NAMESPACE template class Constrained_triang_anim { CGAL::Window_stream _W; public: typedef CGAL::Window_stream VDEVICE; typedef typename GT::GEOMETRY GEOM; typedef typename GT::Base PMDEC; typedef typename PMDEC::Point Point; Constrained_triang_anim() : _W(400,400) { _W.set_show_coordinates(true); _W.init(-120,120,-120,5); _W.display(); } VDEVICE& device() { return _W; } void post_init_animation(GT& gpst) { PM_visualizor V(_W,gpst); V.point(V.target(gpst.e_search)) = Point(-120,0); // to draw we have to embed the virtual search vertex V.draw_skeleton(CGAL::BLUE); _W.read_mouse(); } void pre_event_animation(GT& gpst) { } void post_event_animation(GT& gpst) { PM_visualizor V(_W,gpst); V.draw_ending_bundle(gpst.event,CGAL::GREEN); _W.read_mouse(); } void post_completion_animation(GT& gpst) { _W.clear(); PM_visualizor V(_W,gpst); V.draw_skeleton(CGAL::BLACK); _W.read_mouse(); } }; CGAL_END_NAMESPACE #endif // CGAL_PM_CONSTR_TRIANG_ANIM_H @ \section{A Test of the plane map triangulation} We produce a simple homogeneous kernel, a plane map and use a segment overlay sweep to create an input structure for the constrained triangulation algorithm. <>= #include #include #include #include #include "Affine_geometry.h" #undef CGAL_CFG_NO_TMPL_IN_TMPL_PARAM #include #include #include #include #include #include #include // KERNEL: typedef CGAL::Homogeneous Hom_kernel; typedef CGAL::Affine_geometry Aff_kernel; typedef Aff_kernel::Segment_2 Segment; // HALFEDGE DATA STRUCTURE: struct HDS_traits { typedef Aff_kernel::Point_2 Point; typedef bool Mark; }; typedef CGAL::HalfedgeDS_default HDS; typedef CGAL::PM_decorator< HDS > PM_dec; typedef PM_dec::Halfedge_handle Halfedge_handle; typedef PM_dec::Vertex_handle Vertex_handle; typedef PM_dec::Halfedge_const_handle Halfedge_const_handle; // SEGMENT OVERLAY: template class PM_dec_output : public PMDEC { public: typedef PMDEC Base; typedef typename Base::Plane_map Plane_map; typedef typename Base::Point Point; typedef typename Base::Vertex_handle Vertex_handle; typedef typename Base::Halfedge_handle Halfedge_handle; typedef I ITERATOR; PM_dec_output(HDS& H) : Base(H) {} PM_dec_output(const PM_dec_output& P) : Base(P) {} Vertex_handle new_vertex(const Point& 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 {} }; // PM_dec_output typedef std::list::const_iterator Seg_iterator; typedef CGAL::Segment_overlay_traits< Seg_iterator, PM_dec_output, Aff_kernel> PM_seg_overlay; typedef CGAL::generic_sweep PM_seg_overlay_sweep; // CONSTRAINED TRIANGULATIONS: typedef CGAL::Constrained_triang_traits CTT; typedef CGAL::generic_sweep Constrained_triang_sweep; typedef CGAL::Constrained_triang_anim CTA; typedef CGAL::sweep_observer CTS_observer; // MAIN PROGRAM: int main(int argc, char* argv[]) { // SETDTHREAD(19); CGAL::set_pretty_mode ( cerr ); HDS H; Aff_kernel AK; CTS_observer Obs; Obs.device().message("Insert segments to triangulate."); std::list L; Segment s; if ( argc == 2 ) { std::ifstream log(argv[1]); while ( log >> s ) { L.push_back(s); Obs.device() << s; } } while ( Obs.device() >> s ) L.push_back(s); std::string fname(argv[0]); fname += ".log"; std::ofstream log(fname.c_str()); for (Seg_iterator sit = L.begin(); sit != L.end(); ++sit) log << *sit << " "; log.close(); PM_seg_overlay::OUTPUT D(H); PM_seg_overlay_sweep OV(PM_seg_overlay::INPUT(L.begin(),L.end()),D,AK); CTT::INPUT I; Constrained_triang_sweep CT(I,H,AK); Obs.attach(CT); OV.sweep(); CT.sweep(); return 0; } @ \begin{ignore} <>= // ============================================================================ // // 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$ // <>= // package : Nef_2 // chapter : Nef Polyhedra // // source : nef_2d/Constrained_triang.lw // revision : $Id$ // revision_date : $Date$ // // author(s) : Michael Seel // maintainer : Michael Seel // coordinator : Michael Seel // // implementation: Constrained triangulation of a plane map // ============================================================================ @ \end{ignore} \end{ignoreindiss} %KILLSTART REP \newpage \bibliographystyle{alpha} \bibliography{comp_geo,general,diss} \newpage \section{Appendix} \input manpages/AffineGeometryTraits_2.man \input manpages/GenericSweepTraits.man \input manpages/PMConstDecorator.man \input manpages/PMDecorator.man \input manpages/PM_checker.man \end{document} %KILLEND REP