CDT++
Causal Dynamical Triangulations in C++
Loading...
Searching...
No Matches
Manifold.hpp
Go to the documentation of this file.
1/*******************************************************************************
2 Causal Dynamical Triangulations in C++ using CGAL
3
4 Copyright © 2018 Adam Getchell
5 ******************************************************************************/
6
10
11#ifndef CDT_PLUSPLUS_MANIFOLD_HPP
12#define CDT_PLUSPLUS_MANIFOLD_HPP
13
14#include <cstddef>
15#include <unordered_set>
16
17#include "Geometry.hpp"
18
19namespace manifolds
20{
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>
30 {
31 return foliated_triangulations::make_causal_vertices<dimension>(vertices,
32 timevalues);
33 }
34
37 template <int dimension>
38 class Manifold;
39
41 template <>
42 class [[nodiscard("This contains data!")]] Manifold<3>
43 {
45 using Geometry = Geometry_3;
46
49
52
53 public:
56 static int constexpr dimension = 3;
57
59 static topology_type constexpr topology = topology_type::SPHERICAL;
60
62 ~Manifold() = default;
63
65 Manifold() = default;
66
68 Manifold(Manifold const& other) = default;
69
71 auto operator=(Manifold const& other) -> Manifold& = default;
72
74 Manifold(Manifold&& other) = default;
75
77 auto operator=(Manifold&& other) -> Manifold& = default;
78
83 friend void swap(Manifold& swap_from, Manifold& swap_into) noexcept
84 {
85#ifndef NDEBUG
86 spdlog::debug("{} called.\n", __PRETTY_FUNCTION__);
87#endif
88 using std::swap;
89 swap(swap_from.m_triangulation, swap_into.m_triangulation);
90 swap(swap_from.m_geometry, swap_into.m_geometry);
91 } // swap
92
95 explicit Manifold(Triangulation t_foliated_triangulation)
96 : m_triangulation{std::move(t_foliated_triangulation)}
97 , m_geometry{get_triangulation()}
98 {}
99
105 Manifold(Int_precision const t_desired_simplices,
106 Int_precision const t_desired_timeslices,
107 double const t_initial_radius = INITIAL_RADIUS,
108 double const t_foliation_spacing = FOLIATION_SPACING)
109 : Manifold{
110 Triangulation{t_desired_simplices, t_desired_timeslices,
111 t_initial_radius, t_foliation_spacing}
112 }
113 {}
114
120 explicit Manifold(Causal_vertices_t<3> const& causal_vertices,
121 double const t_initial_radius = INITIAL_RADIUS,
122 double const t_foliation_spacing = FOLIATION_SPACING)
123 : Manifold{
124 Triangulation{causal_vertices, t_initial_radius,
125 t_foliation_spacing}
126 }
127 {}
128
130 void update()
131 try
132 {
133#ifndef NDEBUG
134 spdlog::debug("{} called.\n", __PRETTY_FUNCTION__);
135#endif
136 update_triangulation();
137 update_geometry();
138 }
139 catch (std::system_error const& ex)
140 {
141 spdlog::trace("Exception thrown: {}\n", ex.what());
142 } // update
143
145 [[nodiscard]] auto get_triangulation() const noexcept
146 -> Triangulation const&
147 {
148 return std::cref(m_triangulation);
149 } // get_triangulation
150
152 [[nodiscard]] auto get_delaunay() const noexcept
153 {
154 return get_triangulation().get_delaunay();
155 } // get_delaunay
156
158 [[nodiscard]] auto triangulation() -> Triangulation&
159 {
160 return m_triangulation;
161 } // triangulation
162
164 [[nodiscard]] auto get_geometry() const -> Geometry const&
165 {
166 return m_geometry;
167 } // get_geometry
168
171 [[nodiscard]] auto is_foliated() const -> bool
172 {
173 return m_triangulation.is_foliated();
174 } // is_foliated
175
178 [[nodiscard]] auto is_delaunay() const -> bool
179 {
180 return m_triangulation.is_delaunay();
181 } // is_delaunay
182
185 [[nodiscard]] auto is_valid() const -> bool
186 {
187 return m_triangulation.is_tds_valid();
188 } // is_valid
189
191 [[nodiscard]] auto is_correct() const -> bool
192 {
193 return m_triangulation.is_correct();
194 } // is_correct
195
200 template <typename VertexType>
201 [[nodiscard]] auto is_vertex(VertexType&& t_vertex_candidate) const -> bool
202 {
203 return m_triangulation.get_delaunay().is_vertex(
204 std::forward<VertexType>(t_vertex_candidate));
205 } // is_vertex
206
210 [[nodiscard]] auto is_edge(
211 Edge_handle_t<3> const& t_edge_candidate) const noexcept -> bool
212 {
213 return m_triangulation.get_delaunay().tds().is_edge(
214 t_edge_candidate.first, t_edge_candidate.second,
215 t_edge_candidate.third);
216 } // is_edge
217
219 [[nodiscard]] auto dimensionality() const
220 {
221 return m_triangulation.dimension();
222 }
223
225 [[nodiscard]] auto initial_radius() const
226 {
227 return m_triangulation.initial_radius();
228 }
229
231 [[nodiscard]] auto foliation_spacing() const
232 {
233 return m_triangulation.foliation_spacing();
234 }
235
237 [[nodiscard]] auto N3() const { return m_geometry.N3; }
238
240 [[nodiscard]] auto N3_31() const { return m_geometry.N3_31; }
241
243 [[nodiscard]] auto N3_22() const { return m_geometry.N3_22; }
244
246 [[nodiscard]] auto N3_13() const { return m_geometry.N3_13; }
247
249 [[nodiscard]] auto N3_31_13() const { return m_geometry.N3_31_13; }
250
252 [[nodiscard]] auto simplices() const
253 {
254 return static_cast<Int_precision>(m_triangulation.get_cells().size());
255 } // number_of_simplices
256
258 [[nodiscard]] auto N2() const { return m_geometry.N2; }
259
262 [[nodiscard]] auto N2_SL() const -> auto const&
263 {
264 return m_triangulation.N2_SL();
265 } // N2_SL
266
268 [[nodiscard]] auto faces() const
269 {
270 return static_cast<Int_precision>(
271 m_triangulation.number_of_finite_facets());
272 } // faces
273
275 [[nodiscard]] auto N1() const { return m_geometry.N1; }
276
278 [[nodiscard]] auto N1_SL() const { return m_triangulation.N1_SL(); }
279
281 [[nodiscard]] auto N1_TL() const { return m_triangulation.N1_TL(); }
282
284 [[nodiscard]] auto edges() const
285 {
286 return static_cast<Int_precision>(
287 m_triangulation.number_of_finite_edges());
288 } // edges
289
291 [[nodiscard]] auto N0() const { return m_geometry.N0; }
292
294 [[nodiscard]] auto vertices() const
295 {
296 return static_cast<Int_precision>(m_triangulation.number_of_vertices());
297 } // vertices
298
300 [[nodiscard]] auto min_time() const
301 {
302 return m_triangulation.min_time();
303 } // min_time
304
306 [[nodiscard]] auto max_time() const
307 {
308 return m_triangulation.max_time();
309 } // max_time
310
312 template <typename VertexHandle>
313 [[nodiscard]] auto degree(VertexHandle&& t_vertex) const -> decltype(auto)
314 {
315 return m_triangulation.degree(std::forward<VertexHandle>(t_vertex));
316 } // degree
317
319 template <typename... Ts>
320 [[nodiscard]] auto incident_cells(Ts&&... args) const noexcept
321 -> decltype(auto)
322 {
323 return m_triangulation.incident_cells(std::forward<Ts>(args)...);
324 } // incident_cells
325
327 [[nodiscard]] auto get_timelike_edges() const noexcept -> auto const&
328 {
329 return m_triangulation.get_timelike_edges();
330 } // get_timelike_edges
331
333 [[nodiscard]] auto get_spacelike_edges() const -> auto const&
334 {
335 return m_triangulation.get_spacelike_edges();
336 } // get_spacelike_edges
337
339 [[nodiscard]] auto get_vertices() const noexcept -> auto const&
340 {
341 return m_triangulation.get_vertices();
342 } // get_vertices
343
345 [[nodiscard]] auto get_vertices_span() const noexcept -> auto
346 {
347 return m_triangulation.get_vertices_span();
348 } // get_vertices_span
349
352 [[nodiscard]] auto check_simplices() const -> bool
353 {
354 return this->simplices() == this->N3() &&
355 m_triangulation.check_all_cells();
356 } // check_simplices
357
360 {
361 m_triangulation.print_volume_per_timeslice();
362 } // print_volume_per_timeslice
363
365 void print_vertices() const { m_triangulation.print_vertices(); }
366
369 void print_cells() const { m_triangulation.print_cells(); }
370
372 void print() const
373 try
374 {
375 fmt::print(
376 "Manifold has {} vertices and {} edges and {} faces and {} "
377 "simplices.\n",
378 this->N0(), this->N1(), this->N2(), this->N3());
379 }
380 catch (...)
381 {
382 fmt::print(stderr, "print() went wrong ...\n");
383 throw;
384 } // print
385
387 void print_details() const
388 try
389 {
390 fmt::print(
391 "There are {} (3,1) simplices and {} (2,2) simplices and {} (1,3) "
392 "simplices.\n",
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());
396 }
397 catch (...)
398 {
399 fmt::print(stderr, "print_details() went wrong ...\n");
400 throw;
401 } // print_details
402
406 auto get_vertex(Point_t<3> const& point) const -> Vertex_handle_t<3>
407 {
408 Vertex_handle_t<3> result;
409 m_triangulation.get_delaunay().is_vertex(point, result);
410 return result;
411 } // get_vertex
412
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>
422 {
423 Cell_handle_t<3> result;
424 m_triangulation.get_delaunay().is_cell(vh1, vh2, vh3, vh4, result);
425 return result;
426 } // get_cell_handle
427
428 private:
431 try
432 {
433#ifndef NDEBUG
434 spdlog::debug("{} called.\n", __PRETTY_FUNCTION__);
435#endif
436 // Constructing a new triangulation updates all data structures
437 Triangulation local_triangulation(m_triangulation.get_delaunay());
438 swap(local_triangulation, m_triangulation);
439 }
440 catch (std::system_error const& ex)
441 {
442 fmt::print("Exception thrown: {}\n", ex.what());
443 } // update_triangulation
444
447 void update_geometry() noexcept
448 {
449#ifndef NDEBUG
450 spdlog::debug("{} called.\n", __PRETTY_FUNCTION__);
451#endif
452 // constructing a new geometry updates all data structures
453 Geometry geom(m_triangulation);
454 swap(geom, m_geometry);
455 } // update_geometry
456 };
457
458 using Manifold_3 = Manifold<3>;
459
461 template <>
462 class [[nodiscard("This contains data!")]] Manifold<4>
463 {
466
467 public:
470 static int constexpr dimension = 4;
471
473 static topology_type constexpr topology = topology_type::SPHERICAL;
474 };
475
476 using Manifold_4 = Manifold<4>;
477
478} // namespace manifolds
479
480#endif // CDT_PLUSPLUS_MANIFOLD_HPP
Geometric scalars of the Manifold used to calculate the Regge action.
std::int_fast32_t Int_precision
Definition: Settings.hpp:31
static double constexpr INITIAL_RADIUS
Default foliated triangulation spacings.
Definition: Settings.hpp:47
topology_type
clang-15 does not support std::format
Definition: Utilities.hpp:43
auto is_foliated() const -> bool
Verifies the triangulation is properly foliated.
auto get_vertices_span() const noexcept -> std::span< Vertex_handle 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 &
void print_volume_per_timeslice() const
Print the number of spacelike faces per timeslice.
auto get_timelike_edges() const noexcept -> Edge_container const &
auto get_vertices() const noexcept -> Vertex_container const &
auto incident_cells(VertexHandle &&t_vh) const noexcept -> decltype(auto)
Perfect forwarding to TriangulationDataStructure_3::incident_cells.
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()
void update_triangulation()
Update the triangulation.
Definition: Manifold.hpp:430
Manifold(Triangulation t_foliated_triangulation)
Construct manifold from a Foliated triangulation.
Definition: Manifold.hpp:95
void print() const
Print manifold.
Definition: Manifold.hpp:372
auto is_edge(Edge_handle_t< 3 > const &t_edge_candidate) const noexcept -> bool
Forwarding to FoliatedTriangulation_3.is_edge()
Definition: Manifold.hpp:210
auto get_geometry() const -> Geometry const &
Definition: Manifold.hpp:164
auto degree(VertexHandle &&t_vertex) const -> decltype(auto)
Perfect forwarding to FoliatedTriangulation_3.degree()
Definition: Manifold.hpp:313
auto check_simplices() const -> bool
Definition: Manifold.hpp:352
void print_details() const
Print details of the manifold.
Definition: Manifold.hpp:387
Manifold()=default
Default ctor.
auto get_delaunay() const noexcept
Definition: Manifold.hpp:152
auto N2_SL() const -> auto const &
Definition: Manifold.hpp:262
auto foliation_spacing() const
Radial separation between timeslices.
Definition: Manifold.hpp:231
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()
Definition: Manifold.hpp:369
auto get_triangulation() const noexcept -> Triangulation const &
Definition: Manifold.hpp:145
friend void swap(Manifold &swap_from, Manifold &swap_into) noexcept
Non-member swap function for Manifolds.
Definition: Manifold.hpp:83
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.
Definition: Manifold.hpp:105
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.
Definition: Manifold.hpp:406
void update()
Update the Manifold data structures.
Definition: Manifold.hpp:130
auto incident_cells(Ts &&... args) const noexcept -> decltype(auto)
Perfect forwarding to FoliatedTriangulation_3.incident_cells()
Definition: Manifold.hpp:320
auto initial_radius() const
Initial radius of the first timeslice.
Definition: Manifold.hpp:225
auto simplices() const
Definition: Manifold.hpp:252
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.
Definition: Manifold.hpp:120
auto get_vertices_span() const noexcept -> auto
Call FoliatedTriangulation_3.get_vertices_span()
Definition: Manifold.hpp:345
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.
Definition: Manifold.hpp:419
Geometry m_geometry
The data structure of scalar values for computations.
Definition: Manifold.hpp:51
void print_vertices() const
Print values of a vertex->info()
Definition: Manifold.hpp:365
auto triangulation() -> Triangulation &
Definition: Manifold.hpp:158
auto dimensionality() const
Definition: Manifold.hpp:219
auto operator=(Manifold const &other) -> Manifold &=default
Default copy assignment.
auto get_timelike_edges() const noexcept -> auto const &
Call to triangulation_.get_timelike_edges()
Definition: Manifold.hpp:327
~Manifold()=default
Default dtor.
void print_volume_per_timeslice() const
Print the codimension 1 volume of simplices (faces) per timeslice.
Definition: Manifold.hpp:359
void update_geometry() noexcept
Update geometry data of the manifold when the triangulation has been changed.
Definition: Manifold.hpp:447
auto get_vertices() const noexcept -> auto const &
Call FoliatedTriangulation_3.get_vertices()
Definition: Manifold.hpp:339
auto is_valid() const -> bool
Forwarding to FoliatedTriangulation.is_tds_valid()
Definition: Manifold.hpp:185
auto is_foliated() const -> bool
Forwarding to FoliatedTriangulation_3.is_foliated()
Definition: Manifold.hpp:171
auto get_spacelike_edges() const -> auto const &
Call triangulation.get_spacelike_edges()
Definition: Manifold.hpp:333
Triangulation m_triangulation
The data structure of geometric and combinatorial relationships.
Definition: Manifold.hpp:48
auto is_delaunay() const -> bool
Forwarding to FoliatedTriangulation.is_delaunay()
Definition: Manifold.hpp:178
auto is_correct() const -> bool
Definition: Manifold.hpp:191
auto is_vertex(VertexType &&t_vertex_candidate) const -> bool
Perfect forwarding to FoliatedTriangulation_3.is_vertex()
Definition: Manifold.hpp:201
Geometry_4 m_geometry
The data structure of scalar values for computations.
Definition: Manifold.hpp:465
Definition: group.cpp:53
3D Geometry
Definition: Geometry.hpp:27
Int_precision N1
Number of 1D edges.
Definition: Geometry.hpp:47
Int_precision N3_31
Number of (3,1) simplices.
Definition: Geometry.hpp:32
Int_precision N3_13
Number of (1,3) simplices.
Definition: Geometry.hpp:35
Int_precision N2
Number of 2D faces.
Definition: Geometry.hpp:44
Int_precision N3_22
Number of (2,2) simplices.
Definition: Geometry.hpp:41
Int_precision N3_31_13
Number of (3,1) + (1,3) simplices.
Definition: Geometry.hpp:38
Int_precision N3
Number of 3D simplices.
Definition: Geometry.hpp:29
Int_precision N0
Number of vertices.
Definition: Geometry.hpp:56