CDT++
Causal Dynamical Triangulations in C++
Loading...
Searching...
No Matches
bistellar-flip.cpp
Go to the documentation of this file.
1/*******************************************************************************
2Causal Dynamical Triangulations in C++ using CGAL
3Copyright © 2023 Adam Getchell
4******************************************************************************/
5
14
15#ifdef NDEBUG
16#define DOCTEST_CONFIG_DISABLE
17#include <CGAL/draw_triangulation_3.h>
18#endif
19
20#define DOCTEST_CONFIG_IMPLEMENT
21
22#include <doctest/doctest.h>
23#include <spdlog/spdlog.h>
24
25#include <numbers>
26
27#include "Ergodic_moves_3.hpp"
28
29static inline std::floating_point auto constexpr SQRT_2 =
30 std::numbers::sqrt2_v<double>;
31static inline std::floating_point auto constexpr INV_SQRT_2 = 1 / SQRT_2;
32
33auto bistellar_triangulation_vertices() -> std::vector<Point_t<3>>
34{
35 std::vector vertices{
36 Point_t<3>{ 0, 0, 0},
37 Point_t<3>{ INV_SQRT_2, 0, INV_SQRT_2},
38 Point_t<3>{ 0, INV_SQRT_2, INV_SQRT_2},
39 Point_t<3>{-INV_SQRT_2, 0, INV_SQRT_2},
40 Point_t<3>{ 0, -INV_SQRT_2, INV_SQRT_2},
41 Point_t<3>{ 0, 0, 2}
42 };
43 return vertices;
44}
45
46auto main(int const argc, char* argv[]) -> int
47try
48{
49 // Doctest integration into code
50 doctest::Context context;
51 context.setOption("no-breaks",
52 true); // don't break in debugger when assertions fail
53 context.applyCommandLine(argc, argv);
54
55 int const res = context.run(); // run tests unless --no-run is specified
56 if (context.shouldExit())
57 { // important - query flags (and --exit) rely on the user doing this
58 return res; // propagate the result of the tests
59 }
60
61 context.clearFilters(); // important - otherwise the context filters will be
62 // used during the next evaluation of RUN_ALL_TESTS,
63 // which will lead to wrong results
64
65#ifdef NDEBUG
66 fmt::print("Before bistellar flip.\n");
67 auto vertices = bistellar_triangulation_vertices();
68 ergodic_moves::Delaunay const dt{vertices.begin(), vertices.end()};
69 manifolds::Manifold_3 const manifold{
71 };
72#ifdef ENABLE_VISUALIZATION
73 CGAL::draw(manifold.get_delaunay());
74#endif
75 fmt::print("After bistellar flip.\n");
76 manifold.print_cells();
78#endif
79}
80catch (std::exception const& e)
81{
82 spdlog::critical("Error: {}\n", e.what());
83 return EXIT_FAILURE;
84}
85
86catch (...)
87{
88 spdlog::critical("Something went wrong ... Exiting.\n");
89 return EXIT_FAILURE;
90}
91
92SCENARIO("Perform bistellar flip on Delaunay triangulation" *
93 doctest::test_suite("bistellar"))
94{
95 GIVEN("A triangulation setup for a bistellar flip")
96 {
97 auto vertices = bistellar_triangulation_vertices();
98 ergodic_moves::Delaunay triangulation(vertices.begin(), vertices.end());
99 WHEN("We have a valid triangulation")
100 {
101 CHECK(triangulation.is_valid());
102 THEN("We can perform a bistellar flip")
103 {
104 // Obtain top and bottom vertices by re-inserting, which returns the
105 // Vertex_handle
106 auto top = triangulation.insert(Point_t<3>{0, 0, 2});
107 auto bottom = triangulation.insert(Point_t<3>{0, 0, 0});
108 auto edges = foliated_triangulations::collect_edges<3>(triangulation);
109 auto pivot_edge = ergodic_moves::find_pivot_edge(triangulation, edges);
110 REQUIRE_MESSAGE(pivot_edge, "No pivot edge found.");
111
112 // Check this didn't actually change vertices in the triangulation
113 REQUIRE_EQ(vertices.size(), 6);
114
115 if (pivot_edge)
116 {
117 auto flipped_triangulation = ergodic_moves::bistellar_flip(
118 triangulation, pivot_edge.value(), top, bottom);
119
120 REQUIRE_MESSAGE(flipped_triangulation, "Bistellar flip failed.");
121 if (flipped_triangulation)
122 {
125 WARN(flipped_triangulation->is_valid());
126 }
127 }
128 }
129 }
130 }
131}
Pachner moves on 2+1 dimensional foliated Delaunay triangulations.
void print_delaunay(TriangulationType const &t_triangulation)
Print triangulation statistics.
Definition: Utilities.hpp:144
SCENARIO("Perform bistellar flip on Delaunay triangulation" *doctest::test_suite("bistellar"))