Navigation: Up, Table of Contents, Bibliography, Index, Title Page

Planar Map

Introduction

In this chapter we introduce the planar map which is an embedding of a topological map into the plane. It is derived from the topological map class with additional geometric considerations and functionalities (e.g., point location). In this section we briefly review the geometric concepts added in the planar map, the combinatorial concepts were reviewed in Chapter  reference, ``Topological Map''.

Curve:  We will use the term curve to describe the image of a continuous 1-1 mapping into the plane of any one of the following: the closed unit interval (arc), the open unit interval (unbounded curve), or the unit circle (closed curve). In all cases a curve is non self-intersecting.

x-Monotone curve:  We use a convention that a curve c is x-monotone if c either intersects any vertical line in at most one point, or c is a vertical segment.

Planar subdivision (planar map):  A planar subdivision (or planar map) is an embedding of a planar topological map T into the plane, such that each edge of T is embedded as a bounded x-monotone curve and each vertex is embedded as a planar point. In this embedding no edge intersects another edge's interior.

A face of the subdivision is a maximal connected region of the plane that does not contain any vertex or edge. We consider a face to be open, and its boundary is formed by vertices and halfedges of the subdivision. The halfedges are oriented around a face so that the face they bound is to their left. This means that halfedges on the outer boundary of a face are traversed in counterclockwise order, and halfedges on the inner boundaries (holes) of a face are traversed in clockwise order. Edges around a vertex are also traversed in clockwise order.

Planar Map


begin of advanced section

Point Location Strategies

The planar map users can define which algorithm to use in the point location queries. This is done with a point location class passed to the map in the constructor. The class passed should be derived from the base class Pm_point_location_base which is a pure virtual base class that defines the interface between the algorithm implemented by the users and the planar map. This follows the known strategy pattern  [GHJV95]. The indirection overhead due to the virtual functions is negligible since the optimal point location algorithms (e.g., the one implemented in our default strategy) take (logn) time.

We have derived three concrete classes for point location strategies which are presented at Section reference -  reference.

Default Point Location Strategy

The Pm_default_point_location<Planar_map> class implements the incremental randomized algorithm introduced by Mulmuley  [Mul90] as presented by Seidel  [Sei91], [dBvKOS97]. It uses a trapezoidal map and a search structure to achieve an expected query time of O(logn). Updating of these structures takes expected O(logn) time. There is also a possibility to assure a deterministic worst-case query time of O(logn), with a preprocessing step. The trade-off is in a longer building time.

#include <CGAL/Pm_default_point_location.h>

pl ( bool preprocess = false);
constructs the point location object. If preprocess is true then preprocessing is done to guarantee a deterministic worst case query time of O(logn). However, the preprocessing step can slow down considerably the building time of the structure.

Naive Point Location Strategy

The Pm_naive_point_location<Planar_map> class implements the naive algorithm which goes over all the vertices and halfedges of the planar map. Therefore, the time complexity of a query is linear in the complexity of the planar map. Since it does not use an additional structure, the updating operations are empty and therefore modifying operations such as split_edge and merge_edge take constant time.

#include <CGAL/Pm_naive_point_location.h>

Walk-along-a-Line Point Location Strategy

The Pm_walk_along_line_point_location<Planar_map> class implements a walk-along-a-line algorithm which improves the naive one. Unlike the naive algorithm which goes over all the edges of the planar map, the walk algorithm starts with the unbounded face and ``walks'' along the zone of the vertical ray emanating from the query point. This reduces the number of halfedges that are traversed in a query. Like the naive strategy, the walk strategy does not use an additional structure. The updating operations are empty and therefore modifying operations such as split_edge and merge_edge take constant time.

#include <CGAL/Pm_walk_along_line_point_location.h>

Trade-off Issues

The main trade-off between the two strategies implemented, is between time and storage. Using the naive or walk strategies takes more time but saves storage space.

Another trade-off depends on the need for point location queries compared to the need for other functions. If the users do not need point location queries, but do need other modifying functions (e.g., remove_edge, split_edge and merge_edge) then using the naive or walk strategies is preferable. Note that using the insert function invokes the point location query, therefore when using the naive or walk strategies it is recommended to use the specialized insertion functions : insert_in_face_interior, insert_from_vertex and insert_at_vertices. For example, when using the planar map to represent polygons (e.g., when computing boolean operations on polygons) it might be preferable to use the walk strategy with the specialized insertion functions.

Another trade-off is between the two modes of the default strategy. If preprocessing is not used, the building of the structure is faster. However, for some input sequences the structure might be unbalanced and therefore queries can take longer.


end of advanced section

Requirements for Planar Map Dcel Classes

The planar map class is parameterized with the interface classes Dcel and Traits . The Dcel class is a container class that defines the underlying combinatorial data structure (halfedge data structure) used by the planar map. The Traits class defines the abstract interface between planar maps and the geometric primitives they use.

The requirements for a Dcel class of a planar map are identical to those required for a Dcel class of a topological map. In addition the Vertex and Halfedge classes for the planar map are required to have geometric functionality. The requirements are stated below.

Model for a Planar Map Dcel Class

The Pm_dcel<V,H,F> is a model of the Dcel concept defined above. Its template parameters are the vertex, halfedge and face classes which conform to the requirements given above.

#include <CGAL/Pm_default_dcel.h>

Models for Dcel Vertex, Halfedge and Face

The Pm_vertex_base<Point> is a model for the Dcel vertex concept defined above. The Pm_halfedge_base<Curve> is a model for the Dcel halfedge concept defined above. The Pm_face_base is a model for the Dcel face concept defined above.

#include <CGAL/Pm_default_dcel.h>

The Default Dcel Structure

The Pm_default_dcel<Traits> is a model of the Dcel concept defined above. Its template parameter is the traits class defined below. It is a wrapper for Pm_dcel instantiated with Pm_vertex_base<Traits::Point>, Pm_halfedge_base<Traits::X_curve> and Pm_face_base.

#include <CGAL/Pm_default_dcel.h>

Requirements for Planar Map Traits Classes

The planar map class is parameterized with the interface class Traits which defines the abstract interface between planar maps and the primitives they use. The following requirement catalog lists the primitives, i.e., types, member functions etc., that must be defined for a Traits class that can be used to parameterize planar maps.

Predefined Traits Classes

We supply a default traits class for exact arithmetic - Pm_segment_exact_traits<R> where R is the kernel representation type (Homogeneous or Cartesian). There is also a class for finite precision arithmetic - Pm_segment_epsilon_traits<R>, which is not recommended unless it is known that robustness issues will not arise. Both traits handle finite line segments in the plane and use the CGAL kernel (Point is of type Point_2<R> and X_curve is of type Segment_2<R>). We also supply a traits class that uses LEDA's rational kernel and makes use of its efficient predicates.


begin of advanced section


end of advanced section

Example Program

Figure:  The map generated by the example program

The following program creates a planar map from five curves (see figure  reference). It uses the floating-point arithmetic interface class (Pm_segment_epsilon_traits<R>). The first part of the code initializes five segments. The next part inserts them to the map and checks the validity of the map. In the last part, a vertical ray shooting is performed from the point p.

#include <CGAL/Cartesian.h>
#include <CGAL/Pm_segment_epsilon_traits.h>
#include <CGAL/Pm_default_dcel.h>
#include <CGAL/Planar_map_2.h>


typedef double                                           number_type;
typedef Cartesian<number_type>                      coord_t;
typedef Pm_segment_epsilon_traits<coord_t>          pmtraits;
typedef typename pmtraits::Point                                  point;
typedef typename pmtraits::X_curve                                curve;

typedef Pm_default_dcel<pmtraits>                   pmdcel;


int main()
{
  // creating an instance of Planar_map_2<pmdcel,pmtraits>
  Planar_map_2<pmdcel,pmtraits> pm;

  curve cv[5];
  int i;

  point a1(100, 0), a2(20, 50), a3(180, 50), a4(100, 100);

  // those curves are about to be inserted to pm
  cv[0] = curve(a1, a2);
  cv[1] = curve(a1, a3);
  cv[2] = curve(a2, a3);
  cv[3] = curve(a2, a4);
  cv[4] = curve(a3, a4);
  
  cout << "the curves of the map :" << endl; 
  for (i = 0; i < 5; i++)
    cout << cv[i] << endl;

  cout << endl;

  // insert the five curves to the map
  cout << "inserting the curves to the map..." << endl;
  for (i = 0; i < 5; i++)
  {
    cout << "inserting curve" << i << endl;
    pm.insert(cv[i]);
  }
  
  // check the validity of the map
  cout << "check map validity... ";
  if (pm.is_valid())
    cout << "map valid!" << endl;
  else
    cout << "map invalid!" << endl;
  cout << endl;
  
  // vertical ray shooting upward from p
  point p(95, 30);
  typename Planar_map_2<pmdcel,pmtraits>::Halfedge_handle e;  
  typename Planar_map_2<pmdcel,pmtraits>::Locate_type lt;

  cout << endl << "upward vertical ray shooting from " << p << endl; 
  e=pm.vertical_ray_shoot(p, lt, true);
  if (lt != Planar_map_2<pmdcel,pmtraits>::UNBOUNDED_FACE) {
    cout << "returned the curve :" << e->curve() << endl;
  }
  return 0;  
}

The output of this program is

the curves of the map :
100 0 20 50
100 0 180 50
20 50 180 50
20 50 100 100
180 50 100 100

inserting the curves to the map...
inserting curve0
inserting curve1
inserting curve2
inserting curve3
inserting curve4
check map validity... map valid!

upward vertical ray shooting from 95 30
returned the curve : 20 50 180 50


Navigation: Up, Table of Contents, Bibliography, Index, Title Page
The GALIA project. Jan 18, 2000.