11#ifndef CDT_PLUSPLUS_MANIFOLD_HPP
12#define CDT_PLUSPLUS_MANIFOLD_HPP
15#include <unordered_set>
26 template <
int dimension>
27 auto make_causal_vertices(std::span<Point_t<dimension>
const> vertices,
28 std::span<size_t const>
const timevalues)
29 -> Causal_vertices_t<dimension>
31 return foliated_triangulations::make_causal_vertices<dimension>(vertices,
37 template <
int dimension>
42 class [[nodiscard(
"This contains data!")]]
Manifold<3>
56 static int constexpr dimension = 3;
86 spdlog::debug(
"{} called.\n", __PRETTY_FUNCTION__);
89 swap(swap_from.m_triangulation, swap_into.m_triangulation);
90 swap(swap_from.m_geometry, swap_into.m_geometry);
96 : m_triangulation{std::move(t_foliated_triangulation)}
97 , m_geometry{get_triangulation()}
108 double const t_foliation_spacing = FOLIATION_SPACING)
111 t_initial_radius, t_foliation_spacing}
120 explicit Manifold(Causal_vertices_t<3>
const& causal_vertices,
122 double const t_foliation_spacing = FOLIATION_SPACING)
134 spdlog::debug(
"{} called.\n", __PRETTY_FUNCTION__);
136 update_triangulation();
139 catch (std::system_error
const& ex)
141 spdlog::trace(
"Exception thrown: {}\n", ex.what());
148 return std::cref(m_triangulation);
154 return get_triangulation().get_delaunay();
160 return m_triangulation;
200 template <
typename VertexType>
201 [[nodiscard]]
auto is_vertex(VertexType&& t_vertex_candidate)
const ->
bool
204 std::forward<VertexType>(t_vertex_candidate));
211 Edge_handle_t<3>
const& t_edge_candidate)
const noexcept ->
bool
214 t_edge_candidate.first, t_edge_candidate.second,
215 t_edge_candidate.third);
237 [[nodiscard]]
auto N3()
const {
return m_geometry.
N3; }
240 [[nodiscard]]
auto N3_31()
const {
return m_geometry.
N3_31; }
243 [[nodiscard]]
auto N3_22()
const {
return m_geometry.
N3_22; }
246 [[nodiscard]]
auto N3_13()
const {
return m_geometry.
N3_13; }
258 [[nodiscard]]
auto N2()
const {
return m_geometry.
N2; }
262 [[nodiscard]]
auto N2_SL() const -> auto const&
264 return m_triangulation.
N2_SL();
275 [[nodiscard]]
auto N1()
const {
return m_geometry.
N1; }
278 [[nodiscard]]
auto N1_SL()
const {
return m_triangulation.
N1_SL(); }
281 [[nodiscard]]
auto N1_TL()
const {
return m_triangulation.
N1_TL(); }
291 [[nodiscard]]
auto N0()
const {
return m_geometry.
N0; }
312 template <
typename VertexHandle>
313 [[nodiscard]]
auto degree(VertexHandle&& t_vertex)
const ->
decltype(
auto)
315 return m_triangulation.
degree(std::forward<VertexHandle>(t_vertex));
319 template <
typename... Ts>
354 return this->simplices() == this->N3() &&
376 "Manifold has {} vertices and {} edges and {} faces and {} "
378 this->N0(), this->
N1(), this->N2(), this->N3());
382 fmt::print(stderr,
"print() went wrong ...\n");
391 "There are {} (3,1) simplices and {} (2,2) simplices and {} (1,3) "
393 this->N3_31(), this->N3_22(), this->N3_13());
394 fmt::print(
"There are {} timelike edges and {} spacelike edges.\n",
395 this->N1_TL(), this->N1_SL());
399 fmt::print(stderr,
"print_details() went wrong ...\n");
406 auto get_vertex(Point_t<3>
const& point)
const -> Vertex_handle_t<3>
408 Vertex_handle_t<3> result;
409 m_triangulation.
get_delaunay().is_vertex(point, result);
419 auto get_cell(Vertex_handle_t<3>
const& vh1, Vertex_handle_t<3>
const& vh2,
420 Vertex_handle_t<3>
const& vh3,
421 Vertex_handle_t<3>
const& vh4)
const -> Cell_handle_t<3>
423 Cell_handle_t<3> result;
424 m_triangulation.
get_delaunay().is_cell(vh1, vh2, vh3, vh4, result);
434 spdlog::debug(
"{} called.\n", __PRETTY_FUNCTION__);
438 swap(local_triangulation, m_triangulation);
440 catch (std::system_error
const& ex)
442 fmt::print(
"Exception thrown: {}\n", ex.what());
450 spdlog::debug(
"{} called.\n", __PRETTY_FUNCTION__);
454 swap(geom, m_geometry);
458 using Manifold_3 = Manifold<3>;
462 class [[nodiscard(
"This contains data!")]]
Manifold<4>
470 static int constexpr dimension = 4;
Geometric scalars of the Manifold used to calculate the Regge action.
std::int_fast32_t Int_precision
static double constexpr INITIAL_RADIUS
Default foliated triangulation spacings.
topology_type
clang-15 does not support std::format
3D Foliated triangulation
void print_vertices() const
Print values of a vertex.
auto is_foliated() const -> bool
Verifies the triangulation is properly foliated.
auto get_vertices_span() const noexcept -> std::span< Vertex_handle const >
auto initial_radius() const
auto degree(VertexHandle &&t_vertex) const
Perfect forwarding to TriangulationDataStructure_3::degree.
auto N2_SL() const -> std::multimap< Int_precision, TriangulationTraits< 3 >::Facet > const &
auto get_spacelike_edges() const -> Edge_container const &
auto number_of_finite_edges() const
auto is_correct() const -> bool
void print_volume_per_timeslice() const
Print the number of spacelike faces per timeslice.
auto get_delaunay() const -> Delaunay const &
auto get_timelike_edges() const noexcept -> Edge_container const &
auto foliation_spacing() const
auto get_vertices() const noexcept -> Vertex_container const &
auto get_cells() const -> Cell_container const &
auto incident_cells(VertexHandle &&t_vh) const noexcept -> decltype(auto)
Perfect forwarding to TriangulationDataStructure_3::incident_cells.
auto is_tds_valid() const -> bool
auto check_all_cells() const -> bool
Check that all cells are correctly classified.
void print_cells() const
Print timevalues of each vertex in the cell and the resulting cell->info()
auto number_of_finite_facets() const
auto number_of_vertices() const
auto is_delaunay() const -> bool
void update_triangulation()
Update the triangulation.
Manifold(Triangulation t_foliated_triangulation)
Construct manifold from a Foliated triangulation.
void print() const
Print manifold.
auto is_edge(Edge_handle_t< 3 > const &t_edge_candidate) const noexcept -> bool
Forwarding to FoliatedTriangulation_3.is_edge()
auto get_geometry() const -> Geometry const &
auto degree(VertexHandle &&t_vertex) const -> decltype(auto)
Perfect forwarding to FoliatedTriangulation_3.degree()
auto check_simplices() const -> bool
void print_details() const
Print details of the manifold.
Manifold()=default
Default ctor.
auto get_delaunay() const noexcept
auto N2_SL() const -> auto const &
auto foliation_spacing() const
Radial separation between timeslices.
Manifold(Manifold const &other)=default
Default copy ctor.
void print_cells() const
Print timevalues of each vertex in the cell and the resulting cell->info()
auto get_triangulation() const noexcept -> Triangulation const &
friend void swap(Manifold &swap_from, Manifold &swap_into) noexcept
Non-member swap function for Manifolds.
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.
auto get_vertex(Point_t< 3 > const &point) const -> Vertex_handle_t< 3 >
Obtains a vertex handle from a point.
void update()
Update the Manifold data structures.
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 operator=(Manifold &&other) -> Manifold &=default
Default move assignment.
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.
auto get_vertices_span() const noexcept -> auto
Call FoliatedTriangulation_3.get_vertices_span()
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.
Geometry m_geometry
The data structure of scalar values for computations.
void print_vertices() const
Print values of a vertex->info()
auto triangulation() -> Triangulation &
auto dimensionality() const
auto operator=(Manifold const &other) -> Manifold &=default
Default copy assignment.
auto get_timelike_edges() const noexcept -> auto const &
Call to triangulation_.get_timelike_edges()
~Manifold()=default
Default dtor.
void print_volume_per_timeslice() const
Print the codimension 1 volume of simplices (faces) per timeslice.
void update_geometry() noexcept
Update geometry data of the manifold when the triangulation has been changed.
auto get_vertices() const noexcept -> auto const &
Call FoliatedTriangulation_3.get_vertices()
auto is_valid() const -> bool
Forwarding to FoliatedTriangulation.is_tds_valid()
auto is_foliated() const -> bool
Forwarding to FoliatedTriangulation_3.is_foliated()
auto get_spacelike_edges() const -> auto const &
Call triangulation.get_spacelike_edges()
Triangulation m_triangulation
The data structure of geometric and combinatorial relationships.
auto is_delaunay() const -> bool
Forwarding to FoliatedTriangulation.is_delaunay()
auto is_correct() const -> bool
auto is_vertex(VertexType &&t_vertex_candidate) const -> bool
Perfect forwarding to FoliatedTriangulation_3.is_vertex()
Geometry_4 m_geometry
The data structure of scalar values for computations.
Int_precision N1
Number of 1D edges.
Int_precision N3_31
Number of (3,1) simplices.
Int_precision N3_13
Number of (1,3) simplices.
Int_precision N2
Number of 2D faces.
Int_precision N3_22
Number of (2,2) simplices.
Int_precision N3_31_13
Number of (3,1) + (1,3) simplices.
Int_precision N3
Number of 3D simplices.
Int_precision N0
Number of vertices.