|
CDT++
Causal Dynamical Triangulations in C++
|
#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. | |
Definition at line 42 of file Manifold.hpp.
|
private |
Definition at line 45 of file Manifold.hpp.
|
private |
Definition at line 44 of file Manifold.hpp.
|
inlineexplicit |
Construct manifold from a Foliated triangulation.
| t_foliated_triangulation | Triangulation used to construct manifold |
Definition at line 95 of file Manifold.hpp.
|
inline |
Construct manifold using arguments.
| t_desired_simplices | Number of desired simplices |
| t_desired_timeslices | Number of desired timeslices |
| t_initial_radius | Radius of first timeslice |
| t_foliation_spacing | Radial separation between timeslices |
Definition at line 105 of file Manifold.hpp.
|
inlineexplicit |
Construct manifold from Causal_vertices Pass-by-value-then-move.
| causal_vertices | Causal_vertices to place into the Manifold |
| t_initial_radius | Radius of first timeslice |
| t_foliation_spacing | Radial separation between timeslices |
Definition at line 120 of file Manifold.hpp.
|
inline |
Definition at line 352 of file Manifold.hpp.
References foliated_triangulations::FoliatedTriangulation< 3 >::check_all_cells().
|
inline |
Perfect forwarding to FoliatedTriangulation_3.degree()
Definition at line 313 of file Manifold.hpp.
References foliated_triangulations::FoliatedTriangulation< 3 >::degree().
|
inline |
Definition at line 219 of file Manifold.hpp.
References foliated_triangulations::FoliatedTriangulation< 3 >::dimension().
|
inline |
Definition at line 284 of file Manifold.hpp.
References foliated_triangulations::FoliatedTriangulation< 3 >::number_of_finite_edges().
|
inline |
Definition at line 268 of file Manifold.hpp.
References foliated_triangulations::FoliatedTriangulation< 3 >::number_of_finite_facets().
|
inline |
Radial separation between timeslices.
Definition at line 231 of file Manifold.hpp.
References foliated_triangulations::FoliatedTriangulation< 3 >::foliation_spacing().
|
inline |
Get a cell handle from 4 vertex handles.
| vh1 | The first vertex handle |
| vh2 | The second vertex handle |
| vh3 | The third vertex handle |
| vh4 | The fourth vertex handle |
Definition at line 419 of file Manifold.hpp.
References foliated_triangulations::FoliatedTriangulation< 3 >::get_delaunay().
Referenced by GIVEN().
|
inlinenoexcept |
Definition at line 152 of file Manifold.hpp.
Referenced by GIVEN().
|
inline |
Definition at line 164 of file Manifold.hpp.
|
inline |
Call triangulation.get_spacelike_edges()
Definition at line 333 of file Manifold.hpp.
References foliated_triangulations::FoliatedTriangulation< 3 >::get_spacelike_edges().
|
inlinenoexcept |
Call to triangulation_.get_timelike_edges()
Definition at line 327 of file Manifold.hpp.
References foliated_triangulations::FoliatedTriangulation< 3 >::get_timelike_edges().
|
inlinenoexcept |
Definition at line 145 of file Manifold.hpp.
|
inline |
Obtains a vertex handle from a point.
| point | The point to search for |
Definition at line 406 of file Manifold.hpp.
References foliated_triangulations::FoliatedTriangulation< 3 >::get_delaunay().
Referenced by GIVEN().
|
inlinenoexcept |
Call FoliatedTriangulation_3.get_vertices()
Definition at line 339 of file Manifold.hpp.
References foliated_triangulations::FoliatedTriangulation< 3 >::get_vertices().
|
inlinenoexcept |
Call FoliatedTriangulation_3.get_vertices_span()
Definition at line 345 of file Manifold.hpp.
References foliated_triangulations::FoliatedTriangulation< 3 >::get_vertices_span().
|
inlinenoexcept |
Perfect forwarding to FoliatedTriangulation_3.incident_cells()
Definition at line 320 of file Manifold.hpp.
References foliated_triangulations::FoliatedTriangulation< 3 >::incident_cells().
|
inline |
Initial radius of the first timeslice.
Definition at line 225 of file Manifold.hpp.
References foliated_triangulations::FoliatedTriangulation< 3 >::initial_radius().
|
inline |
Definition at line 191 of file Manifold.hpp.
References foliated_triangulations::FoliatedTriangulation< 3 >::is_correct().
Referenced by GIVEN().
|
inline |
Forwarding to FoliatedTriangulation.is_delaunay()
Definition at line 178 of file Manifold.hpp.
References foliated_triangulations::FoliatedTriangulation< 3 >::is_delaunay().
|
inlinenoexcept |
Forwarding to FoliatedTriangulation_3.is_edge()
| t_edge_candidate | The edge to test |
Definition at line 210 of file Manifold.hpp.
References foliated_triangulations::FoliatedTriangulation< 3 >::get_delaunay().
|
inline |
Forwarding to FoliatedTriangulation_3.is_foliated()
Definition at line 171 of file Manifold.hpp.
References foliated_triangulations::FoliatedTriangulation< 3 >::is_foliated().
|
inline |
Forwarding to FoliatedTriangulation.is_tds_valid()
Definition at line 185 of file Manifold.hpp.
References foliated_triangulations::FoliatedTriangulation< 3 >::is_tds_valid().
|
inline |
Perfect forwarding to FoliatedTriangulation_3.is_vertex()
| VertexType | The vertex type |
| t_vertex_candidate | The vertex to check |
Definition at line 201 of file Manifold.hpp.
References foliated_triangulations::FoliatedTriangulation< 3 >::get_delaunay().
|
inline |
Definition at line 306 of file Manifold.hpp.
References foliated_triangulations::FoliatedTriangulation< 3 >::max_time().
|
inline |
Definition at line 300 of file Manifold.hpp.
References foliated_triangulations::FoliatedTriangulation< 3 >::min_time().
|
inline |
Definition at line 291 of file Manifold.hpp.
References Geometry< 3 >::N0.
|
inline |
Definition at line 275 of file Manifold.hpp.
References Geometry< 3 >::N1.
|
inline |
Definition at line 278 of file Manifold.hpp.
References foliated_triangulations::FoliatedTriangulation< 3 >::N1_SL().
|
inline |
Definition at line 281 of file Manifold.hpp.
References foliated_triangulations::FoliatedTriangulation< 3 >::N1_TL().
|
inline |
Definition at line 258 of file Manifold.hpp.
References Geometry< 3 >::N2.
|
inline |
Definition at line 262 of file Manifold.hpp.
References foliated_triangulations::FoliatedTriangulation< 3 >::N2_SL().
|
inline |
Definition at line 237 of file Manifold.hpp.
References Geometry< 3 >::N3.
|
inline |
Definition at line 246 of file Manifold.hpp.
References Geometry< 3 >::N3_13.
|
inline |
Definition at line 243 of file Manifold.hpp.
References Geometry< 3 >::N3_22.
|
inline |
Definition at line 240 of file Manifold.hpp.
References Geometry< 3 >::N3_31.
|
inline |
Definition at line 249 of file Manifold.hpp.
References Geometry< 3 >::N3_31_13.
|
inline |
|
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().
|
inline |
Print details of the manifold.
Definition at line 387 of file Manifold.hpp.
|
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().
|
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().
|
inline |
Definition at line 252 of file Manifold.hpp.
References foliated_triangulations::FoliatedTriangulation< 3 >::get_cells().
|
inline |
Definition at line 158 of file Manifold.hpp.
|
inline |
Update the Manifold data structures.
Definition at line 130 of file Manifold.hpp.
|
inlineprivatenoexcept |
Update geometry data of the manifold when the triangulation has been changed.
Definition at line 447 of file Manifold.hpp.
|
inlineprivate |
Update the triangulation.
Definition at line 430 of file Manifold.hpp.
References foliated_triangulations::FoliatedTriangulation< 3 >::get_delaunay().
|
inline |
Definition at line 294 of file Manifold.hpp.
References foliated_triangulations::FoliatedTriangulation< 3 >::number_of_vertices().
Non-member swap function for Manifolds.
Used for no-except updates of manifolds after moves.
| swap_from | The value to be swapped from. Assumed to be discarded. |
| swap_into | The value to be swapped into. |
Definition at line 83 of file Manifold.hpp.
|
staticconstexpr |
Dimensionality of the manifold.
Used to determine the manifold dimension at compile-time
Definition at line 56 of file Manifold.hpp.
|
private |
The data structure of scalar values for computations.
Definition at line 51 of file Manifold.hpp.
|
private |
The data structure of geometric and combinatorial relationships.
Definition at line 48 of file Manifold.hpp.
|
staticconstexpr |
Topology of the manifold.
Definition at line 59 of file Manifold.hpp.