Geogram  Version 1.8.9-rc
A programming library of geometric algorithms
GEO::MeshInTriangle Class Reference

Meshes a single triangle with the constraints that come from the intersections with the other triangles. More...

#include <geogram/mesh/mesh_surface_intersection_internal.h>

Inheritance diagram for GEO::MeshInTriangle:
GEO::CDTBase2d

Classes

class  Edge
 An edge of the mesh. More...
 
class  Vertex
 A vertex of the triangulation. More...
 

Public Types

typedef exact::vec3h ExactPoint
 

Public Member Functions

 MeshInTriangle (MeshSurfaceIntersection &EM)
 
const Meshmesh () const
 Gets the readonly initial mesh. More...
 
Meshtarget_mesh ()
 Gets the target mesh. More...
 
void set_dry_run (bool x)
 In dry run mode, the computed local triangulations are not inserted in the global mesh. This is for benchmarking. Default is off.
 
void save_constraints (const std::string &filename)
 For debugging, save constraints to a file. More...
 
void begin_facet (index_t f)
 
index_t add_vertex (index_t f2, TriangleRegion R1, TriangleRegion R2)
 
void add_edge (index_t f2, TriangleRegion AR1, TriangleRegion AR2, TriangleRegion BR1, TriangleRegion BR2)
 
void commit ()
 Creates new vertices and new triangles in target mesh.
 
void clear () override
 
void save (const std::string &filename) const override
 Saves this CDT to a geogram mesh file. 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 fro 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 get_constraints (Mesh &M, bool with_edges=true) const
 For debugging, copies the constraints to a mesh.
 
vec3 mesh_vertex (index_t v) const
 
vec3 mesh_facet_vertex (index_t f, index_t lv) const
 
vec2 mesh_vertex_UV (index_t v) const
 
vec2 mesh_facet_vertex_UV (index_t f, index_t lv) const
 
void log_err () const
 
Sign orient2d (index_t v1, index_t v2, index_t v3) const override
 Tests the orientation of three vertices. More...
 
Sign incircle (index_t v1, index_t v2, index_t v3, index_t v4) const override
 Tests the relative position of a point with respect to the circumscribed circle of a triangle. 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...
 
void get_edge_edge_intersection (index_t e1, index_t e2, ExactPoint &I) const
 Computes the intersection between two edges. More...
 
void get_edge_edge_intersection_2D (index_t e1, index_t e2, ExactPoint &I) const
 Auxilliary function used by get_edge_edge_intersection() for the special case when the two edges are coplanar. More...
 
void begin_insert_transaction () override
 
void commit_insert_transaction () override
 
void rollback_insert_transaction () override
 
- 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.
 

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

Detailed Description

Meshes a single triangle with the constraints that come from the intersections with the other triangles.

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.

Definition at line 70 of file mesh_surface_intersection_internal.h.

Member Function Documentation

◆ clear()

void GEO::MeshInTriangle::clear ( )
overridevirtual
See also
CDT2d::clear()

Reimplemented from GEO::CDTBase2d.

◆ create_intersection()

index_t GEO::MeshInTriangle::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.

◆ get_edge_edge_intersection()

void GEO::MeshInTriangle::get_edge_edge_intersection ( index_t  e1,
index_t  e2,
ExactPoint I 
) const
protected

Computes the intersection between two edges.

Parameters
[in]e1,e2the two edges
[out]Ithe intersection

◆ get_edge_edge_intersection_2D()

void GEO::MeshInTriangle::get_edge_edge_intersection_2D ( index_t  e1,
index_t  e2,
ExactPoint I 
) const
protected

Auxilliary function used by get_edge_edge_intersection() for the special case when the two edges are coplanar.

Parameters
[in]e1,e2the two edges
[out]Ithe intersection

◆ incircle()

Sign GEO::MeshInTriangle::incircle ( index_t  v1,
index_t  v2,
index_t  v3,
index_t  v4 
) const
overrideprotectedvirtual

Tests the relative position of a point with respect to the circumscribed circle of a triangle.

Parameters
[in]v1,v2,v3the three vertices of the triangle oriented anticlockwise
[in]v4the point to be tested
Return values
POSITIVEif the point is inside the circle
ZEROif the point is on the circle
NEGATIVEif the point is outside the circle

Implements GEO::CDTBase2d.

◆ mesh()

const Mesh& GEO::MeshInTriangle::mesh ( ) const
inline

Gets the readonly initial mesh.

Returns
a const reference to a copy of the initial mesh

Definition at line 261 of file mesh_surface_intersection_internal.h.

◆ orient2d()

Sign GEO::MeshInTriangle::orient2d ( index_t  v1,
index_t  v2,
index_t  v3 
) const
overrideprotectedvirtual

Tests the orientation of three vertices.

Parameters
[in]v1,v2,v3the three vertices
Return values
POSITIVEif they are in the trigonometric order
ZEROif they are aligned
NEGATIVEif they are in the anti-trigonometric order

Implements GEO::CDTBase2d.

◆ save()

void GEO::MeshInTriangle::save ( const std::string &  filename) const
overridevirtual

Saves this CDT to a geogram mesh file.

Parameters
[in]filenamewhere to save this CDT

Implements GEO::CDTBase2d.

◆ save_constraints()

void GEO::MeshInTriangle::save_constraints ( const std::string &  filename)
inline

For debugging, save constraints to a file.

Parameters
[in]filenamea mesh filename where to solve the constraints (.obj or .geogram)

Definition at line 287 of file mesh_surface_intersection_internal.h.

◆ target_mesh()

Mesh& GEO::MeshInTriangle::target_mesh ( )
inline

Gets the target mesh.

Returns
a reference to the target mesh

Definition at line 269 of file mesh_surface_intersection_internal.h.


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