BubbleProfiler  0.3.0
by Peter Athron, Csaba Balazs, Michael Bardsley, Andrew Fowlie, Dylan Harries & Graham White
action.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 
32 #include "action.hpp"
33 #include "shooting.hpp"
34 
35 extern "C" double action(const double E,
36  const double alpha,
37  const int dim,
46  using namespace BubbleProfiler;
47  const double false_min = 0.;
48  const double true_min = 1.;
49  const double barrier = 0.75 / alpha - 1.;
50 
51  const auto potential = [E, alpha](double phi) {
52  return -E * (-alpha * phi * phi * phi * phi
53  + phi * phi * phi
54  + 0.5 * (4. * alpha - 3.) * phi * phi);
55  };
56 
57  const auto potential_first = [E, alpha](double phi) {
58  return -E * (-4. * alpha * phi * phi * phi
59  + 3. * phi * phi
60  + (4. * alpha - 3.) * phi);
61  };
62 
63  const auto potential_second = [E, alpha](double phi) {
64  return -E * (-12. * alpha * phi * phi
65  + 6. * phi + (4. * alpha - 3.));
66  };
67 
68  Shooting one_dim;
71  one_dim.set_shooting_abs_tol(settings.shoot_ode_abs);
72  one_dim.set_shooting_rel_tol(settings.shoot_ode_rel);
73  one_dim.set_action_abs_tol(settings.action_ode_abs);
74  one_dim.set_action_rel_tol(settings.action_ode_rel);
75  one_dim.set_drho_frac(settings.drho_frac);
76  one_dim.set_bisect_lambda_max(settings.bisect_lambda_max);
77  one_dim.set_max_iterations(settings.iter_max);
78  one_dim.set_max_periods(settings.periods_max);
79  one_dim.set_f_y_max(settings.f_y_max);
80  one_dim.set_f_y_min(settings.f_y_min);
81  one_dim.set_y_max(settings.y_max);
82 
83  one_dim.solve(potential, potential_first, potential_second,
84  false_min, true_min, barrier,
85  dim, Shooting::Solver_options::Compute_action);
86 
87  return one_dim.get_euclidean_action();
88 }
void set_shooting_rel_tol(double tol)
Definition: shooting.hpp:58
Functions for calculating action for one-dimensional potential parameterized by and ...
void set_max_iterations(int i)
Definition: shooting.hpp:76
double action(const double E, const double alpha, const int dim, BubbleProfiler::Shooting_settings settings)
Definition: action.cpp:35
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_bisect_lambda_max(double l)
Definition: shooting.hpp:74
double get_euclidean_action() const
Definition: shooting.cpp:40
void set_f_y_max(double f)
Definition: shooting.hpp:80
void set_action_rel_tol(double tol)
Definition: shooting.hpp:62
void set_f_y_min(double f)
Definition: shooting.hpp:82
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
void set_y_max(double y)
Definition: shooting.hpp:84
Solve the one-dimensional problem using the shooting method.
Definition: shooting.hpp:44
One-dimensional shooting method.
void set_shooting_abs_tol(double tol)
Definition: shooting.hpp:56
void set_bisection_precision_bits(int b)
Definition: shooting.hpp:52