CDT++
Causal Dynamical Triangulations in C++
Loading...
Searching...
No Matches
Foliated_triangulation_test.cpp
Go to the documentation of this file.
1/*******************************************************************************
2 Causal Dynamical Triangulations in C++ using CGAL
3
4 Copyright © 2018 Adam Getchell
5 ******************************************************************************/
6
12
14
15#include <doctest/doctest.h>
16
17#include <numbers>
18
19using namespace std;
20using namespace foliated_triangulations;
21
22static inline auto constexpr RADIUS_2 = 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 auto constexpr INV_SQRT_2 = 1.0 / SQRT_2;
26
27SCENARIO("FoliatedTriangulation special member and swap properties" *
28 doctest::test_suite("foliated_triangulation"))
29{
30 spdlog::debug("FoliatedTriangulation special member and swap properties.\n");
31 GIVEN("A FoliatedTriangulation_3 class.")
32 {
33 WHEN("It's properties are examined.")
34 {
35 THEN("It is no-throw destructible.")
36 {
37 REQUIRE(is_nothrow_destructible_v<FoliatedTriangulation_3>);
38 spdlog::debug("It is no-throw destructible.\n");
39 }
40 THEN("It is default constructible.")
41 {
42 REQUIRE(is_default_constructible_v<FoliatedTriangulation_3>);
43 spdlog::debug("It is default destructible.\n");
44 }
45 THEN("It is NOT trivially default constructible.")
46 {
47 CHECK_FALSE(
48 is_trivially_default_constructible_v<FoliatedTriangulation_3>);
49 }
50 THEN("Delaunay_triangulation_3 is NOT no-throw default constructible.")
51 {
52 CHECK_FALSE(is_nothrow_default_constructible_v<Delaunay_t<3>>);
53 }
54 THEN(
55 "Therefore FoliatedTriangulation is NOT no-throw default "
56 "constructible.")
57 {
58 CHECK_FALSE(
59 is_nothrow_default_constructible_v<FoliatedTriangulation_3>);
60 }
61 THEN("It is no-throw copy constructible.")
62 {
63 REQUIRE(is_nothrow_copy_constructible_v<FoliatedTriangulation_3>);
64 spdlog::debug("It is no-throw copy constructible.\n");
65 }
66 THEN("It is no-throw copy assignable.")
67 {
68 REQUIRE(is_nothrow_copy_assignable_v<FoliatedTriangulation_3>);
69 spdlog::debug("It is no-throw copy assignable.\n");
70 }
71 THEN("It is no-throw move constructible.")
72 {
73 REQUIRE(is_nothrow_move_constructible_v<FoliatedTriangulation_3>);
74 spdlog::debug("It is no-throw move constructible.\n");
75 }
76 THEN("It is no-throw move assignable.")
77 {
78 REQUIRE(is_nothrow_move_assignable_v<FoliatedTriangulation_3>);
79 spdlog::debug("It is no-throw move assignable.\n");
80 }
81 THEN("It is no-throw swappable.")
82 {
83 REQUIRE(is_nothrow_swappable_v<FoliatedTriangulation_3>);
84 spdlog::debug("It is no-throw swappable.\n");
85 }
86 THEN("It is constructible from a Delaunay triangulation.")
87 {
88 REQUIRE(is_constructible_v<FoliatedTriangulation_3, Delaunay_t<3>>);
89 spdlog::debug("It is constructible from a Delaunay triangulation.\n");
90 }
91 THEN("It is constructible from parameters.")
92 {
93 REQUIRE(is_constructible_v<FoliatedTriangulation_3, Int_precision,
94 Int_precision, double, double>);
95 spdlog::debug("It is constructible from parameters.\n");
96 }
97 THEN("It is constructible from Causal_vertices.")
98 {
99 REQUIRE(
100 is_constructible_v<FoliatedTriangulation_3, Causal_vertices_t<3>>);
101 spdlog::debug("It is constructible from Causal_vertices.\n");
102 }
103 THEN("It is constructible from Causal_vertices and INITIAL_RADIUS.")
104 {
105 REQUIRE(is_constructible_v<FoliatedTriangulation_3,
106 Causal_vertices_t<3>, double>);
107 spdlog::debug(
108 "It is constructible from Causal_vertices and INITIAL_RADIUS.\n");
109 }
110 THEN(
111 "It is constructible from Causal_vertices, INITIAL_RADIUS, and "
112 "RADIAL_SEPARATION.")
113 {
114 REQUIRE(is_constructible_v<FoliatedTriangulation_3,
115 Causal_vertices_t<3>, double, double>);
116 spdlog::debug(
117 "It is constructible from Causal_vertices, INITIAL_RADIUS, and "
118 "RADIAL_SEPARATION.\n");
119 }
120 }
121 }
122}
123
124SCENARIO("FoliatedTriangulation free functions" *
125 doctest::test_suite("foliated_triangulation"))
126{
127 spdlog::debug("foliated_triangulations:: free functions.\n");
128
129 GIVEN("A vector of points and timevalues.")
130 {
131 vector const Vertices{Point_t<3>(1, 0, 0), Point_t<3>(0, 1, 0),
132 Point_t<3>(0, 0, 1),
133 Point_t<3>(RADIUS_2, RADIUS_2, RADIUS_2)};
134 vector<size_t> const Timevalues{1, 1, 1, 2};
135 WHEN("Causal vertices are created.")
136 {
137 auto causal_vertices = make_causal_vertices<3>(Vertices, Timevalues);
138 THEN("They are correct.")
139 {
140 REQUIRE_EQ(causal_vertices.size(), 4);
141 REQUIRE_EQ(causal_vertices[0].first, Point_t<3>(1, 0, 0));
142 REQUIRE_EQ(causal_vertices[0].second, 1);
143 REQUIRE_EQ(causal_vertices[1].first, Point_t<3>(0, 1, 0));
144 REQUIRE_EQ(causal_vertices[1].second, 1);
145 REQUIRE_EQ(causal_vertices[2].first, Point_t<3>(0, 0, 1));
146 REQUIRE_EQ(causal_vertices[2].second, 1);
147 REQUIRE_EQ(causal_vertices[3].first,
148 Point_t<3>(RADIUS_2, RADIUS_2, RADIUS_2));
149 REQUIRE_EQ(causal_vertices[3].second, 2);
150 }
151 }
152 }
153 GIVEN("A mismatched set of points and timevalues.")
154 {
155 vector const Vertices{Point_t<3>(1, 0, 0), Point_t<3>(0, 1, 0),
156 Point_t<3>(0, 0, 1),
157 Point_t<3>(RADIUS_2, RADIUS_2, RADIUS_2)};
158 vector<size_t> const Timevalues{1, 1, 1};
159 WHEN("Causal vertices are created.")
160 {
161 THEN("An exception is thrown.")
162 {
163 REQUIRE_THROWS(make_causal_vertices<3>(Vertices, Timevalues));
164 }
165 }
166 }
167
168 GIVEN("A small foliated 3D triangulation.")
169 {
170 vector Vertices{
171 Point_t<3>{ 1, 0, 0},
172 Point_t<3>{ 0, 1, 0},
173 Point_t<3>{ 0, 0, 1},
174 Point_t<3>{RADIUS_2, RADIUS_2, RADIUS_2}
175 };
176 vector<std::size_t> timevalues{1, 1, 1, 2};
177 auto vertices = make_causal_vertices<3>(Vertices, timevalues);
178 FoliatedTriangulation_3 triangulation(vertices);
179 auto print = [&triangulation](auto& vertex) {
180 fmt::print(
181 "Vertex: ({}) Timevalue: {} is a vertex: {} and is "
182 "infinite: {}\n",
183 utilities::point_to_str(vertex->point()), vertex->info(),
184 triangulation.get_delaunay().tds().is_vertex(vertex),
185 triangulation.is_infinite(vertex));
186 };
187
188 REQUIRE(triangulation.is_initialized());
189 WHEN("check_vertices() is called.")
190 {
191 THEN("The vertices are correct.")
192 {
193 CHECK(foliated_triangulations::check_vertices<3>(
194 triangulation.get_delaunay(), 1.0, 1.0));
195 }
196 }
197 WHEN("check_cells() is called.")
198 {
199 THEN("Cells are correctly classified.")
200 {
201 CHECK(foliated_triangulations::check_cells<3>(
202 triangulation.get_delaunay()));
203 // Human verification
204 triangulation.print_cells();
205 }
206 }
207
208 WHEN("We print a cell in the triangulation.")
209 {
210 THEN("A cell is printed correctly.")
211 {
212 foliated_triangulations::print_cell<3>(triangulation.get_cells().at(0));
213 }
214 }
215
216 WHEN("We ask for a container of vertices given a container of cells.")
217 {
218 auto&& all_vertices =
219 get_vertices_from_cells<3>(triangulation.get_cells());
220 THEN("We get back the correct number of vertices.")
221 {
222 REQUIRE_EQ(all_vertices.size(), 4);
223 // Human verification
224 ranges::for_each(all_vertices, print);
225 }
226 }
227 }
228 GIVEN(
229 "A minimal triangulation with non-default initial radius and radial "
230 "separation.")
231 {
232 auto constexpr desired_simplices = 2;
233 auto constexpr desired_timeslices = 2;
234 auto constexpr initial_radius = 3.0;
235 auto constexpr foliation_spacing = 2.0;
236 FoliatedTriangulation_3 const triangulation(
237 desired_simplices, desired_timeslices, initial_radius,
238 foliation_spacing);
239 THEN("The triangulation is initialized correctly.")
240 {
241 REQUIRE(triangulation.is_initialized());
242 }
243 THEN("The initial radius and radial separation are correct.")
244 {
245 REQUIRE_EQ(triangulation.initial_radius(), initial_radius);
246 REQUIRE_EQ(triangulation.foliation_spacing(), foliation_spacing);
247 // Human verification
248 fmt::print(
249 "The triangulation has an initial radius of {} and a radial "
250 "separation of {}\n",
251 initial_radius, foliation_spacing);
252 }
253 THEN("Each vertex has a valid timevalue.")
254 {
255 for (std::span const checked_vertices(triangulation.get_vertices());
256 auto const& vertex : checked_vertices)
257 {
258 CHECK(triangulation.does_vertex_radius_match_timevalue(vertex));
259 fmt::print(
260 "Vertex ({}) with timevalue of {} has a squared radius of {} and a "
261 "squared expected radius of {} with an expected timevalue of {}.\n",
262 utilities::point_to_str(vertex->point()), vertex->info(),
263 squared_radius<3>(vertex),
264 std::pow(triangulation.expected_radius(vertex), 2),
265 triangulation.expected_timevalue(vertex));
266 }
267 }
268 }
269 GIVEN("A triangulation setup for a (4,4) move")
270 {
271 vector vertices{
272 Point_t<3>{ 0, 0, 0},
273 Point_t<3>{ INV_SQRT_2, 0, INV_SQRT_2},
274 Point_t<3>{ 0, INV_SQRT_2, INV_SQRT_2},
275 Point_t<3>{-INV_SQRT_2, 0, INV_SQRT_2},
276 Point_t<3>{ 0, -INV_SQRT_2, INV_SQRT_2},
277 Point_t<3>{ 0, 0, 2}
278 };
279 vector<size_t> timevalue{1, 2, 2, 2, 2, 3};
280 auto causal_vertices = make_causal_vertices<3>(vertices, timevalue);
281 FoliatedTriangulation_3 const triangulation(causal_vertices, 0, 1);
282 // Verify we have 6 vertices, 13 edges, 12 facets, and 4 cells
283 REQUIRE_EQ(triangulation.number_of_vertices(), 6);
284 REQUIRE_EQ(triangulation.number_of_finite_edges(), 13);
285 REQUIRE_EQ(triangulation.number_of_finite_facets(), 12);
286 REQUIRE(triangulation.number_of_finite_cells() == 4);
287 CHECK_EQ(triangulation.initial_radius(), 0);
288 CHECK_EQ(triangulation.foliation_spacing(), 1);
289 REQUIRE(triangulation.is_delaunay());
290 REQUIRE(triangulation.is_correct());
291 WHEN("We collect edges.")
292 {
293 auto edges = foliated_triangulations::collect_edges<3>(
294 triangulation.get_delaunay());
295 THEN("We have 13 edges.") { REQUIRE_EQ(edges.size(), 13); }
296 }
297 WHEN("We have a point in the triangulation.")
298 {
299 THEN("We can obtain the vertex")
300 {
301 auto vertex = foliated_triangulations::find_vertex<3>(
302 triangulation.get_delaunay(), Point_t<3>{0, 0, 0});
303 REQUIRE_MESSAGE(vertex, "Vertex not found.");
304 if (vertex)
305 {
306 CHECK_EQ(vertex.value()->point(), Point_t<3>{0, 0, 0});
307 CHECK_EQ(vertex.value()->info(), 1);
308 // Human verification
309 fmt::print(
310 "Point(0,0,0) was found as vertex ({}) with a timevalue of {}.\n",
311 utilities::point_to_str(vertex.value()->point()),
312 vertex.value()->info());
313 }
314 }
315 WHEN("We choose a point not in the triangulation.")
316 {
317 THEN("No vertex is found.")
318 {
319 auto vertex = foliated_triangulations::find_vertex<3>(
320 triangulation.get_delaunay(), Point_t<3>{3, 3, 3});
321 REQUIRE_FALSE(vertex);
322 // Human verification
323 fmt::print("Point(3,3,3) was not found.\n");
324 }
325 }
326 WHEN("We check vertices in a cell.")
327 {
328 THEN("The correct vertices yields the correct cell.")
329 {
330 auto v_1 = foliated_triangulations::find_vertex<3>(
331 triangulation.get_delaunay(), Point_t<3>{0, 0, 0});
332 auto v_2 = foliated_triangulations::find_vertex<3>(
333 triangulation.get_delaunay(),
334 Point_t<3>{0, INV_SQRT_2, INV_SQRT_2});
335 auto v_3 = foliated_triangulations::find_vertex<3>(
336 triangulation.get_delaunay(),
337 Point_t<3>{0, -INV_SQRT_2, INV_SQRT_2});
338 auto v_4 = foliated_triangulations::find_vertex<3>(
339 triangulation.get_delaunay(),
340 Point_t<3>{-INV_SQRT_2, 0, INV_SQRT_2});
341 REQUIRE_MESSAGE(v_1, "Vertex v_1 not found.");
342 REQUIRE_MESSAGE(v_2, "Vertex v_2 not found.");
343 REQUIRE_MESSAGE(v_3, "Vertex v_3 not found.");
344 REQUIRE_MESSAGE(v_4, "Vertex v_4 not found.");
345 if (v_1 && v_2 && v_3 && v_4)
346 {
347 auto cell = foliated_triangulations::find_cell<3>(
348 triangulation.get_delaunay(), v_1.value(), v_2.value(),
349 v_3.value(), v_4.value());
350 CHECK(cell);
351 // Human verification
352 triangulation.print_cells();
353 }
354 }
355 THEN("The incorrect vertices do not return a cell.")
356 {
357 auto v_1 = foliated_triangulations::find_vertex<3>(
358 triangulation.get_delaunay(), Point_t<3>{0, 0, 0});
359 auto v_2 = foliated_triangulations::find_vertex<3>(
360 triangulation.get_delaunay(),
361 Point_t<3>{INV_SQRT_2, 0, INV_SQRT_2});
362 auto v_3 = foliated_triangulations::find_vertex<3>(
363 triangulation.get_delaunay(),
364 Point_t<3>{0, INV_SQRT_2, INV_SQRT_2});
365 auto v_4 = foliated_triangulations::find_vertex<3>(
366 triangulation.get_delaunay(), Point_t<3>{0, 0, 2});
367 REQUIRE_MESSAGE(v_1, "Vertex v_1 not found.");
368 REQUIRE_MESSAGE(v_2, "Vertex v_2 not found.");
369 REQUIRE_MESSAGE(v_3, "Vertex v_3 not found.");
370 REQUIRE_MESSAGE(v_4, "Vertex v_4 not found.");
371 if (v_1 && v_2 && v_3 && v_4)
372 {
373 auto cell = foliated_triangulations::find_cell<3>(
374 triangulation.get_delaunay(), v_1.value(), v_2.value(),
375 v_3.value(), v_4.value());
376 REQUIRE_FALSE(cell);
377 }
378 }
379 }
380 }
381 WHEN("A container of cells is printed.")
382 {
383 THEN("The container is printed correctly.")
384 {
385 foliated_triangulations::print_cells<3>(triangulation.get_cells());
386 }
387 }
388 WHEN("We choose a cell in the triangulation.")
389 {
390 auto cell = triangulation.get_cells().at(0);
391 THEN("We can print it's neighbors.")
392 {
393 foliated_triangulations::print_cell<3>(cell);
394 foliated_triangulations::print_neighboring_cells<3>(cell);
395 }
396 }
397 }
398}
399
400SCENARIO("FoliatedTriangulation_3 initialization" *
401 doctest::test_suite("foliated_triangulation"))
402{
403 spdlog::debug("FoliatedTriangulation initialization.\n");
404 GIVEN("A 3D foliated triangulation.")
405 {
406 WHEN("It is default constructed.")
407 {
408 THEN("The default Delaunay triangulation is valid.")
409 {
410 FoliatedTriangulation_3 const triangulation;
411 REQUIRE(triangulation.is_initialized());
412 REQUIRE_EQ(triangulation.max_time(), 0);
413 REQUIRE_EQ(triangulation.min_time(), 0);
414 REQUIRE_EQ(triangulation.initial_radius(), INITIAL_RADIUS);
415 REQUIRE_EQ(triangulation.foliation_spacing(), FOLIATION_SPACING);
416 }
417 }
418 WHEN(
419 "It is constructed from a Delaunay triangulation with 4 causal "
420 "vertices.")
421 {
422 vector Vertices{
423 Point_t<3>{ 1, 0, 0},
424 Point_t<3>{ 0, 1, 0},
425 Point_t<3>{ 0, 0, 1},
426 Point_t<3>{RADIUS_2, RADIUS_2, RADIUS_2}
427 };
428 vector<std::size_t> timevalues{1, 1, 1, 2};
429 auto vertices = make_causal_vertices<3>(Vertices, timevalues);
430 FoliatedTriangulation_3 const triangulation(vertices);
431 THEN("Triangulation is valid and foliated.")
432 {
433 REQUIRE(triangulation.is_initialized());
434 REQUIRE_EQ(triangulation.dimension(), 3);
435 REQUIRE_EQ(triangulation.number_of_vertices(), 4);
436 REQUIRE_EQ(triangulation.number_of_finite_edges(), 6);
437 REQUIRE_EQ(triangulation.number_of_finite_facets(), 4);
438 REQUIRE_EQ(triangulation.number_of_finite_cells(), 1);
439 REQUIRE_EQ(triangulation.max_time(), 2);
440 REQUIRE_EQ(triangulation.min_time(), 1);
441 REQUIRE_EQ(triangulation.initial_radius(), INITIAL_RADIUS);
442 REQUIRE_EQ(triangulation.foliation_spacing(), FOLIATION_SPACING);
443 REQUIRE(triangulation.is_foliated());
444 // Human verification
445 triangulation.print_cells();
446 }
447 }
448 WHEN("Constructing the minimum triangulation.")
449 {
450 auto constexpr desired_simplices = 2;
451 auto constexpr desired_timeslices = 2;
452 FoliatedTriangulation_3 triangulation(desired_simplices,
453 desired_timeslices);
454 THEN("Triangulation is valid and foliated.")
455 {
456 REQUIRE(triangulation.is_initialized());
457 }
458 THEN("The triangulation has sensible values.")
459 {
460 // // We have 1 to 8 vertex_count
461 auto vertex_count{triangulation.number_of_vertices()};
462 CHECK_GE(vertex_count, 1);
463 CHECK_LE(vertex_count, 8);
464 // // We have 1 to 12 simplex_count
465 auto simplex_count{triangulation.number_of_finite_cells()};
466 CHECK_GE(simplex_count, 1);
467 CHECK_LE(simplex_count, 12);
468
469 // Human verification
470 triangulation.print();
471 }
472 THEN("The vertices have correct timevalues.")
473 {
474 auto check = [&triangulation](Vertex_handle_t<3> const& vertex) {
475 CHECK(triangulation.does_vertex_radius_match_timevalue(vertex));
476 };
477 ranges::for_each(triangulation.get_vertices(), check);
478 // Human verification
479 auto print = [&triangulation](Vertex_handle_t<3> const& vertex) {
480 fmt::print(
481 "Vertex: ({}) Timevalue: {} has a squared radius of {} and "
482 "a squared expected radius of {} with an expected timevalue of "
483 "{}.\n",
484 utilities::point_to_str(vertex->point()), vertex->info(),
485 squared_radius<3>(vertex),
486 std::pow(triangulation.expected_radius(vertex), 2),
487 triangulation.expected_timevalue(vertex));
488 };
489 ranges::for_each(triangulation.get_vertices(), print);
490 }
491 }
492 WHEN(
493 "Constructing the minimal triangulation with non-default initial "
494 "radius and separation.")
495 {
496 auto constexpr desired_simplices = 2;
497 auto constexpr desired_timeslices = 2;
498 auto constexpr initial_radius = 3.0;
499 auto constexpr radial_factor = 2.0;
500 FoliatedTriangulation_3 const triangulation(
501 desired_simplices, desired_timeslices, initial_radius, radial_factor);
502 THEN("The triangulation is initialized correctly.")
503 {
504 REQUIRE(triangulation.is_initialized());
505 }
506 THEN("The initial radius and radial separation are correct.")
507 {
508 REQUIRE_EQ(triangulation.initial_radius(), initial_radius);
509 REQUIRE_EQ(triangulation.foliation_spacing(), radial_factor);
510 }
511 }
512 WHEN(
513 "Constructing a small triangulation with fractional initial radius and "
514 "separation.")
515 {
516 auto constexpr desired_simplices = 24;
517 auto constexpr desired_timeslices = 3;
518 auto constexpr initial_radius = 1.5;
519 auto constexpr radial_factor = 1.1;
520 FoliatedTriangulation_3 const triangulation(
521 desired_simplices, desired_timeslices, initial_radius, radial_factor);
522 THEN("The triangulation is initialized correctly.")
523 {
524 REQUIRE(triangulation.is_initialized());
525 }
526 THEN("The initial radius and radial separation are correct.")
527 {
528 REQUIRE_EQ(triangulation.initial_radius(), initial_radius);
529 REQUIRE_EQ(triangulation.foliation_spacing(), radial_factor);
530 }
531 }
532 WHEN("Constructing a medium triangulation.")
533 {
534 auto constexpr desired_simplices = 6400;
535 auto constexpr desired_timeslices = 7;
536 FoliatedTriangulation_3 const triangulation(desired_simplices,
537 desired_timeslices);
538 THEN("Triangulation is valid and foliated.")
539 {
540 REQUIRE(triangulation.is_initialized());
541 }
542 THEN("The triangulation has sensible values.")
543 {
544 REQUIRE_EQ(triangulation.min_time(), 1);
545 // Human verification
546 triangulation.print();
547 }
548 THEN("Data members are correctly populated.")
549 {
550 triangulation.print();
551 // Every cell is classified as (3,1), (2,2), or (1,3)
552 CHECK_EQ(triangulation.get_cells().size(),
553 triangulation.get_three_one().size() +
554 triangulation.get_two_two().size() +
555 triangulation.get_one_three().size());
556 // Every cell is properly labelled
557 CHECK(triangulation.check_all_cells());
558
559 CHECK_FALSE(triangulation.N2_SL().empty());
560
561 CHECK_GT(triangulation.max_time(), 0);
562 CHECK_GT(triangulation.min_time(), 0);
563 CHECK_GT(triangulation.max_time(), triangulation.min_time());
564 auto check_timelike = [](Edge_handle_t<3> const& edge) {
565 CHECK(classify_edge<3>(edge));
566 };
567 ranges::for_each(triangulation.get_timelike_edges(), check_timelike);
568
569 auto check_spacelike = [](Edge_handle_t<3> const& edge) {
570 CHECK(!classify_edge<3>(edge));
571 };
572 ranges::for_each(triangulation.get_spacelike_edges(), check_spacelike);
573 // Human verification
574 fmt::print("There are {} edges.\n",
575 triangulation.number_of_finite_edges());
576 fmt::print("There are {} timelike edges and {} spacelike edges.\n",
577 triangulation.N1_TL(), triangulation.N1_SL());
578 fmt::print(
579 "There are {} vertices with a max timevalue of {} and a min "
580 "timevalue of {}.\n",
581 triangulation.number_of_vertices(), triangulation.max_time(),
582 triangulation.min_time());
583 triangulation.print_volume_per_timeslice();
584 }
585 }
586 }
587}
588
589SCENARIO("FoliatedTriangulation_3 copying" *
590 doctest::test_suite("foliated_triangulation"))
591{
592 spdlog::debug("FoliatedTriangulation_3 copying.\n");
593 GIVEN("A FoliatedTriangulation_3")
594 {
595 auto constexpr desired_simplices = 6400;
596 auto constexpr desired_timeslices = 7;
597 FoliatedTriangulation_3 triangulation(desired_simplices,
598 desired_timeslices);
599 WHEN("It is copied")
600 {
601 auto ft2 = triangulation;
602 THEN("The two objects are distinct.")
603 {
604 auto* ft_ptr = &triangulation;
605 auto* ft2_ptr = &ft2;
606 CHECK_NE(ft_ptr, ft2_ptr);
607 }
608 THEN("The foliated triangulations have identical properties.")
609 {
610 CHECK_EQ(triangulation.is_initialized(), ft2.is_initialized());
611 CHECK_EQ(triangulation.number_of_finite_cells(),
612 ft2.number_of_finite_cells());
613 CHECK_EQ(triangulation.min_time(), ft2.min_time());
614 CHECK_EQ(triangulation.get_cells().size(), ft2.get_cells().size());
615 CHECK_EQ(triangulation.get_three_one().size(),
616 ft2.get_three_one().size());
617 CHECK_EQ(triangulation.get_two_two().size(), ft2.get_two_two().size());
618 CHECK_EQ(triangulation.get_one_three().size(),
619 ft2.get_one_three().size());
620 CHECK_EQ(triangulation.N2_SL().size(), ft2.N2_SL().size());
621 }
622 }
623 }
624}
625
626SCENARIO("Detecting and fixing problems with vertices and cells" *
627 doctest::test_suite("foliated_triangulation"))
628{
629 spdlog::debug("Detecting and fixing problems with vertices and cells.\n");
630 GIVEN("A FoliatedTriangulation_3.")
631 {
632 WHEN("Constructing a triangulation with 4 correct vertices.")
633 {
634 vector Vertices{
635 Point_t<3>{ 1, 0, 0},
636 Point_t<3>{ 0, 1, 0},
637 Point_t<3>{ 0, 0, 1},
638 Point_t<3>{RADIUS_2, RADIUS_2, RADIUS_2}
639 };
640 vector<std::size_t> timevalues{1, 1, 1, 2};
641 auto vertices = make_causal_vertices<3>(Vertices, timevalues);
642 FoliatedTriangulation_3 triangulation(vertices);
643 THEN("No errors in the vertices are detected.")
644 {
645 CHECK(triangulation.check_all_vertices());
646 // Human verification
647 triangulation.print_vertices();
648 }
649 THEN("No errors in the simplex are detected.")
650 {
651 CHECK(triangulation.is_correct());
652 CHECK_FALSE(check_timevalues<3>(triangulation.get_delaunay()));
653 // Human verification
654 triangulation.print_cells();
655 }
656 THEN("No errors in the triangulation foliation are detected")
657 {
658 CHECK_FALSE(foliated_triangulations::fix_timevalues<3>(
659 triangulation.delaunay()));
660 // Human verification
661 utilities::print_delaunay(triangulation.get_delaunay());
662 }
663 AND_WHEN("The vertices are mis-labelled.")
664 {
665 auto break_vertices = [](Vertex_handle_t<3> const& vertex) {
666 vertex->info() = 0;
667 };
668 ranges::for_each(triangulation.get_vertices(), break_vertices);
669
670 THEN("The incorrect vertex labelling is identified.")
671 {
672 CHECK_FALSE(triangulation.check_all_vertices());
673 auto bad_vertices = triangulation.find_incorrect_vertices();
674 CHECK_FALSE(bad_vertices.empty());
675 // Human verification
676 fmt::print("=== Wrong vertex info! ===\n");
677 triangulation.print_vertices();
678 }
679 AND_THEN("The incorrect vertex labelling is fixed.")
680 {
681 CHECK_FALSE(triangulation.check_all_vertices());
682 auto bad_vertices = triangulation.find_incorrect_vertices();
683 CHECK_FALSE(bad_vertices.empty());
684
685 CHECK(triangulation.fix_vertices());
686 CHECK(triangulation.check_all_vertices());
687 fmt::print("=== Corrected vertex info ===\n");
688 triangulation.print_vertices();
689 }
690 }
691 AND_WHEN("The cells are mis-labelled.")
692 {
693 auto break_cells = [](Cell_handle_t<3> const& cell) {
694 cell->info() = 0;
695 };
696 ranges::for_each(triangulation.get_cells(), break_cells);
697 THEN("The incorrect cell labelling is identified.")
698 {
699 CHECK_FALSE(triangulation.check_all_cells());
700 // Human verification
701 fmt::print("=== Wrong cell info! ===\n");
702 triangulation.print_cells();
703 }
704 THEN("The incorrect cell labelling is fixed.")
705 {
706 CHECK_FALSE(triangulation.check_all_cells());
707 CHECK(triangulation.fix_cells());
708 // Human verification
709 fmt::print("=== Corrected cell info ===\n");
710 triangulation.print_cells();
711 CHECK(triangulation.check_all_cells());
712 }
713 }
714 }
715 WHEN(
716 "Constructing a triangulation with an incorrect high timevalue "
717 "vertex.")
718 {
719 vector vertices{
720 Point_t<3>{ 1, 0, 0},
721 Point_t<3>{ 0, 1, 0},
722 Point_t<3>{ 0, 0, 1},
723 Point_t<3>{RADIUS_2, RADIUS_2, RADIUS_2}
724 };
725 vector<std::size_t> timevalues{1, 1, 1, std::numeric_limits<int>::max()};
726 auto causal_vertices = make_causal_vertices<3>(vertices, timevalues);
727 FoliatedTriangulation_3 const triangulation(causal_vertices);
728 THEN("The vertex is fixed on construction.")
729 {
730 CHECK_FALSE(triangulation.fix_vertices());
731 CHECK(triangulation.is_initialized());
732 triangulation.print_cells();
733 }
734 }
735 WHEN("Constructing a triangulation with an incorrect low value vertex.")
736 {
737 vector vertices{
738 Point_t<3>{0, 0, 0},
739 Point_t<3>{0, 1, 0},
740 Point_t<3>{1, 0, 0},
741 Point_t<3>{0, 0, 1}
742 };
743 vector<std::size_t> timevalues{0, 2, 2, 2};
744 auto causal_vertices = make_causal_vertices<3>(vertices, timevalues);
745 FoliatedTriangulation_3 const triangulation(causal_vertices);
746 THEN("The vertex is fixed on construction.")
747 {
748 CHECK_FALSE(triangulation.fix_vertices());
749 CHECK(triangulation.is_initialized());
750 triangulation.print_cells();
751 }
752 }
753 WHEN(
754 "Constructing a triangulation with two incorrect low values and two "
755 "incorrect high values.")
756 {
757 vector vertices{
758 Point_t<3>{0, 0, 0},
759 Point_t<3>{0, 1, 0},
760 Point_t<3>{1, 0, 0},
761 Point_t<3>{0, 0, 1}
762 };
763 vector<std::size_t> timevalues{0, 0, 2, 2};
764 auto causal_vertices = make_causal_vertices<3>(vertices, timevalues);
765 FoliatedTriangulation_3 const triangulation(causal_vertices);
766 THEN("The vertices are fixed on construction.")
767 {
768 CHECK_FALSE(triangulation.fix_vertices());
769 CHECK(triangulation.is_initialized());
770 triangulation.print_cells();
771 }
772 AND_THEN("The cell type is correct.")
773 {
774 CHECK_FALSE(triangulation.fix_vertices());
775 CHECK_FALSE(triangulation.fix_cells());
776 CHECK(triangulation.is_initialized());
777 triangulation.print_cells();
778 }
779 }
780 WHEN(
781 "Constructing a triangulation with all vertices on the same timeslice.")
782 {
783 vector vertices{
784 Point_t<3>{1, 0, 0},
785 Point_t<3>{0, 1, 0},
786 Point_t<3>{0, 0, 1},
787 Point_t<3>{0, 0, -1}
788 };
789 vector<std::size_t> timevalues{1, 1, 1, 1};
790 auto causal_vertices = make_causal_vertices<3>(vertices, timevalues);
791 FoliatedTriangulation_3 const triangulation(causal_vertices);
792 THEN("The vertex error is detected.")
793 {
794 CHECK_FALSE(triangulation.is_initialized());
795 auto cell = triangulation.get_delaunay().finite_cells_begin();
796 CHECK_EQ(expected_cell_type<3>(cell), Cell_type::ACAUSAL);
797 // Human verification
798 triangulation.print_cells();
799 }
800 }
801 WHEN("Constructing a triangulation with an unfixable vertex.")
802 {
803 vector vertices{
804 Point_t<3>{1, 0, 0},
805 Point_t<3>{0, 1, 0},
806 Point_t<3>{0, 0, 1},
807 Point_t<3>{0, 0, 2},
808 Point_t<3>{2, 0, 0},
809 Point_t<3>{0, 3, 0}
810 };
811 vector<std::size_t> timevalues{1, 1, 1, 2, 2, 3};
812 auto causal_vertices = make_causal_vertices<3>(vertices, timevalues);
813 Delaunay_t<3> const delaunay_triangulation{causal_vertices.begin(),
814 causal_vertices.end()};
815 // Passing in a Delaunay triangulation directly allows us to skip the
816 // normal construction process with sanity checks on the triangulation,
817 // which is what we're testing here individually.
818 FoliatedTriangulation_3 triangulation(delaunay_triangulation);
819 THEN("The incorrect cell can be identified.")
820 {
821 auto bad_cells = check_timevalues<3>(delaunay_triangulation);
822 CHECK_MESSAGE(bad_cells.has_value(), "No bad cells found.");
823 if (bad_cells)
824 {
825 fmt::print("Bad cells:\n");
826 print_cells<3>(bad_cells.value());
827 }
828 }
829 AND_THEN("The incorrect vertex can be identified.")
830 {
831 auto bad_cells = check_timevalues<3>(delaunay_triangulation);
832 CHECK_MESSAGE(bad_cells.has_value(), "No bad cells found.");
833 if (bad_cells)
834 {
835 auto bad_vertex = find_bad_vertex<3>(bad_cells->front());
836 fmt::print("Bad vertex ({}) has timevalues {}.\n",
837 utilities::point_to_str(bad_vertex->point()),
838 bad_vertex->info());
839 CHECK_EQ(bad_vertex->info(), 3);
840 }
841 }
842 AND_THEN("The triangulation is fixed.")
843 {
844 fmt::print("Unfixed triangulation:\n");
845 triangulation.print_cells();
846 CHECK(foliated_triangulations::fix_timevalues<3>(
847 triangulation.delaunay()));
848 CHECK(triangulation.is_initialized());
849 fmt::print("Fixed triangulation:\n");
850 print_cells<3>(collect_cells<3>(triangulation.delaunay()));
851 }
852 }
853 }
854}
855
856SCENARIO("FoliatedTriangulation_3 functions from Delaunay3" *
857 doctest::test_suite("foliated_triangulation"))
858{
859 spdlog::debug("FoliatedTriangulation_3 functions from Delaunay3.\n");
860 GIVEN("A FoliatedTriangulation_3.")
861 {
862 WHEN("Constructing a small triangulation.")
863 {
864 vector vertices{
865 Point_t<3>{1, 0, 0},
866 Point_t<3>{0, 1, 0},
867 Point_t<3>{0, 0, 1},
868 Point_t<3>{0, 0, 2},
869 Point_t<3>{2, 0, 0},
870 Point_t<3>{0, 3, 0}
871 };
872 vector<std::size_t> timevalues{1, 1, 1, 2, 2, 3};
873 auto causal_vertices = make_causal_vertices<3>(vertices, timevalues);
874 FoliatedTriangulation_3 triangulation(causal_vertices);
875 THEN("The Foliated triangulation is initially wrong.")
876 {
877 CHECK_FALSE(triangulation.is_initialized());
878 // Human verification
879#ifndef NDEBUG
880 fmt::print("Unfixed triangulation:\n");
881 triangulation.print_cells();
882#endif
883 }
884 THEN("After being fixed, Delaunay3 functions work as expected.")
885 {
886 // Fix the triangulation
887 CHECK(triangulation.is_fixed());
888 CHECK_EQ(triangulation.number_of_finite_cells(), 2);
889 fmt::print("Base Delaunay number of cells: {}\n",
890 triangulation.number_of_finite_cells());
891 CHECK_EQ(triangulation.number_of_finite_facets(), 7);
892 fmt::print("Base Delaunay number of faces: {}\n",
893 triangulation.number_of_finite_facets());
894 triangulation.print_volume_per_timeslice();
895 CHECK_EQ(triangulation.number_of_finite_edges(), 9);
896 fmt::print("Base Delaunay number of edges: {}\n",
897 triangulation.number_of_finite_edges());
898 triangulation.print_edges();
899 CHECK_EQ(triangulation.number_of_vertices(), 5);
900 fmt::print("Base Delaunay number of vertices: {}\n",
901 triangulation.number_of_vertices());
902 CHECK_EQ(triangulation.dimension(), 3);
903 fmt::print("Base Delaunay dimension is: {}\n",
904 triangulation.dimension());
905 // Human verification
906#ifndef NDEBUG
907 utilities::print_delaunay(triangulation.delaunay());
908#endif
909 }
910 }
911 WHEN("Constructing the default triangulation.")
912 {
913 FoliatedTriangulation_3 const triangulation;
914 REQUIRE(triangulation.is_initialized());
915 THEN("is_infinite() identifies a single infinite vertex.")
916 {
917 auto&& vertices = triangulation.get_delaunay().tds().vertices();
918 auto&& vertex = vertices.begin();
919 CHECK_EQ(vertices.size(), 1);
920 CHECK(triangulation.get_delaunay().tds().is_vertex(vertex));
921 CHECK(triangulation.is_infinite(vertex));
922 }
923 }
924 WHEN("Constructing a triangulation with 4 causal vertices.")
925 {
926 vector vertices{
927 Point_t<3>{ 1, 0, 0},
928 Point_t<3>{ 0, 1, 0},
929 Point_t<3>{ 0, 0, 1},
930 Point_t<3>{RADIUS_2, RADIUS_2, RADIUS_2}
931 };
932 vector<std::size_t> timevalues{1, 1, 1, 2};
933 auto causal_vertices = make_causal_vertices<3>(vertices, timevalues);
934 FoliatedTriangulation_3 triangulation(causal_vertices);
935 REQUIRE(triangulation.is_initialized());
936 THEN("The degree of each vertex is 4 (including infinite vertex).")
937 {
938 auto check = [&triangulation](Vertex_handle_t<3> const& vertex) {
939 CHECK_EQ(triangulation.degree(vertex), 4);
940 };
941 ranges::for_each(triangulation.get_vertices(), check);
942 }
943 }
944 }
945}
Create foliated spherical triangulations.
std::int_fast32_t Int_precision
Definition: Settings.hpp:31
static double constexpr INITIAL_RADIUS
Default foliated triangulation spacings.
Definition: Settings.hpp:47
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"))
auto degree(VertexHandle &&t_vertex) const
Perfect forwarding to TriangulationDataStructure_3::degree.
auto get_vertices() const noexcept -> Vertex_container const &