CDT++
Causal Dynamical Triangulations in C++
Loading...
Searching...
No Matches
Apply_move_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 "Apply_move.hpp"
13
14#include <doctest/doctest.h>
15
16#include "Ergodic_moves_3.hpp"
17
18using namespace std;
19
20SCENARIO("Apply an ergodic move to 2+1 manifolds" *
21 doctest::test_suite("apply"))
22{
23 GIVEN("A 2+1 dimensional spherical manifold.")
24 {
25 auto constexpr desired_simplices = 9600;
26 auto constexpr desired_timeslices = 7;
27 manifolds::Manifold_3 manifold(desired_simplices, desired_timeslices);
28 REQUIRE(manifold.is_correct());
29 // Copy of manifold
30 auto manifold_before = manifold;
31 WHEN("A null move is applied to the manifold.")
32 {
33 spdlog::debug("Applying null move to manifold.\n");
34 if (auto result = apply_move(manifold, ergodic_moves::null_move); result)
35 {
36 manifold = result.value();
37 // Update Geometry and Foliated_triangulation with new info
38 manifold.update();
39 }
40 else
41 {
42 spdlog::debug("{}", result.error());
43 REQUIRE(result.has_value());
44 }
45 THEN("The resulting manifold is valid and unchanged.")
46 {
47 CHECK(manifold.is_valid());
48 CHECK_EQ(manifold_before.simplices(), manifold.simplices());
49 CHECK_EQ(manifold_before.faces(), manifold.faces());
50 CHECK_EQ(manifold_before.edges(), manifold.edges());
51 CHECK_EQ(manifold_before.vertices(), manifold.vertices());
52 // Human verification
53 fmt::print("Old manifold.\n");
54 manifold_before.print_details();
55 fmt::print("New manifold after null move:\n");
56 manifold.print_details();
57 }
58 }
59 WHEN("A (2,3) move is applied to the manifold.")
60 {
61 spdlog::debug("Applying (2,3) move to manifold.\n");
62 if (auto result = apply_move(manifold, ergodic_moves::do_23_move); result)
63 {
64 manifold = result.value();
65 // Update Geometry and Foliated_triangulation with new info
66 manifold.update();
67 }
68 else
69 {
70 spdlog::debug("{}", result.error());
71 REQUIRE(result.has_value());
72 }
73 THEN("The resulting manifold has the applied move.")
74 {
75 CHECK(ergodic_moves::check_move(manifold_before, manifold,
76 move_tracker::move_type::TWO_THREE));
77 // Human verification
78 fmt::print("Old manifold.\n");
79 manifold_before.print_details();
80 fmt::print("New manifold after (2,3) move:\n");
81 manifold.print_details();
82 }
83 }
84 WHEN("A (3,2) move is applied to the manifold.")
85 {
86 spdlog::debug("Applying (3,2) move to manifold.\n");
87 if (auto result = apply_move(manifold, ergodic_moves::do_32_move); result)
88 {
89 manifold = result.value();
90 // Update Geometry and Foliated_triangulation with new info
91 manifold.update();
92 }
93 else
94 {
95 spdlog::debug("{}", result.error());
96 // Stop further tests
97 REQUIRE(result.has_value());
98 }
99 THEN("The resulting manifold has the applied move.")
100 {
101 CHECK(ergodic_moves::check_move(manifold_before, manifold,
102 move_tracker::move_type::THREE_TWO));
103 // Human verification
104 fmt::print("Old manifold.\n");
105 manifold_before.print_details();
106 fmt::print("New manifold after (3,2) move:\n");
107 manifold.print_details();
108 }
109 }
110 WHEN("A (2,6) move is applied to the manifold.")
111 {
112 spdlog::debug("Applying (2,6) move to manifold.\n");
113 if (auto result = apply_move(manifold, ergodic_moves::do_26_move); result)
114 {
115 manifold = result.value();
116 // Update Geometry and Foliated_triangulation with new 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 resulting manifold has the applied move.")
126 {
127 CHECK(ergodic_moves::check_move(manifold_before, manifold,
128 move_tracker::move_type::TWO_SIX));
129 // Human verification
130 fmt::print("Old manifold.\n");
131 manifold_before.print_details();
132 fmt::print("New manifold after (2,6) move:\n");
133 manifold.print_details();
134 }
135 }
136 WHEN("A (6,2) move is applied to the manifold.")
137 {
138 spdlog::debug("Applying (6,2) move to manifold.\n");
139 auto result = apply_move(manifold, ergodic_moves::do_62_move);
140 if (result)
141 {
142 manifold = result.value();
143 // Update Geometry and Foliated_triangulation with new info
144 manifold.update();
145 THEN("The resulting manifold has the applied move.")
146 {
147 CHECK(ergodic_moves::check_move(manifold_before, manifold,
148 move_tracker::move_type::SIX_TWO));
149 // Human verification
150 fmt::print("Old manifold.\n");
151 manifold_before.print_details();
152 fmt::print("New manifold after (6,2) move:\n");
153 manifold.print_details();
154 }
155 }
156 else
157 {
158 spdlog::warn("Cannot apply (6,2) move: {}", result.error());
159 doctest::skip("No valid (6,2) move exists in this manifold configuration.");
160 }
161 }
162 WHEN("A (4,4) move is applied to the manifold.")
163 {
164 spdlog::debug("Applying (4,4) move to manifold.\n");
165 if (auto result = apply_move(manifold, ergodic_moves::do_44_move); result)
166 {
167 manifold = result.value();
168 // Update Geometry and Foliated_triangulation with new info
169 manifold.update();
170 }
171 else
172 {
173 spdlog::debug("{}", result.error());
174 // Stop further tests
175 REQUIRE(result.has_value());
176 }
177 THEN("The resulting manifold has the applied move.")
178 {
179 CHECK(ergodic_moves::check_move(manifold_before, manifold,
180 move_tracker::move_type::FOUR_FOUR));
181 // Human verification
182 fmt::print("Old manifold.\n");
183 manifold_before.print_details();
184 fmt::print("New manifold after (4,4) move:\n");
185 manifold.print_details();
186 }
187 }
188 }
189}
Apply Pachner moves to foliated Delaunay triangulations.
auto constexpr apply_move(ManifoldType &&t_manifold, FunctionType t_move) noexcept -> decltype(auto)
An applicative function similar to std::apply on a manifold.
Definition: Apply_move.hpp:32
Pachner moves on 2+1 dimensional foliated Delaunay triangulations.
SCENARIO("Perform bistellar flip on Delaunay triangulation" *doctest::test_suite("bistellar"))