BubbleProfiler  0.3.0
by Peter Athron, Csaba Balazs, Michael Bardsley, Andrew Fowlie, Dylan Harries & Graham White
tabulate.cpp
Go to the documentation of this file.
1 /*
2  * This file is part of BubbleProfiler.
3  *
4  * BubbleProfiler is free software: you can redistribute it and/or modify
5  * it under the terms of the GNU General Public License as published by
6  * the Free Software Foundation, either version 3 of the License, or
7  * (at your option) any later version.
8  *
9  * BubbleProfiler is distributed in the hope that it will be useful,
10  * but WITHOUT ANY WARRANTY; without even the implied warranty of
11  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
12  * GNU General Public License for more details.
13  *
14  * You should have received a copy of the GNU General Public License
15  * along with BubbleProfiler. If not, see <http://www.gnu.org/licenses/>.
16  */
17 
34 #include "error.hpp"
35 #include "gsl_root_finder.hpp"
36 #include "kink_profile_guesser.hpp"
37 #include "logging_manager.hpp"
41 #include "shooting.hpp"
42 
43 #include <chrono>
44 #include <iostream>
45 #include <map>
46 #include <string>
47 #include <tuple>
48 #include <vector>
49 
50 namespace BubbleProfiler {
51 
52 namespace {
53 
54 const std::string shooting_label = "shooting";
55 const std::string perturbative_label = "perturbative";
56 
57 } // anonymous namespace
58 
60  int dim{4};
61  double E{1.};
62  double min_alpha{0.51};
63  double max_alpha{0.74};
64  double alpha_step{0.01};
65  bool use_shooting{true};
66  bool use_perturbative{false};
67  bool verbose{false};
68 };
69 
71  double action{0.};
72  double time_in_ms{0.};
73  int error{0};
74  std::string info{""};
75 };
76 
77 using Profiler_data = std::map<std::string, std::vector<Profiler_result> >;
78 
79 void print_usage()
80 {
81  std::cout <<
82  "Usage: quartic_tabulate [OPTION] <dim> <min_alpha> <max_alpha> <alpha_step>\n\n"
83  "Calculate action for one-dimensional quartic potentials with alpha\n"
84  "varying between the given minimum and maximum values.\n\n"
85  "Example: quartic_tabulate 4 0.51 0.74 0.1\n\n"
86  "Options:\n"
87  " -a, --all-methods print results for all algorithms\n"
88  " -h, --help print this help message\n"
89  " -p, --perturbative use the perturbative algorithm\n"
90  " -v, --verbose produce verbose output"
91  << std::endl;
92 }
93 
94 bool starts_with(const std::string& option, const std::string& prefix)
95 {
96  return !option.compare(0, prefix.size(), prefix);
97 }
98 
99 void parse_positional_args(const std::vector<std::string>& args,
100  Tabulate_options& options)
101 {
102  const auto n_positional_args = args.size();
103 
104  if (n_positional_args == 0) {
105  throw Setup_error(
106  "number of dimensions must be provided");
107  } else if (n_positional_args == 1) {
108  throw Setup_error(
109  "minimum value of alpha must be provided");
110  } else if (n_positional_args == 2) {
111  throw Setup_error(
112  "maximum value of alpha must be provided");
113  } else if (n_positional_args == 3) {
114  throw Setup_error(
115  "alpha increment must be provided");
116  } else if (n_positional_args > 4) {
117  throw Setup_error(
118  "unrecognized command line argument '"
119  + args[4] + "'");
120  }
121 
122  try {
123  options.dim = std::stoi(args[0]);
124  } catch (const std::invalid_argument& e) {
125  throw Setup_error(
126  "invalid number of dimensions '" + args[0] + "'");
127  }
128 
129  try {
130  options.min_alpha = std::stod(args[1]);
131  } catch (const std::invalid_argument& e) {
132  throw Setup_error(
133  "invalid minimum value for alpha, '" + args[1] + "'");
134  }
135 
136  try {
137  options.max_alpha = std::stod(args[2]);
138  } catch (const std::invalid_argument& e) {
139  throw Setup_error(
140  "invalid maximum value for alpha, '" + args[2] + "'");
141  }
142 
143  try {
144  options.alpha_step = std::stod(args[3]);
145  } catch (const std::invalid_argument& e) {
146  throw Setup_error(
147  "invalid alpha increment, '" + args[3] + "'");
148  }
149 }
150 
151 Tabulate_options parse_cmd_line_args(int argc, const char* argv[])
152 {
153  Tabulate_options options;
154 
155  std::vector<std::string> positional_args;
156  int i = 1;
157  while (i < argc) {
158  const std::string opt(argv[i++]);
159 
160  if (!starts_with(opt, "-")) {
161  positional_args.push_back(opt);
162  continue;
163  }
164 
165  if (opt == "-a" || opt == "--all-methods") {
166  options.use_shooting = true;
167  options.use_perturbative = true;
168  continue;
169  }
170 
171  if (opt == "-h" || opt == "--help") {
172  print_usage();
173  exit(EXIT_SUCCESS);
174  }
175 
176  if (opt == "-p" || opt == "--perturbative") {
177  options.use_shooting = false;
178  options.use_perturbative = true;
179  continue;
180  }
181 
182  if (opt == "-v" || opt == "--verbose") {
183  options.verbose = true;
184  continue;
185  }
186 
187  throw Setup_error(
188  "unrecognized command line argument '" + opt + "'");
189  }
190 
191  parse_positional_args(positional_args, options);
192 
193  return options;
194 }
195 
197  const Restricted_quartic_potential& potential, int dim)
198 {
199  using Clock_t = std::chrono::high_resolution_clock;
200 
201  const double false_min = potential.get_local_minimum_location();
202  const double barrier = potential.get_local_maximum_location();
203  const double true_min = potential.get_global_minimum_location();
204 
205  Shooting profiler;
206  profiler.set_bisection_precision_bits(5);
207  profiler.set_action_arrived_rel(1.e-3);
208  profiler.set_shooting_abs_tol(1.e-4);
209  profiler.set_shooting_rel_tol(1.e-4);
210  profiler.set_action_abs_tol(1.e-6);
211  profiler.set_action_rel_tol(1.e-6);
212  profiler.set_drho_frac(1.e-3);
213  profiler.set_bisect_lambda_max(5);
214  profiler.set_max_iterations(100000);
215  profiler.set_max_periods(1.e2);
216  profiler.set_f_y_max(1.e6);
217  profiler.set_f_y_min(1.e-3);
218  profiler.set_y_max(1.e1);
219 
220  Profiler_result result;
221  const auto start_time = Clock_t::now();
222  try {
223  profiler.solve(potential, false_min, true_min, barrier, dim,
224  Shooting::Solver_options::Compute_action);
225  result.action = profiler.get_euclidean_action();
226  } catch (const Error& e) {
227  result.error = 1;
228  result.info = e.what();
229  }
230  const auto stop_time = Clock_t::now();
231  std::chrono::duration<double, std::milli> t_in_ms = stop_time - start_time;
232  result.time_in_ms = t_in_ms.count();
233 
234  return result;
235 }
236 
238  Restricted_quartic_potential potential, int dim)
239 {
240  using Clock_t = std::chrono::high_resolution_clock;
242 
243  RK4_perturbative_profiler profiler;
244 
245  profiler.set_initial_step_size(1.e-3);
246 
247  profiler.set_initial_guesser(std::make_shared<Kink_profile_guesser>());
248 
249  const double action_tol = 1.e-3;
250  const double fields_tol = 1.e-3;
251  profiler.set_convergence_tester(
252  std::make_shared<Relative_convergence_tester>(
253  action_tol, fields_tol));
254 
255  auto root_finder = std::make_shared<Root_finder>();
256  profiler.set_root_finder(root_finder);
257 
258  Eigen::VectorXd true_vacuum(1);
259  true_vacuum << potential.get_global_minimum_location();
260 
261  Eigen::VectorXd false_vacuum(1);
262  false_vacuum << potential.get_local_minimum_location();
263 
264  profiler.set_number_of_dimensions(dim);
265  profiler.set_true_vacuum_loc(true_vacuum);
266  profiler.set_false_vacuum_loc(false_vacuum);
267 
268  Profiler_result result;
269  const auto start_time = Clock_t::now();
270  try {
271  profiler.calculate_bubble_profile(potential);
272  result.action = profiler.get_euclidean_action();
273  } catch (const Error& e) {
274  result.error = 1;
275  result.info = e.what();
276  }
277  const auto stop_time = Clock_t::now();
278  std::chrono::duration<double, std::milli> t_in_ms = stop_time - start_time;
279  result.time_in_ms = t_in_ms.count();
280 
281  return result;
282 }
283 
285  const std::string& profiler,
286  const Restricted_quartic_potential& potential, int dim)
287 {
288  if (profiler == shooting_label) {
289  return run_shooting_profiler(potential, dim);
290  } else if (profiler == perturbative_label) {
291  return run_perturbative_profiler(potential, dim);
292  }
293 
294  throw Setup_error(
295  "unrecognized profiler '" + profiler + "'");
296 }
297 
298 void write_table_header(const std::vector<std::string>& profiler_labels)
299 {
300  std::cout << "# "
301  << std::setw(12) << "dim" << ' '
302  << std::setw(16) << "E" << ' '
303  << std::setw(16) << "alpha";
304  for (const auto& profiler : profiler_labels) {
305  std::cout << ' '
306  << std::setw(16) << profiler + "_action" << ' '
307  << std::setw(16) << profiler + "_time/ms" << ' '
308  << std::setw(12) << profiler + "_error";
309  }
310  std::cout << std::endl;
311 }
312 
313 void write_table(int dim, double E, const std::vector<double>& alpha_data,
314  const Profiler_data& profiler_data)
315 {
316  std::vector<std::string> profilers;
317  for (const auto& x : profiler_data) {
318  profilers.push_back(x.first);
319  }
320 
321  write_table_header(profilers);
322 
323  const std::size_t n_rows = alpha_data.size();
324  for (std::size_t i = 0; i < n_rows; ++i) {
325  std::cout << " "
326  << std::setw(12) << dim << ' '
327  << std::setw(16) << std::setprecision(8)
328  << std::scientific << E << ' '
329  << std::setw(16) << std::setprecision(8)
330  << std::scientific << alpha_data[i];
331 
332  std::string info = "";
333  for (const auto& profiler : profilers) {
334  const auto& data = profiler_data.at(profiler)[i];
335  std::cout << ' '
336  << std::setw(16) << std::setprecision(8)
337  << std::scientific << data.action << ' '
338  << std::setw(16) << std::setprecision(8)
339  << std::scientific << data.time_in_ms << ' '
340  << std::setw(12) << data.error;
341 
342  if (!data.info.empty()) {
343  if (info.empty()) {
344  info = "\t# ";
345  }
346  info += profiler + ": " + data.info + ",";
347  }
348  }
349 
350  std::cout << info << std::endl;
351  }
352 }
353 
354 } // namespace BubbleProfiler
355 
356 int main(int argc, const char* argv[])
357 {
358  using namespace BubbleProfiler;
359 
360  if (argc == 1) {
361  std::cerr << "Usage: quartic_tabulate [OPTION] <dim> <min_alpha> <max_alpha> <alpha_step>\n"
362  << "Try 'quartic_tabulate --help' for more information."
363  << std::endl;
364  return EXIT_FAILURE;
365  }
366 
367  int exit_code = 0;
368  try {
369  Tabulate_options options = parse_cmd_line_args(argc, argv);
370 
371  if (options.verbose) {
374  }
375 
376  if (std::abs(options.alpha_step)
377  < std::numeric_limits<double>::epsilon()) {
378  throw Setup_error(
379  "alpha increment too small");
380  }
381 
382  const bool forward_stepping = options.alpha_step > 0.;
383  const bool start_gt_end = options.min_alpha > options.max_alpha;
384  if ((forward_stepping && start_gt_end) ||
385  (!forward_stepping && !start_gt_end)) {
386  std::swap(options.min_alpha, options.max_alpha);
387  }
388 
389  std::vector<double> alpha_data;
390  Profiler_data profiler_data;
391  if (options.use_shooting) {
392  profiler_data[shooting_label] = std::vector<Profiler_result>();
393  }
394 
395  if (options.use_perturbative) {
396  profiler_data[perturbative_label] = std::vector<Profiler_result>();
397  }
398 
399  double alpha = options.min_alpha;
400  while (alpha <= options.max_alpha) {
401  Restricted_quartic_potential potential(alpha, options.E);
402 
403  alpha_data.push_back(alpha);
404  for (auto& x : profiler_data) {
405  profiler_data[x.first].push_back(
406  run_profiler(x.first, potential, options.dim));
407  }
408 
409  alpha += options.alpha_step;
410  }
411 
412  write_table(options.dim, options.E, alpha_data, profiler_data);
413 
414  } catch (const Error& e) {
415  std::cerr << "Error: " << e.what() << std::endl;
416  exit_code = EXIT_FAILURE;
417  } catch (...) {
418  std::cerr << "Error: unrecognized error occurred."
419  << std::endl;
420  exit_code = EXIT_FAILURE;
421  }
422 
423  return exit_code;
424 }
void run_profiler(const Bubble_profiler_inputs &input)
Calculates the action and bubble profile for the given potential.
void set_shooting_rel_tol(double tol)
Definition: shooting.hpp:58
void parse_positional_args(const std::vector< std::string > &args, Fubini_options &options)
void set_initial_guesser(std::shared_ptr< Profile_guesser > g)
void set_max_iterations(int i)
Definition: shooting.hpp:76
void set_false_vacuum_loc(const Eigen::VectorXd &fv)
std::tuple< double, Field_profiles > run_shooting_profiler(const Bubble_profiler_inputs &input)
Calculates the action and bubble profile using the shooting method.
Bubble_profiler_inputs parse_cmd_line_args(int argc, const char *argv[])
void solve(const std::function< double(double)> &potential_, const std::function< double(double)> &potential_first_, const std::function< double(double)> &potential_second_, double false_min_, double true_min_, double barrier_, int dim_=3, unsigned int options=(Solver_options::Compute_action|Solver_options::Compute_profile))
Definition: shooting.cpp:81
void set_root_finder(std::shared_ptr< Root_finder< Eigen::VectorXd > > rf)
double get_euclidean_action() const
Retrieve the action after running calculate_bubble_profile.
void set_bisect_lambda_max(double l)
Definition: shooting.hpp:74
void write_table(int dim, double E, const std::vector< double > &alpha_data, const Profiler_data &profiler_data)
Definition: tabulate.cpp:313
void set_convergence_tester(std::shared_ptr< Profile_convergence_tester > ct)
auto abs(const Eigen::MatrixBase< Derived > &m) -> decltype(m.cwiseAbs())
static Logging_manager & get_manager()
double get_euclidean_action() const
Definition: shooting.cpp:40
std::map< std::string, std::vector< Profiler_result > > Profiler_data
Definition: tabulate.cpp:77
void set_f_y_max(double f)
Definition: shooting.hpp:80
int main(int argc, const char *argv[])
Definition: tabulate.cpp:356
void set_action_rel_tol(double tol)
Definition: shooting.hpp:62
void set_f_y_min(double f)
Definition: shooting.hpp:82
void write_table_header(const std::vector< std::string > &profiler_labels)
Definition: tabulate.cpp:298
void set_action_arrived_rel(double tol)
Definition: shooting.hpp:54
void set_action_abs_tol(double tol)
Definition: shooting.hpp:60
void set_max_periods(double p)
Definition: shooting.hpp:78
void set_drho_frac(double frac)
Definition: shooting.hpp:64
std::tuple< double, Field_profiles > run_perturbative_profiler(const Bubble_profiler_inputs &input, Profiler &profiler)
Calculates the action and bubble profile using the perturbative method.
void set_y_max(double y)
Definition: shooting.hpp:84
Solve the one-dimensional problem using the shooting method.
Definition: shooting.hpp:44
Exception indicating general setup error.
Definition: error.hpp:32
One-dimensional shooting method.
bool starts_with(const std::string &option, const std::string &prefix)
void set_shooting_abs_tol(double tol)
Definition: shooting.hpp:56
void set_bisection_precision_bits(int b)
Definition: shooting.hpp:52
Bounce solver using perturbative method.