CGAL 5.4.1 - 3D Periodic Triangulations
CGAL::Periodic_3_regular_triangulation_3< PT, TDS > Class Template Reference

#include <CGAL/Periodic_3_regular_triangulation_3.h>

Inherits from

CGAL::Periodic_3_triangulation_3< PT, TDS >.

Definition

Types

typedef TDS::Point_3 Bare_point
 
typedef TDS::Weighted_point_3 Weighted_point
 The type for points p of weighted points \( {p}^{(w)}=(p,w_p)\).
 

Creation

 Periodic_3_regular_triangulation_3 (const Iso_cuboid &domain=Iso_cuboid(0, 0, 0, 1, 1, 1), const Geom_traits &traits=Geom_traits())
 Creates an empty periodic regular triangulation rt, with domain as original domain and possibly specifying a traits class traits. More...
 
 Periodic_3_regular_triangulation_3 (const Periodic_3_regular_triangulation_3 &rt1)
 Copy constructor.
 
template<class InputIterator >
 Periodic_3_regular_triangulation_3 (InputIterator first, InputIterator last, const Iso_cuboid &domain=Iso_cuboid(0, 0, 0, 1, 1, 1), const Geom_traits &traits=Geom_traits(), bool is_large_point_set=false)
 Equivalent to constructing an empty triangulation with the optional domain and traits class arguments and calling insert(first,last). More...
 

Access Functions

size_type number_of_hidden_points () const
 Returns the number of hidden points.
 

Insertion

The following methods insert points in the triangulation ensuring the property that all power spheres are regular.

The inserted weighted points must lie in the original domain (see Section The Flat Torus of the user manual). Note that the insertion of a new point can cause a switch from computing in the 27-sheeted covering space to computing in the 1-sheeted covering space, which invalidates some Vertex_handles and Cell_handles.

Vertex_handle insert (const Weighted_point &p, Cell_handle start=Cell_handle())
 Inserts the point p in the triangulation and returns the corresponding vertex. More...
 
Vertex_handle insert (const Weighted_point &p, Locate_type lt, Cell_handle loc, int li, int lj)
 Inserts point p in the triangulation and returns the corresponding vertex. More...
 

The following method allows one to insert several points.

template<class InputIterator >
std::ptrdiff_t insert (InputIterator first, InputIterator last, bool is_large_point_set=false)
 Inserts the points in the iterator range \( \left[\right.\)first, last \( \left.\right)\). More...
 

Removal

void remove (Vertex_handle v)
 Removes the vertex v from the triangulation. More...
 
template<class InputIterator >
std::ptrdiff_t remove (InputIterator first, InputIterator beyond)
 Removes the vertices specified by the iterator range (first, beyond) of value type Vertex_handle. More...
 

Queries

Bounded_side side_of_power_sphere (Cell_handle c, const Weighted_point &p, const Offset &off=Offset(0, 0, 0)) const
 Returns a value indicating the position of the (weighted point, offset) pair (p,off) with respect to the power sphere of c. More...
 
Vertex_handle nearest_power_vertex (const Bare_point &p, Cell_handle c=Cell_handle()) const
 Returns the vertex of the triangulation which is nearest to \( p\) with respect to the power distance. More...
 

A point-offset pair (p,off) is said to be in conflict with a cell c iff rt.side_of_power_sphere(c, p, off) returns ON_BOUNDED_SIDE.

The set of cells that are in conflict with (p,off) is star-shaped.

template<class OutputIteratorBoundaryFacets , class OutputIteratorCells , class OutputIteratorInternalFacets >
Triple< OutputIteratorBoundaryFacets, OutputIteratorCells, OutputIteratorInternalFacets > find_conflicts (const Weighted_point &p, Cell_handle c, OutputIteratorBoundaryFacets bfit, OutputIteratorCells cit, OutputIteratorInternalFacets ifit) const
 Computes the conflict hole induced by p. More...
 
template<class OutputIterator >
OutputIterator vertices_in_conflict (const Weighted_point &p, Cell_handle c, OutputIterator res) const
 Similar to find_conflicts(), but reports the vertices that are on the boundary of the conflict hole of p, in the output iterator res. More...
 

Power diagram

CGAL offers several functions to display the Voronoi diagram of a set of points in 3D.

Note that a traits class providing exact constructions should be used in order to guarantee the computation of the Voronoi diagram (as opposed to computing the triangulation only, which requires only exact predicates).

Bare_point canonical_dual (Cell_handle c) const
 Returns the representative of the weighted circumcenter of the four vertices of c that lies in the original domain domain.
 
Bare_point dual (Cell_handle c) const
 Returns the circumcenter of the four vertices of c.
 
Periodic_segment_3 dual (Facet f) const
 Returns the dual of the facet f, which is a periodic segment.
 
Periodic_segment_3 dual (Cell_handle c, int i) const
 same as the previous method for the facet (c,i). More...
 
template<class OutputIterator >
OutputIterator dual (Edge e, OutputIterator pts) const
 Returns in the output iterator the points of the dual polygon of the edge e in the same order as the Facet_circulator returns facets incident to the edge e. More...
 
template<class OutputIterator >
OutputIterator dual (Cell_handle c, int i, int j, OutputIterator pts) const
 same as the previous method for the edge (c,i,j). More...
 
template<class OutputIterator >
OutputIterator dual (Vertex_handle v, OutputIterator pts) const
 Returns in the output iterator the points of the dual polyhedron of the vertex v in no particular order. More...
 
template<class Stream >
Stream & draw_dual (Stream &os)
 Sends the set of duals to all the facets of rt into os.
 
Geom_traits::FT dual_volume (Vertex_handle v) const
 Returns the volume of the cell of the power diagram dual to v.
 
Bare_point dual_centroid (Vertex_handle v) const
 Returns the centroid of the cell of the power diagram dual to v.
 

Checking

These methods are mainly a debugging help for the users of advanced features.

bool is_valid (bool verbose=false) const
 Checks the combinatorial validity of the triangulation and the validity of its geometric embedding (see Section Representation). More...
 
bool is_valid (Cell_handle c, bool verbose=false) const
 Checks the combinatorial and geometric validity of the cell (see Section Representation). More...
 

Additional Inherited Members

- Public Types inherited from CGAL::Periodic_3_triangulation_3< PT, TDS >
enum  Locate_type
 The enum Locate_type is defined by Periodic_3_triangulation_3 to specify which case occurs when locating a point in the triangulation. More...
 
enum  Iterator_type
 The enum Iterator_type is defined by Periodic_3_triangulation_3 to specify the behavior of geometric iterators. More...
 
typedef PT Geom_traits
 
typedef TDS Triangulation_data_structure
 
typedef Geom_traits::Periodic_3_offset_3 Offset
 
typedef Geom_traits::Iso_cuboid_3 Iso_cuboid
 A type representing an axis-aligned cuboid. More...
 
typedef array< int, 3 > Covering_sheets
 Integer triple to store the number of sheets in each direction of space.
 
typedef TDS::Vertex::Point Point
 The point type of the triangulation. More...
 
typedef Geom_traits::Point_3 Point_3
 The geometric basic 3D point type.
 
typedef Geom_traits::Segment_3 Segment
 
typedef Geom_traits::Triangle_3 Triangle
 
typedef Geom_traits::Tetrahedron_3 Tetrahedron
 
typedef std::pair< Point, OffsetPeriodic_point
 Represents a point-offset pair. More...
 
typedef std::pair< Point_3, OffsetPeriodic_point_3
 Represents a point-offset pair. More...
 
typedef array< Periodic_point, 2 > Periodic_segment
 A periodic segment. More...
 
typedef array< Periodic_point_3, 2 > Periodic_segment_3
 A periodic segment. More...
 
typedef array< Periodic_point, 3 > Periodic_triangle
 A periodic triangle. More...
 
typedef array< Periodic_point_3, 3 > Periodic_triangle_3
 A periodic triangle. More...
 
typedef array< Periodic_point, 4 > Periodic_tetrahedron
 A periodic tetrahedron. More...
 
typedef array< Periodic_point_3, 4 > Periodic_tetrahedron_3
 A periodic tetrahedron. More...
 
typedef Triangulation_data_structure::Vertex Vertex
 
typedef Triangulation_data_structure::Cell Cell
 
typedef Triangulation_data_structure::Edge Edge
 
typedef Triangulation_data_structure::Facet Facet
 
typedef Triangulation_data_structure::Vertex_handle Vertex_handle
 handle to a vertex
 
typedef Triangulation_data_structure::Cell_handle Cell_handle
 handle to a cell
 
typedef Triangulation_data_structure::size_type size_type
 Size type (an unsigned integral type)
 
typedef Triangulation_data_structure::difference_type difference_type
 Difference type (a signed integral type)
 
typedef Triangulation_data_structure::Cell_iterator Cell_iterator
 iterator over cells
 
typedef Triangulation_data_structure::Facet_iterator Facet_iterator
 iterator over facets
 
typedef Triangulation_data_structure::Edge_iterator Edge_iterator
 iterator over edges
 
typedef Triangulation_data_structure::Vertex_iterator Vertex_iterator
 iterator over vertices
 
typedef unspecified_type Unique_vertex_iterator
 iterator over the vertices whose corresponding points lie in the original domain, i.e. for each set of periodic copies the Unique_vertex_iterator iterates over exactly one representative.
 
typedef Triangulation_data_structure::Cell_circulator Cell_circulator
 circulator over all cells incident to a given edge
 
typedef Triangulation_data_structure::Facet_circulator Facet_circulator
 circulator over all facets incident to a given edge
 
typedef unspecified_type Periodic_tetrahedron_iterator
 iterator over the tetrahedra corresponding to cells of the triangulation.
 
typedef unspecified_type Periodic_triangle_iterator
 iterator over the triangles corresponding to facets of the triangulation.
 
typedef unspecified_type Periodic_segment_iterator
 iterator over the segments corresponding to edges of the triangulation.
 
typedef unspecified_type Periodic_point_iterator
 iterator over the points corresponding to vertices of the triangulation.
 
- Public Member Functions inherited from CGAL::Periodic_3_triangulation_3< PT, TDS >
 Periodic_3_triangulation_3 (const Iso_cuboid &domain=Iso_cuboid(0, 0, 0, 1, 1, 1), const Geom_traits &traits=Geom_traits())
 Introduces an empty triangulation t with domain as original domain. More...
 
 Periodic_3_triangulation_3 (const Periodic_3_triangulation_3 &tr)
 Copy constructor. More...
 
Periodic_3_triangulation_3operator= (const Periodic_3_triangulation_3 &tr)
 The triangulation tr is duplicated, and modifying the copy after the duplication does not modify the original. More...
 
void swap (Periodic_3_triangulation_3 &tr)
 The triangulations tr and t are swapped. More...
 
void clear ()
 Deletes all vertices and all cells of t.
 
const Geom_traitsgeom_traits () const
 Returns a const reference to the geometric traits object.
 
const Triangulation_data_structuretds () const
 Returns a const reference to the triangulation data structure.
 
Iso_cuboid domain () const
 Returns the original domain.
 
Covering_sheets number_of_sheets () const
 This is an advanced function. More...
 
Offset neighbor_offset (Cell_handle ch, int i) const
 Get the offset between the origins of the internal offset coordinate systems of two neighboring cells with respect from ch to its i-th neighbor.
 
Triangulation_data_structuretds ()
 This is an advanced function. More...
 
void set_domain (const Iso_cuboid dom)
 This is an advanced function. More...
 
bool is_extensible_triangulation_in_1_sheet_h1 () const
 The current triangulation remains a triangulation in the 1-sheeted covering space even after adding points if this method returns true. More...
 
bool is_extensible_triangulation_in_1_sheet_h2 () const
 The same as is_extensible_triangulation_in_1_sheet_h1() but with a more precise heuristic, i.e. it might answer true in cases in which is_extensible_triangulation_in_1_sheet_h1() would not. More...
 
bool is_triangulation_in_1_sheet () const
 Returns true if the current triangulation would still be a triangulation in the 1-sheeted covering space, returns false otherwise.
 
void convert_to_1_sheeted_covering () const
 Converts the current triangulation into the same periodic triangulation in the 1-sheeted covering space. More...
 
void convert_to_27_sheeted_covering () const
 Converts the current triangulation into the same periodic triangulation in the 27-sheeted covering space. More...
 
size_type number_of_vertices () const
 Returns the number of vertices. More...
 
size_type number_of_cells () const
 Returns the number of cells. More...
 
size_type number_of_stored_vertices () const
 Returns the number of vertices in the data structure. More...
 
size_type number_of_stored_cells () const
 Returns the number of cells in the data structure. More...
 
size_type number_of_edges () const
 Returns the number of edges. More...
 
size_type number_of_facets () const
 Returns the number of facets. More...
 
size_type number_of_stored_edges () const
 Returns the number of edges in the data structure. More...
 
size_type number_of_stored_facets () const
 Returns the number of facets in the data structure. More...
 
Periodic_point periodic_point (const Vertex_handle v) const
 Returns the periodic point given by vertex v. More...
 
Periodic_point periodic_point (const Cell_handle c, int i) const
 If t is represented in the 1-sheeted covering space, this function returns the periodic point given by the \( i\)-th vertex of cell c, that is the point in the original domain and the offset of the vertex in c. More...
 
Periodic_segment periodic_segment (const Cell_handle c, int i, int j) const
 Returns the periodic segment formed by the two point-offset pairs corresponding to the two vertices of edge (c,i,j). More...
 
Periodic_segment periodic_segment (const Cell_handle c, Offset offset, int i, int j) const
 Returns the periodic segment formed by the two point-offset pairs corresponding to the two vertices of edge (c,i,j). More...
 
Periodic_segment periodic_segment (const Edge &e) const
 Same as the previous method for edge e.
 
Periodic_triangle periodic_triangle (const Cell_handle c, int i) const
 Returns the periodic triangle formed by the three point-offset pairs corresponding to the three vertices of facet (c,i). More...
 
Periodic_triangle periodic_triangle (const Cell_handle c, Offset offset, int i) const
 Returns the periodic triangle formed by the three point-offset pairs corresponding to the three vertices of facet (c,i). More...
 
Periodic_triangle periodic_triangle (const Facet &f) const
 Same as the previous method for facet f.
 
Periodic_tetrahedron periodic_tetrahedron (const Cell_handle c) const
 Returns the periodic tetrahedron formed by the four point-offset pairs corresponding to the four vertices of c.
 
Periodic_tetrahedron periodic_tetrahedron (const Cell_handle c, Offset offset) const
 Returns the periodic tetrahedron formed by the four point-offset pairs corresponding to the four vertices of c. More...
 
Point_3 construct_point (const PP &pp) const
 Converts the periodic point of type PP to a Point_3. More...
 
Point_3 construct_point (const Point &p) const
 Converts the Point p to a Point_3.
 
Point_3 construct_point (const P &p1, const Offset &o1) const
 Same as above, with offsets.
 
Segment construct_segment (const PS &s) const
 Converts the periodic segment of type PS to a Segment. More...
 
Segment construct_segment (const P &p1, const P &p2) const
 Creates a segment from two points. More...
 
Segment construct_segment (const P &p1, const P &p2, const Offset &o1, const Offset &o2) const
 Same as above, with offsets.
 
Triangle construct_triangle (const PT &t) const
 Converts the periodic triangle of type PT to a Triangle. More...
 
Triangle construct_triangle (const P &p1, const P &p2, const P &p3) const
 Creates a triangle from three points. More...
 
Triangle construct_triangle (const P &p1, const P &p2, const P &p3, const Offset &o1, const Offset &o2, const Offset &o3) const
 Same as above, with offsets.
 
Tetrahedron construct_tetrahedron (const PT &t) const
 Converts the periodic tetrahedron of type PT to a Tetrahedron. More...
 
Tetrahedron construct_tetrahedron (const P &p1, const P &p2, const P &p3, const P &p4) const
 Creates a tetrahedron from four points. More...
 
Tetrahedron construct_tetrahedron (const P &p1, const P &p2, const P &p3, const P &p4, const Offset &o1, const Offset &o2, const Offset &o3, const Offset &o4) const
 Same as above, with offsets.
 
bool is_vertex (const Point &p, Vertex_handle &v) const
 Tests whether p is a vertex of t by locating p in the triangulation. More...
 
bool is_vertex (Vertex_handle v) const
 Tests whether v is a vertex of t.
 
bool is_edge (Vertex_handle u, Vertex_handle v, Cell_handle &c, int &i, int &j) const
 Tests whether (u,v) is an edge of t. More...
 
bool is_edge (Vertex_handle u, const Offset &offu, Vertex_handle v, const Offset &offv, Cell_handle &c, int &i, int &j) const
 Tests whether ((u,offu),(v,offu)) is an edge of t. More...
 
bool is_facet (Vertex_handle u, Vertex_handle v, Vertex_handle w, Cell_handle &c, int &i, int &j, int &k) const
 Tests whether (u,v,w) is a facet of t. More...
 
bool is_facet (Vertex_handle u, const Offset &offu, Vertex_handle v, const Offset &offv, Vertex_handle w, const Offset &offw, Cell_handle &c, int &i, int &j, int &k) const
 Tests whether ((u,offu),(v,offv),(w,offw)) is a facet of t. More...
 
bool is_cell (Cell_handle c) const
 Tests whether c is a cell of t.
 
bool is_cell (Vertex_handle u, Vertex_handle v, Vertex_handle w, Vertex_handle x, Cell_handle &c, int &i, int &j, int &k, int &l) const
 Tests whether (u,v,w,x) is a cell of t. More...
 
bool is_cell (Vertex_handle u, Vertex_handle v, Vertex_handle w, Vertex_handle x, Cell_handle &c) const
 Tests whether (u,v,w,x) is a cell of t and computes this cell c. More...
 
bool is_cell (Vertex_handle u, const Offset &offu, Vertex_handle v, const Offset &offv, Vertex_handle w, const Offset &offw, Vertex_handle x, const Offset &offx, Cell_handle &c, int &i, int &j, int &k, int &l) const
 Tests whether ((u,offu),(v,offv),(w,offv),(x,offx)) is a cell of t. More...
 
bool is_cell (Vertex_handle u, const Offset &offu, Vertex_handle v, const Offset &offv, Vertex_handle w, const Offset &offw, Vertex_handle x, const Offset &offx, Cell_handle &c) const
 Tests whether ((u,offu),(v,offv),(w,offv),(x,offx)) is a cell of t and computes this cell c. More...
 
bool has_vertex (const Facet &f, Vertex_handle v, int &j) const
 If v is a vertex of f, then j is the index of v in the cell f.first, and the method returns true.
 
bool has_vertex (Cell_handle c, int i, Vertex_handle v, int &j) const
 Same for facet (c,i). More...
 
bool has_vertex (const Facet &f, Vertex_handle v) const
 
bool has_vertex (Cell_handle c, int i, Vertex_handle v) const
 Same as the first two methods, but these two methods do not return the index of the vertex.
 
bool are_equal (Cell_handle c, int i, Cell_handle n, int j) const
 
bool are_equal (const Facet &f, const Facet &g) const
 
bool are_equal (const Facet &f, Cell_handle n, int j) const
 
Cell_handle locate (const Point &query, Cell_handle start=Cell_handle()) const
 Returns the cell that contains the query in its interior. More...
 
Cell_handle locate (const Point &query, Offset &locate_offset, Cell_handle start=Cell_handle()) const
 Returns the cell that contains the query in its interior. More...
 
Cell_handle locate (const Point &query, Locate_type &lt, int &li, int &lj, Cell_handle start=Cell_handle()) const
 The \( k\)-face that contains query in its interior is returned, by means of the cell returned together with lt, which is set to the locate type of the query (VERTEX, EDGE, FACET, CELL) and two indices li and lj that specify the \( k\)-face of the cell containing query. More...
 
Cell_handle locate (const Point &query, Offset &locate_offset, Locate_type &lt, int &li, int &lj, Cell_handle start=Cell_handle()) const
 The \( k\)-face that contains query in its interior is returned, by means of the cell returned together with lt, which is set to the locate type of the query (VERTEX, EDGE, FACET, CELL) and two indices li and lj that specify the \( k\)-face of the cell containing query. More...
 
Cell_handle inexact_locate (const Point &query, Cell_handle start=Cell_handle()) const
 Same as locate() but uses inexact predicates. More...
 
Bounded_side side_of_cell (const Point &p, Cell_handle c, Locate_type &lt, int &li, int &lj) const
 Returns a value indicating on which side of the oriented boundary of c the point p lies. More...
 
Vertex_iterator vertices_begin () const
 Starts at an arbitrary vertex. More...
 
Vertex_iterator vertices_end () const
 Past-the-end iterator.
 
Edge_iterator edges_begin () const
 Starts at an arbitrary edge. More...
 
Edge_iterator edges_end () const
 Past-the-end iterator.
 
Facet_iterator facets_begin () const
 Starts at an arbitrary facet. More...
 
Facet_iterator facets_end () const
 Past-the-end iterator.
 
Cell_iterator cells_begin () const
 Starts at an arbitrary cell. More...
 
Cell_iterator cells_end () const
 Past-the-end iterator.
 
Unique_vertex_iterator unique_vertices_begin () const
 Starts at an arbitrary vertex. More...
 
Unique_vertex_iterator unique_vertices_end () const
 Past-the-end iterator.
 
Periodic_point_iterator periodic_points_begin (Iterator_type it=STORED) const
 Iterates over the points of the triangulation. More...
 
Periodic_point_iterator periodic_points_end (Iterator_type it=STORED) const
 Past-the-end iterator. More...
 
Periodic_segment_iterator periodic_segments_begin (Iterator_type it=STORED) const
 Iterates over the segments of the triangulation. More...
 
Periodic_segment_iterator periodic_segments_end (Iterator_type it=STORED) const
 Past-the-end iterator. More...
 
Periodic_triangle_iterator periodic_triangles_begin (Iterator_type it=STORED) const
 Iterates over the triangles of the triangulation. More...
 
Periodic_triangle_iterator periodic_triangles_end (Iterator_type it=STORED) const
 Past-the-end iterator. More...
 
Periodic_tetrahedron_iterator periodic_tetrahedra_begin (Iterator_type it=STORED) const
 Iterates over the tetrahedra of the triangulation. More...
 
Periodic_tetrahedron_iterator periodic_tetrahedra_end (Iterator_type it=STORED) const
 Past-the-end iterator. More...
 
Cell_circulator incident_cells (Edge e) const
 Starts at an arbitrary cell incident to e.
 
Cell_circulator incident_cells (Cell_handle c, int i, int j) const
 As above for edge (i,j) of c.
 
Cell_circulator incident_cells (Edge e, Cell_handle start) const
 Starts at cell start. More...
 
Cell_circulator incident_cells (Cell_handle c, int i, int j, Cell_handle start) const
 As above for edge (i,j) of c.
 
Facet_circulator incident_facets (Edge e) const
 Starts at an arbitrary facet incident to e.
 
Facet_circulator incident_facets (Cell_handle c, int i, int j) const
 As above for edge (i,j) of c.
 
Facet_circulator incident_facets (Edge e, Facet start) const
 Starts at facet start. More...
 
Facet_circulator incident_facets (Edge e, Cell_handle start, int f) const
 Starts at facet of index f in start.
 
Facet_circulator incident_facets (Cell_handle c, int i, int j, Facet start) const
 As above for edge (i,j) of c.
 
Facet_circulator incident_facets (Cell_handle c, int i, int j, Cell_handle start, int f) const
 As above for edge (i,j) of c and facet (start,f).
 
OutputIterator incident_cells (Vertex_handle v, OutputIterator cells) const
 Copies the Cell_handles of all cells incident to v to the output iterator cells. More...
 
OutputIterator incident_facets (Vertex_handle v, OutputIterator facets) const
 Copies the Facets incident to v to the output iterator facets. More...
 
OutputIterator incident_edges (Vertex_handle v, OutputIterator edges) const
 Copies the Edges incident to v to the output iterator edges. More...
 
OutputIterator adjacent_vertices (Vertex_handle v, OutputIterator vertices) const
 Copies the Vertex_handles of all vertices adjacent to v to the output iterator vertices. More...
 
size_type degree (Vertex_handle v) const
 Returns the degree of a vertex, that is, the number of adjacent vertices. More...
 
int mirror_index (Cell_handle c, int i) const
 Returns the index of c in its \( i^{th}\) neighbor. More...
 
Vertex_handle mirror_vertex (Cell_handle c, int i) const
 Returns the vertex of the \( i^{th}\) neighbor of c that is opposite to c. More...
 
Facet mirror_facet (Facet f) const
 Returns the same facet viewed from the other adjacent cell.
 
bool is_valid (bool verbose=false) const
 Checks the combinatorial validity of the triangulation. More...
 
bool is_valid (Cell_handle c, bool verbose=false) const
 Checks the combinatorial validity of the cell by calling the is_valid method of the Triangulation_data_structure cell class. More...
 

Constructor & Destructor Documentation

◆ Periodic_3_regular_triangulation_3() [1/2]

template<typename PT , typename TDS >
CGAL::Periodic_3_regular_triangulation_3< PT, TDS >::Periodic_3_regular_triangulation_3 ( const Iso_cuboid domain = Iso_cuboid(0, 0, 0, 1, 1, 1),
const Geom_traits traits = Geom_traits() 
)

Creates an empty periodic regular triangulation rt, with domain as original domain and possibly specifying a traits class traits.

Precondition
domain is a cube.

◆ Periodic_3_regular_triangulation_3() [2/2]

template<typename PT , typename TDS >
template<class InputIterator >
CGAL::Periodic_3_regular_triangulation_3< PT, TDS >::Periodic_3_regular_triangulation_3 ( InputIterator  first,
InputIterator  last,
const Iso_cuboid domain = Iso_cuboid(0, 0, 0, 1, 1, 1),
const Geom_traits traits = Geom_traits(),
bool  is_large_point_set = false 
)

Equivalent to constructing an empty triangulation with the optional domain and traits class arguments and calling insert(first,last).

Precondition
The value_type of first and last are Weighted_points lying inside domain and domain is a cube. Their weights are non-negative and smaller than 1/64 times the squared cube edge length. If the fourth argument is_large_point_set is set to true a heuristic for optimizing the insertion of large point sets is applied.

Member Function Documentation

◆ dual() [1/4]

template<typename PT , typename TDS >
Periodic_segment_3 CGAL::Periodic_3_regular_triangulation_3< PT, TDS >::dual ( Cell_handle  c,
int  i 
) const

same as the previous method for the facet (c,i).

Precondition
\( i\in\{0,1,2,3\}\)

◆ dual() [2/4]

template<typename PT , typename TDS >
template<class OutputIterator >
OutputIterator CGAL::Periodic_3_regular_triangulation_3< PT, TDS >::dual ( Edge  e,
OutputIterator  pts 
) const

Returns in the output iterator the points of the dual polygon of the edge e in the same order as the Facet_circulator returns facets incident to the edge e.

The points form the dual polygon in \( \mathbb R^3\), so they do not necessarily all lie inside the original domain.

◆ dual() [3/4]

template<typename PT , typename TDS >
template<class OutputIterator >
OutputIterator CGAL::Periodic_3_regular_triangulation_3< PT, TDS >::dual ( Cell_handle  c,
int  i,
int  j,
OutputIterator  pts 
) const

same as the previous method for the edge (c,i,j).

Precondition
\( i,j\in\{0,1,2,3\}, i\neq j\)

◆ dual() [4/4]

template<typename PT , typename TDS >
template<class OutputIterator >
OutputIterator CGAL::Periodic_3_regular_triangulation_3< PT, TDS >::dual ( Vertex_handle  v,
OutputIterator  pts 
) const

Returns in the output iterator the points of the dual polyhedron of the vertex v in no particular order.

The points form the dual polyhedron in \( \mathbb R^3\), so they do not necessarily lie all inside the original domain.

◆ find_conflicts()

template<typename PT , typename TDS >
template<class OutputIteratorBoundaryFacets , class OutputIteratorCells , class OutputIteratorInternalFacets >
Triple<OutputIteratorBoundaryFacets,OutputIteratorCells,OutputIteratorInternalFacets> CGAL::Periodic_3_regular_triangulation_3< PT, TDS >::find_conflicts ( const Weighted_point p,
Cell_handle  c,
OutputIteratorBoundaryFacets  bfit,
OutputIteratorCells  cit,
OutputIteratorInternalFacets  ifit 
) const

Computes the conflict hole induced by p.

The starting cell c must be in conflict. Then this function returns respectively in the output iterators:

  • cit: the cells in conflict.
  • bfit: the facets on the boundary, that is, the facets (t, i) where the cell t is in conflict, but t->neighbor(i) is not.
  • ifit: the facets inside the hole, that is, delimiting two cells in conflict.

Returns the Triple composed of the resulting output iterators.

Precondition
c is in conflict with p and p lies in the original domain domain. Its weight is non-negative and smaller than 1/64 times the squared cube edge length.

◆ insert() [1/3]

template<typename PT , typename TDS >
Vertex_handle CGAL::Periodic_3_regular_triangulation_3< PT, TDS >::insert ( const Weighted_point p,
Cell_handle  start = Cell_handle() 
)

Inserts the point p in the triangulation and returns the corresponding vertex.

The optional argument start is used as a starting place for the point location.

If this insertion creates a vertex, this vertex is returned.

If p coincides with an existing vertex and has a greater weight, then the existing weighted point becomes hidden (see Periodic_3RegularTriangulationCellBase_3) and p replaces it as vertex of the triangulation.

If p coincides with an already existing vertex (both point and weights being equal), then this vertex is returned and the triangulation remains unchanged.

Otherwise, if p does not appear as a vertex of the triangulation, then it is stored as a hidden point and this method returns the default constructed handle.

Precondition
p lies in the original domain. Its weight is non-negative and smaller than 1/64 times the squared cube edge length.

◆ insert() [2/3]

template<typename PT , typename TDS >
Vertex_handle CGAL::Periodic_3_regular_triangulation_3< PT, TDS >::insert ( const Weighted_point p,
Locate_type  lt,
Cell_handle  loc,
int  li,
int  lj 
)

Inserts point p in the triangulation and returns the corresponding vertex.

Similar to the above insert() function, but takes as additional parameter the return values of a previous location query. See description of Periodic_3_triangulation_3::locate().

Precondition
p lies in the original domain. Its weight is non-negative and smaller than 1/64 times the squared cube edge length.

◆ insert() [3/3]

template<typename PT , typename TDS >
template<class InputIterator >
std::ptrdiff_t CGAL::Periodic_3_regular_triangulation_3< PT, TDS >::insert ( InputIterator  first,
InputIterator  last,
bool  is_large_point_set = false 
)

Inserts the points in the iterator range \( \left[\right.\)first, last \( \left.\right)\).

Returns the number of inserted points. This function uses spatial sorting (cf. chapter Spatial Sorting) and therefore is not guaranteed to insert the points following the order of InputIterator. If the third argument is_large_point_set is set to true, a heuristic for optimizing the insertion of large point sets is applied.

Precondition
The value_type of first and last are Weighted_points lying inside the original domain. Their weights are non-negative and smaller than 1/64 times the squared cube edge length.

◆ is_valid() [1/2]

template<typename PT , typename TDS >
bool CGAL::Periodic_3_regular_triangulation_3< PT, TDS >::is_valid ( bool  verbose = false) const

Checks the combinatorial validity of the triangulation and the validity of its geometric embedding (see Section Representation).

Also checks that the power spheres of all cells are regular.

When verbose is set to true, messages describing the first invalidity encountered (if any) are printed.

◆ is_valid() [2/2]

template<typename PT , typename TDS >
bool CGAL::Periodic_3_regular_triangulation_3< PT, TDS >::is_valid ( Cell_handle  c,
bool  verbose = false 
) const

Checks the combinatorial and geometric validity of the cell (see Section Representation).

Also checks that its power sphere is regular.

When verbose is set to true, messages are printed to give a precise indication of the kind of invalidity encountered (if any).

◆ nearest_power_vertex()

template<typename PT , typename TDS >
Vertex_handle CGAL::Periodic_3_regular_triangulation_3< PT, TDS >::nearest_power_vertex ( const Bare_point p,
Cell_handle  c = Cell_handle() 
) const

Returns the vertex of the triangulation which is nearest to \( p\) with respect to the power distance.

This means that the power of the query point p with respect to the weighted point in the returned vertex is smaller than the power of p with respect to the weighted point in any other vertex. Ties are broken arbitrarily. It always returns a vertex corresponding to a point inside domain even if computing in a multiply sheeted covering space. The default constructed handle is returned if the triangulation is empty. The optional argument c is a hint specifying where to start the search.

Precondition
c is a cell of rt and p lies in the original domain domain.

◆ remove() [1/2]

template<typename PT , typename TDS >
void CGAL::Periodic_3_regular_triangulation_3< PT, TDS >::remove ( Vertex_handle  v)

Removes the vertex v from the triangulation.

When computing in the 27-sheeted covering space it removes all 27 copies of v.

◆ remove() [2/2]

template<typename PT , typename TDS >
template<class InputIterator >
std::ptrdiff_t CGAL::Periodic_3_regular_triangulation_3< PT, TDS >::remove ( InputIterator  first,
InputIterator  beyond 
)

Removes the vertices specified by the iterator range (first, beyond) of value type Vertex_handle.

remove() is called for each element of the range. The number of vertices removed is returned; this number does not account for periodic copies of removed vertices.

Precondition
The iterator must not iterate over several periodic copies of the same vertex, use e.g. the Unique_vertex_iterator.

◆ side_of_power_sphere()

template<typename PT , typename TDS >
Bounded_side CGAL::Periodic_3_regular_triangulation_3< PT, TDS >::side_of_power_sphere ( Cell_handle  c,
const Weighted_point p,
const Offset off = Offset(0, 0, 0) 
) const

Returns a value indicating the position of the (weighted point, offset) pair (p,off) with respect to the power sphere of c.

More precisely, it returns:

  • ON_BOUNDED_SIDE if the angle between the weighted point (p,off) and the power sphere of c is larger than \( \pi/2\) or if the ball of (p,off) is included in the power sphere of c.
  • ON_BOUNDARY if the ball (p,off) is orthogonal to the power sphere of c.
  • ON_UNBOUNDED_SIDE if the angle between the weighted point (p,off) and the power sphere of c is less than \( \pi/2\) or if the two balls do not intersect.

◆ vertices_in_conflict()

template<typename PT , typename TDS >
template<class OutputIterator >
OutputIterator CGAL::Periodic_3_regular_triangulation_3< PT, TDS >::vertices_in_conflict ( const Weighted_point p,
Cell_handle  c,
OutputIterator  res 
) const

Similar to find_conflicts(), but reports the vertices that are on the boundary of the conflict hole of p, in the output iterator res.

Returns the resulting output iterator.

Precondition
c is in conflict with p and p lies in the original domain domain. Its weight is non-negative and smaller than 1/64 times the squared cube edge length.