Graphite  Version 3
An experimental 3D geometry processing program
GEO::ExactCDT2d Class Reference

Constrained Delaunay Triangulation with vertices that are exact points. Can be used to implement 2D CSG. More...

#include <geogram/delaunay/CDT_2d.h>

Inheritance diagram for GEO::ExactCDT2d:
GEO::CDTBase2d

Public Types

typedef exact::vec2h ExactPoint
 

Public Member Functions

 ExactCDT2d ()
 ExactCDT2d constructor.
 
 ~ExactCDT2d () override
 ExactCDT2d destructor.
 
void clear () override
 Removes everything from this triangulation. More...
 
index_t insert (const ExactPoint &p, index_t id=0, index_t hint=index_t(-1))
 Inserts a point. More...
 
void insert_constraint (index_t v1, index_t v2, index_t operand_bits=0)
 Inserts a constraint. More...
 
void create_enclosing_quad (const ExactPoint &p1, const ExactPoint &p2, const ExactPoint &p3, const ExactPoint &p4)
 Creates a first large enclosing quad. More...
 
void create_enclosing_rectangle (double x1, double y1, double x2, double y2)
 Creates a first large enclosing rectangle. More...
 
const ExactPointpoint (index_t v) const
 Gets a point by vertex index. More...
 
index_t vertex_id (index_t v) const
 Gets a vertex id by vertex index. More...
 
void set_vertex_id (index_t v, index_t id)
 Sets a vertex id by vertex index. More...
 
void classify_triangles (const std::string &boolean_expression, bool mark_only=false)
 Used by 2D CSG operations, discards triangles according to a boolean operation. More...
 
void save (const std::string &filename) const override
 
- Public Member Functions inherited from GEO::CDTBase2d
 CDTBase2d ()
 CDTBase2d constructor.
 
virtual ~CDTBase2d ()
 CDTBase2d destructor.
 
void insert_constraint (index_t i, index_t j)
 Inserts a constraint. More...
 
void remove_external_triangles (bool remove_internal_holes=false)
 Recursively removes all the triangles adjacent to the border, and keeps what's surrounded by constraints. More...
 
void set_delaunay (bool delaunay)
 Specifies whether a constrained Delaunay triangulation should be constructed, or just a plain constrained triangulation. More...
 
index_t nT () const
 Gets the number of triangles.
 
index_t nv () const
 Gets the number of vertices.
 
index_t ncnstr () const
 Gets the number of constraints.
 
index_t Tv (index_t t, index_t lv) const
 Gets a vertex of a triangle. More...
 
index_t Tv_find (index_t t, index_t v) const
 Finds the local index of a vertex in a triangle. More...
 
index_t Tadj (index_t t, index_t le) const
 Gets a triangle adjacent to a triangle. More...
 
index_t Tadj_find (index_t t1, index_t t2) const
 Finds the edge accross which a triangle is adjacent to another one. More...
 
index_t vT (index_t v) const
 Gets a triangle incident to a given vertex. More...
 
index_t Tedge_cnstr_first (index_t t, index_t le) const
 Gets the constraint associated with an edge. More...
 
index_t edge_cnstr_next (index_t ecit) const
 Gets the successor of an edge constraint iterator. More...
 
index_t edge_cnstr (index_t ecit) const
 Gets an edge constraint from an edge constraint iterator. More...
 
index_t Tedge_cnstr_nb (index_t t, index_t le) const
 Gets the number of constraints associated with a triange edge. More...
 
bool Tedge_is_Delaunay (index_t t, index_t le) const
 Tests whether a triangle edge is Delaunay. More...
 
void check_consistency () const
 Checks both combinatorics and geometry, aborts on unconsistency.
 

Protected Member Functions

void add_point (const ExactPoint &p, index_t id=index_t(-1))
 
void begin_insert_transaction () override
 
void commit_insert_transaction () override
 
void rollback_insert_transaction () override
 
Sign orient2d (index_t i, index_t j, index_t k) const override
 Computes the orientation predicate in 2d. More...
 
Sign incircle (index_t i, index_t j, index_t k, index_t l) const override
 Incircle predicate. More...
 
index_t create_intersection (index_t E1, index_t i, index_t j, index_t E2, index_t k, index_t l) override
 Given two segments that have an intersection, create the intersection. More...
 
- Protected Member Functions inherited from GEO::CDTBase2d
index_t insert (index_t v, index_t hint=index_t(-1))
 Inserts a new point. More...
 
void create_enclosing_triangle (index_t v1, index_t v2, index_t v3)
 Creates the combinatorics for a first large enclosing triangle. More...
 
void create_enclosing_quad (index_t v1, index_t v2, index_t v3, index_t v4)
 Creates the combinatorics for a first large enclosing quad. More...
 
void Tset_flag (index_t t, index_t flag)
 Sets a triangle flag. More...
 
void Treset_flag (index_t t, index_t flag)
 Resets a triangle flag. More...
 
bool Tflag_is_set (index_t t, index_t flag)
 Tests a triangle flag. More...
 
bool Tis_in_list (index_t t) const
 Tests whether a triangle is in a DList. More...
 
void insert_vertex_in_edge (index_t v, index_t t, index_t le, DList &S)
 Inserts a vertex in an edge. More...
 
void insert_vertex_in_edge (index_t v, index_t t, index_t le)
 Inserts a vertex in an edge. More...
 
void insert_vertex_in_triangle (index_t v, index_t t, DList &S)
 Inserts a vertex in a triangle. More...
 
index_t find_intersected_edges (index_t i, index_t j, DList &Q)
 Finds the edges intersected by a constraint. More...
 
void walk_constraint_v (CDT2d_ConstraintWalker &W)
 Used by find_intersected_edges()
 
void walk_constraint_t (CDT2d_ConstraintWalker &W, DList &Q)
 Used by find_intersected_edges()
 
void constrain_edges (index_t i, index_t j, DList &Q, DList &N)
 Constrains an edge by iteratively flipping the intersected edges. More...
 
void Delaunayize_vertex_neighbors (index_t from_v)
 Restores Delaunay condition starting from the triangles incident to a given vertex. More...
 
void Delaunayize_vertex_neighbors (index_t v, DList &S)
 Restores Delaunay condition starting from the triangles incident to a given vertex. More...
 
void Delaunayize_new_edges (DList &N)
 Restores Delaunay condition for a set of edges after inserting a constrained edge. More...
 
void Tset (index_t t, index_t v1, index_t v2, index_t v3, index_t adj1, index_t adj2, index_t adj3, index_t e1cnstr=index_t(-1), index_t e2cnstr=index_t(-1), index_t e3cnstr=index_t(-1))
 Sets all the combinatorial information of a triangle and edge flags. More...
 
void Trot (index_t t, index_t lv)
 Rotates indices in triangle t in such a way that a given vertex becomes vertex 0. More...
 
void swap_edge (index_t t1, bool swap_t1_t2=false)
 Swaps an edge. More...
 
void Tadj_set (index_t t, index_t le, index_t adj)
 Sets a triangle adjacency relation. More...
 
index_t Topp (index_t t, index_t e=0) const
 Gets the neighboring triangle vertex opposite to a given vertex.
 
void Tadj_back_connect (index_t t1, index_t le1, index_t prev_t2_adj_e2)
 After having changed connections from triangle to a neighbor, creates connections from neighbor to triangle. More...
 
index_t Tnew ()
 Creates a new triangle. More...
 
void Tset_edge_cnstr_first (index_t t, index_t le, index_t ecit)
 Sets the constraints list associated with an edge. More...
 
void Tadd_edge_cnstr (index_t t, index_t le, index_t cnstr_id)
 Adds a constraint to a triangle edge. More...
 
void Tadd_edge_cnstr_with_neighbor (index_t t, index_t le, index_t cnstr_id)
 Adds a constraint to a triangle edge and to the neighboring edge if it exists. More...
 
bool Tedge_is_constrained (index_t t, index_t le) const
 Tests whether an edge is constrained. More...
 
void for_each_T_around_v (index_t v, std::function< bool(index_t t, index_t lv)> doit)
 Calls a user-defined function for each triangle around a vertex. More...
 
index_t locate (index_t v, index_t hint=index_t(-1), Sign *orient=nullptr) const
 Locates a vertex. More...
 
bool is_convex_quad (index_t t) const
 Tests whether triange t and its neighbor accross edge 0 form a strictly convex quad. More...
 
void remove_marked_triangles ()
 Removes all the triangles that have the flag T_MARKED_FLAG set. More...
 
void Tcheck (index_t t) const
 Consistency check for a triangle. More...
 
void debug_Tcheck (index_t t) const
 Consistency check for a triangle in debug mode, ignored in release mode. More...
 
void check_combinatorics () const
 Consistency combinatorial check for all the triangles. More...
 
void debug_check_combinatorics () const
 Consistency combinatorial check for all the triangles in debug mode, ignored in release mode. More...
 
virtual void check_geometry () const
 Consistency geometrical check for all the triangles. More...
 
void debug_check_geometry () const
 Consistency geometrical check for all the triangles in debug mode, ignored in release mode. More...
 
void debug_check_consistency () const
 Checks both combinatorics and geometry in debug mode, ignored in release mode, aborts on unconsistency.
 
bool segment_segment_intersect (index_t u1, index_t u2, index_t v1, index_t v2) const
 Tests whether two segments have a frank intersection. More...
 
bool segment_edge_intersect (index_t v1, index_t v2, index_t t, index_t le) const
 Tests whether an edge triangle and a segment have a frank intersection. More...
 
void check_edge_intersections (index_t v1, index_t v2, const DList &Q)
 Checks that the edges stored in a DList exactly correspond to all edge intersections between a segment and the triangle edges. More...
 
index_t eT (Edge E)
 Gets a triangle incident a a given edge. More...
 
index_t locate_naive (index_t v, index_t hint=index_t(-1), Sign *orient=nullptr) const
 Simpler version of locate() kept for reference. More...
 
void constrain_edges_naive (index_t i, index_t j, DList &Q, vector< Edge > &N)
 Simpler version of constrain_edges() kept for reference. More...
 
void Delaunayize_new_edges_naive (vector< Edge > &N)
 Simpler version of Delaunayize_new_edges() that uses a vector instead of a DList, kept for reference.
 

Protected Attributes

vector< ExactPointpoint_
 
vector< double > length_
 
vector< index_tid_
 
vector< index_tcnstr_operand_bits_
 
vector< index_tfacet_inclusion_bits_
 
std::map< trindex, Signpred_cache_
 
bool use_pred_cache_insert_buffer_
 
std::vector< std::pair< trindex, Sign > > pred_cache_insert_buffer_
 
vector< bindexconstraints_
 
- Protected Attributes inherited from GEO::CDTBase2d
index_t nv_
 
index_t ncnstr_
 
vector< index_tT_
 
vector< index_tTadj_
 
vector< index_tv2T_
 
vector< uint8_t > Tflags_
 
vector< index_tTecnstr_first_
 
vector< index_tecnstr_val_
 
vector< index_tecnstr_next_
 
vector< index_tTnext_
 
vector< index_tTprev_
 
bool delaunay_
 
Sign orient_012_
 
bool exact_incircle_
 
bool exact_intersections_
 

Additional Inherited Members

- Protected Types inherited from GEO::CDTBase2d
enum  { DLIST_S_ID =0 , DLIST_Q_ID =1 , DLIST_N_ID =2 , DLIST_NB =3 }
 Constants for list_id.
 
enum  { T_MARKED_FLAG = DLIST_NB , T_VISITED_FLAG = DLIST_NB+1 }
 Constants for triangle flags.
 
typedef std::pair< index_t, index_tEdge
 
- Static Protected Member Functions inherited from GEO::CDTBase2d
static index_t find_3 (const index_t *T, index_t v)
 Finds the index of an integer in an array of three integers. More...
 

Detailed Description

Constrained Delaunay Triangulation with vertices that are exact points. Can be used to implement 2D CSG.

Points are represented using exact 2d homogeneous coordinates. Unlike CDT2d, this ensures exact representation of constraints intersections with guaranteed behavior. Under the hood, it inherits CDTBase2d (constrained Delaunay triangulation), and redefines orient2d(), incircle2d() and create_intersection() using vectors with homogeneous coordinates stored as arithmetic expansions (vec2HE) or arbitrary-precision floating point numbers (vec2HEx) if compiled with Tessael's geogramplus extension package.

See also
CDT2d

Definition at line 1537 of file CDT_2d.h.

Member Function Documentation

◆ classify_triangles()

void GEO::ExactCDT2d::classify_triangles ( const std::string &  boolean_expression,
bool  mark_only = false 
)

Used by 2D CSG operations, discards triangles according to a boolean operation.

Discards all the triangles that are outside the object defined by the boolean expression. It uses the operand bits associated with the constraints.

Parameters
[in]boolean_expressiona string with the boolean expression, as defined by BooleanExpression constructor. Each variable corresponds to an operand bit associated with the constraints. There can be up to 32 operands.
mark_onlyif set, triangles to be discarded are marked (but not discarded)
See also
BooleanExpression, insert_constraint()

◆ clear()

void GEO::ExactCDT2d::clear ( )
overridevirtual

Removes everything from this triangulation.

Reimplemented from GEO::CDTBase2d.

◆ create_enclosing_quad()

void GEO::ExactCDT2d::create_enclosing_quad ( const ExactPoint p1,
const ExactPoint p2,
const ExactPoint p3,
const ExactPoint p4 
)

Creates a first large enclosing quad.

Parameters
[in]p1,p2,p3,p4the four vertices of the quad

The quad needs to be convex. create_enclosing_rectangle() or create_enclosing_quad() need to be called before anything else

◆ create_enclosing_rectangle()

void GEO::ExactCDT2d::create_enclosing_rectangle ( double  x1,
double  y1,
double  x2,
double  y2 
)
inline

Creates a first large enclosing rectangle.

Parameters
[in]x1,y1,x2,y2rectangle bounds

create_enclosing_triangle(), create_enclosing_rectangle() or create_enclosing_quad() need to be called before anything else

Definition at line 1604 of file CDT_2d.h.

◆ create_intersection()

index_t GEO::ExactCDT2d::create_intersection ( index_t  E1,
index_t  i,
index_t  j,
index_t  E2,
index_t  k,
index_t  l 
)
overrideprotectedvirtual

Given two segments that have an intersection, create the intersection.

The intersection is given both as the indices of segment extremities (i,j) and (k,l), that one can use to retreive the points in derived classes, and constraint indices E1 and E2, that derived classes may use to retreive symbolic information attached to the constraint

Parameters
[in]E1the index of the first edge, corresponding to the value of ncnstr() when insert_constraint() was called for that edge
[in]i,jthe vertices of the first segment
[in]E2the index of the second edge, corresponding to the value of ncnstr() when insert_constraint() was called for that edge
[in]k,lthe vertices of the second segment
Returns
the index of a newly created vertex that corresponds to the intersection between [i , j] and [k , l]

Implements GEO::CDTBase2d.

◆ incircle()

Sign GEO::ExactCDT2d::incircle ( index_t  i,
index_t  j,
index_t  k,
index_t  l 
) const
overrideprotectedvirtual

Incircle predicate.

Parameters
[in]i,j,kthe three vertices of a triangle
[in]lanother vertex
Return values
POSITIVEif l is inside the circumscribed circle of the triangle
ZEROif l is on the circumscribed circle of the triangle
NEGATIVEif l is outside the circumscribed circle of the triangle

Implements GEO::CDTBase2d.

◆ insert()

index_t GEO::ExactCDT2d::insert ( const ExactPoint p,
index_t  id = 0,
index_t  hint = index_t(-1) 
)

Inserts a point.

Parameters
[in]pthe point to be inserted
[in]hinta triangle not too far away from the point to be inserted
[in]idan opaque identifier attached to the vertex that can be used by client code for instance to keep relations with a mesh. It can be queried using the vertex_id() function.
Returns
the index of the created point. Duplicated points are detected (and then the index of the existing point is returned)

◆ insert_constraint()

void GEO::ExactCDT2d::insert_constraint ( index_t  v1,
index_t  v2,
index_t  operand_bits = 0 
)
inline

Inserts a constraint.

Parameters
[in]v1,v2the two extremities of the constraint, as returned by insert()
[in]operand_bitsoptional bitfield used by 2D CSG, indicating on which primitive boundaries the constraint is. Each bit set corresponds to a primitive. It is used by the classify() function.
See also
CDTBase::insert_constraint() and classify()

Definition at line 1580 of file CDT_2d.h.

◆ orient2d()

Sign GEO::ExactCDT2d::orient2d ( index_t  i,
index_t  j,
index_t  k 
) const
overrideprotectedvirtual

Computes the orientation predicate in 2d.

Computes the sign of the signed area of the triangle p0, p1, p2.

Parameters
[in]p0,p1,p2vertices of the triangle as 2d vectors with homogeneous coordinates stored as expansion_nt (arbitrary precision).
Return values
POSITIVEif the triangle is oriented counter-clockwise
ZEROif the triangle is flat
NEGATIVEif the triangle is oriented clockwise

Implements GEO::CDTBase2d.

◆ point()

const ExactPoint& GEO::ExactCDT2d::point ( index_t  v) const
inline

Gets a point by vertex index.

Parameters
[in]vvertex index
Returns
the point at index v

Definition at line 1620 of file CDT_2d.h.

◆ save()

void GEO::ExactCDT2d::save ( const std::string &  filename) const
overridevirtual
See also
CDTBase2d::save()

Implements GEO::CDTBase2d.

◆ set_vertex_id()

void GEO::ExactCDT2d::set_vertex_id ( index_t  v,
index_t  id 
)
inline

Sets a vertex id by vertex index.

Parameters
[in]vvertex index
[in]idvertex id

Definition at line 1640 of file CDT_2d.h.

◆ vertex_id()

index_t GEO::ExactCDT2d::vertex_id ( index_t  v) const
inline

Gets a vertex id by vertex index.

Parameters
[in]vvertex index
Returns
the point at index v

Definition at line 1630 of file CDT_2d.h.


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