BubbleProfiler  0.3.0
by Peter Athron, Csaba Balazs, Michael Bardsley, Andrew Fowlie, Dylan Harries & Graham White
Public Types | Public Member Functions | Private Types | Private Member Functions | Private Attributes | List of all members
BubbleProfiler::Shooting Class Reference

Solve the one-dimensional problem using the shooting method. More...

#include <shooting.hpp>

Collaboration diagram for BubbleProfiler::Shooting:
Collaboration graph
[legend]

Public Types

enum  Solver_options : unsigned int { Compute_action = 0x01, Compute_profile = 0x02 }
 

Public Member Functions

void set_bisection_precision_bits (int b)
 
void set_action_arrived_rel (double tol)
 
void set_shooting_abs_tol (double tol)
 
void set_shooting_rel_tol (double tol)
 
void set_action_abs_tol (double tol)
 
void set_action_rel_tol (double tol)
 
void set_drho_frac (double frac)
 
void set_evolve_change_rel (double frac)
 
void set_bisect_lambda_max (double l)
 
void set_max_iterations (int i)
 
void set_max_periods (double p)
 
void set_f_y_max (double f)
 
void set_f_y_min (double f)
 
void set_y_max (double y)
 
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))
 
void solve (const Potential &potential, double false_min_, double true_min_, double barrier_, int dim_=3, unsigned int options=(Solver_options::Compute_action|Solver_options::Compute_profile))
 
const Field_profilesget_bubble_profile () const
 
double get_euclidean_action () const
 
double quadratic_potential_solution (double) const
 
double quadratic_potential_solution (double, double) const
 

Private Types

using evolve_type = boost::array< double, 3 >
 
using state_type = boost::array< double, 2 >
 
using error_stepper_type = boost::numeric::odeint::runge_kutta_dopri5< state_type >
 
using controlled_stepper_type = boost::numeric::odeint::controlled_runge_kutta< error_stepper_type >
 

Private Member Functions

void initialize_potential_parameters (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_)
 
void ode (const state_type &x, state_type &dxdt, double rho)
 ODE for bounce solution. More...
 
double calculate_roll_time (double lambda, double dVdphi, double d2Vdphi2) const
 Find the (dimensionless) time at which the field changes by a small amount, \(m \hat\rho\). More...
 
double calculate_roll_velocity (double lambda, double y, double dVdphi, double d2Vdphi2) const
 Find velocity \(\hat\phi(\hat\rho)\). More...
 
evolve_type evolve (double lambda) const
 Uses an approximation solution to evolve the system forward until the field changes by a required amount. More...
 
double bubble_scale () const
 Characteristic scale of bubble. More...
 
double energy (double phi, double dot_phi) const
 
double shoot (double lambda)
 Performs a single shot of our shooting method. More...
 
double unmap (double lambda) const
 
double shooting ()
 Solve for \(\phi_0\) by bisection. More...
 
double integrand (double dot_phi, double rho) const
 Action integrand from kinetic term without pre-factor. More...
 
double action (double lambda)
 
double interior (double lambda)
 
Field_profiles calculate_bubble_profile (double lambda)
 

Private Attributes

int dim {3}
 
double lambda_sol {0.}
 
double euclidean_action {0.}
 
Field_profiles profiles {}
 
bool action_computed {false}
 
bool profile_computed {false}
 
controlled_stepper_type shoot_stepper {}
 
controlled_stepper_type action_stepper {}
 
double drho_guess {0.}
 
double scale {0.}
 
double barrier {0.}
 
double false_min {0.}
 
double true_min {1.}
 
double sign_min {1.}
 
double p_barrier {0.}
 
double p_false_min {0.}
 
double p_true_min {0.}
 
std::function< double(double)> potential {nullptr}
 
std::function< double(double)> potential_first {nullptr}
 
std::function< double(double)> potential_second {nullptr}
 
int shoot_bisect_bits {5}
 
double action_arrived_rel {1.e-3}
 
double shoot_ode_abs {1.e-4}
 
double shoot_ode_rel {1.e-4}
 
double action_ode_abs {1.e-6}
 
double action_ode_rel {1.e-6}
 
double drho_frac {1.e-3}
 
double evolve_change_rel {1.e-2}
 
double bisect_lambda_max {5.}
 
int iter_max {100000}
 
double periods_max {1.e2}
 
double f_y_max {1.e6}
 
double f_y_min {1.e-3}
 
double y_max {1.e1}
 

Detailed Description

Solve the one-dimensional problem using the shooting method.

Examples:
action.cpp, and scale.cpp.

Definition at line 44 of file shooting.hpp.

Member Typedef Documentation

using BubbleProfiler::Shooting::controlled_stepper_type = boost::numeric::odeint::controlled_runge_kutta<error_stepper_type>
private

Definition at line 121 of file shooting.hpp.

using BubbleProfiler::Shooting::error_stepper_type = boost::numeric::odeint::runge_kutta_dopri5<state_type>
private

Definition at line 119 of file shooting.hpp.

using BubbleProfiler::Shooting::evolve_type = boost::array<double, 3>
private

Definition at line 116 of file shooting.hpp.

using BubbleProfiler::Shooting::state_type = boost::array<double, 2>
private

Definition at line 117 of file shooting.hpp.

Member Enumeration Documentation

Enumerator
Compute_action 
Compute_profile 

Definition at line 46 of file shooting.hpp.

Member Function Documentation

double BubbleProfiler::Shooting::action ( double  lambda)
private
double BubbleProfiler::Shooting::bubble_scale ( ) const
private

Characteristic scale of bubble.

Found from period of oscillation about barrier. We approximate the potential at the barrier by a quadratic such that

\[ T = \frac{2\pi}{\sqrt{-V^{\prime\prime}(\phi_b)}}. \]

Definition at line 277 of file shooting.cpp.

References barrier, potential_second, and BubbleProfiler::Sqrt().

Referenced by initialize_potential_parameters().

Field_profiles BubbleProfiler::Shooting::calculate_bubble_profile ( double  lambda)
private
Returns
The bubble profile, \(\phi(\rho)\)

Definition at line 490 of file shooting.cpp.

References BubbleProfiler::Abs(), action_arrived_rel, action_ode_abs, action_ode_rel, drho_guess, evolve(), false_min, iter_max, ode(), periods_max, quadratic_potential_solution(), scale, sign_min, and unmap().

Referenced by solve().

double BubbleProfiler::Shooting::calculate_roll_time ( double  lambda,
double  dVdphi,
double  d2Vdphi2 
) const
private

Find the (dimensionless) time at which the field changes by a small amount, \(m \hat\rho\).

Returns
\(m \hat\rho\)
Parameters
lambda\(\lambda\), related to \(\phi_0\)
dVdphi\(\left.\frac{dV}{d\phi}\right|_\lambda\)
d2Vdphi2\(\left.\frac{d^2V}{d\phi^2}\right|_\lambda\)

Definition at line 164 of file shooting.cpp.

References BubbleProfiler::Abs(), BubbleProfiler::approx_root_neg_4(), BubbleProfiler::approx_root_pos_4(), BubbleProfiler::asinc(), BubbleProfiler::asinch(), barrier, dim, evolve_change_rel, f_y_max, false_min, BubbleProfiler::Log(), true_min, and unmap().

Referenced by evolve().

double BubbleProfiler::Shooting::calculate_roll_velocity ( double  lambda,
double  y,
double  dVdphi,
double  d2Vdphi2 
) const
private

Find velocity \(\hat\phi(\hat\rho)\).

Returns
\(\hat\phi(\hat\rho)\)
Parameters
lambda\(\lambda\), related to \(\phi_0\)
y\(y = m \hat\rho\)
dVdphi\(\left.\frac{dV}{d\phi}\right|_\lambda\)
d2Vdphi2\(\left.\frac{d^2V}{d\phi^2}\right|_\lambda\)

Definition at line 200 of file shooting.cpp.

References BubbleProfiler::Abs(), BubbleProfiler::BesselI(), BubbleProfiler::BesselJ(), BubbleProfiler::Cos(), BubbleProfiler::Cosh(), dim, evolve_change_rel, f_y_max, f_y_min, false_min, sign_min, BubbleProfiler::Sqrt(), and unmap().

Referenced by evolve().

double BubbleProfiler::Shooting::energy ( double  phi,
double  dot_phi 
) const
private
Returns
Energy at \((\phi, \dot\phi)\)

Definition at line 283 of file shooting.cpp.

References p_false_min, and potential.

Referenced by shoot().

Shooting::evolve_type BubbleProfiler::Shooting::evolve ( double  lambda) const
private

Uses an approximation solution to evolve the system forward until the field changes by a required amount.

The required amount is specified by Shooting::rchange using:

\[ |\phi(\rho) - \phi_0| = Shooting::rchange |\phi_f - \phi_0|. \]

This method is used to avoid solving the ODE during the time period where friction dominates, which can be very slow, and avoids a singularity at \(\rho =0\).

Parameters
lambda\(\lambda\)
Returns
\((\hat\rho, \phi(\hat\rho), \dot\phi(\hat\rho))\)

Definition at line 251 of file shooting.cpp.

References BubbleProfiler::Abs(), calculate_roll_time(), calculate_roll_velocity(), evolve_change_rel, false_min, potential_first, potential_second, sign_min, BubbleProfiler::Sqrt(), and unmap().

Referenced by action(), calculate_bubble_profile(), interior(), and shoot().

const Field_profiles & BubbleProfiler::Shooting::get_bubble_profile ( ) const
double BubbleProfiler::Shooting::get_euclidean_action ( ) const

Euclidean action, \(S\)

Examples:
action.cpp, and scale.cpp.

Definition at line 40 of file shooting.cpp.

References action_computed, and euclidean_action.

Referenced by action(), main(), BubbleProfiler::run_shooting_profiler(), set_y_max(), and BubbleProfiler::solve_shooting().

void BubbleProfiler::Shooting::initialize_potential_parameters ( 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_ 
)
private
double BubbleProfiler::Shooting::integrand ( double  dot_phi,
double  rho 
) const
private

Action integrand from kinetic term without pre-factor.

We use a trick that

\[ I = \int d\rho \rho^{d-1} \left[ \frac12 \dot\phi^2 + V(\phi) - V(\phi_f)\right] \equiv I_1 + I_2 \]

but by the fact that action is stationary under \(\rho \to a \rho\)

\[ d I_2 = (2 - d) I_1 \]

such that

\[ I = I_1 + I_2 = \frac{2}{d} I_1 = \frac{1}{d} \int d\rho \dot\phi^2 \rho^{d-1} \]

Parameters
dot_phi\(\dot\phi\)
rho\(\rho\)
Returns
Action without prefactor

Definition at line 409 of file shooting.cpp.

References dim.

Referenced by action().

double BubbleProfiler::Shooting::interior ( double  lambda)
private
void BubbleProfiler::Shooting::ode ( const state_type x,
state_type dxdt,
double  rho 
)
private

ODE for bounce solution.

Rewrite second-order as two coupled first order equations,

\[ \dot\phi = \chi \]

and

\[ \dot\chi = -\frac{2}{\rho} \phi + V^\prime(\phi) \]

Parameters
x\((\phi, \chi)\)
dxdt\((\dot\phi, \dot\chi)\)
rho\(\rho\)

Definition at line 158 of file shooting.cpp.

References dim, and potential_first.

Referenced by action(), calculate_bubble_profile(), and shoot().

double BubbleProfiler::Shooting::quadratic_potential_solution ( double  rho) const

This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.

Definition at line 572 of file shooting.cpp.

References lambda_sol.

Referenced by calculate_bubble_profile(), BubbleProfiler::Shooting_profile_guesser::calculate_field_profiles(), and set_y_max().

double BubbleProfiler::Shooting::quadratic_potential_solution ( double  rho,
double  lambda 
) const
void BubbleProfiler::Shooting::set_action_abs_tol ( double  tol)
inline

Absolute tolerance in ODE solver for action integral

Examples:
action.cpp.

Definition at line 60 of file shooting.hpp.

References action_ode_abs.

Referenced by action(), main(), and BubbleProfiler::run_shooting_profiler().

void BubbleProfiler::Shooting::set_action_arrived_rel ( double  tol)
inline

Relative tolerance for judging when \(\phi \approx \phi_f\)

Examples:
action.cpp.

Definition at line 54 of file shooting.hpp.

References action_arrived_rel.

Referenced by action(), main(), and BubbleProfiler::run_shooting_profiler().

void BubbleProfiler::Shooting::set_action_rel_tol ( double  tol)
inline

Relative tolerance in ODE solver for action integral

Examples:
action.cpp.

Definition at line 62 of file shooting.hpp.

References action_ode_rel.

Referenced by action(), main(), and BubbleProfiler::run_shooting_profiler().

void BubbleProfiler::Shooting::set_bisect_lambda_max ( double  l)
inline

Maximum value of \(\lambda\)

Examples:
action.cpp.

Definition at line 74 of file shooting.hpp.

References bisect_lambda_max.

Referenced by action(), main(), and BubbleProfiler::run_shooting_profiler().

void BubbleProfiler::Shooting::set_bisection_precision_bits ( int  b)
inline

Bits of precision in bisection of \(\phi_0\)

Examples:
action.cpp.

Definition at line 52 of file shooting.hpp.

References shoot_bisect_bits.

Referenced by action(), main(), BubbleProfiler::run_shooting_profiler(), and BubbleProfiler::solve_shooting().

void BubbleProfiler::Shooting::set_drho_frac ( double  frac)
inline

The initial ODE step size, \(d\rho\), is this fraction of the bubble scale

Examples:
action.cpp.

Definition at line 64 of file shooting.hpp.

References drho_frac.

Referenced by action(), main(), and BubbleProfiler::run_shooting_profiler().

void BubbleProfiler::Shooting::set_evolve_change_rel ( double  frac)
inline

The EOM is integrated analytically until

\[ \frac{\left|\phi(\hat\rho) - \phi_0 \right|}{\left|\phi_f - \phi_0\right|}, \]

changes by this amount.

Definition at line 72 of file shooting.hpp.

References evolve_change_rel.

Referenced by BubbleProfiler::run_shooting_profiler().

void BubbleProfiler::Shooting::set_f_y_max ( double  f)
inline

Threshold for asymptotic approximation in analytic evolution

Examples:
action.cpp.

Definition at line 80 of file shooting.hpp.

References f_y_max.

Referenced by action(), main(), and BubbleProfiler::run_shooting_profiler().

void BubbleProfiler::Shooting::set_f_y_min ( double  f)
inline

Threshold for asymptotic approximation in analytic evolution

Examples:
action.cpp.

Definition at line 82 of file shooting.hpp.

References f_y_min.

Referenced by action(), main(), and BubbleProfiler::run_shooting_profiler().

void BubbleProfiler::Shooting::set_max_iterations ( int  i)
inline

Maximum number of iterations of shooting method

Examples:
action.cpp.

Definition at line 76 of file shooting.hpp.

References iter_max.

Referenced by action(), main(), and BubbleProfiler::run_shooting_profiler().

void BubbleProfiler::Shooting::set_max_periods ( double  p)
inline

Maximum number periods before the ODE integration aborts

Examples:
action.cpp.

Definition at line 78 of file shooting.hpp.

References periods_max.

Referenced by action(), main(), and BubbleProfiler::run_shooting_profiler().

void BubbleProfiler::Shooting::set_shooting_abs_tol ( double  tol)
inline

Absolute tolerance in ODE solver for a shot

Examples:
action.cpp.

Definition at line 56 of file shooting.hpp.

References shoot_ode_abs.

Referenced by action(), main(), and BubbleProfiler::run_shooting_profiler().

void BubbleProfiler::Shooting::set_shooting_rel_tol ( double  tol)
inline

Relative tolerance in ODE solver for a shot

Examples:
action.cpp.

Definition at line 58 of file shooting.hpp.

References shoot_ode_rel.

Referenced by action(), main(), and BubbleProfiler::run_shooting_profiler().

void BubbleProfiler::Shooting::set_y_max ( double  y)
inline

Threshold for asymptotic approximation in analytic evolution

Examples:
action.cpp.

Definition at line 84 of file shooting.hpp.

References get_bubble_profile(), get_euclidean_action(), potential, quadratic_potential_solution(), solve(), and y_max.

Referenced by action(), main(), and BubbleProfiler::run_shooting_profiler().

double BubbleProfiler::Shooting::shoot ( double  lambda)
private

Performs a single shot of our shooting method.

Evolves field from an initial guess for \(\phi_0\), given by \(\lambda\) and tests whether it under or over shoots the false vacuum.

Parameters
lambda\(\lambda\)
Returns
Whether overshot or undershot

Definition at line 288 of file shooting.cpp.

References BubbleProfiler::Abs(), drho_guess, energy(), evolve(), false_min, iter_max, ode(), periods_max, scale, shoot_stepper, sign_min, and unmap().

Referenced by shooting().

double BubbleProfiler::Shooting::shooting ( )
private

Solve for \(\phi_0\) by bisection.

Bisect between over and under shots, beginning from the interval \((\phi_b, \phi_t)\).

Returns
\(\lambda\) corresponding to \(\phi_0\)

Definition at line 380 of file shooting.cpp.

References bisect_lambda_max, iter_max, shoot(), and shoot_bisect_bits.

Referenced by solve().

void BubbleProfiler::Shooting::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) 
)

This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts. Solver that accepts potential and first and second derivatives

Examples:
action.cpp, and scale.cpp.

Definition at line 81 of file shooting.cpp.

References action(), action_computed, action_ode_abs, action_ode_rel, action_stepper, calculate_bubble_profile(), euclidean_action, initialize_potential_parameters(), lambda_sol, profile_computed, profiles, shoot_ode_abs, shoot_ode_rel, shoot_stepper, and shooting().

Referenced by action(), main(), BubbleProfiler::run_shooting_profiler(), set_y_max(), solve(), BubbleProfiler::Shooting_profile_guesser::solve_one_dimensional(), and BubbleProfiler::solve_shooting().

void BubbleProfiler::Shooting::solve ( const Potential potential,
double  false_min_,
double  true_min_,
double  barrier_,
int  dim_ = 3,
unsigned int  options = (Solver_options::Compute_action |                                      Solver_options::Compute_profile) 
)

This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts. Solver that accepts a Potential object

Definition at line 128 of file shooting.cpp.

References BubbleProfiler::Potential::get_number_of_fields(), BubbleProfiler::Potential::partial(), and solve().

double BubbleProfiler::Shooting::unmap ( double  lambda) const
private

Member Data Documentation

double BubbleProfiler::Shooting::action_arrived_rel {1.e-3}
private

Definition at line 167 of file shooting.hpp.

Referenced by action(), calculate_bubble_profile(), and set_action_arrived_rel().

bool BubbleProfiler::Shooting::action_computed {false}
private

Definition at line 131 of file shooting.hpp.

Referenced by get_euclidean_action(), and solve().

double BubbleProfiler::Shooting::action_ode_abs {1.e-6}
private

Definition at line 170 of file shooting.hpp.

Referenced by calculate_bubble_profile(), set_action_abs_tol(), and solve().

double BubbleProfiler::Shooting::action_ode_rel {1.e-6}
private

Definition at line 171 of file shooting.hpp.

Referenced by calculate_bubble_profile(), set_action_rel_tol(), and solve().

controlled_stepper_type BubbleProfiler::Shooting::action_stepper {}
private

Definition at line 136 of file shooting.hpp.

Referenced by action(), and solve().

double BubbleProfiler::Shooting::barrier {0.}
private

Location of the barrier, \(\phi_b\)

Definition at line 144 of file shooting.hpp.

Referenced by bubble_scale(), calculate_roll_time(), initialize_potential_parameters(), and unmap().

double BubbleProfiler::Shooting::bisect_lambda_max {5.}
private

Definition at line 174 of file shooting.hpp.

Referenced by set_bisect_lambda_max(), and shooting().

int BubbleProfiler::Shooting::dim {3}
private
double BubbleProfiler::Shooting::drho_frac {1.e-3}
private

Definition at line 172 of file shooting.hpp.

Referenced by initialize_potential_parameters(), and set_drho_frac().

double BubbleProfiler::Shooting::drho_guess {0.}
private

Guess of step size, \(d\rho\)

Definition at line 139 of file shooting.hpp.

Referenced by action(), calculate_bubble_profile(), initialize_potential_parameters(), and shoot().

double BubbleProfiler::Shooting::euclidean_action {0.}
private

Solution for Euclidean action, \(S\)

Definition at line 129 of file shooting.hpp.

Referenced by get_euclidean_action(), and solve().

double BubbleProfiler::Shooting::evolve_change_rel {1.e-2}
private
double BubbleProfiler::Shooting::f_y_max {1.e6}
private

Definition at line 177 of file shooting.hpp.

Referenced by calculate_roll_time(), calculate_roll_velocity(), and set_f_y_max().

double BubbleProfiler::Shooting::f_y_min {1.e-3}
private

Definition at line 178 of file shooting.hpp.

Referenced by calculate_roll_velocity(), and set_f_y_min().

double BubbleProfiler::Shooting::false_min {0.}
private

Location of the false minimum, \(\phi_f\)

Definition at line 146 of file shooting.hpp.

Referenced by action(), calculate_bubble_profile(), calculate_roll_time(), calculate_roll_velocity(), evolve(), initialize_potential_parameters(), interior(), and shoot().

int BubbleProfiler::Shooting::iter_max {100000}
private

Definition at line 175 of file shooting.hpp.

Referenced by action(), calculate_bubble_profile(), set_max_iterations(), shoot(), and shooting().

double BubbleProfiler::Shooting::lambda_sol {0.}
private

Solution for \(\lambda\)

Definition at line 127 of file shooting.hpp.

Referenced by quadratic_potential_solution(), and solve().

double BubbleProfiler::Shooting::p_barrier {0.}
private

Potential at barrier, \(V(\phi_b)\)

Definition at line 152 of file shooting.hpp.

Referenced by initialize_potential_parameters().

double BubbleProfiler::Shooting::p_false_min {0.}
private

Potential at false minimum, \(V(\phi_f)\)

Definition at line 154 of file shooting.hpp.

Referenced by energy(), and initialize_potential_parameters().

double BubbleProfiler::Shooting::p_true_min {0.}
private

Potential at true minimum, \(V(\phi_t)\)

Definition at line 156 of file shooting.hpp.

Referenced by initialize_potential_parameters().

double BubbleProfiler::Shooting::periods_max {1.e2}
private

Definition at line 176 of file shooting.hpp.

Referenced by action(), calculate_bubble_profile(), set_max_periods(), and shoot().

std::function<double(double)> BubbleProfiler::Shooting::potential {nullptr}
private

Potential, \(V(\phi)\)

Definition at line 159 of file shooting.hpp.

Referenced by energy(), initialize_potential_parameters(), and set_y_max().

std::function<double(double)> BubbleProfiler::Shooting::potential_first {nullptr}
private

Derivative of potential, \(V^{\prime}(\phi)\)

Definition at line 161 of file shooting.hpp.

Referenced by evolve(), initialize_potential_parameters(), interior(), ode(), and quadratic_potential_solution().

std::function<double(double)> BubbleProfiler::Shooting::potential_second {nullptr}
private

Second derivative of potential, \(V^{\prime\prime}(\phi)\)

Definition at line 163 of file shooting.hpp.

Referenced by bubble_scale(), evolve(), initialize_potential_parameters(), interior(), and quadratic_potential_solution().

bool BubbleProfiler::Shooting::profile_computed {false}
private

Definition at line 132 of file shooting.hpp.

Referenced by get_bubble_profile(), and solve().

Field_profiles BubbleProfiler::Shooting::profiles {}
private

Definition at line 130 of file shooting.hpp.

Referenced by get_bubble_profile(), and solve().

double BubbleProfiler::Shooting::scale {0.}
private

Estimated scale of bubble

Definition at line 141 of file shooting.hpp.

Referenced by action(), calculate_bubble_profile(), initialize_potential_parameters(), and shoot().

int BubbleProfiler::Shooting::shoot_bisect_bits {5}
private

Definition at line 166 of file shooting.hpp.

Referenced by set_bisection_precision_bits(), and shooting().

double BubbleProfiler::Shooting::shoot_ode_abs {1.e-4}
private

Definition at line 168 of file shooting.hpp.

Referenced by set_shooting_abs_tol(), and solve().

double BubbleProfiler::Shooting::shoot_ode_rel {1.e-4}
private

Definition at line 169 of file shooting.hpp.

Referenced by set_shooting_rel_tol(), and solve().

controlled_stepper_type BubbleProfiler::Shooting::shoot_stepper {}
private

Definition at line 135 of file shooting.hpp.

Referenced by shoot(), and solve().

double BubbleProfiler::Shooting::sign_min {1.}
private

\(\textrm{sign}\,\phi_f - \phi_t\)

Definition at line 150 of file shooting.hpp.

Referenced by action(), calculate_bubble_profile(), calculate_roll_velocity(), evolve(), initialize_potential_parameters(), and shoot().

double BubbleProfiler::Shooting::true_min {1.}
private

Location of the true minimum, \(\phi_t\)

Definition at line 148 of file shooting.hpp.

Referenced by calculate_roll_time(), initialize_potential_parameters(), quadratic_potential_solution(), and unmap().

double BubbleProfiler::Shooting::y_max {1.e1}
private

Definition at line 179 of file shooting.hpp.

Referenced by interior(), and set_y_max().


The documentation for this class was generated from the following files: