CGAL 5.3.2 - 3D Periodic Triangulations
CGAL::Periodic_3_triangulation_3< Traits, TDS > Class Template Reference

#include <CGAL/Periodic_3_triangulation_3.h>

Definition

The class Periodic_3_triangulation_3 represents a 3-dimensional triangulation of a point set in \( \mathbb T_c^3\).

Template Parameters
Traitsmust be a model of the concept Periodic_3TriangulationTraits_3.
TDSmust be a model of the concept TriangulationDataStructure_3 with vertex and cell are models of Periodic_3TriangulationDSVertexBase_3 and Periodic_3TriangulationDSCellBase_3, Its default value is:
See also
CGAL::Periodic_3_Delaunay_triangulation_3<Traits, TDS>
CGAL::Periodic_3_regular_triangulation_3<Traits, TDS>

Related Functions

(Note that these are not member functions.)

istream & operator>> (istream &is, Periodic_3_triangulation_3 &t)
 Reads a triangulation from is and stores it in t. More...
 
ostream & operator<< (ostream &os, const Periodic_3_triangulation_3 &t)
 Writes the triangulation t into os. More...
 

Types

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

Only vertices ( \( 0\)-faces) and cells ( \( 3\)-faces) are stored.

Edges ( \( 1\)-faces) and facets ( \( 2\)-faces) are not explicitly represented and thus there are no corresponding classes (see Section Representation).

typedef Triangulation_data_structure::Vertex Vertex
 
typedef Triangulation_data_structure::Cell Cell
 
typedef Triangulation_data_structure::Edge Edge
 
typedef Triangulation_data_structure::Facet Facet
 

The vertices and faces of the triangulations are accessed through handles, iterators and circulators.

A handle is a type which supports the two dereference operators and operator->. The Handle concept is documented in the support library. Iterators and circulators are bidirectional and non-mutable. The edges and facets of the triangulation can also be visited through iterators and circulators which are bidirectional and non-mutable.

Iterators and circulators are convertible to the corresponding handles, thus the user can pass them directly as arguments to the functions.

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
 

Geometric Iterators:

The following iterators have value type Periodic_segment, Periodic_triangle, and Periodic_tetrahedron, which have inner type Point.

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.
 

Enums

enum  Locate_type {
  VERTEX =0, EDGE, FACET, CELL,
  EMPTY
}
 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 { STORED =0, UNIQUE, STORED_COVER_DOMAIN, UNIQUE_COVER_DOMAIN }
 The enum Iterator_type is defined by Periodic_3_triangulation_3 to specify the behavior of geometric iterators. More...
 

Creation

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

Assignment

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.
 
template<class Traits , class TDS1 , class TDS2 >
bool operator== (const Periodic_3_triangulation_3< Traits, TDS1 > &t1, const Periodic_3_triangulation_3< Traits, TDS2 > &t2)
 Equality operator. More...
 
template<class Traits , class TDS1 , class TDS2 >
bool operator!= (const Periodic_3_triangulation_3< Traits, TDS1 > &t1, const Periodic_3_triangulation_3< Traits, TDS2 > &t2)
 The opposite of operator==.
 

Access Functions

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.
 

Non const access

void set_domain (const Iso_cuboid dom)
 This is an advanced function. More...
 
Triangulation_data_structuretds ()
 This is an advanced function. More...
 

Non-constant-time queries and conversions

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

Constant-time access functions

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

Non-constant-time access functions

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

Geometric Access Functions

The following functions return object of types Periodic_segment, Periodic_triangle, and Periodic_tetrahedron, which have inner type Point.

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

Warning
The following functions were renamed with CGAL 4.11 to clarify that they return geometric objects with inner type Point_3.

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

template<typename PP >
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.
 
template<typename P >
Point_3 construct_point (const P &p1, const Offset &o1) const
 Same as above, with offsets.
 
template<typename PS >
Segment construct_segment (const PS &s) const
 Converts the periodic segment of type PS to a Segment. More...
 
template<typename P >
Segment construct_segment (const P &p1, const P &p2) const
 Creates a segment from two points. More...
 
template<typename P >
Segment construct_segment (const P &p1, const P &p2, const Offset &o1, const Offset &o2) const
 Same as above, with offsets.
 
template<typename PT >
Triangle construct_triangle (const PT &t) const
 Converts the periodic triangle of type PT to a Triangle. More...
 
template<typename P >
Triangle construct_triangle (const P &p1, const P &p2, const P &p3) const
 Creates a triangle from three points. More...
 
template<typename P >
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.
 
template<typename PT >
Tetrahedron construct_tetrahedron (const PT &t) const
 Converts the periodic tetrahedron of type PT to a Tetrahedron. More...
 
template<typename P >
Tetrahedron construct_tetrahedron (const P &p1, const P &p2, const P &p3, const P &p4) const
 Creates a tetrahedron from four points. More...
 
template<typename P >
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.
 

Queries

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

There is a method has_vertex in the cell class.

The analogous methods for facets are defined here.

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.
 

The following three methods test whether two facets have the same vertices.

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
 

Point Location

The class Periodic_3_triangulation_3 provides three functions to locate a given point with respect to a triangulation.

It provides also functions to test if a given point is inside a face or not. Note that the class Periodic_3_Delaunay_triangulation_3 also provides a nearest_vertex() function.

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 inexact_locate (const Point &query, Cell_handle start=Cell_handle()) const
 Same as locate() but uses inexact predicates. 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...
 
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...
 

Cell, Face, Edge and Vertex Iterators

The following iterators allow the user to visit cells, facets, edges and vertices of the stored triangulation, i.e. in case of computing in a multiply sheeted covering space all stored periodic copies of each item are returned.

These iterators are non-mutable, bidirectional and their value types are respectively Cell, Facet, Edge and Vertex. They are all invalidated by any change in the triangulation.

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.
 

Geometric Iterators

The following iterators allow the user to obtain geometric primitives corresponding to cells, facets, edges, and vertices of the triangulation.

These iterators are non-mutable, bidirectional and their value types are respectively Periodic_point, Periodic_segment, Periodic_triangle, and Periodic_tetrahedron. They are all invalidated by any change in the triangulation. If the periodic triangulation is not computed in the 1-sheeted covering space, these iterators can be used to retain only the geometric primitives in the original domain. This can be controlled using the enum Iterator_type.

STORED
STORED_COVER_DOMAIN
UNIQUE
UNIQUE_COVER_DOMAIN
The four different modes of the geometric iterators: STORED, STORED_COVER_DOMAIN, UNIQUE, UNIQUE_COVER_DOMAIN. Note that in case of computing in the 1-sheeted covering space, STORED and UNIQUE give the same result.
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 and Facet Circulators

The following circulators respectively visit all cells or all facets incident to a given edge.

They are non-mutable and bidirectional. They are invalidated by any modification of one of the cells traversed.

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

Traversal of the Incident Cells and Facets, and the Adjacent Vertices of a Given Vertex

template<class OutputIterator >
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...
 
template<class OutputIterator >
OutputIterator incident_facets (Vertex_handle v, OutputIterator facets) const
 Copies the Facets incident to v to the output iterator facets. More...
 
template<class OutputIterator >
OutputIterator incident_edges (Vertex_handle v, OutputIterator edges) const
 Copies the Edges incident to v to the output iterator edges. More...
 
template<class OutputIterator >
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...
 

Traversal Between Adjacent Cells

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.
 

Checking

Advanced

The responsibility of keeping a valid triangulation belongs to the user when using advanced operations allowing a direct manipulation of cells and vertices.

We provide the user with the following methods to help debugging.

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

Member Typedef Documentation

◆ Iso_cuboid

template<typename Traits, typename TDS>
typedef Geom_traits::Iso_cuboid_3 CGAL::Periodic_3_triangulation_3< Traits, TDS >::Iso_cuboid

A type representing an axis-aligned cuboid.

It must be a model of Traits::Iso_cuboid_3. Used to represent the original domain.

◆ Periodic_point

template<typename Traits, typename TDS>
typedef std::pair< Point, Offset > CGAL::Periodic_3_triangulation_3< Traits, TDS >::Periodic_point

Represents a point-offset pair.

The point in the pair lies in the original domain. Note that the inner type is Point.

◆ Periodic_point_3

template<typename Traits, typename TDS>
typedef std::pair< Point_3, Offset > CGAL::Periodic_3_triangulation_3< Traits, TDS >::Periodic_point_3

Represents a point-offset pair.

The point in the pair lies in the original domain. Note that the inner type is Point_3.

◆ Periodic_segment

template<typename Traits, typename TDS>
typedef array< Periodic_point, 2> CGAL::Periodic_3_triangulation_3< Traits, TDS >::Periodic_segment

A periodic segment.

Note that the inner type is Periodic_point.

◆ Periodic_segment_3

template<typename Traits, typename TDS>
typedef array< Periodic_point_3, 2> CGAL::Periodic_3_triangulation_3< Traits, TDS >::Periodic_segment_3

A periodic segment.

Note that the inner type is Periodic_point_3.

◆ Periodic_tetrahedron

template<typename Traits, typename TDS>
typedef array< Periodic_point, 4> CGAL::Periodic_3_triangulation_3< Traits, TDS >::Periodic_tetrahedron

A periodic tetrahedron.

Note that the inner type is Periodic_point.

◆ Periodic_tetrahedron_3

template<typename Traits, typename TDS>
typedef array< Periodic_point_3, 4> CGAL::Periodic_3_triangulation_3< Traits, TDS >::Periodic_tetrahedron_3

A periodic tetrahedron.

Note that the inner type is Periodic_point_3.

◆ Periodic_triangle

template<typename Traits, typename TDS>
typedef array< Periodic_point, 3> CGAL::Periodic_3_triangulation_3< Traits, TDS >::Periodic_triangle

A periodic triangle.

Note that the inner type is Periodic_point.

◆ Periodic_triangle_3

template<typename Traits, typename TDS>
typedef array< Periodic_point_3, 3> CGAL::Periodic_3_triangulation_3< Traits, TDS >::Periodic_triangle_3

A periodic triangle.

Note that the inner type is Periodic_point_3.

◆ Point

template<typename Traits, typename TDS>
typedef TDS::Vertex::Point CGAL::Periodic_3_triangulation_3< Traits, TDS >::Point

The point type of the triangulation.

Note
This type is equal to Geom_traits::Point_3 when considering periodic Delaunay triangulations and to Geom_traits::Weighted_point_3 when considering periodic weighted Delaunay triangulations.

Member Enumeration Documentation

◆ Iterator_type

template<typename Traits, typename TDS>
enum CGAL::Periodic_3_triangulation_3::Iterator_type

The enum Iterator_type is defined by Periodic_3_triangulation_3 to specify the behavior of geometric iterators.

The elements of the enum have the following meaning:

See also
CGAL::Periodic_3_triangulation_3
Enumerator
STORED 

Return all geometric primitives as they are stored internally in Triangulation_data_structure_3.

UNIQUE 

Return only one representative of each geometric primitive even if the triangulation is computed in a multiply sheeted covering space.

Choose the representative whose maximum offset is minimal but non-negative in each direction of space.

STORED_COVER_DOMAIN 

Same as STORED but return additionally all primitives whose intersection with the original domain of the current covering space is non-empty.

UNIQUE_COVER_DOMAIN 

Same as UNIQUE but return additionally all primitives whose intersection with the original domain is non-empty.

◆ Locate_type

template<typename Traits, typename TDS>
enum CGAL::Periodic_3_triangulation_3::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.

If the triangulation does not contain any points EMPTY is returned.

See also
CGAL::Periodic_3_triangulation_3
Enumerator
VERTEX 
EDGE 
FACET 
CELL 
EMPTY 

Constructor & Destructor Documentation

◆ Periodic_3_triangulation_3() [1/2]

template<typename Traits, typename TDS>
CGAL::Periodic_3_triangulation_3< Traits, 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.

Precondition
domain is a cube.

◆ Periodic_3_triangulation_3() [2/2]

template<typename Traits, typename TDS>
CGAL::Periodic_3_triangulation_3< Traits, TDS >::Periodic_3_triangulation_3 ( const Periodic_3_triangulation_3< Traits, TDS > &  tr)

Copy constructor.

All vertices and faces are duplicated.

Member Function Documentation

◆ adjacent_vertices()

template<typename Traits, typename TDS>
template<class OutputIterator >
OutputIterator CGAL::Periodic_3_triangulation_3< Traits, TDS >::adjacent_vertices ( Vertex_handle  v,
OutputIterator  vertices 
) const

Copies the Vertex_handles of all vertices adjacent to v to the output iterator vertices.

Returns the resulting output iterator.

Precondition
v \( \neq\) Vertex_handle(), t.is_vertex(v).

◆ cells_begin()

template<typename Traits, typename TDS>
Cell_iterator CGAL::Periodic_3_triangulation_3< Traits, TDS >::cells_begin ( ) const

Starts at an arbitrary cell.

Iterates over all cells. Returns cells_end() if t.number_of_vertices() \( =0\).

◆ construct_point()

template<typename Traits, typename TDS>
template<typename PP >
Point_3 CGAL::Periodic_3_triangulation_3< Traits, TDS >::construct_point ( const PP &  pp) const

Converts the periodic point of type PP to a Point_3.

The type PP can be either Periodic_point or Periodic_point_3.

◆ construct_segment() [1/2]

template<typename Traits, typename TDS>
template<typename PS >
Segment CGAL::Periodic_3_triangulation_3< Traits, TDS >::construct_segment ( const PS &  s) const

Converts the periodic segment of type PS to a Segment.

The type PS can be either Periodic_segment or Periodic_segment_3.

◆ construct_segment() [2/2]

template<typename Traits, typename TDS>
template<typename P >
Segment CGAL::Periodic_3_triangulation_3< Traits, TDS >::construct_segment ( const P &  p1,
const P &  p2 
) const

Creates a segment from two points.

The type P can be either Point or Point_3.

◆ construct_tetrahedron() [1/2]

template<typename Traits, typename TDS>
template<typename PT >
Tetrahedron CGAL::Periodic_3_triangulation_3< Traits, TDS >::construct_tetrahedron ( const PT &  t) const

Converts the periodic tetrahedron of type PT to a Tetrahedron.

The type PT can be either Periodic_tetrahedron or Periodic_tetrahedron_3.

◆ construct_tetrahedron() [2/2]

template<typename Traits, typename TDS>
template<typename P >
Tetrahedron CGAL::Periodic_3_triangulation_3< Traits, TDS >::construct_tetrahedron ( const P &  p1,
const P &  p2,
const P &  p3,
const P &  p4 
) const

Creates a tetrahedron from four points.

The type P can be either Point or Point_3.

◆ construct_triangle() [1/2]

template<typename Traits, typename TDS>
template<typename PT >
Triangle CGAL::Periodic_3_triangulation_3< Traits, TDS >::construct_triangle ( const PT &  t) const

Converts the periodic triangle of type PT to a Triangle.

The type PT can be either Periodic_triangle or Periodic_triangle_3.

◆ construct_triangle() [2/2]

template<typename Traits, typename TDS>
template<typename P >
Triangle CGAL::Periodic_3_triangulation_3< Traits, TDS >::construct_triangle ( const P &  p1,
const P &  p2,
const P &  p3 
) const

Creates a triangle from three points.

The type P can be either Point or Point_3.

◆ convert_to_1_sheeted_covering()

template<typename Traits, typename TDS>
void CGAL::Periodic_3_triangulation_3< Traits, TDS >::convert_to_1_sheeted_covering ( ) const

Converts the current triangulation into the same periodic triangulation in the 1-sheeted covering space.

Attention
It is not recommended to interfere with the built-in covering management. Especially a premature conversion to the 1-sheeted covering space might lead to problems when modifying the triangulation later.

◆ convert_to_27_sheeted_covering()

template<typename Traits, typename TDS>
void CGAL::Periodic_3_triangulation_3< Traits, TDS >::convert_to_27_sheeted_covering ( ) const

Converts the current triangulation into the same periodic triangulation in the 27-sheeted covering space.

Attention
It is not recommended to interfere with the built-in covering management. Especially a premature conversion to the 1-sheeted covering space might lead to problems when modifying the triangulation later.

◆ degree()

template<typename Traits, typename TDS>
size_type CGAL::Periodic_3_triangulation_3< Traits, TDS >::degree ( Vertex_handle  v) const

Returns the degree of a vertex, that is, the number of adjacent vertices.

Precondition
v \( \neq\) Vertex_handle(), t.is_vertex(v).

◆ edges_begin()

template<typename Traits, typename TDS>
Edge_iterator CGAL::Periodic_3_triangulation_3< Traits, TDS >::edges_begin ( ) const

Starts at an arbitrary edge.

Iterates over all edges. Returns edges_end() if t.number_of_vertices() \( =0\).

◆ facets_begin()

template<typename Traits, typename TDS>
Facet_iterator CGAL::Periodic_3_triangulation_3< Traits, TDS >::facets_begin ( ) const

Starts at an arbitrary facet.

Iterates over all facets. Returns facets_end() if t.number_of_vertices() \( =0\).

◆ has_vertex()

template<typename Traits, typename TDS>
bool CGAL::Periodic_3_triangulation_3< Traits, TDS >::has_vertex ( Cell_handle  c,
int  i,
Vertex_handle  v,
int &  j 
) const

Same for facet (c,i).

Computes the index j of v in c.

◆ incident_cells() [1/2]

template<typename Traits, typename TDS>
Cell_circulator CGAL::Periodic_3_triangulation_3< Traits, TDS >::incident_cells ( Edge  e,
Cell_handle  start 
) const

Starts at cell start.

Precondition
start is incident to e.

◆ incident_cells() [2/2]

template<typename Traits, typename TDS>
template<class OutputIterator >
OutputIterator CGAL::Periodic_3_triangulation_3< Traits, TDS >::incident_cells ( Vertex_handle  v,
OutputIterator  cells 
) const

Copies the Cell_handles of all cells incident to v to the output iterator cells.

Returns the resulting output iterator.

Precondition
v \( \neq\) Vertex_handle(), t.is_vertex(v).

◆ incident_edges()

template<typename Traits, typename TDS>
template<class OutputIterator >
OutputIterator CGAL::Periodic_3_triangulation_3< Traits, TDS >::incident_edges ( Vertex_handle  v,
OutputIterator  edges 
) const

Copies the Edges incident to v to the output iterator edges.

Returns the resulting output iterator.

Precondition
v \( \neq\) Vertex_handle(), t.is_vertex(v).

◆ incident_facets() [1/2]

template<typename Traits, typename TDS>
Facet_circulator CGAL::Periodic_3_triangulation_3< Traits, TDS >::incident_facets ( Edge  e,
Facet  start 
) const

Starts at facet start.

Precondition
start is incident to e.

◆ incident_facets() [2/2]

template<typename Traits, typename TDS>
template<class OutputIterator >
OutputIterator CGAL::Periodic_3_triangulation_3< Traits, TDS >::incident_facets ( Vertex_handle  v,
OutputIterator  facets 
) const

Copies the Facets incident to v to the output iterator facets.

Returns the resulting output iterator.

Precondition
v \( \neq\) Vertex_handle(), t.is_vertex(v).

◆ inexact_locate()

template<typename Traits, typename TDS>
Cell_handle CGAL::Periodic_3_triangulation_3< Traits, TDS >::inexact_locate ( const Point query,
Cell_handle  start = Cell_handle() 
) const

Same as locate() but uses inexact predicates.

This function returns a handle on a cell that is a good approximation of the exact location of query, while being faster. Note that it may return a handle on a cell whose interior does not contain query.

Note that this function is available only if the Cartesian coordinates of query are accessible with functions x(), y() and z().

◆ is_cell() [1/4]

template<typename Traits, typename TDS>
bool CGAL::Periodic_3_triangulation_3< Traits, TDS >::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.

If the cell c is found, the method computes the indices i, j, k and l of the vertices u, v, w and x in c, in this order.

Precondition
u, v, w and x are vertices of t.

◆ is_cell() [2/4]

template<typename Traits, typename TDS>
bool CGAL::Periodic_3_triangulation_3< Traits, TDS >::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.

Precondition
u, v, w and x are vertices of t.

◆ is_cell() [3/4]

template<typename Traits, typename TDS>
bool CGAL::Periodic_3_triangulation_3< Traits, TDS >::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.

If the cell c is found, the method computes the indices i, j, k and l of the vertices u, v, w and x in c, in this order.

Precondition
u, v, w and x are vertices of t.

◆ is_cell() [4/4]

template<typename Traits, typename TDS>
bool CGAL::Periodic_3_triangulation_3< Traits, TDS >::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.

Precondition
u, v, w and x are vertices of t.

◆ is_edge() [1/2]

template<typename Traits, typename TDS>
bool CGAL::Periodic_3_triangulation_3< Traits, TDS >::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.

If the edge is found, it gives a cell c having this edge and the indices i and j of the vertices u and v in c, in this order.

Precondition
u and v are vertices of t.

◆ is_edge() [2/2]

template<typename Traits, typename TDS>
bool CGAL::Periodic_3_triangulation_3< Traits, TDS >::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.

If the edge is found, it gives a cell c having this edge and the indices i and j of the vertices u and v in c, in this order.

Precondition
u and v are vertices of t.

◆ is_extensible_triangulation_in_1_sheet_h1()

template<typename Traits, typename TDS>
bool CGAL::Periodic_3_triangulation_3< Traits, TDS >::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.

This test relies on a heuristic, i.e. if it answers false nothing is known. This test runs in constant-time when not computing in the 1-sheeted covering space. (This test uses the length of the longest edge in the triangulation as a criterion [1].)

◆ is_extensible_triangulation_in_1_sheet_h2()

template<typename Traits, typename TDS>
bool CGAL::Periodic_3_triangulation_3< Traits, TDS >::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.

However, it is much less time efficient when not computing in the 1-sheeted covering space. (This test uses the diameter of the largest empty ball in the input point set as a criterion [1].)

◆ is_facet() [1/2]

template<typename Traits, typename TDS>
bool CGAL::Periodic_3_triangulation_3< Traits, TDS >::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.

If the facet is found, it computes a cell c having this facet and the indices i, j and k of the vertices u, v and w in c, in this order.

Precondition
u, v and w are vertices of t.

◆ is_facet() [2/2]

template<typename Traits, typename TDS>
bool CGAL::Periodic_3_triangulation_3< Traits, TDS >::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.

If the facet is found, it computes a cell c having this facet and the indices i, j and k of the vertices u, v and w in c, in this order.

Precondition
u, v and w are vertices of t.

◆ is_valid() [1/2]

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

Checks the combinatorial validity of the triangulation.

Checks also the validity of its geometric embedding (see Section Representation). When verbose is set to true, messages describing the first invalidity encountered are printed.

◆ is_valid() [2/2]

template<typename Traits, typename TDS>
bool CGAL::Periodic_3_triangulation_3< Traits, TDS >::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.

Also checks the geometric validity of c, if c is finite. (See Section Representation.)

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

◆ is_vertex()

template<typename Traits, typename TDS>
bool CGAL::Periodic_3_triangulation_3< Traits, TDS >::is_vertex ( const Point p,
Vertex_handle v 
) const

Tests whether p is a vertex of t by locating p in the triangulation.

If p is found, the associated vertex v is given.

Precondition
p lies in the original domain domain.

◆ locate() [1/4]

template<typename Traits, typename TDS>
Cell_handle CGAL::Periodic_3_triangulation_3< Traits, TDS >::locate ( const Point query,
Cell_handle  start = Cell_handle() 
) const

Returns the cell that contains the query in its interior.

If query lies on a facet, an edge or on a vertex, one of the cells having query on its boundary is returned.

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

Precondition
query lies in the original domain domain.

◆ locate() [2/4]

template<typename Traits, typename TDS>
Cell_handle CGAL::Periodic_3_triangulation_3< Traits, TDS >::locate ( const Point query,
Offset locate_offset,
Cell_handle  start = Cell_handle() 
) const

Returns the cell that contains the query in its interior.

If query lies on a facet, an edge or on a vertex, one of the cells having query on its boundary is returned.

locate_offset is the offset that must be used by the function periodic_tetrahedron() together with the returned cell, so that the constructed Periodic_tetrahedron contains the query point.

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

Precondition
query lies in the original domain domain.

◆ locate() [3/4]

template<typename Traits, typename TDS>
Cell_handle CGAL::Periodic_3_triangulation_3< Traits, TDS >::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.

If the \( k\)-face is a cell, li and lj have no meaning; if it is a facet (resp. vertex), li gives the index of the facet (resp. vertex) and lj has no meaning; if it is an edge, li and lj give the indices of its vertices.

If there is no vertex in the triangulation yet, lt is set to EMPTY and locate returns the default constructed handle.

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

Precondition
query lies in the original domain domain.

◆ locate() [4/4]

template<typename Traits, typename TDS>
Cell_handle CGAL::Periodic_3_triangulation_3< Traits, TDS >::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.

If the \( k\)-face is a cell, li and lj have no meaning; if it is a facet (resp. vertex), li gives the index of the facet (resp. vertex) and lj has no meaning; if it is an edge, li and lj give the indices of its vertices.

If there is no vertex in the triangulation yet, lt is set to EMPTY and locate returns the default constructed handle.

locate_offset is the offset that must be used by the function periodic_tetrahedron() together with the returned cell, so that the constructed Periodic_tetrahedron contains the query point.

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

Precondition
query lies in the original domain domain.

◆ mirror_index()

template<typename Traits, typename TDS>
int CGAL::Periodic_3_triangulation_3< Traits, TDS >::mirror_index ( Cell_handle  c,
int  i 
) const

Returns the index of c in its \( i^{th}\) neighbor.

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

◆ mirror_vertex()

template<typename Traits, typename TDS>
Vertex_handle CGAL::Periodic_3_triangulation_3< Traits, TDS >::mirror_vertex ( Cell_handle  c,
int  i 
) const

Returns the vertex of the \( i^{th}\) neighbor of c that is opposite to c.

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

◆ number_of_cells()

template<typename Traits, typename TDS>
size_type CGAL::Periodic_3_triangulation_3< Traits, TDS >::number_of_cells ( ) const

Returns the number of cells.

Counts all cells that are representatives of the same tetrahedron in \( \mathbb T_c^3\) as one cell.

◆ number_of_edges()

template<typename Traits, typename TDS>
size_type CGAL::Periodic_3_triangulation_3< Traits, TDS >::number_of_edges ( ) const

Returns the number of edges.

Counts all edges that are representatives of the same segment in \( \mathbb T_c^3\) as one edge.

◆ number_of_facets()

template<typename Traits, typename TDS>
size_type CGAL::Periodic_3_triangulation_3< Traits, TDS >::number_of_facets ( ) const

Returns the number of facets.

Counts all facets that are representatives of the same triangle in \( \mathbb T_c^3\) as one facet.

◆ number_of_sheets()

template<typename Traits, typename TDS>
Covering_sheets CGAL::Periodic_3_triangulation_3< Traits, TDS >::number_of_sheets ( ) const

This is an advanced function.

Advanced

Returns the number of sheets of the covering space the triangulation is currently computed in.

◆ number_of_stored_cells()

template<typename Traits, typename TDS>
size_type CGAL::Periodic_3_triangulation_3< Traits, TDS >::number_of_stored_cells ( ) const

Returns the number of cells in the data structure.

This is the same as the number of sheets times number_of_cells().

◆ number_of_stored_edges()

template<typename Traits, typename TDS>
size_type CGAL::Periodic_3_triangulation_3< Traits, TDS >::number_of_stored_edges ( ) const

Returns the number of edges in the data structure.

This is the same as the number of sheets times number_of_edges().

◆ number_of_stored_facets()

template<typename Traits, typename TDS>
size_type CGAL::Periodic_3_triangulation_3< Traits, TDS >::number_of_stored_facets ( ) const

Returns the number of facets in the data structure.

This is the same as the number of sheets times number_of_facets().

◆ number_of_stored_vertices()

template<typename Traits, typename TDS>
size_type CGAL::Periodic_3_triangulation_3< Traits, TDS >::number_of_stored_vertices ( ) const

Returns the number of vertices in the data structure.

This is the same as the number of sheets times number_of_vertices().

◆ number_of_vertices()

template<typename Traits, typename TDS>
size_type CGAL::Periodic_3_triangulation_3< Traits, TDS >::number_of_vertices ( ) const

Returns the number of vertices.

Counts all vertices that are representatives of the same point in \( \mathbb T_c^3\) as one vertex.

◆ operator=()

template<typename Traits, typename TDS>
Periodic_3_triangulation_3& CGAL::Periodic_3_triangulation_3< Traits, TDS >::operator= ( const Periodic_3_triangulation_3< Traits, TDS > &  tr)

The triangulation tr is duplicated, and modifying the copy after the duplication does not modify the original.

The previous triangulation held by t is deleted.

◆ periodic_point() [1/2]

template<typename Traits, typename TDS>
Periodic_point CGAL::Periodic_3_triangulation_3< Traits, TDS >::periodic_point ( const Vertex_handle  v) const

Returns the periodic point given by vertex v.

If t is represented in the 1-sheeted covering space, the offset is always zero. Otherwise v can correspond to a periodic copy outside domain of an input point.

◆ periodic_point() [2/2]

template<typename Traits, typename TDS>
Periodic_point CGAL::Periodic_3_triangulation_3< Traits, TDS >::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.

If t is represented in the 27-sheeted covering space, this offset is possibly added to another offset determining the periodic copy.

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

◆ periodic_points_begin()

template<typename Traits, typename TDS>
Periodic_point_iterator CGAL::Periodic_3_triangulation_3< Traits, TDS >::periodic_points_begin ( Iterator_type  it = STORED) const

Iterates over the points of the triangulation.

Its behavior is defined by the Iterator_type it as described on CGAL::Periodic_3_triangulation_3::Iterator_type.

◆ periodic_points_end()

template<typename Traits, typename TDS>
Periodic_point_iterator CGAL::Periodic_3_triangulation_3< Traits, TDS >::periodic_points_end ( Iterator_type  it = STORED) const

Past-the-end iterator.

Note that to match another Periodic_point_iterator both must have the same Iterator_type it.

◆ periodic_segment() [1/2]

template<typename Traits, typename TDS>
Periodic_segment CGAL::Periodic_3_triangulation_3< Traits, TDS >::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).

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

◆ periodic_segment() [2/2]

template<typename Traits, typename TDS>
Periodic_segment CGAL::Periodic_3_triangulation_3< Traits, TDS >::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).

A translation in accordance with offset is applied on the point-offet pairs.

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

◆ periodic_segments_begin()

template<typename Traits, typename TDS>
Periodic_segment_iterator CGAL::Periodic_3_triangulation_3< Traits, TDS >::periodic_segments_begin ( Iterator_type  it = STORED) const

Iterates over the segments of the triangulation.

Its behavior is defined by the Iterator_type it as described on CGAL::Periodic_3_triangulation_3::Iterator_type.

◆ periodic_segments_end()

template<typename Traits, typename TDS>
Periodic_segment_iterator CGAL::Periodic_3_triangulation_3< Traits, TDS >::periodic_segments_end ( Iterator_type  it = STORED) const

Past-the-end iterator.

Note that to match another Periodic_segment_iterator both must have the same Iterator_type it.

◆ periodic_tetrahedra_begin()

template<typename Traits, typename TDS>
Periodic_tetrahedron_iterator CGAL::Periodic_3_triangulation_3< Traits, TDS >::periodic_tetrahedra_begin ( Iterator_type  it = STORED) const

Iterates over the tetrahedra of the triangulation.

Its behavior is defined by the Iterator_type it as described on CGAL::Periodic_3_triangulation_3::Iterator_type.

◆ periodic_tetrahedra_end()

template<typename Traits, typename TDS>
Periodic_tetrahedron_iterator CGAL::Periodic_3_triangulation_3< Traits, TDS >::periodic_tetrahedra_end ( Iterator_type  it = STORED) const

Past-the-end iterator.

Note that to match another Periodic_tetrahedron_iterator both must have the same Iterator_type it.

◆ periodic_tetrahedron()

template<typename Traits, typename TDS>
Periodic_tetrahedron CGAL::Periodic_3_triangulation_3< Traits, TDS >::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.

A translation in accordance with offset is applied on the point-offet pairs.

◆ periodic_triangle() [1/2]

template<typename Traits, typename TDS>
Periodic_triangle CGAL::Periodic_3_triangulation_3< Traits, TDS >::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).

The triangle is oriented so that its normal points to the inside of cell c.

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

◆ periodic_triangle() [2/2]

template<typename Traits, typename TDS>
Periodic_triangle CGAL::Periodic_3_triangulation_3< Traits, TDS >::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).

A translation in accordance with offset is applied on the point-offet pairs.

The triangle is oriented so that its normal points to the inside of cell c.

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

◆ periodic_triangles_begin()

template<typename Traits, typename TDS>
Periodic_triangle_iterator CGAL::Periodic_3_triangulation_3< Traits, TDS >::periodic_triangles_begin ( Iterator_type  it = STORED) const

Iterates over the triangles of the triangulation.

Its behavior is defined by the Iterator_type it as described on CGAL::Periodic_3_triangulation_3::Iterator_type.

◆ periodic_triangles_end()

template<typename Traits, typename TDS>
Periodic_triangle_iterator CGAL::Periodic_3_triangulation_3< Traits, TDS >::periodic_triangles_end ( Iterator_type  it = STORED) const

Past-the-end iterator.

Note that to match another Periodic_triangle_iterator both must have the same Iterator_type it.

◆ set_domain()

template<typename Traits, typename TDS>
void CGAL::Periodic_3_triangulation_3< Traits, TDS >::set_domain ( const Iso_cuboid  dom)

This is an advanced function.

Advanced

Changes the domain. Note that this function calls clear(), i.e., it erases the existing triangulation.

◆ side_of_cell()

template<typename Traits, typename TDS>
Bounded_side CGAL::Periodic_3_triangulation_3< Traits, TDS >::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 precisely, it returns:

  • ON_BOUNDED_SIDE if p is inside the cell.
  • ON_BOUNDARY if p on the boundary of the cell. Then lt together with li and lj give the precise location on the boundary. (See the descriptions of the locate methods.)
  • ON_UNBOUNDED_SIDE if p lies outside the cell.
    Precondition
    query lies in the original domain domain.

◆ swap()

template<typename Traits, typename TDS>
void CGAL::Periodic_3_triangulation_3< Traits, TDS >::swap ( Periodic_3_triangulation_3< Traits, TDS > &  tr)

The triangulations tr and t are swapped.

t.swap(tr) should be preferred to t = tr or to t(tr) if tr is deleted after that. Indeed, there is no copy of cells and vertices, thus this method runs in constant time.

◆ tds()

template<typename Traits, typename TDS>
Triangulation_data_structure& CGAL::Periodic_3_triangulation_3< Traits, TDS >::tds ( )

This is an advanced function.

Advanced

Returns a reference to the triangulation data structure.

The responsibility of keeping a valid triangulation belongs to the user when using advanced operations allowing a direct manipulation of the tds. This method is mainly a help for users implementing their own triangulation algorithms.

◆ unique_vertices_begin()

template<typename Traits, typename TDS>
Unique_vertex_iterator CGAL::Periodic_3_triangulation_3< Traits, TDS >::unique_vertices_begin ( ) const

Starts at an arbitrary vertex.

Iterates over all 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. Returns unique_vertices_end() if t.number_of_vertices() \( =0\).

◆ vertices_begin()

template<typename Traits, typename TDS>
Vertex_iterator CGAL::Periodic_3_triangulation_3< Traits, TDS >::vertices_begin ( ) const

Starts at an arbitrary vertex.

Iterates over all vertices. Returns vertices_end() if t.number_of_vertices() \( =0\).

Friends And Related Function Documentation

◆ operator<<()

template<typename Traits, typename TDS>
ostream & operator<< ( ostream &  os,
const Periodic_3_triangulation_3< Traits, TDS > &  t 
)
related

Writes the triangulation t into os.

The information in the iostream is:

  • the original domain
  • the number of sheets of the covering space as in number_of_sheets()
  • the number of vertices
  • the non-combinatorial information of vertices (point resp. point-offset pairs, etc.)
  • the number of cells
  • the indices of the vertices of each cell
  • the indices of the neighbors of each cell, where the index corresponds to the preceding list of cells
  • the offsets corresponding to the vertices of the cells
  • the non-combinatorial information of each cell

◆ operator==()

template<typename Traits, typename TDS>
template<class Traits , class TDS1 , class TDS2 >
bool operator== ( const Periodic_3_triangulation_3< Traits, TDS1 > &  t1,
const Periodic_3_triangulation_3< Traits, TDS2 > &  t2 
)
related

Equality operator.

Returns true iff there exist a bijection between the vertices of t1 and those of t2 and a bijection between the cells of t1 and those of t2, which preserve the geometry of the triangulation, that is, the points of each corresponding pair of vertices are equal, and the tetrahedra corresponding to each pair of cells are equal (up to a permutation of their vertices).

◆ operator>>()

template<typename Traits, typename TDS>
istream & operator>> ( istream &  is,
Periodic_3_triangulation_3< Traits, TDS > &  t 
)
related

Reads a triangulation from is and stores it in t.

Precondition
is has the below described format.

The information in the iostream is:

  • the original domain
  • the number of sheets of the covering space as in number_of_sheets()
  • the number of vertices
  • the non-combinatorial information of vertices (point resp. point-offset pairs, etc.)
  • the number of cells
  • the indices of the vertices of each cell
  • the indices of the neighbors of each cell, where the index corresponds to the preceding list of cells
  • the offsets corresponding to the vertices of the cells
  • the non-combinatorial information of each cell