Geogram
Version 1.9.1-rc
A programming library of geometric algorithms
|
Base class for constrained Delaunay triangulation. More...
#include <geogram/delaunay/CDT_2d.h>
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_t > | Edge |
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_t > | T_ |
vector< index_t > | Tadj_ |
vector< index_t > | v2T_ |
vector< uint8_t > | Tflags_ |
vector< index_t > | Tecnstr_first_ |
vector< index_t > | ecnstr_val_ |
vector< index_t > | ecnstr_next_ |
vector< index_t > | Tnext_ |
vector< index_t > | Tprev_ |
bool | delaunay_ |
Sign | orient_012_ |
bool | exact_incircle_ |
bool | exact_intersections_ |
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:
|
inlineprotected |
Checks that the edges stored in a DList exactly correspond to all edge intersections between a segment and the triangle edges.
[in] | v1,v2 | the two vertices of the constrained segment |
[in] | Q | a list of triangle. For each triangle in Q, edge 0 is supposed to have an intersection with v1 , v2 |
|
protectedvirtual |
Consistency geometrical check for all the triangles.
aborts if inconsistency is detected
Constrains an edge by iteratively flipping the intersected edges.
[in] | i,j | the extremities of the edge |
[in] | Q | the list of intersected edges, computed by find_intersected_edges() |
[out] | N | optional DList with the new edges that need to be re-Delaunized by find_intersected_edges(), ignored if uninitialized |
|
protected |
Simpler version of constrain_edges() kept for reference.
|
protected |
Creates the combinatorics for a first large enclosing quad.
[in] | v1,v2,v3,v4 | the four vertices of the quad, in 0,1,2,3 |
create_enclosing_triangle() or create_enclosing_quad() need to be called before anything else
Creates the combinatorics for a first large enclosing triangle.
[in] | v1,v2,v3 | the three vertices of the first triangle, in 0,1,2 |
create_enclosing_triangle() or create_enclosing_quad() need to be called before anything else
|
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
[in] | E1 | the index of the first edge, corresponding to the value of ncnstr() when insert_constraint() was called for that edge |
[in] | i,j | the vertices of the first segment |
[in] | E2 | the index of the second edge, corresponding to the value of ncnstr() when insert_constraint() was called for that edge |
[in] | k,l | the vertices of the second segment |
i
, j
] and [k
, l
] Implemented in GEO::MeshInTriangle, GEO::ExactCDT2d, and GEO::CDT2d.
|
inlineprotected |
|
inlineprotected |
|
inlineprotected |
|
protected |
|
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.
[in] | from_v | the vertex. Cannot be a vertex incident to the border. |
Restores Delaunay condition starting from the triangles incident to a given vertex.
[in] | v | the vertex |
[in] | S | a 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.
Gets an edge constraint from an edge constraint iterator.
[in] | ecit | the edge constraint iterator. Needs to be a valid iterator, different from index_t(-1). |
Gets the successor of an edge constraint iterator.
[in] | ecit | the edge constraint iterator |
ecit
or index_t(-1) if ecit
is the last of the list.
|
inlineprotected |
Finds the index of an integer in an array of three integers.
[in] | T | a const pointer to an array of three integers |
[in] | v | the integer to retrieve in T |
v
in T
T
are different and one of them is equal to v
. Finds the edges intersected by a constraint.
[in] | i,j | the two vertices of the constraint |
[out] | Q | for 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.
|
inlineprotected |
Calls a user-defined function for each triangle around a vertex.
[in] | v | the vertex |
[in] | doit | the 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. |
|
protectedpure virtual |
Incircle predicate.
[in] | i,j,k | the three vertices of a triangle |
[in] | l | another vertex |
POSITIVE | if l is inside the circumscribed circle of the triangle |
ZERO | if l is on the circumscribed circle of the triangle |
NEGATIVE | if l is outside the circumscribed circle of the triangle |
Implemented in GEO::MeshInTriangle, GEO::ExactCDT2d, and GEO::CDT2d.
Inserts a new point.
[in] | v | the index of the new point, supposed to be equal to nv() |
[in] | hint | an optional triangle, not too far away from the point to be inserted |
Inserts a constraint.
[in] | i,j | the indices of the two vertices of the constrained segment |
|
protected |
Inserts a vertex in an edge.
[in] | v | the vertex to be inserted |
[in] | t | a triangle incident to the edge |
[in] | le | the local index of the edge in t |
[out] | S | DList of created triangles, ignored if uninitialized |
Inserts a vertex in a triangle.
[in] | v | the vertex to be inserted |
[in] | t | the triangle |
[out] | S | optional DList of created triangles |
|
protected |
Tests whether triange t and its neighbor accross edge 0 form a strictly convex quad.
true | if triange t and its neighbor accross edge 0 form a strictly convex quad |
false | otherwise |
|
protected |
Locates a vertex.
[in] | v | the vertex index |
[in] | hint | an optional triangle, not too far away from the point to be inserted |
[out] | orient | a 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. |
v
|
protectedpure virtual |
Orientation predicate.
[in] | i,j,k | three vertices |
Implemented in GEO::MeshInTriangle, GEO::ExactCDT2d, and GEO::CDT2d.
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.
[in] | remove_internal_holes | if 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).
|
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.
|
pure virtual |
Saves this CDT to a geogram mesh file.
[in] | filename | where to save this CDT |
Implemented in GEO::MeshInTriangle, GEO::ExactCDT2d, and GEO::CDT2d.
|
inlineprotected |
Tests whether an edge triangle and a segment have a frank intersection.
[in] | v1,v2 | the two extremities of the segment |
[in] | t | a triangle |
[in] | le | local edge index in 0,1,2 |
true | if edge le of t has a frank intersection with edge v1 , v2 |
false | otherwise |
|
inline |
|
protected |
Swaps an edge.
Swaps edge 0 of t1
. Vertex 0 of t1
is vertex 0 of the two new triangles.
[in] | t1 | a triangle index. Its edge opposite to vertex 0 is swapped |
[in] | swap_t1_t2 | if 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). |
|
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
[in] | t1 | a triangle |
[in] | le1 | a local edge of t1 , in 0,1,2 |
[in] | prev_t2_adj_e2 | the triangle adjacent to t2 that t1 will replace, where t2 = Tadj(t1,le1) |
|
inlineprotected |
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:
where 'cnstr' corresponds to the value of ncnstr() when insert_constraint() was called for that constraint.
[in] | t | a triangle |
[in] | le | local edge index, in 0,1,2 |
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.
[in] | t | the triangle |
[in] | le | the local index of the edge (0,1,2) in the triangle |
Tests whether a triangle edge is Delaunay.
returns true also for constrained edges and edges on borders
|
inlineprotected |
|
inlineprotected |
|
protected |
|
protected |
true if incircle() is 100% exact
|
protected |
|
protected |
|
protected |