Geogram  Version 1.9.1
A programming library of geometric algorithms
GEO::CDTBase2d Class Referenceabstract

Base class for constrained Delaunay triangulation. More...

#include <geogram/delaunay/CDT_2d.h>

Inheritance diagram for GEO::CDTBase2d:
GEO::CDT2d GEO::ExactCDT2d GEO::MeshInTriangle

Classes

struct  DList
 Doubly connected triangle list. More...
 

Public Member Functions

 CDTBase2d ()
 CDTBase2d constructor.
 
virtual ~CDTBase2d ()
 CDTBase2d destructor.
 
virtual void clear ()
 Removes everything from this triangulation.
 
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...
 
virtual void save (const std::string &filename) const =0
 Saves this CDT to a geogram mesh file. 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 Types

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
 

Protected Member Functions

virtual void begin_insert_transaction ()
 
virtual void commit_insert_transaction ()
 
virtual void rollback_insert_transaction ()
 
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...
 
virtual Sign orient2d (index_t i, index_t j, index_t k) const =0
 Orientation predicate. More...
 
virtual Sign incircle (index_t i, index_t j, index_t k, index_t l) const =0
 Incircle predicate. More...
 
virtual index_t create_intersection (index_t E1, index_t i, index_t j, index_t E2, index_t k, index_t l)=0
 Given two segments that have an intersection, create the intersection. 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.
 

Static Protected Member Functions

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...
 

Protected Attributes

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_
 

Detailed Description

Base class for constrained Delaunay triangulation.

Manages the combinatorics of the constrained Delaunay triangulation. The points need to be stored elsewhere, and manipulated through indices, with two predicates:

  • orient2d(i,j,k)
  • incircle(i,j,k,l) and one construction:
  • create_intersection(i,j,k,l) See CDT2d for an example of implementation

Definition at line 75 of file CDT_2d.h.

Member Function Documentation

◆ check_combinatorics()

void GEO::CDTBase2d::check_combinatorics ( ) const
inlineprotected

Consistency combinatorial check for all the triangles.

aborts if inconsistency is detected

Definition at line 1151 of file CDT_2d.h.

◆ check_edge_intersections()

void GEO::CDTBase2d::check_edge_intersections ( index_t  v1,
index_t  v2,
const DList Q 
)
protected

Checks that the edges stored in a DList exactly correspond to all edge intersections between a segment and the triangle edges.

Parameters
[in]v1,v2the two vertices of the constrained segment
[in]Qa list of triangle. For each triangle in Q, edge 0 is supposed to have an intersection with v1 , v2

◆ check_geometry()

virtual void GEO::CDTBase2d::check_geometry ( ) const
protectedvirtual

Consistency geometrical check for all the triangles.

aborts if inconsistency is detected

◆ constrain_edges()

void GEO::CDTBase2d::constrain_edges ( index_t  i,
index_t  j,
DList Q,
DList N 
)
protected

Constrains an edge by iteratively flipping the intersected edges.

Parameters
[in]i,jthe extremities of the edge
[in]Qthe list of intersected edges, computed by find_intersected_edges()
[out]Noptional DList with the new edges that need to be re-Delaunized by find_intersected_edges(), ignored if uninitialized

◆ constrain_edges_naive()

void GEO::CDTBase2d::constrain_edges_naive ( index_t  i,
index_t  j,
DList Q,
vector< Edge > &  N 
)
protected

Simpler version of constrain_edges() kept for reference.

See also
constrain_edges()

◆ create_enclosing_quad()

void GEO::CDTBase2d::create_enclosing_quad ( index_t  v1,
index_t  v2,
index_t  v3,
index_t  v4 
)
protected

Creates the combinatorics for a first large enclosing quad.

Parameters
[in]v1,v2,v3,v4the four vertices of the quad, in 0,1,2,3

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

◆ create_enclosing_triangle()

void GEO::CDTBase2d::create_enclosing_triangle ( index_t  v1,
index_t  v2,
index_t  v3 
)
protected

Creates the combinatorics for a first large enclosing triangle.

Parameters
[in]v1,v2,v3the three vertices of the first triangle, in 0,1,2

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

◆ create_intersection()

virtual index_t GEO::CDTBase2d::create_intersection ( index_t  E1,
index_t  i,
index_t  j,
index_t  E2,
index_t  k,
index_t  l 
)
protectedpure virtual

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]

Implemented in GEO::MeshInTriangle, GEO::ExactCDT2d, and GEO::CDT2d.

◆ debug_check_combinatorics()

void GEO::CDTBase2d::debug_check_combinatorics ( ) const
inlineprotected

Consistency combinatorial check for all the triangles in debug mode, ignored in release mode.

aborts if inconsistency is detected

Definition at line 1162 of file CDT_2d.h.

◆ debug_check_geometry()

void GEO::CDTBase2d::debug_check_geometry ( ) const
inlineprotected

Consistency geometrical check for all the triangles in debug mode, ignored in release mode.

aborts if inconsistency is detected

Definition at line 1179 of file CDT_2d.h.

◆ debug_Tcheck()

void GEO::CDTBase2d::debug_Tcheck ( index_t  t) const
inlineprotected

Consistency check for a triangle in debug mode, ignored in release mode.

aborts if inconsistency is detected

Parameters
[in]tthe triangle to be tested

Definition at line 1139 of file CDT_2d.h.

◆ Delaunayize_new_edges()

void GEO::CDTBase2d::Delaunayize_new_edges ( DList N)
protected

Restores Delaunay condition for a set of edges after inserting a constrained edge.

Parameters
[in]Nthe edges for which Delaunay condition should be restored.

◆ Delaunayize_vertex_neighbors() [1/2]

void GEO::CDTBase2d::Delaunayize_vertex_neighbors ( index_t  from_v)
protected

Restores Delaunay condition starting from the triangles incident to a given vertex.

This version uses internally a stack, initialized with the triangles incident to the vertex.

Parameters
[in]from_vthe vertex. Cannot be a vertex incident to the border.

◆ Delaunayize_vertex_neighbors() [2/2]

void GEO::CDTBase2d::Delaunayize_vertex_neighbors ( index_t  v,
DList S 
)
protected

Restores Delaunay condition starting from the triangles incident to a given vertex.

Parameters
[in]vthe vertex
[in]Sa stack of triangles, initialized with the triangles incident to the vertex. Each triangle t is Trot()-ed in such a way that the vertex v corresponds to Vt(t,0)

Each time a triangle edge is swapped, the two new neighbors are recursively examined.

◆ edge_cnstr()

index_t GEO::CDTBase2d::edge_cnstr ( index_t  ecit) const
inline

Gets an edge constraint from an edge constraint iterator.

Parameters
[in]ecitthe edge constraint iterator. Needs to be a valid iterator, different from index_t(-1).
Returns
the edge constraint associated with the iterator.
See also
Tedge_cnstr_first(), edge_cnstr_next()

Definition at line 261 of file CDT_2d.h.

◆ edge_cnstr_next()

index_t GEO::CDTBase2d::edge_cnstr_next ( index_t  ecit) const
inline

Gets the successor of an edge constraint iterator.

Parameters
[in]ecitthe edge constraint iterator
Returns
the edge constraint iterator to the successor of ecit or index_t(-1) if ecit is the last of the list.
See also
Tedge_cnstr_first(), edge_cnstr()

Definition at line 250 of file CDT_2d.h.

◆ eT()

index_t GEO::CDTBase2d::eT ( Edge  E)
inlineprotected

Gets a triangle incident a a given edge.

Parameters
[in]Ethe edge
Returns
a triangle with E as its edge 0

Definition at line 1260 of file CDT_2d.h.

◆ find_3()

static index_t GEO::CDTBase2d::find_3 ( const index_t T,
index_t  v 
)
inlinestaticprotected

Finds the index of an integer in an array of three integers.

Parameters
[in]Ta const pointer to an array of three integers
[in]vthe integer to retrieve in T
Returns
the index (0,1 or 2) of v in T
Precondition
The three entries of T are different and one of them is equal to v.

Definition at line 1084 of file CDT_2d.h.

◆ find_intersected_edges()

index_t GEO::CDTBase2d::find_intersected_edges ( index_t  i,
index_t  j,
DList Q 
)
protected

Finds the edges intersected by a constraint.

Parameters
[in]i,jthe two vertices of the constraint
[out]Qfor each intersected edge, a triangle t will be pushed-back to Q, such that vT(t,1) and vT(t,2) are the extremities of the intersected edge. In addition, each triangle t is marked.

If a vertex k that is exactly on the constraint is found, then traversal stops there and k is returned. One can find the remaining intersections by continuing to call the function with (k,j) until j is returned.

Returns
the first vertex on [i,j] encountered when traversing the segment [i,j].

◆ for_each_T_around_v()

void GEO::CDTBase2d::for_each_T_around_v ( index_t  v,
std::function< bool(index_t t, index_t lv)>  doit 
)
inlineprotected

Calls a user-defined function for each triangle around a vertex.

Parameters
[in]vthe vertex
[in]doitthe function, that takes as argument the current triangle t and the local index lv of v in t. The function returns true if iteration is finished and can be exited, false otherwise.

Definition at line 975 of file CDT_2d.h.

◆ incircle()

virtual Sign GEO::CDTBase2d::incircle ( index_t  i,
index_t  j,
index_t  k,
index_t  l 
) const
protectedpure virtual

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

Implemented in GEO::MeshInTriangle, GEO::ExactCDT2d, and GEO::CDT2d.

◆ insert()

index_t GEO::CDTBase2d::insert ( index_t  v,
index_t  hint = index_t(-1) 
)
protected

Inserts a new point.

Parameters
[in]vthe index of the new point, supposed to be equal to nv()
[in]hintan optional triangle, not too far away from the point to be inserted
Returns
the index of the created point. May be different from v if the point already existed in the triangulation

◆ insert_constraint()

void GEO::CDTBase2d::insert_constraint ( index_t  i,
index_t  j 
)

Inserts a constraint.

Parameters
[in]i,jthe indices of the two vertices of the constrained segment

◆ insert_vertex_in_edge() [1/2]

void GEO::CDTBase2d::insert_vertex_in_edge ( index_t  v,
index_t  t,
index_t  le 
)
inlineprotected

Inserts a vertex in an edge.

Parameters
[in]vthe vertex to be inserted
[in]ta triangle incident to the edge
[in]lethe local index of the edge in t

Definition at line 647 of file CDT_2d.h.

◆ insert_vertex_in_edge() [2/2]

void GEO::CDTBase2d::insert_vertex_in_edge ( index_t  v,
index_t  t,
index_t  le,
DList S 
)
protected

Inserts a vertex in an edge.

Parameters
[in]vthe vertex to be inserted
[in]ta triangle incident to the edge
[in]lethe local index of the edge in t
[out]SDList of created triangles, ignored if uninitialized

◆ insert_vertex_in_triangle()

void GEO::CDTBase2d::insert_vertex_in_triangle ( index_t  v,
index_t  t,
DList S 
)
protected

Inserts a vertex in a triangle.

Parameters
[in]vthe vertex to be inserted
[in]tthe triangle
[out]Soptional DList of created triangles

◆ is_convex_quad()

bool GEO::CDTBase2d::is_convex_quad ( index_t  t) const
protected

Tests whether triange t and its neighbor accross edge 0 form a strictly convex quad.

Return values
trueif triange t and its neighbor accross edge 0 form a strictly convex quad
falseotherwise

◆ locate()

index_t GEO::CDTBase2d::locate ( index_t  v,
index_t  hint = index_t(-1),
Sign orient = nullptr 
) const
protected

Locates a vertex.

Parameters
[in]vthe vertex index
[in]hintan optional triangle, not too far away from the point to be inserted
[out]orienta pointer to the three orientations in the triangle. If one of them is zero, the point is on an edge, and if two of them are zero, it is on a vertex.
Returns
a triangle that contains v

◆ locate_naive()

index_t GEO::CDTBase2d::locate_naive ( index_t  v,
index_t  hint = index_t(-1),
Sign orient = nullptr 
) const
protected

Simpler version of locate() kept for reference.

See also
locate()

◆ orient2d()

virtual Sign GEO::CDTBase2d::orient2d ( index_t  i,
index_t  j,
index_t  k 
) const
protectedpure virtual

Orientation predicate.

Parameters
[in]i,j,kthree vertices
Returns
the sign of det(pj-pi,pk-pi)

Implemented in GEO::MeshInTriangle, GEO::ExactCDT2d, and GEO::CDT2d.

◆ remove_external_triangles()

void GEO::CDTBase2d::remove_external_triangles ( bool  remove_internal_holes = false)

Recursively removes all the triangles adjacent to the border, and keeps what's surrounded by constraints.

Parameters
[in]remove_internal_holesif set, triangles inside the internal closed loops of constrained edges are removed as well.

If remove_internal_holes is set, closed loops inside holes are considered as "matter" (and kept), and so on and so forth. This also works if there are overlapping constraints (what counts is the number of constraints associated with each triangle edge). Note that this does not work if there is a chain of constrained internal edges (as opposed to a loop).

◆ remove_marked_triangles()

void GEO::CDTBase2d::remove_marked_triangles ( )
protected

Removes all the triangles that have the flag T_MARKED_FLAG set.

This compresses the triangle array, and updates triangle adjacencies, as well as the vertex to triangle array. Note that it uses Tnext_'s storage for internal bookkeeping, so this function should not be used if there exists a non-empty DList.

◆ save()

virtual void GEO::CDTBase2d::save ( const std::string &  filename) const
pure virtual

Saves this CDT to a geogram mesh file.

Parameters
[in]filenamewhere to save this CDT

Implemented in GEO::MeshInTriangle, GEO::ExactCDT2d, and GEO::CDT2d.

◆ segment_edge_intersect()

bool GEO::CDTBase2d::segment_edge_intersect ( index_t  v1,
index_t  v2,
index_t  t,
index_t  le 
) const
inlineprotected

Tests whether an edge triangle and a segment have a frank intersection.

Parameters
[in]v1,v2the two extremities of the segment
[in]ta triangle
[in]lelocal edge index in 0,1,2
Return values
trueif edge le of t has a frank intersection with edge v1 , v2
falseotherwise

Definition at line 1234 of file CDT_2d.h.

◆ segment_segment_intersect()

bool GEO::CDTBase2d::segment_segment_intersect ( index_t  u1,
index_t  u2,
index_t  v1,
index_t  v2 
) const
inlineprotected

Tests whether two segments have a frank intersection.

Parameters
[in]u1,u2the two extremities of the first segment
[in]v1,v2the two extremities of the second segment
Return values
trueif u1 , u2 has a frank intersection with v1 , v2
falseotherwise

Definition at line 1215 of file CDT_2d.h.

◆ set_delaunay()

void GEO::CDTBase2d::set_delaunay ( bool  delaunay)
inline

Specifies whether a constrained Delaunay triangulation should be constructed, or just a plain constrained triangulation.

Parameters
[in]delaunaytrue if a Delaunay triangulation should be constructed, false otherwise.

Definition at line 124 of file CDT_2d.h.

◆ swap_edge()

void GEO::CDTBase2d::swap_edge ( index_t  t1,
bool  swap_t1_t2 = false 
)
protected

Swaps an edge.

Swaps edge 0 of t1. Vertex 0 of t1 is vertex 0 of the two new triangles.

Parameters
[in]t1a triangle index. Its edge opposite to vertex 0 is swapped
[in]swap_t1_t2if set, swap which triangle will be t1 and which triangle will be Tadj(t1,0) in the new pair of triange (needed for two configurations of the optimized constraint enforcement algorithm).

◆ Tadd_edge_cnstr()

void GEO::CDTBase2d::Tadd_edge_cnstr ( index_t  t,
index_t  le,
index_t  cnstr_id 
)
inlineprotected

Adds a constraint to a triangle edge.

Parameters
[in]ta triangle
[in]lelocal edge index, in 0,1,2
[in]cnstr_idthe constraint

Definition at line 903 of file CDT_2d.h.

◆ Tadd_edge_cnstr_with_neighbor()

void GEO::CDTBase2d::Tadd_edge_cnstr_with_neighbor ( index_t  t,
index_t  le,
index_t  cnstr_id 
)
inlineprotected

Adds a constraint to a triangle edge and to the neighboring edge if it exists.

Parameters
[in]ta triangle
[in]lelocal edge index, in 0,1,2
[in]cnstr_idthe constraint

Definition at line 936 of file CDT_2d.h.

◆ Tadj()

index_t GEO::CDTBase2d::Tadj ( index_t  t,
index_t  le 
) const
inline

Gets a triangle adjacent to a triangle.

Parameters
[in]tthe triangle
[in]lethe local edge index, in 0,1,2
Returns
the triangle adjacent to t accross le, or index_t(-1) if there is no such triangle

Definition at line 180 of file CDT_2d.h.

◆ Tadj_back_connect()

void GEO::CDTBase2d::Tadj_back_connect ( index_t  t1,
index_t  le1,
index_t  prev_t2_adj_e2 
)
inlineprotected

After having changed connections from triangle to a neighbor, creates connections from neighbor to triangle.

edge flags are copied from the neighbor to t1. If there is no triangle accross le1, then nothing is done

Parameters
[in]t1a triangle
[in]le1a local edge of t1, in 0,1,2
[in]prev_t2_adj_e2the triangle adjacent to t2 that t1 will replace, where t2 = Tadj(t1,le1)

Definition at line 852 of file CDT_2d.h.

◆ Tadj_find()

index_t GEO::CDTBase2d::Tadj_find ( index_t  t1,
index_t  t2 
) const
inline

Finds the edge accross which a triangle is adjacent to another one.

Parameters
[in]t1,t2the two triangles
Returns
the local edge index le such that Tadj(t1,le) = t2

Definition at line 193 of file CDT_2d.h.

◆ Tadj_set()

void GEO::CDTBase2d::Tadj_set ( index_t  t,
index_t  le,
index_t  adj 
)
inlineprotected

Sets a triangle adjacency relation.

Parameters
[in]ta triangle
[in]lelocal edge index, in 0,1,2
[in]adjthe triangle adjacent to t accross le

Definition at line 820 of file CDT_2d.h.

◆ Tcheck()

void GEO::CDTBase2d::Tcheck ( index_t  t) const
inlineprotected

Consistency check for a triangle.

aborts if inconsistency is detected

Parameters
[in]tthe triangle to be tested

Definition at line 1117 of file CDT_2d.h.

◆ Tedge_cnstr_first()

index_t GEO::CDTBase2d::Tedge_cnstr_first ( index_t  t,
index_t  le 
) const
inline

Gets the constraint associated with an edge.

When constraining segments on a CDT2d by calling insert_constraint(), some segments in the triangulation may be included in several constraints (it the constraints are co-linear and overlapping). One iterates on the constraints associated with an edge as follows:

for(
index_t ecit = Tedge_cnstr_first(t,le);
ecit != index_t(-1);
ecit = edge_cnstr_next(ecit)
) {
index_t cnstr = edge_cnstr(ecit);
... // do something with cnstr
}
index_t Tedge_cnstr_first(index_t t, index_t le) const
Gets the constraint associated with an edge.
Definition: CDT_2d.h:237
index_t edge_cnstr_next(index_t ecit) const
Gets the successor of an edge constraint iterator.
Definition: CDT_2d.h:250
index_t edge_cnstr(index_t ecit) const
Gets an edge constraint from an edge constraint iterator.
Definition: CDT_2d.h:261
geo_index_t index_t
The type for storing and manipulating indices.
Definition: numeric.h:329

where 'cnstr' corresponds to the value of ncnstr() when insert_constraint() was called for that constraint.

Parameters
[in]ta triangle
[in]lelocal edge index, in 0,1,2
Returns
the edge constraints iterator associated with this edge or index_t(-1) if the edge is not constrained.
See also
edge_cnstr_next(), edge_cnstr()

Definition at line 237 of file CDT_2d.h.

◆ Tedge_cnstr_nb()

index_t GEO::CDTBase2d::Tedge_cnstr_nb ( index_t  t,
index_t  le 
) const
inline

Gets the number of constraints associated with a triange edge.

There can be several constraints associated with the same edge, whenever there are overlapping constraints. For instance, this function is useful to test the parity of the number of constraints when classifying inside/outside triangles in a CSG operation.

Parameters
[in]tthe triangle
[in]lethe local index of the edge (0,1,2) in the triangle
Returns
the number of constraints associated with the edge

Definition at line 276 of file CDT_2d.h.

◆ Tedge_is_constrained()

bool GEO::CDTBase2d::Tedge_is_constrained ( index_t  t,
index_t  le 
) const
inlineprotected

Tests whether an edge is constrained.

Parameters
[in]ta triangle
[in]lelocal edge index, in 0,1,2
Return values
trueif the edge is constrained
falseotherwise

Definition at line 962 of file CDT_2d.h.

◆ Tedge_is_Delaunay()

bool GEO::CDTBase2d::Tedge_is_Delaunay ( index_t  t,
index_t  le 
) const

Tests whether a triangle edge is Delaunay.

returns true also for constrained edges and edges on borders

◆ Tflag_is_set()

bool GEO::CDTBase2d::Tflag_is_set ( index_t  t,
index_t  flag 
)
inlineprotected

Tests a triangle flag.

Parameters
[in]tthe triangle
[in]flagthe flag, in 0..7
Return values
trueif the flag is set
falseotherwise

Definition at line 368 of file CDT_2d.h.

◆ Tis_in_list()

bool GEO::CDTBase2d::Tis_in_list ( index_t  t) const
inlineprotected

Tests whether a triangle is in a DList.

Parameters
[in]tthe triangle
Return values
trueif the triangle is in a list
falseotherwise

Definition at line 398 of file CDT_2d.h.

◆ Tnew()

index_t GEO::CDTBase2d::Tnew ( )
inlineprotected

Creates a new triangle.

Returns
the index of the new triange

Definition at line 870 of file CDT_2d.h.

◆ Treset_flag()

void GEO::CDTBase2d::Treset_flag ( index_t  t,
index_t  flag 
)
inlineprotected

Resets a triangle flag.

Parameters
[in]tthe triangle
[in]flagthe flag, in 0..7

Definition at line 355 of file CDT_2d.h.

◆ Trot()

void GEO::CDTBase2d::Trot ( index_t  t,
index_t  lv 
)
inlineprotected

Rotates indices in triangle t in such a way that a given vertex becomes vertex 0.

On exit, vertex lv of t becomes vertex 0

Parameters
[in]ta triangle index
[in]lvlocal vertex index in 0,1,2

Definition at line 783 of file CDT_2d.h.

◆ Tset()

void GEO::CDTBase2d::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) 
)
inlineprotected

Sets all the combinatorial information of a triangle and edge flags.

Parameters
[in]tthe triangle
[in]v1,v2,v3the three vertices
[in]adj1,adj2,adj3the three triangles adjacent to t
[in]e1cnstr,e2cnstr,e3cnstroptional edge constraints

Definition at line 740 of file CDT_2d.h.

◆ Tset_edge_cnstr_first()

void GEO::CDTBase2d::Tset_edge_cnstr_first ( index_t  t,
index_t  le,
index_t  ecit 
)
inlineprotected

Sets the constraints list associated with an edge.

Parameters
[in]ta triangle
[in]lelocal edge index, in 0,1,2
[in]ecitthe edge constraint iterator that points to the first constraint associated with the edge

Definition at line 889 of file CDT_2d.h.

◆ Tset_flag()

void GEO::CDTBase2d::Tset_flag ( index_t  t,
index_t  flag 
)
inlineprotected

Sets a triangle flag.

Parameters
[in]tthe triangle
[in]flagthe flag, in 0..7

Definition at line 344 of file CDT_2d.h.

◆ Tv()

index_t GEO::CDTBase2d::Tv ( index_t  t,
index_t  lv 
) const
inline

Gets a vertex of a triangle.

Parameters
[in]tthe triangle
[in]lvthe local index of the vertex, in 0,1,2
Returns
the global index of the vertex

Definition at line 155 of file CDT_2d.h.

◆ Tv_find()

index_t GEO::CDTBase2d::Tv_find ( index_t  t,
index_t  v 
) const
inline

Finds the local index of a vertex in a triangle.

Parameters
[in]tthe triangle
[in]vthe vertex
Returns
lv such that Tv(t,lv) = v

Definition at line 167 of file CDT_2d.h.

◆ vT()

index_t GEO::CDTBase2d::vT ( index_t  v) const
inline

Gets a triangle incident to a given vertex.

Parameters
[in]va vertex
Returns
a triangle t such that there exists lv in 0,1,2 such that Tv(t,lv) = v

Definition at line 205 of file CDT_2d.h.

Member Data Documentation

◆ delaunay_

bool GEO::CDTBase2d::delaunay_
protected

if set, compute a CDT, else just a CT

Definition at line 1324 of file CDT_2d.h.

◆ ecnstr_next_

vector<index_t> GEO::CDTBase2d::ecnstr_next_
protected

edge constraints list, links

Definition at line 1321 of file CDT_2d.h.

◆ ecnstr_val_

vector<index_t> GEO::CDTBase2d::ecnstr_val_
protected

edge constraints list, indices

Definition at line 1320 of file CDT_2d.h.

◆ exact_incircle_

bool GEO::CDTBase2d::exact_incircle_
protected

true if incircle() is 100% exact

Definition at line 1326 of file CDT_2d.h.

◆ exact_intersections_

bool GEO::CDTBase2d::exact_intersections_
protected

true if intersections are 100% exact

Definition at line 1327 of file CDT_2d.h.

◆ orient_012_

Sign GEO::CDTBase2d::orient_012_
protected

global triangles orientation

Definition at line 1325 of file CDT_2d.h.

◆ T_

vector<index_t> GEO::CDTBase2d::T_
protected

triangles vertices array

Definition at line 1315 of file CDT_2d.h.

◆ Tadj_

vector<index_t> GEO::CDTBase2d::Tadj_
protected

triangles adjacency array

Definition at line 1316 of file CDT_2d.h.

◆ Tecnstr_first_

vector<index_t> GEO::CDTBase2d::Tecnstr_first_
protected

index in edge constraints list

Definition at line 1319 of file CDT_2d.h.

◆ Tflags_

vector<uint8_t> GEO::CDTBase2d::Tflags_
protected

triangle flags

Definition at line 1318 of file CDT_2d.h.

◆ Tnext_

vector<index_t> GEO::CDTBase2d::Tnext_
protected

doubly connected triangle list

Definition at line 1322 of file CDT_2d.h.

◆ Tprev_

vector<index_t> GEO::CDTBase2d::Tprev_
protected

doubly connected triangle list

Definition at line 1323 of file CDT_2d.h.

◆ v2T_

vector<index_t> GEO::CDTBase2d::v2T_
protected

vertex to triangle back pointer

Definition at line 1317 of file CDT_2d.h.


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