CDT++
Causal Dynamical Triangulations in C++
Loading...
Searching...
No Matches
Ergodic_moves_3_test.cpp
Go to the documentation of this file.
1/*******************************************************************************
2 Causal Dynamical Triangulations in C++ using CGAL
3
4 Copyright © 2019 Adam Getchell
5 ******************************************************************************/
6
11
12#include "Ergodic_moves_3.hpp"
13
14#include <doctest/doctest.h>
15
16#include <numbers>
17
18using namespace std;
19using namespace manifolds;
20
21static inline std::floating_point auto constexpr RADIUS_2 =
22 2.0 * std::numbers::inv_sqrt3_v<double>;
23static inline std::floating_point auto constexpr SQRT_2 =
24 std::numbers::sqrt2_v<double>;
25static inline std::floating_point auto constexpr INV_SQRT_2 = 1 / SQRT_2;
26
27SCENARIO("Use check_move to validate successful move" *
28 doctest::test_suite("ergodic"))
29{
30 spdlog::debug("Use check_move to validate successful move.\n");
31 GIVEN("A triangulation setup for a (2,3) move")
32 {
33 vector vertices{
34 Point_t<3>{ 1, 0, 0},
35 Point_t<3>{ 0, 1, 0},
36 Point_t<3>{ 0, 0, 1},
37 Point_t<3>{RADIUS_2, RADIUS_2, RADIUS_2},
38 Point_t<3>{ SQRT_2, SQRT_2, 0}
39 };
40 vector<size_t> timevalues{1, 1, 1, 2, 2};
41 auto causal_vertices = make_causal_vertices<3>(vertices, timevalues);
42 Manifold_3 manifold(causal_vertices);
43 WHEN("A correct (2,3) move is performed.")
44 {
45 spdlog::debug("When a correct (2,3) move is performed.\n");
46 // Copy manifold
47 auto manifold_before = manifold;
48 if (auto result = ergodic_moves::do_23_move(manifold); result)
49 {
50 manifold = result.value();
51 // Update geometry with new triangulation
52 manifold.update();
53 }
54 else
55 {
56 spdlog::debug("{}", result.error());
57 REQUIRE(result.has_value());
58 }
59 // Manual check
60 REQUIRE(manifold.is_correct());
61 CHECK_EQ(manifold.vertices(), 5);
62 CHECK_EQ(manifold.edges(), 10); // +1 timelike edge
63 CHECK_EQ(manifold.faces(), 9); // +2 faces
64 CHECK_EQ(manifold.simplices(), 3); // +1 (2,2) simplex
65 CHECK_EQ(manifold.N3_31(), 1);
66 CHECK_EQ(manifold.N3_22(), 2);
67 CHECK_EQ(manifold.N1_SL(), 4);
68 CHECK_EQ(manifold.N1_TL(), 6);
69 // CHECK(manifold.is_delaunay());
70 THEN("check_move returns true")
71 {
72 CHECK(ergodic_moves::check_move(manifold_before, manifold,
73 move_tracker::move_type::TWO_THREE));
74 }
75 }
76 }
77}
78
80 "Perform ergodic moves on the minimal manifold necessary for (2,3) and (3,2) moves" *
81 doctest::test_suite("ergodic"))
82{
83 spdlog::debug(
84 "Perform ergodic moves on the minimal manifold necessary for (2,3) and (3,2) moves.\n");
85 GIVEN("A triangulation setup for (2,3) moves")
86 {
87 vector vertices{
88 Point_t<3>{ 1, 0, 0},
89 Point_t<3>{ 0, 1, 0},
90 Point_t<3>{ 0, 0, 1},
91 Point_t<3>{RADIUS_2, RADIUS_2, RADIUS_2},
92 Point_t<3>{ SQRT_2, SQRT_2, 0}
93 };
94 vector<size_t> timevalues{1, 1, 1, 2, 2};
95 auto causal_vertices = make_causal_vertices<3>(vertices, timevalues);
96 Manifold_3 manifold(causal_vertices);
97
98 REQUIRE(manifold.is_correct());
99 REQUIRE_EQ(manifold.vertices(), 5);
100 REQUIRE_EQ(manifold.edges(), 9);
101 REQUIRE_EQ(manifold.faces(), 7);
102 REQUIRE_EQ(manifold.simplices(), 2);
103 REQUIRE_EQ(manifold.N3_31(), 1);
104 REQUIRE_EQ(manifold.N3_22(), 1);
105 REQUIRE_EQ(manifold.N1_SL(), 4);
106 REQUIRE_EQ(manifold.N1_TL(), 5);
107 REQUIRE(manifold.is_delaunay());
108 WHEN("A (2,3) move is performed")
109 {
110 spdlog::debug("When a (2,3) move is performed.\n");
111 // Copy manifold
112 auto manifold_before = manifold;
113 if (auto result = ergodic_moves::do_23_move(manifold); result)
114 {
115 manifold = result.value();
116 // Update geometry with new triangulation info
117 manifold.update();
118 }
119 else
120 {
121 spdlog::debug("{}", result.error());
122 // Stop further tests
123 REQUIRE(result.has_value());
124 }
125 THEN("The move is correct and the manifold invariants are maintained")
126 {
127 CHECK(ergodic_moves::check_move(manifold_before, manifold,
128 move_tracker::move_type::TWO_THREE));
129 // Manual check
130 REQUIRE(manifold.is_correct());
131 CHECK_EQ(manifold.vertices(), 5);
132 CHECK_EQ(manifold.edges(), 10); // +1 timelike edge
133 CHECK_EQ(manifold.faces(), 9); // +2 faces
134 CHECK_EQ(manifold.simplices(), 3); // +1 (2,2) simplex
135 CHECK_EQ(manifold.N3_31(), 1);
136 CHECK_EQ(manifold.N3_22(), 2);
137 CHECK_EQ(manifold.N1_SL(), 4);
138 CHECK_EQ(manifold.N1_TL(), 6);
139 // CHECK(manifold.is_delaunay());
140 // Human-readable output
141 manifold.print_details();
142 manifold.print_cells();
143 }
144 }
145 WHEN("A (3,2) move is performed")
146 {
147 spdlog::debug("When a (3,2) move is performed.\n");
148 // First, do a (2,3) move to set up the manifold
149 if (auto start = ergodic_moves::do_23_move(manifold); start)
150 {
151 manifold = start.value();
152 // Update geometry with new triangulation info
153 manifold.update();
154 }
155 else
156 {
157 spdlog::debug(
158 "The (2,3) move to set up the manifold for the (3,2) move "
159 "failed.\n");
160 // Stop further tests
161 REQUIRE(start.has_value());
162 }
163 // Verify we have 1 (3,1) simplex and 2 (2,2) simplices, etc.
164 REQUIRE_EQ(manifold.vertices(), 5);
165 REQUIRE_EQ(manifold.edges(), 10);
166 REQUIRE_EQ(manifold.faces(), 9);
167 REQUIRE_EQ(manifold.simplices(), 3);
168 REQUIRE_EQ(manifold.N3_31(), 1);
169 REQUIRE_EQ(manifold.N3_22(), 2);
170 REQUIRE_EQ(manifold.N1_SL(), 4);
171 REQUIRE_EQ(manifold.N1_TL(), 6);
172
173 // Copy manifold
174 auto manifold_before = manifold;
175 // Do move
176 if (auto result = ergodic_moves::do_32_move(manifold); result)
177 {
178 manifold = result.value();
179 // Update geometry with new triangulation info
180 manifold.update();
181 }
182 else
183 {
184 spdlog::debug("{}", result.error());
185 // Stop further tests
186 REQUIRE(result.has_value());
187 }
188 THEN("The move is correct and the manifold invariants are maintained")
189 {
190 CHECK(ergodic_moves::check_move(manifold_before, manifold,
191 move_tracker::move_type::THREE_TWO));
192 // Manual check
193 REQUIRE(manifold.is_correct());
194 CHECK_EQ(manifold.vertices(), 5);
195 CHECK_EQ(manifold.edges(), 9);
196 CHECK_EQ(manifold.faces(), 7);
197 CHECK_EQ(manifold.simplices(), 2);
198 CHECK_EQ(manifold.N3_31(), 1);
199 CHECK_EQ(manifold.N3_22(), 1);
200 CHECK_EQ(manifold.N1_SL(), 4);
201 CHECK_EQ(manifold.N1_TL(), 5);
202 CHECK(manifold.is_delaunay());
203 // Human-readable output
204 manifold.print_details();
205 manifold.print_cells();
206 }
207 }
208 WHEN("An improperly prepared (3,2) move is performed")
209 {
210 auto result = ergodic_moves::do_32_move(manifold);
211 THEN("The move is not performed")
212 {
213 CHECK_FALSE(result);
214 CHECK_EQ(result.error(), "No (3,2) move possible.\n");
215 }
216 }
217 }
218}
220 "Perform ergodic moves on the minimal manifold necessary (2,6) and (6,2) moves" *
221 doctest::test_suite("ergodic"))
222{
223 spdlog::debug(
224 "Perform ergodic moves on the minimal manifold necessary (2,6) and (6,2) moves.\n");
225 GIVEN("A triangulation setup for a (2,6) move")
226 {
227 vector vertices{
228 Point_t<3>{ 0, 0, 0},
229 Point_t<3>{ 1, 0, 0},
230 Point_t<3>{ 0, 1, 0},
231 Point_t<3>{ 0, 0, 1},
232 Point_t<3>{RADIUS_2, RADIUS_2, RADIUS_2}
233 };
234 vector<size_t> timevalues{0, 1, 1, 1, 2};
235 auto causal_vertices = make_causal_vertices<3>(vertices, timevalues);
236 Manifold_3 manifold(causal_vertices);
237
238 REQUIRE(manifold.is_correct());
239 REQUIRE_EQ(manifold.vertices(), 5);
240 REQUIRE_EQ(manifold.edges(), 9);
241 REQUIRE_EQ(manifold.faces(), 7);
242 REQUIRE_EQ(manifold.simplices(), 2);
243 REQUIRE_EQ(manifold.N3_31(), 1);
244 REQUIRE_EQ(manifold.N3_22(), 0);
245 REQUIRE_EQ(manifold.N3_13(), 1);
246 REQUIRE_EQ(manifold.N3_31_13(), 2);
247 REQUIRE_EQ(manifold.N1_SL(), 3);
248 REQUIRE_EQ(manifold.N1_TL(), 6);
249 REQUIRE(manifold.is_delaunay());
250 WHEN("A (2,6) move is performed")
251 {
252 spdlog::debug("When a (2,6) move is performed.\n");
253 // Copy manifold
254 auto manifold_before = manifold;
255 // Do move and check results
256 if (auto result = ergodic_moves::do_26_move(manifold); result)
257 {
258 manifold = result.value();
259 // Update geometry with new triangulation info
260 manifold.update();
261 }
262 else
263 {
264 spdlog::debug("The (2,6) move failed.\n");
265 // Stop further tests
266 REQUIRE(result.has_value());
267 }
268 THEN("The move is correct and the manifold invariants are maintained")
269 {
270 CHECK(ergodic_moves::check_move(manifold_before, manifold,
271 move_tracker::move_type::TWO_SIX));
272 // Manual check
273 REQUIRE(manifold.is_correct());
274 CHECK_EQ(manifold.vertices(), 6); // +1 vertex
275 CHECK_EQ(manifold.edges(), 14); // +3 spacelike and +2 timelike edges
276 CHECK_EQ(manifold.faces(), 15); // +8 faces
277 CHECK_EQ(manifold.simplices(), 6); // +2 (3,1) and +2 (1,3) simplices
278 CHECK_EQ(manifold.N3_31(), 3);
279 CHECK_EQ(manifold.N3_22(), 0);
280 CHECK_EQ(manifold.N3_13(), 3);
281 CHECK_EQ(manifold.N3_31_13(), 6);
282 CHECK_EQ(manifold.N1_SL(), 6); // +3 spacelike edges
283 CHECK_EQ(manifold.N1_TL(), 8); // +2 timelike edges
284 CHECK(manifold.is_delaunay());
285 // Human-readable output
286 fmt::print("Manifold before (2,6):\n");
287 manifold_before.print_details();
288 manifold_before.print_cells();
289 fmt::print("Manifold after (2,6):\n");
290 manifold.print_details();
291 manifold.print_cells();
292 }
293 }
294 WHEN("A (6,2) move is performed")
295 {
296 spdlog::debug("When a (6,2) move is performed.\n");
297 // First, do a (2,6) move to set up the manifold
298 if (auto start = ergodic_moves::do_26_move(manifold); start)
299 {
300 manifold = start.value();
301 // Update geometry with new triangulation info
302 manifold.update();
303 }
304 else
305 {
306 spdlog::debug(
307 "The (2,6) move to set up the manifold for the (6,2) move "
308 "failed.\n");
309 // Stop further tests
310 REQUIRE(start.has_value());
311 }
312 // Verify we have 3 (3,1) simplices and 3 (1,3) simplices, etc.
313 REQUIRE_EQ(manifold.vertices(), 6);
314 REQUIRE_EQ(manifold.edges(), 14);
315 REQUIRE_EQ(manifold.faces(), 15);
316 REQUIRE_EQ(manifold.simplices(), 6);
317 REQUIRE_EQ(manifold.N3_31(), 3);
318 REQUIRE_EQ(manifold.N3_22(), 0);
319 REQUIRE_EQ(manifold.N3_13(), 3);
320 REQUIRE_EQ(manifold.N3_31_13(), 6);
321 REQUIRE_EQ(manifold.N1_SL(), 6);
322 REQUIRE_EQ(manifold.N1_TL(), 8);
323 REQUIRE(manifold.is_delaunay());
324
325 // Copy manifold
326 auto manifold_before = manifold;
327 // Do move and check results
328 if (auto result = ergodic_moves::do_62_move(manifold); result)
329 {
330 manifold.update();
331 }
332 else { spdlog::info("The (6,2) move failed.\n"); }
333 THEN("The move is correct and the manifold invariants are maintained")
334 {
335 // Check the move
336 CHECK(ergodic_moves::check_move(manifold_before, manifold,
337 move_tracker::move_type::SIX_TWO));
338 // Manual check
339 REQUIRE(manifold.is_correct());
340 CHECK(manifold.get_triangulation().is_foliated());
341 CHECK(manifold.get_triangulation().is_tds_valid());
342 CHECK(manifold.get_triangulation().check_all_cells());
343 CHECK_EQ(manifold.vertices(), 5);
344 CHECK_EQ(manifold.edges(), 9);
345 CHECK_EQ(manifold.faces(), 7);
346 CHECK_EQ(manifold.simplices(), 2);
347 CHECK_EQ(manifold.N3_31(), 1);
348 CHECK_EQ(manifold.N3_22(), 0);
349 CHECK_EQ(manifold.N3_13(), 1);
350 CHECK_EQ(manifold.N3_31_13(), 2);
351 CHECK_EQ(manifold.N1_SL(), 3);
352 CHECK_EQ(manifold.N1_TL(), 6);
353 CHECK(manifold.is_delaunay());
354 // Human-readable output
355 fmt::print("Manifold before (6,2):\n");
356 manifold_before.print_details();
357 manifold_before.print_cells();
358 fmt::print("Manifold after (6,2):\n");
359 manifold.print_details();
360 manifold.print_cells();
361 }
362 }
363 WHEN("An improperly prepared (6,2) move is performed")
364 {
365 auto result = ergodic_moves::do_62_move(manifold);
366 THEN("The move is not performed")
367 {
368 CHECK_FALSE(result);
369 CHECK_EQ(result.error(), "No (6,2) move possible.\n");
370 }
371 }
372 }
373}
374SCENARIO("Perform ergodic moves on the minimal manifold necessary (4,4) moves" *
375 doctest::test_suite("ergodic") * doctest::skip())
376{
377 spdlog::debug(
378 "Perform ergodic moves on the minimal manifold necessary (4,4) moves.\n");
379 GIVEN("A triangulation setup for a (4,4) move")
380 {
381 vector vertices{
382 Point_t<3>{ 0, 0, 0},
383 Point_t<3>{ INV_SQRT_2, 0, INV_SQRT_2},
384 Point_t<3>{ 0, INV_SQRT_2, INV_SQRT_2},
385 Point_t<3>{-INV_SQRT_2, 0, INV_SQRT_2},
386 Point_t<3>{ 0, -INV_SQRT_2, INV_SQRT_2},
387 Point_t<3>{ 0, 0, 2}
388 };
389 vector<size_t> timevalues{0, 1, 1, 1, 1, 2};
390 auto causal_vertices = make_causal_vertices<3>(vertices, timevalues);
391 Manifold_3 manifold(causal_vertices, 0, 1);
392 // Verify we have 4 vertices, 4 edges, 4 faces, and 4 simplices
393 REQUIRE_EQ(manifold.vertices(), 6);
394 REQUIRE_EQ(manifold.edges(), 13);
395 REQUIRE_EQ(manifold.faces(), 12);
396 REQUIRE_EQ(manifold.simplices(), 4);
397 REQUIRE_EQ(manifold.N3_31(), 2);
398 REQUIRE_EQ(manifold.N3_22(), 0);
399 REQUIRE_EQ(manifold.N3_13(), 2);
400 REQUIRE_EQ(manifold.N3_31_13(), 4);
401 REQUIRE_EQ(manifold.N1_SL(), 5);
402 REQUIRE_EQ(manifold.N1_TL(), 8);
403 CHECK_EQ(manifold.initial_radius(), 0);
404 CHECK_EQ(manifold.foliation_spacing(), 1);
405 REQUIRE(manifold.is_delaunay());
406 REQUIRE(manifold.is_correct());
407
408 WHEN("A (4,4) move is performed")
409 {
410 spdlog::debug("When a (4,4) move is performed.\n");
411 // Copy manifold
412 auto manifold_before = manifold;
413 // Human verification
414 fmt::print("Manifold before (4,4):\n");
415 manifold_before.print_details();
416 manifold_before.print_cells();
417 // Do move and check results
418 if (auto result = ergodic_moves::do_44_move(manifold); result)
419 {
420 manifold.update();
421 }
422 else { spdlog::info("The (4,4) move failed.\n"); }
423 fmt::print("Manifold after (4,4):\n");
424 manifold.print_details();
425 manifold.print_cells();
426 THEN("The move is correct and the manifold invariants are maintained")
427 {
428 // Check the move
429 CHECK(ergodic_moves::check_move(manifold_before, manifold,
430 move_tracker::move_type::FOUR_FOUR));
431 }
432 }
433 }
434}
435
436SCENARIO("Test convenience functions needed for bistellar flip" *
437 doctest::test_suite("ergodic"))
438{
439 GIVEN("A triangulation setup for a bistellar flip")
440 {
441 vector vertices{
442 Point_t<3>{ 0, 0, 0},
443 Point_t<3>{ INV_SQRT_2, 0, INV_SQRT_2},
444 Point_t<3>{ 0, INV_SQRT_2, INV_SQRT_2},
445 Point_t<3>{-INV_SQRT_2, 0, INV_SQRT_2},
446 Point_t<3>{ 0, -INV_SQRT_2, INV_SQRT_2},
447 Point_t<3>{ 0, 0, 2}
448 };
449 ergodic_moves::Delaunay triangulation(vertices.begin(), vertices.end());
450 CHECK(triangulation.is_valid());
451 auto edges = foliated_triangulations::collect_edges<3>(triangulation);
452 WHEN("We get all the finite cells in the triangulation")
453 {
454 auto cells = foliated_triangulations::collect_cells<3>(triangulation);
455 THEN("We have 4 cells") { REQUIRE_EQ(cells.size(), 4); }
456 }
457 WHEN("We get all finite edges in the triangulation")
458 {
459 THEN("We have 13 edges") { REQUIRE_EQ(edges.size(), 13); }
460 }
461 WHEN("We find the pivot edge in the triangulation")
462 {
463 auto pivot_edge = ergodic_moves::find_pivot_edge(triangulation, edges);
464 REQUIRE_MESSAGE(pivot_edge, "No pivot edge found.");
465
466 auto Contains = [&vertices](Point_t<3> point) {
467 return ranges::any_of(vertices, [&point](Point_t<3> const& test) {
468 return test == point;
469 });
470 };
471
472 if (pivot_edge)
473 {
474 auto incident_cells = ergodic_moves::get_incident_cells(
475 triangulation, pivot_edge.value());
476 REQUIRE_MESSAGE(incident_cells, "No incident cells found.");
477 THEN("We have a pivot edge")
478 {
479 CHECK_MESSAGE(pivot_edge, "Pivot edge found");
480 REQUIRE(triangulation.tds().is_edge(
481 pivot_edge->first, pivot_edge->second, pivot_edge->third));
482 auto pivot_from_1 =
483 pivot_edge->first->vertex(pivot_edge->second)->point();
484 auto pivot_from_2 =
485 pivot_edge->first->vertex(pivot_edge->third)->point();
486 // Verify Contains
487 REQUIRE_FALSE(Contains(Point_t<3>{0, 0, 1}));
488 REQUIRE(Contains(pivot_from_1));
489 REQUIRE(Contains(pivot_from_2));
490 }
491 if (incident_cells)
492 {
493 THEN("We can obtain the cells incident to that edge")
494 {
495 REQUIRE_EQ(incident_cells->size(), 4);
496 }
497 AND_THEN("We can obtain the vertices from the incident cells")
498 {
499 auto incident_vertices =
500 ergodic_moves::get_vertices(incident_cells.value());
501 REQUIRE_EQ(incident_vertices.size(), 6);
502 }
503 }
504 }
505 }
506 WHEN("We get all finite vertices in the triangulation")
507 {
508 THEN("We have 6 vertices")
509 {
510 auto all_finite_vertices =
511 foliated_triangulations::collect_vertices<3>(triangulation);
512 REQUIRE_EQ(all_finite_vertices.size(), 6);
513 }
514 }
515 }
516}
517
518SCENARIO("Perform bistellar flip on Delaunay triangulation" *
519 doctest::test_suite("ergodic"))
520{
521 GIVEN("A triangulation setup for a bistellar flip")
522 {
523 vector vertices{
524 Point_t<3>{ 0, 0, 0},
525 Point_t<3>{ INV_SQRT_2, 0, INV_SQRT_2},
526 Point_t<3>{ 0, INV_SQRT_2, INV_SQRT_2},
527 Point_t<3>{-INV_SQRT_2, 0, INV_SQRT_2},
528 Point_t<3>{ 0, -INV_SQRT_2, INV_SQRT_2},
529 Point_t<3>{ 0, 0, 2}
530 };
531 ergodic_moves::Delaunay triangulation(vertices.begin(), vertices.end());
532 WHEN("We have a valid triangulation")
533 {
534 CHECK(triangulation.is_valid());
535 THEN("We can perform a bistellar flip")
536 {
537 // Obtain top and bottom vertices by re-inserting, which returns the
538 // Vertex_handle
539 auto top = triangulation.insert(Point_t<3>{0, 0, 2});
540 auto bottom = triangulation.insert(Point_t<3>{0, 0, 0});
541 auto edges = foliated_triangulations::collect_edges<3>(triangulation);
542 auto pivot_edge = ergodic_moves::find_pivot_edge(triangulation, edges);
543 REQUIRE_MESSAGE(pivot_edge, "No pivot edge found.");
544
545 // Check this didn't actually change vertices in the triangulation
546 REQUIRE_EQ(vertices.size(), 6);
547
548 if (pivot_edge)
549 {
550 auto flipped_triangulation = ergodic_moves::bistellar_flip(
551 triangulation, pivot_edge.value(), top, bottom);
552
553 REQUIRE_MESSAGE(flipped_triangulation, "Bistellar flip failed.");
554 if (flipped_triangulation)
555 {
558 WARN(flipped_triangulation->is_valid());
559 }
560 }
561 }
562 }
563 }
564}
Pachner moves on 2+1 dimensional foliated Delaunay triangulations.
auto do_32_move(Manifold &t_manifold) -> Expected
Perform a (3,2) move.
auto get_incident_cells(Delaunay const &triangulation, Edge_handle const edge) -> std::optional< Cell_container >
Return a container of cells incident to an edge.
auto do_62_move(Manifold &t_manifold) -> Expected
Perform a (6,2) move.
auto get_vertices(Cell_container const &cells)
Return a container of all vertices in a container of cells.
auto find_pivot_edge(Delaunay const &triangulation, Edge_container const &edges) -> std::optional< Edge_handle >
SCENARIO("Perform bistellar flip on Delaunay triangulation" *doctest::test_suite("bistellar"))