CDT++
Causal Dynamical Triangulations in C++
Loading...
Searching...
No Matches
Classes | Concepts | Typedefs | Enumerations | Functions | Variables
Foliated_triangulation.hpp File Reference

Create foliated spherical triangulations. More...

#include "Triangulation_traits.hpp"
#include "Utilities.hpp"
+ Include dependency graph for Foliated_triangulation.hpp:
+ This graph shows which files directly or indirectly include this file:

Go to the source code of this file.

Classes

class  foliated_triangulations::FoliatedTriangulation< 3 >
 3D Foliated triangulation More...
 
class  foliated_triangulations::FoliatedTriangulation< 4 >
 4D Triangulation More...
 

Concepts

concept  ContainerType
 This is std::movable from <concepts>
 

Typedefs

template<int dimension>
using Causal_vertices_t = std::vector< std::pair< Point_t< dimension >, Int_precision > >
 
template<int dimension>
using Cell_handle_t = typename TriangulationTraits< dimension >::Cell_handle
 
template<int dimension>
using Delaunay_t = typename TriangulationTraits< dimension >::Delaunay
 
template<int dimension>
using Edge_handle_t = typename TriangulationTraits< dimension >::Edge_handle
 
template<int dimension>
using Face_handle_t = typename TriangulationTraits< dimension >::Face_handle
 
template<int dimension>
using Facet_t = typename TriangulationTraits< dimension >::Facet
 
using foliated_triangulations::FoliatedTriangulation_3 = FoliatedTriangulation< 3 >
 
using foliated_triangulations::FoliatedTriangulation_4 = FoliatedTriangulation< 4 >
 
template<int dimension>
using Point_t = typename TriangulationTraits< dimension >::Point
 
template<int dimension>
using Spherical_points_generator_t = typename TriangulationTraits< dimension >::Spherical_points_generator
 
template<int dimension>
using Vertex_handle_t = typename TriangulationTraits< dimension >::Vertex_handle
 

Enumerations

enum class  Cell_type {
  THREE_ONE = 31 , TWO_TWO = 22 , ONE_THREE = 13 , FOUR_ONE = 41 ,
  THREE_TWO = 32 , TWO_THREE = 23 , ONE_FOUR = 14 , ACAUSAL = 99 ,
  UNCLASSIFIED = 0
}
 (n,m) is number of vertices on (lower, higher) timeslice More...
 

Functions

template<int dimension>
auto foliated_triangulations::check_cells (Delaunay_t< dimension > const &t_triangulation) -> bool
 Check all finite cells in the Delaunay triangulation.
 
template<int dimension>
auto foliated_triangulations::check_timevalues (Delaunay_t< dimension > const &t_triangulation) -> std::optional< std::vector< Cell_handle_t< dimension > > >
 Check cells for correct foliation.
 
template<int dimension>
auto foliated_triangulations::check_vertices (Delaunay_t< dimension > const &t_triangulation, double t_initial_radius, double t_foliation_spacing)
 Check if vertices have the correct timevalues.
 
template<int dimension>
auto foliated_triangulations::classify_edge (Edge_handle_t< dimension > const &t_edge) -> bool
 Predicate to classify edge as timelike or spacelike.
 
template<int dimension>
auto foliated_triangulations::collect_cells (Delaunay_t< dimension > const &t_triangulation) -> std::vector< Cell_handle_t< dimension > >
 Obtain all finite cells in the Delaunay triangulation.
 
template<int dimension>
auto foliated_triangulations::collect_edges (Delaunay_t< dimension > const &delaunay)
 Returns a container of all the finite edges in the triangulation.
 
template<int dimension>
auto foliated_triangulations::collect_vertices (Delaunay_t< dimension > const &t_triangulation)
 Obtain all finite vertices in the Delaunay triangulation.
 
template<int dimension, ContainerType Container>
void foliated_triangulations::debug_print_cells (Container &&t_cells)
 Write to debug log timevalues of each vertex in the cell and the resulting cell->info.
 
template<int dimension>
auto foliated_triangulations::expected_cell_type (Cell_handle_t< dimension > const &t_cell)
 Classifies cells by their timevalues.
 
template<int dimension>
auto foliated_triangulations::expected_timevalue (Vertex_handle_t< dimension > const &t_vertex, double t_initial_radius, double t_foliation_spacing) -> Int_precision
 Find the expected timevalue for a vertex.
 
template<int dimension>
auto foliated_triangulations::filter_cells (std::vector< Cell_handle_t< dimension > > const &t_cells, Cell_type const &t_cell_type) -> std::vector< Cell_handle_t< dimension > >
 
template<int dimension>
auto foliated_triangulations::filter_edges (std::vector< Edge_handle_t< dimension > > const &t_edges, bool t_is_Timelike_pred) -> std::vector< Edge_handle_t< dimension > >
 
template<int dimension>
auto foliated_triangulations::find_bad_vertex (Cell_handle_t< dimension > const &cell) -> Vertex_handle_t< dimension >
 Find the vertex that is causing a cell's foliation to be invalid.
 
template<int dimension>
auto foliated_triangulations::find_cell (Delaunay_t< dimension > const &delaunay, Vertex_handle_t< dimension > const &vh1, Vertex_handle_t< dimension > const &vh2, Vertex_handle_t< dimension > const &vh3, Vertex_handle_t< dimension > const &vh4) -> std::optional< Cell_handle_t< dimension > >
 Returns the cell containing the vertices.
 
template<int dimension>
auto foliated_triangulations::find_incorrect_cells (Delaunay_t< dimension > const &t_triangulation)
 Check all finite cells in the Delaunay triangulation.
 
template<int dimension>
auto foliated_triangulations::find_incorrect_vertices (Delaunay_t< dimension > const &t_triangulation, double t_initial_radius, double t_foliation_spacing)
 Obtain vertices with incorrect timevalues.
 
template<int dimension>
auto foliated_triangulations::find_incorrect_vertices (std::vector< Cell_handle_t< dimension > > const &t_cells, double t_initial_radius, double t_foliation_spacing)
 Obtain vertices with incorrect timevalues.
 
template<int dimension, ContainerType Container>
auto foliated_triangulations::find_max_timevalue (Container &&t_vertices) -> Int_precision
 
template<int dimension, ContainerType Container>
auto foliated_triangulations::find_min_timevalue (Container &&t_vertices) -> Int_precision
 
template<int dimension>
auto foliated_triangulations::find_vertex (Delaunay_t< dimension > const &delaunay, Point_t< dimension > const &point) -> std::optional< Vertex_handle_t< dimension > >
 Returns the vertex containing the given point.
 
template<int dimension>
auto foliated_triangulations::fix_cells (Delaunay_t< dimension > const &t_triangulation) -> bool
 Fix simplices with the wrong type.
 
template<int dimension>
auto foliated_triangulations::fix_timevalues (Delaunay_t< dimension > &t_triangulation) -> bool
 Fix the vertices of a cell to be consistent with the foliation.
 
template<int dimension>
auto foliated_triangulations::fix_vertices (Delaunay_t< dimension > const &t_triangulation, double const t_initial_radius, double const t_foliation_spacing) -> bool
 Fix vertices with incorrect timevalues.
 
template<int dimension>
auto foliated_triangulations::fix_vertices (std::vector< Cell_handle_t< dimension > > const &t_cells, double t_initial_radius, double t_foliation_spacing)
 Fix vertices with incorrect timevalues.
 
template<int dimension>
auto foliated_triangulations::get_vertices_from_cells (std::vector< Cell_handle_t< dimension > > const &t_cells)
 Extracts vertices from cells.
 
template<int dimension>
auto foliated_triangulations::is_cell_type_correct (Cell_handle_t< dimension > const &t_cell) -> bool
 Checks if a cell is classified correctly.
 
template<int dimension>
auto foliated_triangulations::is_vertex_timevalue_correct (Vertex_handle_t< dimension > const &t_vertex, double const t_initial_radius, double const t_foliation_spacing) -> bool
 Checks if vertex timevalue is correct.
 
template<int dimension>
auto foliated_triangulations::make_causal_vertices (std::span< Point_t< dimension > const > vertices, std::span< size_t const > timevalues) -> Causal_vertices_t< dimension >
 Create causal vertices from vertices and timevalues.
 
template<int dimension>
auto foliated_triangulations::make_foliated_ball (Int_precision const t_simplices, Int_precision const t_timeslices, double const initial_radius=INITIAL_RADIUS, double const foliation_spacing=FOLIATION_SPACING)
 Make foliated ball.
 
template<int dimension>
auto foliated_triangulations::make_triangulation (Int_precision const t_simplices, Int_precision t_timeslices, double const initial_radius=INITIAL_RADIUS, double const foliation_spacing=FOLIATION_SPACING) -> Delaunay_t< dimension >
 Make a Delaunay triangulation.
 
template<int dimension>
void foliated_triangulations::print_cell (Cell_handle_t< dimension > cell)
 Print a cell in the triangulation.
 
template<int dimension, typename Container >
void foliated_triangulations::print_cells (Container &&t_cells)
 Print timevalues of each vertex in the cell and the resulting cell->info()
 
template<int dimension>
void foliated_triangulations::print_edge (Edge_handle_t< dimension > const &t_edge)
 Print edge.
 
template<int dimension>
void foliated_triangulations::print_neighboring_cells (Cell_handle_t< dimension > cell)
 Print neighboring cells.
 
template<int dimension>
auto foliated_triangulations::squared_radius (Vertex_handle_t< dimension > const &t_vertex) -> double
 Calculate the squared radius from the origin.
 
template<int dimension, ContainerType Container>
auto foliated_triangulations::volume_per_timeslice (Container &&t_facets) -> std::multimap< Int_precision, Facet_t< 3 > >
 Collect spacelike facets into a container indexed by time value.
 

Variables

template<int dimension>
auto constexpr foliated_triangulations::compare_v_info
 
static int constexpr foliated_triangulations::MAX_FIX_PASSES = 50
 

Detailed Description

Create foliated spherical triangulations.

Author
Adam Getchell

Extends CGAL's Delaunay_triangulation_3 and Triangulation_3 classes to create foliated spherical triangulations of a given dimension.

The dimensionality, number of desired simplices, and number of desired timeslices is given. Successive spheres are created with increasing radii, parameterized by INITIAL_RADIUS and RADIAL_FACTOR. Each vertex at a given radius is assigned a timeslice so that the entire triangulation will have a defined foliation of time.

See also
https://doc.cgal.org/latest/Triangulation_3/index.html#Chapter_3D_Triangulations

Definition in file Foliated_triangulation.hpp.

Typedef Documentation

◆ Causal_vertices_t

template<int dimension>
using Causal_vertices_t = std::vector<std::pair<Point_t<dimension>, Int_precision> >

Definition at line 34 of file Foliated_triangulation.hpp.

◆ Cell_handle_t

template<int dimension>
using Cell_handle_t = typename TriangulationTraits<dimension>::Cell_handle

Definition at line 38 of file Foliated_triangulation.hpp.

◆ Delaunay_t

template<int dimension>
using Delaunay_t = typename TriangulationTraits<dimension>::Delaunay

Definition at line 28 of file Foliated_triangulation.hpp.

◆ Edge_handle_t

template<int dimension>
using Edge_handle_t = typename TriangulationTraits<dimension>::Edge_handle

Definition at line 47 of file Foliated_triangulation.hpp.

◆ Face_handle_t

template<int dimension>
using Face_handle_t = typename TriangulationTraits<dimension>::Face_handle

Definition at line 41 of file Foliated_triangulation.hpp.

◆ Facet_t

template<int dimension>
using Facet_t = typename TriangulationTraits<dimension>::Facet

Definition at line 44 of file Foliated_triangulation.hpp.

◆ FoliatedTriangulation_3

using foliated_triangulations::FoliatedTriangulation_3 = typedef FoliatedTriangulation<3>

Definition at line 1577 of file Foliated_triangulation.hpp.

◆ FoliatedTriangulation_4

using foliated_triangulations::FoliatedTriangulation_4 = typedef FoliatedTriangulation<4>

Definition at line 1584 of file Foliated_triangulation.hpp.

◆ Point_t

template<int dimension>
using Point_t = typename TriangulationTraits<dimension>::Point

Definition at line 31 of file Foliated_triangulation.hpp.

◆ Spherical_points_generator_t

template<int dimension>
using Spherical_points_generator_t = typename TriangulationTraits<dimension>::Spherical_points_generator

Definition at line 53 of file Foliated_triangulation.hpp.

◆ Vertex_handle_t

template<int dimension>
using Vertex_handle_t = typename TriangulationTraits<dimension>::Vertex_handle

Definition at line 50 of file Foliated_triangulation.hpp.

Enumeration Type Documentation

◆ Cell_type

enum class Cell_type
strong

(n,m) is number of vertices on (lower, higher) timeslice

Definition at line 64 of file Foliated_triangulation.hpp.

Function Documentation

◆ check_cells()

template<int dimension>
auto foliated_triangulations::check_cells ( Delaunay_t< dimension > const &  t_triangulation) -> bool

Check all finite cells in the Delaunay triangulation.

Template Parameters
dimensionDimensionality of the Delaunay triangulation
Parameters
t_triangulationThe Delaunay triangulation
Returns
True if there are no finite cells, or all finite cells are correctly classified

Definition at line 584 of file Foliated_triangulation.hpp.

References foliated_triangulations::check_cells().

Referenced by foliated_triangulations::check_cells().

◆ check_timevalues()

template<int dimension>
auto foliated_triangulations::check_timevalues ( Delaunay_t< dimension > const &  t_triangulation) -> std::optional<std::vector<Cell_handle_t<dimension>>>

Check cells for correct foliation.

The timevalues of the vertices of a cell differ by at most one and cannot all be the same. The first case would correspond to the cell (simplex) spanning more than one timeslice; the second would correspond to the cell being purely spacelike. Both of these cases are causally inconsistent. Note that this takes a Delaunay triangulation as input, as it is expected to be called while the Foliated triangulation is still being constructed. Note also that to guard against numeric errors causing invalid cells, fix_vertices() should be called before this function.

Template Parameters
dimensionThe dimensionality of the cells and triangulation
Parameters
t_triangulationThe Delaunay triangulation
Returns
An optional container of invalid cells

Definition at line 787 of file Foliated_triangulation.hpp.

References foliated_triangulations::check_timevalues().

Referenced by foliated_triangulations::check_timevalues().

◆ check_vertices()

template<int dimension>
auto foliated_triangulations::check_vertices ( Delaunay_t< dimension > const &  t_triangulation,
double  t_initial_radius,
double  t_foliation_spacing 
)

Check if vertices have the correct timevalues.

Template Parameters
dimensionDimensionality of the vertices and Delaunay triangulation
Parameters
t_triangulationThe Delaunay triangulation
t_initial_radiusThe initial radius of the radial foliation
t_foliation_spacingThe spacing between successive leaves
Returns
True if all vertices have correct timevalues

Definition at line 367 of file Foliated_triangulation.hpp.

References foliated_triangulations::check_vertices().

Referenced by foliated_triangulations::check_vertices().

◆ classify_edge()

template<int dimension>
auto foliated_triangulations::classify_edge ( Edge_handle_t< dimension > const &  t_edge) -> bool

Predicate to classify edge as timelike or spacelike.

Template Parameters
dimensionThe dimensionality of the simplices
Parameters
t_edgeThe Edge_handle to classify
Returns
True if timelike and false if spacelike

Definition at line 227 of file Foliated_triangulation.hpp.

References foliated_triangulations::classify_edge().

Referenced by foliated_triangulations::classify_edge().

◆ collect_cells()

template<int dimension>
auto foliated_triangulations::collect_cells ( Delaunay_t< dimension > const &  t_triangulation) -> std::vector<Cell_handle_t<dimension>>

Obtain all finite cells in the Delaunay triangulation.

Template Parameters
dimensionDimensionality of the Delaunay triangulation
Parameters
t_triangulationThe triangulation
Returns
A container of finite cells

Definition at line 384 of file Foliated_triangulation.hpp.

References foliated_triangulations::collect_cells().

Referenced by foliated_triangulations::FoliatedTriangulation< 3 >::FoliatedTriangulation(), and foliated_triangulations::collect_cells().

◆ collect_edges()

template<int dimension>
auto foliated_triangulations::collect_edges ( Delaunay_t< dimension > const &  delaunay)

Returns a container of all the finite edges in the triangulation.

Regardless of the dimensionality of the triangulation, the edges are 1-d simplices connecting 0-d vertices.

Template Parameters
dimensionThe dimensionality of the triangulation
Parameters
delaunayThe triangulation
Returns
Container of all the finite edges in the triangulation

Definition at line 114 of file Foliated_triangulation.hpp.

References foliated_triangulations::collect_edges().

Referenced by foliated_triangulations::FoliatedTriangulation< 3 >::FoliatedTriangulation(), and foliated_triangulations::collect_edges().

◆ collect_vertices()

template<int dimension>
auto foliated_triangulations::collect_vertices ( Delaunay_t< dimension > const &  t_triangulation)

Obtain all finite vertices in the Delaunay triangulation.

Template Parameters
dimensionDimensionality of the Delaunay triangulation
Parameters
t_triangulationThe Delaunay triangulation
Returns
A container of finite vertices

Definition at line 345 of file Foliated_triangulation.hpp.

References foliated_triangulations::collect_vertices().

Referenced by foliated_triangulations::FoliatedTriangulation< 3 >::FoliatedTriangulation(), and foliated_triangulations::collect_vertices().

◆ debug_print_cells()

template<int dimension, ContainerType Container>
void foliated_triangulations::debug_print_cells ( Container &&  t_cells)

Write to debug log timevalues of each vertex in the cell and the resulting cell->info.

Template Parameters
dimensionThe dimensionality of the simplices
ContainerThe type of container
Parameters
t_cellsThe cells to write to debug log

Definition at line 667 of file Foliated_triangulation.hpp.

References foliated_triangulations::debug_print_cells().

Referenced by foliated_triangulations::debug_print_cells().

◆ expected_cell_type()

template<int dimension>
auto foliated_triangulations::expected_cell_type ( Cell_handle_t< dimension > const &  t_cell)

Classifies cells by their timevalues.

Template Parameters
dimensionThe dimensionality of the simplices
Parameters
t_cellThe simplex to check
Returns
The type of the simplex

Definition at line 503 of file Foliated_triangulation.hpp.

References foliated_triangulations::expected_cell_type().

Referenced by foliated_triangulations::expected_cell_type().

◆ expected_timevalue()

template<int dimension>
auto foliated_triangulations::expected_timevalue ( Vertex_handle_t< dimension > const &  t_vertex,
double  t_initial_radius,
double  t_foliation_spacing 
) -> Int_precision

Find the expected timevalue for a vertex.

The formula for the expected timevalue is:

\[t=\frac{R-I+S}{S}\]

Where R is radius, I is INITIAL_RADIUS, and S is RADIAL_SEPARATION

Template Parameters
dimensionDimensionality of the vertex
Parameters
t_vertexThe vertex
t_initial_radiusThe initial radius of the radial foliation
t_foliation_spacingThe spacing between successive leaves
Returns
The effective radius of the vertex

Definition at line 309 of file Foliated_triangulation.hpp.

References foliated_triangulations::expected_timevalue().

Referenced by foliated_triangulations::expected_timevalue(), and foliated_triangulations::FoliatedTriangulation< 3 >::print_vertices().

◆ filter_cells()

template<int dimension>
auto foliated_triangulations::filter_cells ( std::vector< Cell_handle_t< dimension > > const &  t_cells,
Cell_type const &  t_cell_type 
) -> std::vector<Cell_handle_t<dimension>>
Template Parameters
dimensionThe dimensionality of the simplices
Parameters
t_cellsThe container of simplices
t_cell_typeThe type of simplex to filter by
Returns
A container of simplices filtered by type

Definition at line 269 of file Foliated_triangulation.hpp.

References foliated_triangulations::filter_cells().

Referenced by foliated_triangulations::FoliatedTriangulation< 3 >::FoliatedTriangulation(), and foliated_triangulations::filter_cells().

◆ filter_edges()

template<int dimension>
auto foliated_triangulations::filter_edges ( std::vector< Edge_handle_t< dimension > > const &  t_edges,
bool  t_is_Timelike_pred 
) -> std::vector<Edge_handle_t<dimension>>
Template Parameters
dimensionThe dimensionality of the simplices
Parameters
t_edgesThe container of edges to filter
t_is_Timelike_predThe predicate to filter by
Returns
A container of is_Timelike edges

Definition at line 250 of file Foliated_triangulation.hpp.

References foliated_triangulations::filter_edges().

Referenced by foliated_triangulations::FoliatedTriangulation< 3 >::FoliatedTriangulation(), and foliated_triangulations::filter_edges().

◆ find_bad_vertex()

template<int dimension>
auto foliated_triangulations::find_bad_vertex ( Cell_handle_t< dimension > const &  cell) -> Vertex_handle_t<dimension>

Find the vertex that is causing a cell's foliation to be invalid.

Template Parameters
dimensionDimensionality of the cell
Parameters
cellThe cell to check
Returns
The offending vertex

Definition at line 809 of file Foliated_triangulation.hpp.

References foliated_triangulations::find_bad_vertex().

Referenced by foliated_triangulations::find_bad_vertex().

◆ find_cell()

template<int dimension>
auto foliated_triangulations::find_cell ( Delaunay_t< dimension > const &  delaunay,
Vertex_handle_t< dimension > const &  vh1,
Vertex_handle_t< dimension > const &  vh2,
Vertex_handle_t< dimension > const &  vh3,
Vertex_handle_t< dimension > const &  vh4 
) -> std::optional<Cell_handle_t<dimension>>

Returns the cell containing the vertices.

Template Parameters
dimensionThe dimensionality of the triangulation
Parameters
delaunayThe triangulation
vh1The first vertex
vh2The second vertex
vh3The third vertex
vh4The fourth vertex
Returns
The cell containing the vertices
See also
https://doc.cgal.org/latest/Triangulation_3/classCGAL_1_1Triangulation__3.html#a8766c9a0c2a84203be31537e5e015646

Definition at line 165 of file Foliated_triangulation.hpp.

References foliated_triangulations::find_cell().

Referenced by foliated_triangulations::find_cell().

◆ find_incorrect_cells()

template<int dimension>
auto foliated_triangulations::find_incorrect_cells ( Delaunay_t< dimension > const &  t_triangulation)

Check all finite cells in the Delaunay triangulation.

Template Parameters
dimensionDimensionality of the Delaunay triangulation
Parameters
t_triangulationThe Delaunay triangulation
Returns
A container of cells that are not classified correctly

Definition at line 600 of file Foliated_triangulation.hpp.

References foliated_triangulations::find_incorrect_cells().

Referenced by foliated_triangulations::find_incorrect_cells().

◆ find_incorrect_vertices() [1/2]

template<int dimension>
auto foliated_triangulations::find_incorrect_vertices ( Delaunay_t< dimension > const &  t_triangulation,
double  t_initial_radius,
double  t_foliation_spacing 
)

Obtain vertices with incorrect timevalues.

Template Parameters
dimensionDimensionality of the vertices and Delaunay triangulation
Parameters
t_triangulationThe Delaunay triangulation
t_initial_radiusThe initial radius of the radial foliation
t_foliation_spacingThe spacing between successive leaves
Returns
A container of vertices with incorrect timevalues

Definition at line 449 of file Foliated_triangulation.hpp.

References foliated_triangulations::find_incorrect_vertices().

◆ find_incorrect_vertices() [2/2]

template<int dimension>
auto foliated_triangulations::find_incorrect_vertices ( std::vector< Cell_handle_t< dimension > > const &  t_cells,
double  t_initial_radius,
double  t_foliation_spacing 
)

Obtain vertices with incorrect timevalues.

Template Parameters
dimensionDimensionality of vertices and cells
Parameters
t_cellsContainer of cells to check
t_initial_radiusThe initial radius of the radial foliation
t_foliation_spacingThe spacing between successive leaves
Returns
A container of vertices with incorrect timevalues

Definition at line 425 of file Foliated_triangulation.hpp.

References foliated_triangulations::find_incorrect_vertices().

Referenced by foliated_triangulations::find_incorrect_vertices().

◆ find_max_timevalue()

template<int dimension, ContainerType Container>
auto foliated_triangulations::find_max_timevalue ( Container &&  t_vertices) -> Int_precision
Template Parameters
dimensionThe dimensionality of the simplices
Parameters
t_verticesThe container of vertices
Returns
The maximum timevalue in the container

Definition at line 192 of file Foliated_triangulation.hpp.

References foliated_triangulations::find_max_timevalue().

Referenced by foliated_triangulations::FoliatedTriangulation< 3 >::FoliatedTriangulation(), and foliated_triangulations::find_max_timevalue().

◆ find_min_timevalue()

template<int dimension, ContainerType Container>
auto foliated_triangulations::find_min_timevalue ( Container &&  t_vertices) -> Int_precision
Template Parameters
dimensionThe dimensionality of the simplices
Parameters
t_verticesThe container of vertices
Returns
The minimum timevalue in the container

Definition at line 210 of file Foliated_triangulation.hpp.

References foliated_triangulations::find_min_timevalue().

Referenced by foliated_triangulations::FoliatedTriangulation< 3 >::FoliatedTriangulation(), and foliated_triangulations::find_min_timevalue().

◆ find_vertex()

template<int dimension>
auto foliated_triangulations::find_vertex ( Delaunay_t< dimension > const &  delaunay,
Point_t< dimension > const &  point 
) -> std::optional<Vertex_handle_t<dimension>>

Returns the vertex containing the given point.

Template Parameters
dimensionThe dimensionality of the triangulation
Parameters
delaunayThe triangulation
pointThe point to find the vertex for
Returns
The vertex containing the given point
See also
https://doc.cgal.org/latest/Triangulation_3/classCGAL_1_1Triangulation__3.html#a5b45572c663e5d2c10f26e7be421e140

Definition at line 142 of file Foliated_triangulation.hpp.

References foliated_triangulations::find_vertex().

Referenced by foliated_triangulations::find_vertex().

◆ fix_cells()

template<int dimension>
auto foliated_triangulations::fix_cells ( Delaunay_t< dimension > const &  t_triangulation) -> bool

Fix simplices with the wrong type.

Template Parameters
dimensionThe dimensionality of the simplices
Parameters
t_triangulationThe Delaunay triangulation
Returns
True if cells->info() was fixed

Definition at line 618 of file Foliated_triangulation.hpp.

References foliated_triangulations::fix_cells().

Referenced by foliated_triangulations::fix_cells().

◆ fix_timevalues()

template<int dimension>
auto foliated_triangulations::fix_timevalues ( Delaunay_t< dimension > &  t_triangulation) -> bool

Fix the vertices of a cell to be consistent with the foliation.

Template Parameters
dimensionDimensionality of the triangulation
Parameters
t_triangulationThe Delaunay triangulation
Returns
True if incorrectly foliated simplices were fixed

Definition at line 843 of file Foliated_triangulation.hpp.

References foliated_triangulations::fix_timevalues().

Referenced by foliated_triangulations::fix_timevalues().

◆ fix_vertices() [1/2]

template<int dimension>
auto foliated_triangulations::fix_vertices ( Delaunay_t< dimension > const &  t_triangulation,
double const  t_initial_radius,
double const  t_foliation_spacing 
) -> bool

Fix vertices with incorrect timevalues.

Changes vertex->info() to the correct timevalue

Template Parameters
dimensionDimensionality of the vertices and Delaunay triangulation
Parameters
t_triangulationThe triangulation
t_initial_radiusThe initial radius of the radial foliation
t_foliation_spacingThe spacing between successive leaves
Returns
True if any vertex->info() was fixed

Definition at line 490 of file Foliated_triangulation.hpp.

References foliated_triangulations::fix_vertices().

◆ fix_vertices() [2/2]

template<int dimension>
auto foliated_triangulations::fix_vertices ( std::vector< Cell_handle_t< dimension > > const &  t_cells,
double  t_initial_radius,
double  t_foliation_spacing 
)

Fix vertices with incorrect timevalues.

Changes vertex->info() to the correct timevalue using foliated_triangulations::expected_timevalue

Template Parameters
dimensionDimensionality of vertices and cells
Parameters
t_cellsContainer of cells to check
t_initial_radiusThe initial radius of the radial foliation
t_foliation_spacing
Returns
True if any vertex->info() was fixed

Definition at line 467 of file Foliated_triangulation.hpp.

References foliated_triangulations::fix_vertices().

Referenced by foliated_triangulations::fix_vertices().

◆ get_vertices_from_cells()

template<int dimension>
auto foliated_triangulations::get_vertices_from_cells ( std::vector< Cell_handle_t< dimension > > const &  t_cells)

Extracts vertices from cells.

Parameters
t_cellsThe cells from which to extract vertices
Returns
All of the vertices contained in the cells

Definition at line 402 of file Foliated_triangulation.hpp.

References foliated_triangulations::get_vertices_from_cells().

Referenced by foliated_triangulations::get_vertices_from_cells().

◆ is_cell_type_correct()

template<int dimension>
auto foliated_triangulations::is_cell_type_correct ( Cell_handle_t< dimension > const &  t_cell) -> bool

Checks if a cell is classified correctly.

Template Parameters
dimensionThe dimensionality of the simplices
Parameters
t_cellThe simplex to check
Returns
True if the cell_info matches expected cell_info

Definition at line 569 of file Foliated_triangulation.hpp.

References foliated_triangulations::is_cell_type_correct().

Referenced by foliated_triangulations::is_cell_type_correct().

◆ is_vertex_timevalue_correct()

template<int dimension>
auto foliated_triangulations::is_vertex_timevalue_correct ( Vertex_handle_t< dimension > const &  t_vertex,
double const  t_initial_radius,
double const  t_foliation_spacing 
) -> bool

Checks if vertex timevalue is correct.

Template Parameters
dimensionDimensionality of the vertex
Parameters
t_vertexThe vertex
t_initial_radiusThe initial radius of the radial foliation
t_foliation_spacingThe spacing between successive leaves
Returns
True if the timevalue of the vertex matches its effective radius

Definition at line 326 of file Foliated_triangulation.hpp.

References foliated_triangulations::is_vertex_timevalue_correct().

Referenced by foliated_triangulations::is_vertex_timevalue_correct().

◆ make_causal_vertices()

template<int dimension>
auto foliated_triangulations::make_causal_vertices ( std::span< Point_t< dimension > const >  vertices,
std::span< size_t const >  timevalues 
) -> Causal_vertices_t<dimension>

Create causal vertices from vertices and timevalues.

Template Parameters
dimensionDimensionality of the manifold
Parameters
verticesThe vertices of the manifold
timevaluesThe timevalue of each vertex
Returns
A container of vertices that have an associated timevalue

Definition at line 89 of file Foliated_triangulation.hpp.

References foliated_triangulations::make_causal_vertices().

Referenced by foliated_triangulations::make_causal_vertices().

◆ make_foliated_ball()

template<int dimension>
auto foliated_triangulations::make_foliated_ball ( Int_precision const  t_simplices,
Int_precision const  t_timeslices,
double const  initial_radius = INITIAL_RADIUS,
double const  foliation_spacing = FOLIATION_SPACING 
)

Make foliated ball.

Makes a solid ball of successive layers of spheres at a given radius.

Template Parameters
dimensionThe dimensionality of the simplices
Parameters
t_simplicesThe desired number of simplices in the triangulation
t_timeslicesThe desired number of timeslices in the triangulation
initial_radiusThe radius of the first time slice
foliation_spacingThe distance between successive time slices
Returns
A container of (vertex, timevalue) pairs

Definition at line 882 of file Foliated_triangulation.hpp.

References foliated_triangulations::make_foliated_ball().

Referenced by foliated_triangulations::make_foliated_ball().

◆ make_triangulation()

template<int dimension>
auto foliated_triangulations::make_triangulation ( Int_precision const  t_simplices,
Int_precision  t_timeslices,
double const  initial_radius = INITIAL_RADIUS,
double const  foliation_spacing = FOLIATION_SPACING 
) -> Delaunay_t<dimension>

Make a Delaunay triangulation.

Template Parameters
dimensionDimensionality of the Delaunay triangulation
Parameters
t_simplicesNumber of desired simplices
t_timeslicesNumber of desired timeslices
initial_radiusRadius of first timeslice
foliation_spacingRadial separation between timeslices
Returns
A Delaunay triangulation with a timevalue for each vertex

Definition at line 916 of file Foliated_triangulation.hpp.

References foliated_triangulations::make_triangulation().

Referenced by foliated_triangulations::FoliatedTriangulation< 3 >::FoliatedTriangulation(), and foliated_triangulations::make_triangulation().

◆ print_cell()

template<int dimension>
void foliated_triangulations::print_cell ( Cell_handle_t< dimension >  cell)

Print a cell in the triangulation.

Template Parameters
dimensionThe dimensionality of the triangulation
Parameters
cellThe cell to print

Definition at line 634 of file Foliated_triangulation.hpp.

References foliated_triangulations::print_cell().

Referenced by foliated_triangulations::print_cell().

◆ print_cells()

template<int dimension, typename Container >
void foliated_triangulations::print_cells ( Container &&  t_cells)

Print timevalues of each vertex in the cell and the resulting cell->info()

Template Parameters
dimensionThe dimensionality of the simplices
Parameters
t_cellsThe cells to print

Definition at line 652 of file Foliated_triangulation.hpp.

References foliated_triangulations::print_cells().

Referenced by foliated_triangulations::print_cells().

◆ print_edge()

template<int dimension>
void foliated_triangulations::print_edge ( Edge_handle_t< dimension > const &  t_edge)

Print edge.

An edge is represented by a cell and two indices which refer to the vertices of the cell connected by the edge. The data is stored in a CGAL::Triple, which is an extension of std::pair.

Template Parameters
dimensionThe dimensionality of the simplices
Parameters
t_edgeThe edge to print
See also
https://doc.cgal.org/latest/TDS_3/index.html#fig__TDS3figrepres
https://doc.cgal.org/latest/STL_Extension/classCGAL_1_1Triple.html

Definition at line 705 of file Foliated_triangulation.hpp.

References foliated_triangulations::print_edge().

Referenced by foliated_triangulations::print_edge().

◆ print_neighboring_cells()

template<int dimension>
void foliated_triangulations::print_neighboring_cells ( Cell_handle_t< dimension >  cell)

Print neighboring cells.

Template Parameters
dimensionThe dimensionality of the simplices
Parameters
cellThe cell to print neighbors of

Definition at line 687 of file Foliated_triangulation.hpp.

References foliated_triangulations::print_neighboring_cells().

Referenced by foliated_triangulations::print_neighboring_cells().

◆ squared_radius()

template<int dimension>
auto foliated_triangulations::squared_radius ( Vertex_handle_t< dimension > const &  t_vertex) -> double

Calculate the squared radius from the origin.

Template Parameters
dimensionThe dimensionality of the simplices
Parameters
t_vertexThe vertex to check
Returns
The squared radial distance of the vertex from the origin

Definition at line 289 of file Foliated_triangulation.hpp.

References foliated_triangulations::squared_radius().

Referenced by foliated_triangulations::squared_radius().

◆ volume_per_timeslice()

template<int dimension, ContainerType Container>
auto foliated_triangulations::volume_per_timeslice ( Container &&  t_facets) -> std::multimap<Int_precision, Facet_t<3>>

Collect spacelike facets into a container indexed by time value.

Warning! Turning on debugging info will generate gigabytes of logs.

Template Parameters
dimensionThe dimensionality of the simplices
Parameters
t_facetsA container of facets
Returns
Container with spacelike facets per timeslice

Definition at line 724 of file Foliated_triangulation.hpp.

References foliated_triangulations::volume_per_timeslice().

Referenced by foliated_triangulations::FoliatedTriangulation< 3 >::FoliatedTriangulation(), and foliated_triangulations::volume_per_timeslice().

Variable Documentation

◆ compare_v_info

template<int dimension>
auto constexpr foliated_triangulations::compare_v_info
constexpr
Initial value:
= [](Vertex_handle_t<dimension> const& lhs,
Vertex_handle_t<dimension> const& rhs) {
return lhs->info() < rhs->info();
}
Template Parameters
dimensionThe dimensionality of the simplices
Returns
True if timevalue of lhs is less than rhs

Definition at line 183 of file Foliated_triangulation.hpp.

◆ MAX_FIX_PASSES

int constexpr foliated_triangulations::MAX_FIX_PASSES = 50
staticconstexpr

Definition at line 81 of file Foliated_triangulation.hpp.