LiTe
 All Classes Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
Classes | Public Types | Public Member Functions | List of all members
LineTes Class Reference

Polygonal tessellation of a bounded rectangular domain. More...

#include <ttessel.h>

Inheritance diagram for LineTes:
TTessel

Classes

class  Seg
 Segment in a line tessellation. More...
 
class  Seg_list
 List of Seg objects. More...
 
class  Seg_sublist
 Sublist of Seg objects. More...
 

Public Types

typedef SegSeg_handle
 Pointer to a Seg object.
 
typedef Seg_list::iterator Seg_list_iterator
 Iterator for accessing segments in a list.
 
typedef Seg_sublist::iterator Seg_sublist_iterator
 Iterator for accessing segments in a sublist of segments.
 

Public Member Functions

 LineTes ()
 Empty LineTes constructor. More...
 
virtual void insert_window (Rectangle)
 Define the rectangular domain delimiting the tessellation. More...
 
virtual void insert_window (Polygon &)
 Define the polygonal domain delimiting the tessellation. More...
 
Halfedge_handle split_from_edge (Halfedge_handle, Point2, Halfedge_handle, Point2)
 Split a tessellation cell from an edge to another one. More...
 
Halfedge_handle split_from_vertex (Vertex_handle, Halfedge_handle, Point2)
 Split a tessellation cell from a vertex to an edge. More...
 
Face_handle suppress_edge (Halfedge_handle)
 Suppress an edge from the tessellation. More...
 
void clear ()
 Clear a tessellation. More...
 
Seg_list_iterator segments_begin ()
 Return an iterator pointing to the beginning of the list of tessellation segments.
 
Seg_list_iterator segments_end ()
 Return an iterator pointing to the end of the list of tessellation segments.
 
Size number_of_segments ()
 Return the total number of segments. More...
 
bool is_valid (bool verbose=false)
 Check whether the LineTes object is valid. More...
 
bool check_all_segments (bool verbose=false)
 Test whether all segments of a tessellation are valid. More...
 
bool check_segment (Seg_handle s, bool verbose=false)
 Test whether a segment of a tessellation is valid. More...
 
bool check_all_halfedges (bool verbose=false)
 Check all halfedges of a line tessellation. More...
 
bool check_halfedge (LineTes::Halfedge_handle e, bool verbose=false)
 Test whether a halfedge is valid. More...
 
bool halfedge_exists (LineTes::Halfedge_handle)
 Test whether a halfedge handle is indeed a halfedge handle of a tessellation.
 
bool check_halfedge_neighbours (LineTes::Halfedge_handle e, bool verbose=false)
 Check the neighbours of a halfedge along its segment. More...
 
bool check_halfedge_prev_hf (LineTes::Halfedge_handle e, bool verbose=false)
 Check the previous neighbour of a halfedge along its segment. More...
 
bool check_halfedge_next_hf (LineTes::Halfedge_handle e, bool verbose=false)
 Check the next neighbour of a halfedge along its segment. More...
 
bool check_halfedge_segment (LineTes::Halfedge_handle e, bool verbose=false)
 Check the segment attribute of a halfedge. More...
 
bool check_halfedge_dir (LineTes::Halfedge_handle e, bool verbose=false)
 Check the direction attribute of a halfedge. More...
 
void print_all_segments (bool endsOnly=false)
 Print on standard output the segment points. More...
 
Polygon get_window ()
 Return the domain that is tessellated.
 
Size number_of_window_edges ()
 Return the number of sides of the tessellated domain.
 
double get_window_perimeter ()
 Return the perimeter of the tessellated domain.
 
int is_on_boundary (Halfedge_handle)
 Test whether a halfedge lies along the domain boundary. More...
 
bool is_on_boundary (Seg_handle s)
 Test whether a segment lies on the domain boundary. More...
 
Seg_handle insert_segment (Segment)
 Insert a new line segment in a line tessellation. More...
 
void remove_ivertices (Size imax=0, bool verbose=false)
 Remove all I-vertices from a line tessellation. More...
 
void remove_lvertices (Size imax=0, bool verbose=false)
 Remove all L-vertices from a line tessellation. More...
 
bool is_a_T_tessellation (bool verbose=false, std::ostream &out=std::clog)
 Test whether a line tessellation is a T-tessellation. More...
 
void write (std::ostream &)
 Write a LineTes object to an output stream. More...
 
void read (std::istream &)
 Read a LineTes object from an input stream. More...
 

Detailed Description

Polygonal tessellation of a bounded rectangular domain.

This class inherits from the Arrangement_2 class provided by CGAL. Its main feature concerns the so-called segments. A segment is defined as a maximal subset of edges that are aligned and connected. The LineTes class is of interest for tessellations with segments that consist of more than one edge.

Constructor & Destructor Documentation

LineTes::LineTes ( )

Empty LineTes constructor.

Equivalent to the empty CGAL::Arrangement_2 constructor

Member Function Documentation

bool LineTes::check_all_halfedges ( bool  verbose = false)

Check all halfedges of a line tessellation.

Parameters
verbose: if true, when a invalid halfedge is found, print some data on the standard output stream for logging.

All halfedges of the line tessellation are inspected. Check is based on the LineTes::check_halfedge method.

bool LineTes::check_all_segments ( bool  verbose = false)

Test whether all segments of a tessellation are valid.

Returns
true if all segments are valid, false otherwise.

All segments are tested using the check_segment method.

bool LineTes::check_halfedge ( LineTes::Halfedge_handle  e,
bool  verbose = false 
)

Test whether a halfedge is valid.

Parameters
e: a handle of the halfedge to be checked.
verbose: if true and halfedge is found to be non valid, some information is sent to the standard output for logging.
Returns
true if the halfedge is valid, false otherwise.

Several tests are performed: existence of the halfedge, validity of its neighbours along its segment, consistency of the segment and direction attributes.

bool LineTes::check_halfedge_dir ( LineTes::Halfedge_handle  e,
bool  verbose = false 
)

Check the direction attribute of a halfedge.

Parameters
e: handle of the halfedge to be checked.
verbose: if true, details are sent to std::clog. Default to false.
Returns
true if the halfedge and its segment directions are consistent, false otherwise.
bool LineTes::check_halfedge_neighbours ( LineTes::Halfedge_handle  e,
bool  verbose = false 
)

Check the neighbours of a halfedge along its segment.

Parameters
e: a handle of the halfedge to be checked.
verbose: if true, details are sent to std::clog. Default to false.
Returns
true if neighbours are valid, false otherwise.
bool LineTes::check_halfedge_next_hf ( LineTes::Halfedge_handle  e,
bool  verbose = false 
)

Check the next neighbour of a halfedge along its segment.

Parameters
e: a handle of the halfedge to be checked.
verbose: if true, details are sent to std::clog. Default to false.
Returns
true if the next neighbour is valid, false otherwise.
bool LineTes::check_halfedge_prev_hf ( LineTes::Halfedge_handle  e,
bool  verbose = false 
)

Check the previous neighbour of a halfedge along its segment.

Parameters
e: a handle of the halfedge to be checked.
verbose: if true, details are sent to std::clog. Default to false.
Returns
true if the previous neighbour is valid, false otherwise.
bool LineTes::check_halfedge_segment ( LineTes::Halfedge_handle  e,
bool  verbose = false 
)

Check the segment attribute of a halfedge.

Parameters
e: handle of the halfedge to be checked.
verbose: if true, details are sent to std::clog. Default to false.
Returns
true if the segment is valid, false otherwise.

Features checked: alignement of the halfedge and its segment, occurence of the segment in the tessellation segment list.

bool LineTes::check_segment ( Seg_handle  s,
bool  verbose = false 
)

Test whether a segment of a tessellation is valid.

Parameters
s: a handle of the segment to be tested.
verbose: if true and the segment is non valid, some information is sent to the standard output for logging.
Returns
true if the segment is valid, false otherwise.

Just check that the halfedge identifying the segment is indeed a tessellation halfedge.

void LineTes::clear ( )

Clear a tessellation.

All segments are removed including the domain sides.

LineTes::Seg_handle LineTes::insert_segment ( Segment  iseg)

Insert a new line segment in a line tessellation.

Parameters
iseg: line segment to be inserted in the tessellation.
Returns
handle to the newly created Seg object.

Eventually the segment is clipped by the tessellation domain before insertion.

void LineTes::insert_window ( Rectangle  r)
virtual

Define the rectangular domain delimiting the tessellation.

The LineTes class is designed for representing tessellation of bounded rectangular domains. The insert_window method should be used as a first step where one defines the domain to be tessellated.

Reimplemented in TTessel.

void LineTes::insert_window ( Polygon p)
virtual

Define the polygonal domain delimiting the tessellation.

Do not use this method at the moment as non-rectangular domains are not fully supported by the LineTes class.

Reimplemented in TTessel.

bool LineTes::is_a_T_tessellation ( bool  verbose = false,
std::ostream &  out = std::clog 
)

Test whether a line tessellation is a T-tessellation.

Parameters
verbose: if true, details are sent to std::clog. Default to false.
out: output stream to be written if verbose is true.
Returns
true if is a T-tessellation, false otherwise.

Check whether all internal vertices are T-vertices and whether all vertices on the domain boundary are either of degree 2 or T-vertices.

int LineTes::is_on_boundary ( Halfedge_handle  e)

Test whether a halfedge lies along the domain boundary.

Test if an edge lies on the window boundary. Return

  • 0 if the edge is internal.
  • 1 if the halfedge bounds an internal (bounded) face.
  • 2 if the halfedge bounds the external (unbounded) face.
bool LineTes::is_on_boundary ( Seg_handle  s)
inline

Test whether a segment lies on the domain boundary.

Parameters
s: a handle to the segment to be tested
bool LineTes::is_valid ( bool  verbose = false)

Check whether the LineTes object is valid.

Parameters
verbose: if true, details are sent to std::clog. Default to false. return true if valid, false otherwise

The following checks are performed:

  • Validity as an arrangement.
  • Validity of all halfedges (method check_all_halfedges).
  • Validity of all segments (method check_all_segments).
Size LineTes::number_of_segments ( )
inline

Return the total number of segments.

Boundary and internal segments are included

void LineTes::print_all_segments ( bool  endsOnly = false)

Print on standard output the segment points.

Parameters
endsOnly: see the documentation of the print method of the Seg class.
void LineTes::read ( std::istream &  is)

Read a LineTes object from an input stream.

Parameters
is: input stream.

Data provided by the input stream should be formatted as described in LineTes::write documentation.

void LineTes::remove_ivertices ( Size  imax = 0,
bool  verbose = false 
)

Remove all I-vertices from a line tessellation.

Parameters
imax: maximal number of I-vertices to be removed. If zero, (default value) all I-vertices are removed.
verbose: If true, information about I-vertex processing is sent to the standard output. Default to false.

An I-vertex is a tessellation vertex of degree 1. Removal is performed either by shortening or lengthening of the incident edge. Shortening is preferred if the removed length is less than the added length.

The sequence of removals is ordered by length variation. At each iteration, the I-vertex to be removed is the one inducing the smallest length variation.

void LineTes::remove_lvertices ( Size  imax = 0,
bool  verbose = false 
)

Remove all L-vertices from a line tessellation.

Parameters
imax: maximal number of I-vertices to be removed. If zero, (default value) all I-vertices are removed.
verbose: if true, information about I-vertex processing is sent to the standard output. Default to false.

An L-vertex is a tessellation vertex of degree 2 with both incident edges which are not aligned. Removal is performed by lengthening one of the incident edges. The edge to be lengthened is the one that induces the shortest length variation.

The sequence of removals is ordered by length variation. At each iteration, the L-vertex to be removed is the one inducing the smallest length variation.

LineTes::Halfedge_handle LineTes::split_from_edge ( Halfedge_handle  e1,
Point2  p1,
Halfedge_handle  e2,
Point2  p2 
)

Split a tessellation cell from an edge to another one.

Parameters
e1: handle of the halfedge where the splitting segment starts.
p1: point on the edge where the splitting segment starts.
e2: handle of the halfedge where the splitting segment ends.
p2: point on the edge where the splitting segment ends.
Returns
handle of one of the halfedges along the splitting segment.
Precondition
e1 should not bound the external face. If so, return NULL_HALFEDGE_HANDLE.
LineTes::Halfedge_handle LineTes::split_from_vertex ( Vertex_handle  v,
Halfedge_handle  e,
Point2  p 
)

Split a tessellation cell from a vertex to an edge.

Parameters
v: handle of the vertex where the splitting segment starts.
e: handle of the edge where the splitting segment ends.
ppoint where the splitting segment ends.
Returns
handle of one of the halfedges along the splitting segment.
LineTes::Face_handle LineTes::suppress_edge ( Halfedge_handle  e)

Suppress an edge from the tessellation.

Parameters
e: handle of the halfedge to be suppressed.
Returns
a handle to the new face created by the edge removal.
Precondition
e should not lie on the domain boundary and it should be a the end of an existing segment. If the condition is not fulfilled, the tessellation is not modified and NULL_HALFEDGE_HANDLE is returned.

The vertices of the suppressed edge are removed if they are of degree 2 and if the connect two consecutive edges on a segment.

void LineTes::write ( std::ostream &  os)

Write a LineTes object to an output stream.

Parameters
os: output stream.

Output is done based on a format without loss of numerical precision. Therefore it is possible to recover a LineTes object from the output. The output format is as follows:

  • A coordinate is output as two integers (numerator and denominator).
  • A point is output as its two Cartesian coordinates.
  • A line segment is output as its two end points.
  • First, vertices of the tessellated domain are sent to the output stream (single line).
  • Subsequently, internal segments are sent to the output stream (one per line).

The documentation for this class was generated from the following files: