14#include <doctest/doctest.h>
19using namespace manifolds;
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;
27SCENARIO(
"Use check_move to validate successful move" *
28 doctest::test_suite(
"ergodic"))
30 spdlog::debug(
"Use check_move to validate successful move.\n");
31 GIVEN(
"A triangulation setup for a (2,3) move")
37 Point_t<3>{RADIUS_2, RADIUS_2, RADIUS_2},
38 Point_t<3>{ SQRT_2, SQRT_2, 0}
40 vector<size_t> timevalues{1, 1, 1, 2, 2};
41 auto causal_vertices = make_causal_vertices<3>(vertices, timevalues);
43 WHEN(
"A correct (2,3) move is performed.")
45 spdlog::debug(
"When a correct (2,3) move is performed.\n");
47 auto manifold_before = manifold;
48 if (
auto result = ergodic_moves::do_23_move(manifold); result)
50 manifold = result.value();
56 spdlog::debug(
"{}", result.error());
57 REQUIRE(result.has_value());
60 REQUIRE(manifold.is_correct());
61 CHECK_EQ(manifold.vertices(), 5);
62 CHECK_EQ(manifold.edges(), 10);
63 CHECK_EQ(manifold.faces(), 9);
64 CHECK_EQ(manifold.simplices(), 3);
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);
70 THEN(
"check_move returns true")
72 CHECK(ergodic_moves::check_move(manifold_before, manifold,
73 move_tracker::move_type::TWO_THREE));
80 "Perform ergodic moves on the minimal manifold necessary for (2,3) and (3,2) moves" *
81 doctest::test_suite(
"ergodic"))
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")
91 Point_t<3>{RADIUS_2, RADIUS_2, RADIUS_2},
92 Point_t<3>{ SQRT_2, SQRT_2, 0}
94 vector<size_t> timevalues{1, 1, 1, 2, 2};
95 auto causal_vertices = make_causal_vertices<3>(vertices, timevalues);
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")
110 spdlog::debug(
"When a (2,3) move is performed.\n");
112 auto manifold_before = manifold;
113 if (
auto result = ergodic_moves::do_23_move(manifold); result)
115 manifold = result.value();
121 spdlog::debug(
"{}", result.error());
123 REQUIRE(result.has_value());
125 THEN(
"The move is correct and the manifold invariants are maintained")
127 CHECK(ergodic_moves::check_move(manifold_before, manifold,
128 move_tracker::move_type::TWO_THREE));
130 REQUIRE(manifold.is_correct());
131 CHECK_EQ(manifold.vertices(), 5);
132 CHECK_EQ(manifold.edges(), 10);
133 CHECK_EQ(manifold.faces(), 9);
134 CHECK_EQ(manifold.simplices(), 3);
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);
141 manifold.print_details();
142 manifold.print_cells();
145 WHEN(
"A (3,2) move is performed")
147 spdlog::debug(
"When a (3,2) move is performed.\n");
149 if (
auto start = ergodic_moves::do_23_move(manifold); start)
151 manifold = start.value();
158 "The (2,3) move to set up the manifold for the (3,2) move "
161 REQUIRE(start.has_value());
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);
174 auto manifold_before = manifold;
176 if (
auto result = ergodic_moves::do_32_move(manifold); result)
178 manifold = result.value();
184 spdlog::debug(
"{}", result.error());
186 REQUIRE(result.has_value());
188 THEN(
"The move is correct and the manifold invariants are maintained")
190 CHECK(ergodic_moves::check_move(manifold_before, manifold,
191 move_tracker::move_type::THREE_TWO));
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());
204 manifold.print_details();
205 manifold.print_cells();
208 WHEN(
"An improperly prepared (3,2) move is performed")
211 THEN(
"The move is not performed")
214 CHECK_EQ(result.error(),
"No (3,2) move possible.\n");
220 "Perform ergodic moves on the minimal manifold necessary (2,6) and (6,2) moves" *
221 doctest::test_suite(
"ergodic"))
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")
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}
234 vector<size_t> timevalues{0, 1, 1, 1, 2};
235 auto causal_vertices = make_causal_vertices<3>(vertices, timevalues);
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")
252 spdlog::debug(
"When a (2,6) move is performed.\n");
254 auto manifold_before = manifold;
256 if (
auto result = ergodic_moves::do_26_move(manifold); result)
258 manifold = result.value();
264 spdlog::debug(
"The (2,6) move failed.\n");
266 REQUIRE(result.has_value());
268 THEN(
"The move is correct and the manifold invariants are maintained")
270 CHECK(ergodic_moves::check_move(manifold_before, manifold,
271 move_tracker::move_type::TWO_SIX));
273 REQUIRE(manifold.is_correct());
274 CHECK_EQ(manifold.vertices(), 6);
275 CHECK_EQ(manifold.edges(), 14);
276 CHECK_EQ(manifold.faces(), 15);
277 CHECK_EQ(manifold.simplices(), 6);
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);
283 CHECK_EQ(manifold.N1_TL(), 8);
284 CHECK(manifold.is_delaunay());
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();
294 WHEN(
"A (6,2) move is performed")
296 spdlog::debug(
"When a (6,2) move is performed.\n");
298 if (
auto start = ergodic_moves::do_26_move(manifold); start)
300 manifold = start.value();
307 "The (2,6) move to set up the manifold for the (6,2) move "
310 REQUIRE(start.has_value());
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());
326 auto manifold_before = manifold;
328 if (
auto result = ergodic_moves::do_62_move(manifold); result)
332 else { spdlog::info(
"The (6,2) move failed.\n"); }
333 THEN(
"The move is correct and the manifold invariants are maintained")
336 CHECK(ergodic_moves::check_move(manifold_before, manifold,
337 move_tracker::move_type::SIX_TWO));
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());
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();
363 WHEN(
"An improperly prepared (6,2) move is performed")
366 THEN(
"The move is not performed")
369 CHECK_EQ(result.error(),
"No (6,2) move possible.\n");
374SCENARIO(
"Perform ergodic moves on the minimal manifold necessary (4,4) moves" *
375 doctest::test_suite(
"ergodic") * doctest::skip())
378 "Perform ergodic moves on the minimal manifold necessary (4,4) moves.\n");
379 GIVEN(
"A triangulation setup for a (4,4) move")
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},
389 vector<size_t> timevalues{0, 1, 1, 1, 1, 2};
390 auto causal_vertices = make_causal_vertices<3>(vertices, timevalues);
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());
408 WHEN(
"A (4,4) move is performed")
410 spdlog::debug(
"When a (4,4) move is performed.\n");
412 auto manifold_before = manifold;
414 fmt::print(
"Manifold before (4,4):\n");
415 manifold_before.print_details();
416 manifold_before.print_cells();
418 if (
auto result = ergodic_moves::do_44_move(manifold); result)
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")
429 CHECK(ergodic_moves::check_move(manifold_before, manifold,
430 move_tracker::move_type::FOUR_FOUR));
436SCENARIO(
"Test convenience functions needed for bistellar flip" *
437 doctest::test_suite(
"ergodic"))
439 GIVEN(
"A triangulation setup for a bistellar flip")
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},
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")
454 auto cells = foliated_triangulations::collect_cells<3>(triangulation);
455 THEN(
"We have 4 cells") { REQUIRE_EQ(cells.size(), 4); }
457 WHEN(
"We get all finite edges in the triangulation")
459 THEN(
"We have 13 edges") { REQUIRE_EQ(edges.size(), 13); }
461 WHEN(
"We find the pivot edge in the triangulation")
464 REQUIRE_MESSAGE(pivot_edge,
"No pivot edge found.");
466 auto Contains = [&vertices](Point_t<3> point) {
467 return ranges::any_of(vertices, [&point](Point_t<3>
const& test) {
468 return test == point;
475 triangulation, pivot_edge.value());
476 REQUIRE_MESSAGE(incident_cells,
"No incident cells found.");
477 THEN(
"We have a pivot edge")
479 CHECK_MESSAGE(pivot_edge,
"Pivot edge found");
480 REQUIRE(triangulation.tds().is_edge(
481 pivot_edge->first, pivot_edge->second, pivot_edge->third));
483 pivot_edge->first->vertex(pivot_edge->second)->point();
485 pivot_edge->first->vertex(pivot_edge->third)->point();
487 REQUIRE_FALSE(Contains(Point_t<3>{0, 0, 1}));
488 REQUIRE(Contains(pivot_from_1));
489 REQUIRE(Contains(pivot_from_2));
493 THEN(
"We can obtain the cells incident to that edge")
495 REQUIRE_EQ(incident_cells->size(), 4);
497 AND_THEN(
"We can obtain the vertices from the incident cells")
499 auto incident_vertices =
501 REQUIRE_EQ(incident_vertices.size(), 6);
506 WHEN(
"We get all finite vertices in the triangulation")
508 THEN(
"We have 6 vertices")
510 auto all_finite_vertices =
511 foliated_triangulations::collect_vertices<3>(triangulation);
512 REQUIRE_EQ(all_finite_vertices.size(), 6);
518SCENARIO(
"Perform bistellar flip on Delaunay triangulation" *
519 doctest::test_suite(
"ergodic"))
521 GIVEN(
"A triangulation setup for a bistellar flip")
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},
531 ergodic_moves::Delaunay triangulation(vertices.begin(), vertices.end());
532 WHEN(
"We have a valid triangulation")
534 CHECK(triangulation.is_valid());
535 THEN(
"We can perform a bistellar flip")
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.");
546 REQUIRE_EQ(vertices.size(), 6);
550 auto flipped_triangulation = ergodic_moves::bistellar_flip(
551 triangulation, pivot_edge.value(), top, bottom);
553 REQUIRE_MESSAGE(flipped_triangulation,
"Bistellar flip failed.");
554 if (flipped_triangulation)
558 WARN(flipped_triangulation->is_valid());
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"))