CDT++
Causal Dynamical Triangulations in C++
Loading...
Searching...
No Matches
initialize.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
10
11#include <boost/program_options.hpp>
12
13#include "Manifold.hpp"
14
15using namespace std;
16namespace po = boost::program_options;
17
18static string_view constexpr USAGE{
19 R"(Causal Dynamical Triangulations in C++ using CGAL.
20
21Copyright (c) 2014 Adam Getchell
22
23A program that generates d-dimensional triangulated spacetimes
24with a defined causal structure. Specify the topology of the triangulation
25(spherical or toroidal), the desired number of simplices, and the
26desired number of timeslices.
27
28Usage:./initialize (--spherical | --toroidal) -n SIMPLICES -t TIMESLICES
29 [-d DIM]
30 [--init INITIAL RADIUS]
31 [--foliate FOLIATION SPACING]
32 [--output]
33
34Optional arguments are in square brackets.
35
36Examples:
37./initialize --spherical --simplices 32000 --timeslices 11 --init 1.0 --foliate 1.0 --output
38./initialize -s -n32000 -t11 -i1.0 -f1.0 -o
39
40Options)"};
41
42auto main(int const argc, char* const argv[]) -> int
43try
44{
45 std::string const intro{USAGE};
46 // Parsed arguments
47 topology_type topology;
48 long long simplices;
49 long long timeslices;
50 long long dimensions;
51 double initial_radius;
52 double foliation_spacing;
53 bool save_file;
54
55 po::options_description description(intro);
56 description.add_options()("help,h", "Show this message")(
57 "version,v", "Show program version")("spherical,s", "Spherical topology")(
58 "toroidal,e", "Toroidal topology")("simplices,n",
59 po::value<long long>(&simplices),
60 "Approximate number of simplices")(
61 "timeslices,t", po::value<long long>(&timeslices),
62 "Number of timeslices")(
63 "dimensions,d", po::value<long long>(&dimensions)->default_value(3),
64 "Dimensionality")("init,i",
65 po::value<double>(&initial_radius)->default_value(1.0),
66 "Initial radius")(
67 "foliate,f", po::value<double>(&foliation_spacing)->default_value(1.0),
68 "Foliation spacing")("output,o", "Save triangulation into OFF file");
69
70 po::variables_map args;
71 po::store(po::parse_command_line(argc, argv, description), args);
72 po::notify(args);
73
74 if (args.count("help"))
75 {
76 cout << description << "\n";
77 return EXIT_SUCCESS;
78 }
79
80 if (args.count("version"))
81 {
82 fmt::print("CDT initializer version 1.0\n");
83 return EXIT_SUCCESS;
84 }
85
86 if (args.count("spherical")) { topology = topology_type::SPHERICAL; }
87 else if (args.count("toroidal")) { topology = topology_type::TOROIDAL; }
88 else
89 {
90 fmt::print("Topology not specified.\n");
91 return EXIT_FAILURE;
92 }
93
94 if (args.count("simplices"))
95 {
96 simplices = args["simplices"].as<long long>();
97 }
98 else
99 {
100 fmt::print("Number of simplices not specified.\n");
101 return EXIT_FAILURE;
102 }
103
104 if (args.count("timeslices"))
105 {
106 timeslices = args["timeslices"].as<long long>();
107 }
108 else
109 {
110 fmt::print("Number of timeslices not specified.\n");
111 return EXIT_FAILURE;
112 }
113
114 if (args.count("output")) { save_file = true; }
115 else { save_file = false; }
116
117 // Initialize triangulation
118 manifolds::Manifold_3 universe;
119
120 // Display job parameters
121 fmt::print("Topology is {}\n", utilities::topology_to_str(topology));
122 fmt::print("Number of dimensions = {}\n", dimensions);
123 fmt::print("Number of desired simplices = {}\n", simplices);
124 fmt::print("Number of desired timeslices = {}\n", timeslices);
125 fmt::print("Initial radius = {}\n", initial_radius);
126 fmt::print("Foliation spacing = {}\n", foliation_spacing);
127
128 if (save_file) { fmt::print("Output will be saved.\n"); }
129
130 if (simplices < 2 || timeslices < 2)
131 {
132 throw invalid_argument(
133 "Simplices and timeslices should be greater or equal to 2.");
134 }
135
136 switch (topology)
137 {
138 case topology_type::SPHERICAL:
139 if (dimensions == 3)
140 {
141 // Start your run
142 manifolds::Manifold_3 populated_universe(
143 static_cast<Int_precision>(simplices),
144 static_cast<Int_precision>(timeslices), initial_radius,
145 foliation_spacing);
146 swap(populated_universe, universe);
147 }
148 else { throw invalid_argument("Currently, dimensions cannot be >3."); }
149 break;
150 case topology_type::TOROIDAL:
151 throw invalid_argument("Toroidal triangulations not yet supported.");
152 }
153 universe.print();
155 fmt::print("Final number of simplices: {}\n", universe.N3());
156 if (save_file) { utilities::write_file(universe); }
157 return EXIT_SUCCESS;
158}
159catch (invalid_argument const& InvalidArgument)
160{
161 spdlog::critical("{}\n", InvalidArgument.what());
162 spdlog::critical("Invalid parameter ... Exiting.\n");
163 return EXIT_FAILURE;
164}
165catch (...)
166{
167 spdlog::critical("Something went wrong ... Exiting.\n");
168 return EXIT_FAILURE;
169}
Data structures for manifolds.
std::int_fast32_t Int_precision
Definition: Settings.hpp:31
void write_file(std::filesystem::path const &filename, TriangulationType triangulation)
Write triangulation to file.
Definition: Utilities.hpp:163
topology_type
clang-15 does not support std::format
Definition: Utilities.hpp:43
void print() const
Print manifold.
Definition: Manifold.hpp:372
void print_volume_per_timeslice() const
Print the codimension 1 volume of simplices (faces) per timeslice.
Definition: Manifold.hpp:359