CDT++
Causal Dynamical Triangulations in C++
Loading...
Searching...
No Matches
Public Member Functions | Static Public Attributes | Private Types | Private Member Functions | Private Attributes | Friends | List of all members
manifolds::Manifold< 3 > Class Reference

3D Manifold More...

#include <Manifold.hpp>

+ Collaboration diagram for manifolds::Manifold< 3 >:

Public Member Functions

 Manifold ()=default
 Default ctor.
 
 Manifold (Causal_vertices_t< 3 > const &causal_vertices, double const t_initial_radius=INITIAL_RADIUS, double const t_foliation_spacing=FOLIATION_SPACING)
 Construct manifold from Causal_vertices Pass-by-value-then-move.
 
 Manifold (Int_precision const t_desired_simplices, Int_precision const t_desired_timeslices, double const t_initial_radius=INITIAL_RADIUS, double const t_foliation_spacing=FOLIATION_SPACING)
 Construct manifold using arguments.
 
 Manifold (Manifold &&other)=default
 Default move ctor.
 
 Manifold (Manifold const &other)=default
 Default copy ctor.
 
 Manifold (Triangulation t_foliated_triangulation)
 Construct manifold from a Foliated triangulation.
 
 ~Manifold ()=default
 Default dtor.
 
auto check_simplices () const -> bool
 
template<typename VertexHandle >
auto degree (VertexHandle &&t_vertex) const -> decltype(auto)
 Perfect forwarding to FoliatedTriangulation_3.degree()
 
auto dimensionality () const
 
auto edges () const
 
auto faces () const
 
auto foliation_spacing () const
 Radial separation between timeslices.
 
auto get_cell (Vertex_handle_t< 3 > const &vh1, Vertex_handle_t< 3 > const &vh2, Vertex_handle_t< 3 > const &vh3, Vertex_handle_t< 3 > const &vh4) const -> Cell_handle_t< 3 >
 Get a cell handle from 4 vertex handles.
 
auto get_delaunay () const noexcept
 
auto get_geometry () const -> Geometry const &
 
auto get_spacelike_edges () const -> auto const &
 Call triangulation.get_spacelike_edges()
 
auto get_timelike_edges () const noexcept -> auto const &
 Call to triangulation_.get_timelike_edges()
 
auto get_triangulation () const noexcept -> Triangulation const &
 
auto get_vertex (Point_t< 3 > const &point) const -> Vertex_handle_t< 3 >
 Obtains a vertex handle from a point.
 
auto get_vertices () const noexcept -> auto const &
 Call FoliatedTriangulation_3.get_vertices()
 
auto get_vertices_span () const noexcept -> auto
 Call FoliatedTriangulation_3.get_vertices_span()
 
template<typename... Ts>
auto incident_cells (Ts &&... args) const noexcept -> decltype(auto)
 Perfect forwarding to FoliatedTriangulation_3.incident_cells()
 
auto initial_radius () const
 Initial radius of the first timeslice.
 
auto is_correct () const -> bool
 
auto is_delaunay () const -> bool
 Forwarding to FoliatedTriangulation.is_delaunay()
 
auto is_edge (Edge_handle_t< 3 > const &t_edge_candidate) const noexcept -> bool
 Forwarding to FoliatedTriangulation_3.is_edge()
 
auto is_foliated () const -> bool
 Forwarding to FoliatedTriangulation_3.is_foliated()
 
auto is_valid () const -> bool
 Forwarding to FoliatedTriangulation.is_tds_valid()
 
template<typename VertexType >
auto is_vertex (VertexType &&t_vertex_candidate) const -> bool
 Perfect forwarding to FoliatedTriangulation_3.is_vertex()
 
auto max_time () const
 
auto min_time () const
 
auto N0 () const
 
auto N1 () const
 
auto N1_SL () const
 
auto N1_TL () const
 
auto N2 () const
 
auto N2_SL () const -> auto const &
 
auto N3 () const
 
auto N3_13 () const
 
auto N3_22 () const
 
auto N3_31 () const
 
auto N3_31_13 () const
 
auto operator= (Manifold &&other) -> Manifold &=default
 Default move assignment.
 
auto operator= (Manifold const &other) -> Manifold &=default
 Default copy assignment.
 
void print () const
 Print manifold.
 
void print_cells () const
 Print timevalues of each vertex in the cell and the resulting cell->info()
 
void print_details () const
 Print details of the manifold.
 
void print_vertices () const
 Print values of a vertex->info()
 
void print_volume_per_timeslice () const
 Print the codimension 1 volume of simplices (faces) per timeslice.
 
auto simplices () const
 
auto triangulation () -> Triangulation &
 
void update ()
 Update the Manifold data structures.
 
auto vertices () const
 

Static Public Attributes

static int constexpr dimension = 3
 Dimensionality of the manifold.
 
static topology_type constexpr topology = topology_type::SPHERICAL
 Topology of the manifold.
 

Private Types

using Geometry = Geometry_3
 
using Triangulation = foliated_triangulations::FoliatedTriangulation_3
 

Private Member Functions

void update_geometry () noexcept
 Update geometry data of the manifold when the triangulation has been changed.
 
void update_triangulation ()
 Update the triangulation.
 

Private Attributes

Geometry m_geometry
 The data structure of scalar values for computations.
 
Triangulation m_triangulation
 The data structure of geometric and combinatorial relationships.
 

Friends

void swap (Manifold &swap_from, Manifold &swap_into) noexcept
 Non-member swap function for Manifolds.
 

Detailed Description

3D Manifold

Definition at line 42 of file Manifold.hpp.

Member Typedef Documentation

◆ Geometry

Definition at line 45 of file Manifold.hpp.

◆ Triangulation

Definition at line 44 of file Manifold.hpp.

Constructor & Destructor Documentation

◆ Manifold() [1/3]

manifolds::Manifold< 3 >::Manifold ( Triangulation  t_foliated_triangulation)
inlineexplicit

Construct manifold from a Foliated triangulation.

Parameters
t_foliated_triangulationTriangulation used to construct manifold

Definition at line 95 of file Manifold.hpp.

◆ Manifold() [2/3]

manifolds::Manifold< 3 >::Manifold ( Int_precision const  t_desired_simplices,
Int_precision const  t_desired_timeslices,
double const  t_initial_radius = INITIAL_RADIUS,
double const  t_foliation_spacing = FOLIATION_SPACING 
)
inline

Construct manifold using arguments.

Parameters
t_desired_simplicesNumber of desired simplices
t_desired_timeslicesNumber of desired timeslices
t_initial_radiusRadius of first timeslice
t_foliation_spacingRadial separation between timeslices

Definition at line 105 of file Manifold.hpp.

◆ Manifold() [3/3]

manifolds::Manifold< 3 >::Manifold ( Causal_vertices_t< 3 > const &  causal_vertices,
double const  t_initial_radius = INITIAL_RADIUS,
double const  t_foliation_spacing = FOLIATION_SPACING 
)
inlineexplicit

Construct manifold from Causal_vertices Pass-by-value-then-move.

Parameters
causal_verticesCausal_vertices to place into the Manifold
t_initial_radiusRadius of first timeslice
t_foliation_spacingRadial separation between timeslices

Definition at line 120 of file Manifold.hpp.

Member Function Documentation

◆ check_simplices()

auto manifolds::Manifold< 3 >::check_simplices ( ) const -> bool
inline
Returns
True if all cells in triangulation are classified and match number in geometry

Definition at line 352 of file Manifold.hpp.

References foliated_triangulations::FoliatedTriangulation< 3 >::check_all_cells().

◆ degree()

template<typename VertexHandle >
auto manifolds::Manifold< 3 >::degree ( VertexHandle &&  t_vertex) const -> decltype(auto)
inline

Perfect forwarding to FoliatedTriangulation_3.degree()

Definition at line 313 of file Manifold.hpp.

References foliated_triangulations::FoliatedTriangulation< 3 >::degree().

◆ dimensionality()

auto manifolds::Manifold< 3 >::dimensionality ( ) const
inline
Returns
Run-time dimensionality of the triangulation data structure

Definition at line 219 of file Manifold.hpp.

References foliated_triangulations::FoliatedTriangulation< 3 >::dimension().

◆ edges()

auto manifolds::Manifold< 3 >::edges ( ) const
inline
Returns
Number of 1D edges in triangulation data structure

Definition at line 284 of file Manifold.hpp.

References foliated_triangulations::FoliatedTriangulation< 3 >::number_of_finite_edges().

◆ faces()

auto manifolds::Manifold< 3 >::faces ( ) const
inline
Returns
Number of 2D faces in triangulation data structure

Definition at line 268 of file Manifold.hpp.

References foliated_triangulations::FoliatedTriangulation< 3 >::number_of_finite_facets().

◆ foliation_spacing()

auto manifolds::Manifold< 3 >::foliation_spacing ( ) const
inline

Radial separation between timeslices.

Definition at line 231 of file Manifold.hpp.

References foliated_triangulations::FoliatedTriangulation< 3 >::foliation_spacing().

◆ get_cell()

auto manifolds::Manifold< 3 >::get_cell ( Vertex_handle_t< 3 > const &  vh1,
Vertex_handle_t< 3 > const &  vh2,
Vertex_handle_t< 3 > const &  vh3,
Vertex_handle_t< 3 > const &  vh4 
) const -> Cell_handle_t<3>
inline

Get a cell handle from 4 vertex handles.

Parameters
vh1The first vertex handle
vh2The second vertex handle
vh3The third vertex handle
vh4The fourth vertex handle
Returns
The cell handle containing the vertex handles

Definition at line 419 of file Manifold.hpp.

References foliated_triangulations::FoliatedTriangulation< 3 >::get_delaunay().

Referenced by GIVEN().

◆ get_delaunay()

auto manifolds::Manifold< 3 >::get_delaunay ( ) const
inlinenoexcept
Returns
A read-only reference to the Delaunay triangulation

Definition at line 152 of file Manifold.hpp.

Referenced by GIVEN().

◆ get_geometry()

auto manifolds::Manifold< 3 >::get_geometry ( ) const -> Geometry const&
inline
Returns
A read-only reference to the Geometry

Definition at line 164 of file Manifold.hpp.

◆ get_spacelike_edges()

auto manifolds::Manifold< 3 >::get_spacelike_edges ( ) const -> auto const&
inline

Call triangulation.get_spacelike_edges()

Definition at line 333 of file Manifold.hpp.

References foliated_triangulations::FoliatedTriangulation< 3 >::get_spacelike_edges().

◆ get_timelike_edges()

auto manifolds::Manifold< 3 >::get_timelike_edges ( ) const -> auto const&
inlinenoexcept

Call to triangulation_.get_timelike_edges()

Definition at line 327 of file Manifold.hpp.

References foliated_triangulations::FoliatedTriangulation< 3 >::get_timelike_edges().

◆ get_triangulation()

auto manifolds::Manifold< 3 >::get_triangulation ( ) const -> Triangulation const&
inlinenoexcept
Returns
A read-only reference to the triangulation

Definition at line 145 of file Manifold.hpp.

◆ get_vertex()

auto manifolds::Manifold< 3 >::get_vertex ( Point_t< 3 > const &  point) const -> Vertex_handle_t<3>
inline

Obtains a vertex handle from a point.

Parameters
pointThe point to search for
Returns
The vertex handle to the vertex containing the point

Definition at line 406 of file Manifold.hpp.

References foliated_triangulations::FoliatedTriangulation< 3 >::get_delaunay().

Referenced by GIVEN().

◆ get_vertices()

auto manifolds::Manifold< 3 >::get_vertices ( ) const -> auto const&
inlinenoexcept

Call FoliatedTriangulation_3.get_vertices()

Definition at line 339 of file Manifold.hpp.

References foliated_triangulations::FoliatedTriangulation< 3 >::get_vertices().

◆ get_vertices_span()

auto manifolds::Manifold< 3 >::get_vertices_span ( ) const -> auto
inlinenoexcept

Call FoliatedTriangulation_3.get_vertices_span()

Definition at line 345 of file Manifold.hpp.

References foliated_triangulations::FoliatedTriangulation< 3 >::get_vertices_span().

◆ incident_cells()

template<typename... Ts>
auto manifolds::Manifold< 3 >::incident_cells ( Ts &&...  args) const -> decltype(auto)
inlinenoexcept

Perfect forwarding to FoliatedTriangulation_3.incident_cells()

Definition at line 320 of file Manifold.hpp.

References foliated_triangulations::FoliatedTriangulation< 3 >::incident_cells().

◆ initial_radius()

auto manifolds::Manifold< 3 >::initial_radius ( ) const
inline

Initial radius of the first timeslice.

Definition at line 225 of file Manifold.hpp.

References foliated_triangulations::FoliatedTriangulation< 3 >::initial_radius().

◆ is_correct()

auto manifolds::Manifold< 3 >::is_correct ( ) const -> bool
inline
Returns
If base data structures are correct

Definition at line 191 of file Manifold.hpp.

References foliated_triangulations::FoliatedTriangulation< 3 >::is_correct().

Referenced by GIVEN().

◆ is_delaunay()

auto manifolds::Manifold< 3 >::is_delaunay ( ) const -> bool
inline

Forwarding to FoliatedTriangulation.is_delaunay()

Returns
True if the Manifold triangulation is Delaunay

Definition at line 178 of file Manifold.hpp.

References foliated_triangulations::FoliatedTriangulation< 3 >::is_delaunay().

◆ is_edge()

auto manifolds::Manifold< 3 >::is_edge ( Edge_handle_t< 3 > const &  t_edge_candidate) const -> bool
inlinenoexcept

Forwarding to FoliatedTriangulation_3.is_edge()

Parameters
t_edge_candidateThe edge to test
Returns
True if the candidate is an edge

Definition at line 210 of file Manifold.hpp.

References foliated_triangulations::FoliatedTriangulation< 3 >::get_delaunay().

◆ is_foliated()

auto manifolds::Manifold< 3 >::is_foliated ( ) const -> bool
inline

Forwarding to FoliatedTriangulation_3.is_foliated()

Returns
True if the Manifold triangulation is foliated

Definition at line 171 of file Manifold.hpp.

References foliated_triangulations::FoliatedTriangulation< 3 >::is_foliated().

◆ is_valid()

auto manifolds::Manifold< 3 >::is_valid ( ) const -> bool
inline

Forwarding to FoliatedTriangulation.is_tds_valid()

Returns
True if the TriangulationDataStructure is valid

Definition at line 185 of file Manifold.hpp.

References foliated_triangulations::FoliatedTriangulation< 3 >::is_tds_valid().

◆ is_vertex()

template<typename VertexType >
auto manifolds::Manifold< 3 >::is_vertex ( VertexType &&  t_vertex_candidate) const -> bool
inline

Perfect forwarding to FoliatedTriangulation_3.is_vertex()

Template Parameters
VertexTypeThe vertex type
Parameters
t_vertex_candidateThe vertex to check
Returns
True if the vertex candidate is a vertex

Definition at line 201 of file Manifold.hpp.

References foliated_triangulations::FoliatedTriangulation< 3 >::get_delaunay().

◆ max_time()

auto manifolds::Manifold< 3 >::max_time ( ) const
inline
Returns
Maximum timeslice value in triangulation data structure

Definition at line 306 of file Manifold.hpp.

References foliated_triangulations::FoliatedTriangulation< 3 >::max_time().

◆ min_time()

auto manifolds::Manifold< 3 >::min_time ( ) const
inline
Returns
Minimum timeslice value in triangulation data structure

Definition at line 300 of file Manifold.hpp.

References foliated_triangulations::FoliatedTriangulation< 3 >::min_time().

◆ N0()

auto manifolds::Manifold< 3 >::N0 ( ) const
inline
Returns
Number of vertices in geometry data structure

Definition at line 291 of file Manifold.hpp.

References Geometry< 3 >::N0.

◆ N1()

auto manifolds::Manifold< 3 >::N1 ( ) const
inline
Returns
Number of 1D edges in geometry data structure

Definition at line 275 of file Manifold.hpp.

References Geometry< 3 >::N1.

◆ N1_SL()

auto manifolds::Manifold< 3 >::N1_SL ( ) const
inline
Returns
Number of spacelike edges in triangulation data structure

Definition at line 278 of file Manifold.hpp.

References foliated_triangulations::FoliatedTriangulation< 3 >::N1_SL().

◆ N1_TL()

auto manifolds::Manifold< 3 >::N1_TL ( ) const
inline
Returns
Number of timelike edges in triangulation data structure

Definition at line 281 of file Manifold.hpp.

References foliated_triangulations::FoliatedTriangulation< 3 >::N1_TL().

◆ N2()

auto manifolds::Manifold< 3 >::N2 ( ) const
inline
Returns
Number of 2D faces in geometry data structure

Definition at line 258 of file Manifold.hpp.

References Geometry< 3 >::N2.

◆ N2_SL()

auto manifolds::Manifold< 3 >::N2_SL ( ) const -> auto const&
inline
Returns
An associative container of spacelike faces indexed by timevalue

Definition at line 262 of file Manifold.hpp.

References foliated_triangulations::FoliatedTriangulation< 3 >::N2_SL().

◆ N3()

auto manifolds::Manifold< 3 >::N3 ( ) const
inline
Returns
Number of 3D simplices in geometry data structure

Definition at line 237 of file Manifold.hpp.

References Geometry< 3 >::N3.

◆ N3_13()

auto manifolds::Manifold< 3 >::N3_13 ( ) const
inline
Returns
Number of (1,3) simplices in geometry data structure

Definition at line 246 of file Manifold.hpp.

References Geometry< 3 >::N3_13.

◆ N3_22()

auto manifolds::Manifold< 3 >::N3_22 ( ) const
inline
Returns
Number of (2,2) simplices in geometry data structure

Definition at line 243 of file Manifold.hpp.

References Geometry< 3 >::N3_22.

◆ N3_31()

auto manifolds::Manifold< 3 >::N3_31 ( ) const
inline
Returns
Number of (3,1) simplices in geometry data structure

Definition at line 240 of file Manifold.hpp.

References Geometry< 3 >::N3_31.

◆ N3_31_13()

auto manifolds::Manifold< 3 >::N3_31_13 ( ) const
inline
Returns
Number of (3,1) and (1,3) simplices in geometry data structure

Definition at line 249 of file Manifold.hpp.

References Geometry< 3 >::N3_31_13.

◆ print()

void manifolds::Manifold< 3 >::print ( ) const
inline

Print manifold.

Definition at line 372 of file Manifold.hpp.

Referenced by GIVEN(), and main().

◆ print_cells()

void manifolds::Manifold< 3 >::print_cells ( ) const
inline

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

Definition at line 369 of file Manifold.hpp.

References foliated_triangulations::FoliatedTriangulation< 3 >::print_cells().

◆ print_details()

void manifolds::Manifold< 3 >::print_details ( ) const
inline

Print details of the manifold.

Definition at line 387 of file Manifold.hpp.

Referenced by GIVEN(), and main().

◆ print_vertices()

void manifolds::Manifold< 3 >::print_vertices ( ) const
inline

Print values of a vertex->info()

Definition at line 365 of file Manifold.hpp.

References foliated_triangulations::FoliatedTriangulation< 3 >::print_vertices().

Referenced by GIVEN().

◆ print_volume_per_timeslice()

void manifolds::Manifold< 3 >::print_volume_per_timeslice ( ) const
inline

Print the codimension 1 volume of simplices (faces) per timeslice.

Definition at line 359 of file Manifold.hpp.

References foliated_triangulations::FoliatedTriangulation< 3 >::print_volume_per_timeslice().

Referenced by main().

◆ simplices()

auto manifolds::Manifold< 3 >::simplices ( ) const
inline
Returns
Number of 3D simplices in triangulation data structure

Definition at line 252 of file Manifold.hpp.

References foliated_triangulations::FoliatedTriangulation< 3 >::get_cells().

◆ triangulation()

auto manifolds::Manifold< 3 >::triangulation ( ) -> Triangulation&
inline
Returns
A mutable reference to the triangulation

Definition at line 158 of file Manifold.hpp.

◆ update()

void manifolds::Manifold< 3 >::update ( )
inline

Update the Manifold data structures.

Definition at line 130 of file Manifold.hpp.

◆ update_geometry()

void manifolds::Manifold< 3 >::update_geometry ( )
inlineprivatenoexcept

Update geometry data of the manifold when the triangulation has been changed.

Definition at line 447 of file Manifold.hpp.

◆ update_triangulation()

void manifolds::Manifold< 3 >::update_triangulation ( )
inlineprivate

Update the triangulation.

Definition at line 430 of file Manifold.hpp.

References foliated_triangulations::FoliatedTriangulation< 3 >::get_delaunay().

◆ vertices()

auto manifolds::Manifold< 3 >::vertices ( ) const
inline
Returns
Number of vertices in triangulation data structure

Definition at line 294 of file Manifold.hpp.

References foliated_triangulations::FoliatedTriangulation< 3 >::number_of_vertices().

Friends And Related Function Documentation

◆ swap

void swap ( Manifold< 3 > &  swap_from,
Manifold< 3 > &  swap_into 
)
friend

Non-member swap function for Manifolds.

Used for no-except updates of manifolds after moves.

Parameters
swap_fromThe value to be swapped from. Assumed to be discarded.
swap_intoThe value to be swapped into.

Definition at line 83 of file Manifold.hpp.

Member Data Documentation

◆ dimension

int constexpr manifolds::Manifold< 3 >::dimension = 3
staticconstexpr

Dimensionality of the manifold.

Used to determine the manifold dimension at compile-time

Definition at line 56 of file Manifold.hpp.

◆ m_geometry

Geometry manifolds::Manifold< 3 >::m_geometry
private

The data structure of scalar values for computations.

Definition at line 51 of file Manifold.hpp.

◆ m_triangulation

Triangulation manifolds::Manifold< 3 >::m_triangulation
private

The data structure of geometric and combinatorial relationships.

Definition at line 48 of file Manifold.hpp.

◆ topology

topology_type constexpr manifolds::Manifold< 3 >::topology = topology_type::SPHERICAL
staticconstexpr

Topology of the manifold.

Definition at line 59 of file Manifold.hpp.


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