Graphite
Version 3
An experimental 3D geometry processing program
|
Constrained Delaunay triangulation. More...
#include <geogram/delaunay/CDT_2d.h>
Public Member Functions | |
void | clear () override |
Removes everything from this triangulation. More... | |
void | create_enclosing_triangle (const vec2 &p1, const vec2 &p2, const vec2 &p3) |
Creates a first large enclosing triangle. More... | |
void | create_enclosing_quad (const vec2 &p1, const vec2 &p2, const vec2 &p3, const vec2 &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... | |
index_t | insert (const vec2 &p, index_t hint=index_t(-1)) |
Inserts a point. More... | |
void | insert (index_t nb_points, const double *points, index_t *indices=nullptr, bool remove_unreferenced_vertices=false) |
Batch-inserts a set of point. More... | |
void | save (const std::string &filename) const override |
Saves this CDT to a geogram mesh file. More... | |
const vec2 | point (index_t v) const |
Gets a point by index. More... | |
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 | |
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 | |
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... | |
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< vec2 > | point_ |
Protected Attributes inherited from GEO::CDTBase2d | |
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_ |
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_t > | Edge |
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... | |
Constrained Delaunay triangulation.
Example:
If some constraints are intersecting, new vertices are generated. They can be accessed using the function vec2 CDT::point(index_t v). Vertices coming from an intersection are between indices nv1 and CDT::nv(), where nv1 is the value of CDT::nv() before inserting the constraints (nv1 corresponds to the number of times CDT::insert() was called plus the number of points in the enclosing polygon). Note that like input points, constraint intersections are represented using double-precision floating point numbers, which is not always sufficient to ensure robustness. If the input has intersecting constraints and bullet-proof guarantees are needed, one can use ExactCDT2d instead.
If you want only a constrained triangulation (not Delaunay), you can call CDT::set_Delaunay(false) before inserting the points.
If you have many points to insert, you can use the function:
It is much much faster than inserting the points one by one. Internally it uses Amenta et.al's BRIO method (multi-resolution spatial sort).
|
overridevirtual |
Removes everything from this triangulation.
Reimplemented from GEO::CDTBase2d.
void GEO::CDT2d::create_enclosing_quad | ( | const vec2 & | p1, |
const vec2 & | p2, | ||
const vec2 & | p3, | ||
const vec2 & | p4 | ||
) |
Creates a first large enclosing quad.
[in] | p1,p2,p3,p4 | the four vertices of the quad |
The quad needs to be convex. create_enclosing_triangle(), create_enclosing_rectangle() or create_enclosing_quad() need to be called before anything else
|
inline |
Creates a first large enclosing rectangle.
[in] | x1,y1,x2,y2 | rectangle bounds |
create_enclosing_triangle(), create_enclosing_rectangle() or create_enclosing_quad() need to be called before anything else
Creates a first large enclosing triangle.
[in] | p1,p2,p3 | the three vertices of the first triangle |
create_enclosing_triangle(), create_enclosing_rectangle() or create_enclosing_quad() need to be called before anything else
|
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
[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
] Implements GEO::CDTBase2d.
|
overrideprotectedvirtual |
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 |
Implements GEO::CDTBase2d.
void GEO::CDT2d::insert | ( | index_t | nb_points, |
const double * | points, | ||
index_t * | indices = nullptr , |
||
bool | remove_unreferenced_vertices = false |
||
) |
Batch-inserts a set of point.
In general, it is much faster than calling insert() multiple times. Internally it uses a spatial sort (Amenta et.al's BRIO method). On exit, the optional indices
array contains the index mapping. indices[i] may be different from i if there were duplicated points. If there may be duplicated points and if one wants to insert constraint using CDTBase2d::insert_constraint(), one needs to call insert_constraint(indices[i],indices[j]) to get the correct translation of the indices
[in] | points | a contiguous array of all point coordinates |
[in] | nb_points | number of points |
[out] | indices | an optional pointer to an array of size nb_points of indices. On exit, indices[i] contains the index of the mesh vertex that corresponds to the i-th point. If there are duplicated points or if remove_unreferenced_vertices is set, then indices[i] may be different from i |
[in] | remove_unreferenced_vertices | if set, then duplicated vertices are not stored in the vertices array. Internally, this systematically changes the order of the points. For this reason, if this flag is set, then one needs to do index mapping with indices , even when there is no duplicated point |
Computes the orientation predicate in 2d.
Computes the sign of the signed area of the triangle p0, p1, p2.
[in] | p0,p1,p2 | vertices of the triangle as 2d vectors with homogeneous coordinates stored as expansion_nt (arbitrary precision). |
POSITIVE | if the triangle is oriented counter-clockwise |
ZERO | if the triangle is flat |
NEGATIVE | if the triangle is oriented clockwise |
Implements GEO::CDTBase2d.
|
overridevirtual |
Saves this CDT to a geogram mesh file.
[in] | filename | where to save this CDT |
Implements GEO::CDTBase2d.