BubbleProfiler  0.3.0
by Peter Athron, Csaba Balazs, Michael Bardsley, Andrew Fowlie, Dylan Harries & Graham White
gaussian_potential_test.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 
18 #include "gaussian_potential.hpp"
19 #include "profile_guesser.hpp"
20 
21 #define private public
22 #include "kink_profile_guesser.hpp"
23 #undef private
24 
25 #include "field_profiles.hpp"
27 
28 #include <boost/program_options.hpp>
29 #include <iostream>
30 #include <Eigen/Core>
31 
32 namespace po = boost::program_options;
33 
35  double gamma;
36  double lambda;
37  int n_fields;
38  std::string output_path{""};
39 };
40 
41 Gaussian_inputs parse_cmd_line_args(int argc, const char* argv[])
42 {
43  Gaussian_inputs inputs;
44 
45  po::options_description desc("gaussian_potential_test options");
46 
47  desc.add_options()
48  ("gamma", po::value<double>()->required(), "scale of secondary gaussian relative to primary")
49  ("lambda", po::value<double>()->required(), "geometric distance between gaussians")
50  ("n_fields", po::value<int>()->required(), "number of fields in potential")
51  ("output", po::value<std::string>(), "output path for field profiles")
52  ("help", "print help message");
53 
54  po::variables_map vm;
55  po::store(
56  po::parse_command_line(
57  argc, argv, desc,
58  po::command_line_style::unix_style ^ po::command_line_style::allow_short),
59  vm);
60 
61  if (vm.count("help")) {
62  std::cout << desc << std::endl;
63  exit(EXIT_SUCCESS);
64  }
65 
66  po::notify(vm);
67 
68  inputs.gamma = vm["gamma"].as<double>();
69  inputs.lambda = vm["lambda"].as<double>();
70  inputs.n_fields = vm["n_fields"].as<int>();
71 
72  if (vm.count("output")) {
73  inputs.output_path = vm["output"].as<std::string>();
74  }
75 
76  return inputs;
77 }
78 
79 int main(int argc, const char* argv[])
80 {
81  using namespace BubbleProfiler;
82 
83  Gaussian_inputs input = parse_cmd_line_args(argc, argv);
84 
85  Gaussian_potential potential(input.gamma, input.lambda, input.n_fields);
86 
87  std::shared_ptr<Kink_profile_guesser> kink_guesser
88  = std::make_shared<Kink_profile_guesser>();
89  std::shared_ptr<Profile_guesser> guesser(kink_guesser);
90 
91  Eigen::VectorXd false_vacuum(Eigen::VectorXd::Zero(input.n_fields));
92  Eigen::VectorXd true_vacuum(Eigen::VectorXd::Zero(input.n_fields));
93 
94  for (int i = 0; i < input.n_fields; ++i) {
95  true_vacuum(i) = input.lambda / std::sqrt(input.n_fields);
96  }
97 
98  Field_profiles ansatz = guesser->get_profile_guess(
99  potential, true_vacuum, input.n_fields, -1.0, -1.0, 1.e-3, 100);
100 
101  double alpha = kink_guesser->alpha;
102  std::cout << "Alpha: " << alpha << '\n';
103 
104  RK4_perturbative_profiler profiler;
105 
106  profiler.set_true_vacuum_loc(true_vacuum);
107  profiler.set_false_vacuum_loc(false_vacuum);
108  profiler.set_initial_guesser(guesser);
109 
110  auto convergence_tester = std::make_shared<Relative_convergence_tester>();
111  convergence_tester->set_max_iterations(30);
112  profiler.set_convergence_tester(convergence_tester);
113 
114  // Build string of field names
115  std::vector<std::string> field_names;
116  for (int i = 0; i < input.n_fields; ++i) {
117  field_names.push_back(std::to_string(i));
118  }
119 
120  if (!input.output_path.empty()) {
121  Profile_observer observer(field_names, input.output_path, potential);
122  profiler.calculate_bubble_profile(potential, observer);
123  }
124  else {
125  profiler.calculate_bubble_profile(potential);
126  }
127 
128  std::cout << "Action: " << profiler.get_euclidean_action() << '\n';
129 
130  return 0;
131 }
contains the definition of the Field_profiles clas
contains definitions of Gaussian_potential class
void set_initial_guesser(std::shared_ptr< Profile_guesser > g)
void set_false_vacuum_loc(const Eigen::VectorXd &fv)
double get_euclidean_action() const
Retrieve the action after running calculate_bubble_profile.
void set_convergence_tester(std::shared_ptr< Profile_convergence_tester > ct)
int main(int argc, const char *argv[])
Potential composed of two multivariate Gaussian functions.
Observes the profile during iteration and writes out profile, perturbations and action for each step ...
Definition: observers.hpp:50
Gaussian_inputs parse_cmd_line_args(int argc, const char *argv[])
Discretized set of field profiles.
Bounce solver using perturbative method.