BubbleProfiler  0.3.0
by Peter Athron, Csaba Balazs, Michael Bardsley, Andrew Fowlie, Dylan Harries & Graham White
Namespaces | Classes | Typedefs | Enumerations | Functions | Variables
BubbleProfiler Namespace Reference

Namespaces

 detail
 
 logging
 

Classes

struct  Action_options
 
class  Algebraic_potential
 Implements a multi dimensional scalar field potential, which may be expressed as almost any algerbraic fucntion. More...
 
struct  Bubble_profiler_inputs
 
class  BVP_solver_error
 Exception indicating failure to integrate BVP. More...
 
class  Controlled_step_size_integrator
 
class  Domain_error
 Exception indicating function evaluation out of allowed domain. More...
 
class  Dummy_observer
 A default observer that does not perform any actions when called. More...
 
class  Error
 
class  Field_profiles
 Discretized set of field profiles. More...
 
class  Fixed_step_size_integrator
 
struct  Fubini_options
 
class  Gaussian_potential
 Potential composed of two multivariate Gaussian functions. More...
 
class  Generalized_fubini_observer
 
class  Generalized_fubini_potential
 
class  GSL_interpolator
 Provides routines for carrying out 1D interpolation. More...
 
class  GSL_root_finder
 
class  Instream_profile_guesser
 Provides guess for field profiles read from an input stream. More...
 
class  IO_error
 Exception indicating generic input/output error. More...
 
class  Kink_profile_guesser
 Computes an initial guess for the bubble profile using a kink ansatz. More...
 
class  Logarithmic_observer
 
struct  Logarithmic_options
 
class  Logarithmic_potential
 
class  NLopt_optimizer
 
class  No_convergence_error
 Exception indicating failure of iterative procedure to converge. More...
 
class  Not_implemented_error
 
class  Numerical_error
 Exception indicating generic numerical error. More...
 
class  Optimizer_error
 Exception indicating that the optimizer failed to converge. More...
 
class  Out_of_bounds_error
 Exception indicating out of bounds access error. More...
 
class  Out_of_memory_error
 Exception indicating that memory allocation failed. More...
 
class  Perturbations_ODE_system
 Class implementing the system of ODEs obeyed by profile perturbations. More...
 
class  Perturbative_profiler
 Bounce solver using perturbative method. More...
 
class  Potential
 Abstract base class for a generic potential. More...
 
class  Profile_convergence_tester
 
class  Profile_guesser
 Abstract class to represent ansatz generators. More...
 
class  Profile_observer
 Observes the profile during iteration and writes out profile, perturbations and action for each step to file. More...
 
struct  Profiler_result
 
class  RAII_guard
 Carries out provided clean-up actions at destruction. More...
 
class  Relative_convergence_tester
 
class  Restricted_quartic_potential
 
class  Root_finder
 
struct  Scale_options
 
class  Setup_error
 Exception indicating general setup error. More...
 
class  Shooting
 Solve the one-dimensional problem using the shooting method. More...
 
class  Shooting_profile_guesser
 Constructs a solution ansatz using the overshoot/undershoot method. More...
 
struct  Shooting_settings
 
struct  State_algebra_dispatcher
 
struct  Tabulate_options
 
class  Thin_wall_error
 Exception indicating a thin-wall ansatz which we can't solve (yet!) More...
 
class  Thin_wall_observer
 
struct  Thin_wall_options
 
class  Thin_wall_potential
 Example thin wall potential as given by Coleman. More...
 

Typedefs

using Constant_step_size_RK4 = Fixed_step_size_integrator< double, Eigen::VectorXd >
 
using Constant_step_size_euler = Fixed_step_size_integrator< double, Eigen::VectorXd, boost::numeric::odeint::euler< Eigen::VectorXd, double, Eigen::VectorXd, double, State_algebra_dispatcher< Eigen::VectorXd >::algebra_type > >
 
using Controlled_step_size_RK4 = Controlled_step_size_integrator< double, Eigen::VectorXd, boost::numeric::odeint::runge_kutta_cash_karp54< Eigen::VectorXd, double, Eigen::VectorXd, double, State_algebra_dispatcher< Eigen::VectorXd >::algebra_type > >
 
using Controlled_step_size_RKD5 = Controlled_step_size_integrator< double, Eigen::VectorXd, boost::numeric::odeint::runge_kutta_dopri5< Eigen::VectorXd, double, Eigen::VectorXd, double, State_algebra_dispatcher< Eigen::VectorXd >::algebra_type > >
 
using Euler_perturbative_profiler = Perturbative_profiler< Constant_step_size_euler >
 
using RK4_perturbative_profiler = Perturbative_profiler< Constant_step_size_RK4 >
 
using CRK4_perturbative_profiler = Perturbative_profiler< Controlled_step_size_RK4 >
 
using CRKD5_perturbative_profiler = Perturbative_profiler< Controlled_step_size_RKD5 >
 
using Profiler_data = std::map< std::string, std::vector< Profiler_result > >
 

Enumerations

enum  Integration_rule {
  Integration_rule::GK15, Integration_rule::GK21, Integration_rule::GK31, Integration_rule::GK41,
  Integration_rule::GK51, Integration_rule::GK61
}
 
enum  Root_finder_status { Root_finder_status::SUCCESS, Root_finder_status::FAIL }
 
enum  Integration_algorithm : int { Integration_algorithm::Runge_kutta_4 = 0, Integration_algorithm::Euler }
 

Functions

double calculate_kinetic_action (const Potential &potential, const Field_profiles &profiles, std::size_t max_intervals=1000, double rel_tol=1.e-4, double abs_tol=1.e-4, Integration_rule rule=Integration_rule::GK31)
 Calculates the action using only the kinetic term contributions. More...
 
double calculate_potential_action (const Potential &potential, const Field_profiles &profiles, std::size_t max_intervals=1000, double rel_tol=1.e-4, double abs_tol=1.e-4, Integration_rule rule=Integration_rule::GK31)
 Calculates the action using only the potential term contributions. More...
 
double calculate_full_action (const Potential &potential, const Field_profiles &profiles, std::size_t max_intervals=1000, double rel_tol=1.e-4, double abs_tol=1.e-4, Integration_rule rule=Integration_rule::GK31)
 Calculates the action using all terms in the action integrand. More...
 
double calculate_action (const Potential &potential, const Field_profiles &profiles, std::size_t max_intervals=1000, double rel_tol=1.e-4, double abs_tol=1.e-4, Integration_rule rule=Integration_rule::GK31, bool use_kinetic=true)
 Calculates the action using either the kinetic or potential terms. More...
 
Eigen::VectorXd interpolate_f_at (const Eigen::VectorXd &x_old, const Eigen::VectorXd &y_old, const Eigen::VectorXd &x_new)
 Computes interpolated function values at a new set of grid points. More...
 
Eigen::VectorXd interpolate_f_at (const Eigen::VectorXd &x_old, const Eigen::VectorXd &y_old, int n_old, const Eigen::VectorXd &x_new, int n_new)
 Computes interpolated function values at a new set of grid points. More...
 
int get_integration_rule_key (Integration_rule)
 
template<typename F >
std::tuple< int, double, double > integrate_gsl_qag (F f, double a, double b, double abs_tol, double rel_tol, std::size_t max_intervals, Integration_rule rule)
 
template<typename F >
std::tuple< int, double, double > integrate_gsl_qags (F f, double a, double b, double abs_tol, double rel_tol, std::size_t max_intervals)
 
template<typename F >
std::tuple< int, double, double, int > integrate_gsl_cquad (F f, double a, double b, double abs_tol, double rel_tol, std::size_t n_intervals)
 
template<typename T >
Abs (T x) noexcept
 
template<typename T >
Abs (const std::complex< T > &z) noexcept
 
template<typename Order , typename Argument >
auto BesselI (Order v, Argument z) -> decltype(boost::math::cyl_bessel_i(v, z))
 
template<typename Order , typename Argument >
auto BesselJ (Order v, Argument z) -> decltype(boost::math::cyl_bessel_j(v, z))
 
template<typename Order , typename Argument >
auto BesselK (Order v, Argument z) -> decltype(boost::math::cyl_bessel_k(v, z))
 
template<typename Order , typename Argument >
auto BesselY (Order v, Argument z) -> decltype(boost::math::cyl_neumann(v, z))
 
template<typename T , typename = typename std::enable_if< std::is_floating_point<T>::value,T>::type>
Exp (T x) noexcept
 
template<typename T , typename = typename std::enable_if< std::is_integral<T>::value,T>::type>
double Exp (T n) noexcept
 
template<typename T >
std::complex< T > Exp (const std::complex< T > &z) noexcept
 
template<typename T , typename = typename std::enable_if< std::is_floating_point<T>::value,T>::type>
Log (T x) noexcept
 
template<typename T , typename = typename std::enable_if< std::is_integral<T>::value,T>::type>
double Log (T n) noexcept
 
template<typename T >
std::complex< T > Log (const std::complex< T > &z) noexcept
 
template<typename T >
int Sign (const T &x)
 
template<typename T , typename = typename std::enable_if< std::is_floating_point<T>::value,T>::type>
Cos (T x) noexcept
 
template<typename T , typename = typename std::enable_if< std::is_integral<T>::value,T>::type>
double Cos (T n) noexcept
 
template<typename T >
std::complex< T > Cos (const std::complex< T > &z) noexcept
 
template<typename T , typename = typename std::enable_if< std::is_floating_point<T>::value,T>::type>
Sin (T x) noexcept
 
template<typename T , typename = typename std::enable_if< std::is_integral<T>::value,T>::type>
double Sin (T n) noexcept
 
template<typename T >
std::complex< T > Sin (const std::complex< T > &z) noexcept
 
template<typename T , typename = typename std::enable_if< std::is_floating_point<T>::value,T>::type>
Tan (T x) noexcept
 
template<typename T , typename = typename std::enable_if< std::is_integral<T>::value,T>::type>
double Tan (T n) noexcept
 
template<typename T >
std::complex< T > Tan (const std::complex< T > &z) noexcept
 
template<typename T >
auto Cot (T x) noexcept-> decltype(Tan(x))
 
template<typename T >
auto Csc (T x) noexcept-> decltype(Sin(x))
 
template<typename T >
auto Sec (T x) noexcept-> decltype(Cos(x))
 
template<typename T , typename = typename std::enable_if< std::is_floating_point<T>::value,T>::type>
ArcCos (T x) noexcept
 
template<typename T , typename = typename std::enable_if< std::is_integral<T>::value,T>::type>
double ArcCos (T n) noexcept
 
template<typename T >
std::complex< T > ArcCos (const std::complex< T > &z) noexcept
 
template<typename T , typename = typename std::enable_if< std::is_floating_point<T>::value,T>::type>
ArcSin (T x) noexcept
 
template<typename T , typename = typename std::enable_if< std::is_integral<T>::value,T>::type>
double ArcSin (T n) noexcept
 
template<typename T >
std::complex< T > ArcSin (const std::complex< T > &z) noexcept
 
template<typename T , typename = typename std::enable_if< std::is_floating_point<T>::value,T>::type>
ArcTan (T x) noexcept
 
template<typename T , typename = typename std::enable_if< std::is_integral<T>::value,T>::type>
double ArcTan (T n) noexcept
 
template<typename T >
std::complex< T > ArcTan (const std::complex< T > &z) noexcept
 
template<typename T , typename = typename std::enable_if< std::is_floating_point<T>::value,T>::type>
Cosh (T x) noexcept
 
template<typename T , typename = typename std::enable_if< std::is_integral<T>::value,T>::type>
double Cosh (T n) noexcept
 
template<typename T >
std::complex< T > Cosh (const std::complex< T > &z) noexcept
 
template<typename T , typename = typename std::enable_if< std::is_floating_point<T>::value,T>::type>
Sinh (T x) noexcept
 
template<typename T , typename = typename std::enable_if< std::is_integral<T>::value,T>::type>
double Sinh (T n) noexcept
 
template<typename T >
std::complex< T > Sinh (const std::complex< T > &z) noexcept
 
template<typename T , typename = typename std::enable_if< std::is_floating_point<T>::value,T>::type>
Tanh (T x) noexcept
 
template<typename T , typename = typename std::enable_if< std::is_integral<T>::value,T>::type>
double Tanh (T n) noexcept
 
template<typename T >
std::complex< T > Tanh (const std::complex< T > &z) noexcept
 
template<typename T >
auto Coth (T x) noexcept-> decltype(Tanh(x))
 
template<typename T >
auto Csch (T x) noexcept-> decltype(Sinh(x))
 
template<typename T >
auto Sech (T x) noexcept-> decltype(Cosh(x))
 
template<typename T , typename = typename std::enable_if< std::is_floating_point<T>::value,T>::type>
ArcCosh (T x) noexcept
 
template<typename T , typename = typename std::enable_if< std::is_integral<T>::value,T>::type>
double ArcCosh (T n) noexcept
 
template<typename T >
std::complex< T > ArcCosh (const std::complex< T > &z) noexcept
 
template<typename T , typename = typename std::enable_if< std::is_floating_point<T>::value,T>::type>
ArcSinh (T x) noexcept
 
template<typename T , typename = typename std::enable_if< std::is_integral<T>::value,T>::type>
double ArcSinh (T n) noexcept
 
template<typename T >
std::complex< T > ArcSinh (const std::complex< T > &z) noexcept
 
template<typename T , typename = typename std::enable_if< std::is_floating_point<T>::value,T>::type>
ArcTanh (T x) noexcept
 
template<typename T , typename = typename std::enable_if< std::is_integral<T>::value,T>::type>
double ArcTanh (T n) noexcept
 
template<typename T >
std::complex< T > ArcTanh (const std::complex< T > &z) noexcept
 
template<typename T , typename = typename std::enable_if< std::is_arithmetic<T>::value,T>::type>
Im (T) noexcept
 
template<typename T >
Im (const std::complex< T > &z) noexcept
 
template<typename T , typename = typename std::enable_if< std::is_arithmetic<T>::value,T>::type>
Re (T x) noexcept
 
template<typename T >
Re (const std::complex< T > &z) noexcept
 
template<typename T , typename = typename std::enable_if< std::is_floating_point<T>::value,T>::type>
Sqrt (T x) noexcept
 
template<typename T , typename = typename std::enable_if< std::is_integral<T>::value,T>::type>
double Sqrt (T n) noexcept
 
bool optimization_succeeded (nlopt::result)
 
double area_n_sphere (int n)
 
double asinch (double a)
 
double asinc (double a)
 
double approx_root_pos_4 (double a)
 
double approx_root_neg_4 (double a)
 
template<typename F >
constexpr RAII_guard< F > make_raii_guard (F f)
 
Eigen::MatrixXd calculate_rotation_to_target (const Eigen::VectorXd &target)
 Calculate the rotation matrix for a rotation to the target vector. More...
 
template<typename Scalar , typename Data >
Scalar linear_lagrange_interpolation_at_point (Scalar x, const Data &xdata, const Data &ydata)
 
template<typename Scalar , typename Data >
Scalar quadratic_lagrange_interpolation_at_point (Scalar x, const Data &xdata, const Data &ydata)
 
double Wm1 (double a)
 
double series (double a)
 
void validate (boost::any &v, const std::vector< std::string > &values, Integration_algorithm *, int)
 
Bubble_profiler_inputs parse_cmd_line_args (int argc, const char *argv[])
 
void write_action (double action, std::ostream &ostr)
 Write the value of the action to requested output stream. More...
 
void write_profiles (const std::vector< std::string > &fields, const Field_profiles &profiles, std::ostream &ostr)
 Write the given field profiles to requested output stream. More...
 
Eigen::VectorXd find_global_min (const Potential &potential, double opt_timeout)
 Find the global minimum of the given potential. More...
 
Eigen::VectorXd find_one_dimensional_barrier (const Potential &potential, const Eigen::VectorXd &true_vacuum_loc, const Eigen::VectorXd &false_vacuum_loc, double opt_timeout)
 Find barrier directly between two minima of a 1D potential. More...
 
void initialize_extrema (const Potential &potential, const Bubble_profiler_inputs &input, Eigen::VectorXd &true_vacuum_loc, Eigen::VectorXd &false_vacuum_loc)
 Locate the requested critical points if not already provided. More...
 
void initialize_extrema (const Potential &potential, const Bubble_profiler_inputs &input, Eigen::VectorXd &true_vacuum_loc, Eigen::VectorXd &false_vacuum_loc, Eigen::VectorXd &barrier_loc)
 Locate the necessary critical points if not already provided. More...
 
std::tuple< double, Field_profilesrun_shooting_profiler (const Bubble_profiler_inputs &input)
 Calculates the action and bubble profile using the shooting method. More...
 
template<class Profiler >
std::tuple< double, Field_profilesrun_perturbative_profiler (const Bubble_profiler_inputs &input, Profiler &profiler)
 Calculates the action and bubble profile using the perturbative method. More...
 
std::tuple< double, Field_profilesrun_perturbative_profiler (const Bubble_profiler_inputs &input)
 Calculates the action and bubble profile using the perturbative method. More...
 
void run_profiler (const Bubble_profiler_inputs &input)
 Calculates the action and bubble profile for the given potential. More...
 
void configure_logging (const Bubble_profiler_inputs &input)
 
void print_usage ()
 
bool starts_with (const std::string &option, const std::string &prefix)
 
std::string get_option_value (const std::string &option, const std::string &sep="=")
 
void parse_positional_args (const std::vector< std::string > &args, Fubini_options &options)
 
double parse_global_min (const std::string &value)
 
std::tuple< double, Field_profilessolve_shooting (const Generalized_fubini_potential &potential, double true_min, bool calculate_profile)
 
std::tuple< double, Field_profilessolve_perturbative (Generalized_fubini_potential &potential, double true_min, const std::string &observer_file)
 
void write_profile (std::ostream &ostr, const Field_profiles &numerical_profile, const Field_profiles &exact_profile)
 
void write_action (std::ostream &ostr, double numerical_action, double exact_action)
 
void parse_positional_args (const std::vector< std::string > &args, Logarithmic_options &options)
 
std::tuple< double, Field_profilessolve_shooting (const Logarithmic_potential &potential, double true_min, bool calculate_profile)
 
std::tuple< double, Field_profilessolve_perturbative (Logarithmic_potential &potential, double true_min, const std::string &observer_file)
 
void parse_positional_args (const std::vector< std::string > &args, Action_options &options)
 
void parse_positional_args (const std::vector< std::string > &args, Tabulate_options &options)
 
Profiler_result run_shooting_profiler (const Restricted_quartic_potential &potential, int dim)
 
Profiler_result run_perturbative_profiler (Restricted_quartic_potential potential, int dim)
 
Profiler_result run_profiler (const std::string &profiler, const Restricted_quartic_potential &potential, int dim)
 
void write_table_header (const std::vector< std::string > &profiler_labels)
 
void write_table (int dim, double E, const std::vector< double > &alpha_data, const Profiler_data &profiler_data)
 
void parse_positional_args (const std::vector< std::string > &args, Scale_options &options)
 
void parse_positional_args (const std::vector< std::string > &args, Thin_wall_options &options)
 
std::tuple< double, Field_profilessolve_shooting (const Thin_wall_potential &potential, int dim, bool calculate_profile)
 
std::tuple< double, Field_profilessolve_perturbative (Thin_wall_potential &potential, int dim, const std::string &observer_file)
 
void write_profile (std::ostream &ostr, const Field_profiles &numerical_profile)
 

Variables

const std::string default_ansatz_file = ""
 
const double default_domain_start = -1.
 
const double default_domain_end = -1.
 
const double default_initial_step_size = 1.e-2
 
const int default_interp_fraction = 1.0
 
const Integration_algorithm default_algorithm
 
const double default_rtol_action = 1.e-3
 
const double default_rtol_fields = 1.e-3
 
const std::string default_output_file = "-"
 
const bool default_write_perturbations = true
 
const bool default_force_output = false
 
const double default_opt_timeout = 10.
 
const int default_max_iterations = -1
 
const bool default_assume_origin_false_vacuum = false
 
const bool default_use_perturbative = false
 
const bool default_write_profiles = false
 
const bool default_shooting_ansatz = false
 
const bool default_verbose = false
 
const int default_shoot_bisect_bits = 5
 
const double default_action_arrived_rel = 1.e-3
 
const double default_shoot_ode_abs = 1.e-4
 
const double default_shoot_ode_rel = 1.e-4
 
const double default_action_ode_abs = 1.e-6
 
const double default_action_ode_rel = 1.e-6
 
const double default_drho_frac = 1.e-3
 
const double default_evolve_change_rel = 1.e-2
 
const double default_bisect_lambda_max = 5.0
 
const int default_iter_max = 100000
 
const double default_periods_max = 1.e2
 
const double default_f_y_max = 1.e6
 
const double default_f_y_min = 1.e-3
 
const double default_y_max = 1.e1
 

Typedef Documentation

using BubbleProfiler::Constant_step_size_euler = typedef Fixed_step_size_integrator< double, Eigen::VectorXd, boost::numeric::odeint::euler< Eigen::VectorXd, double, Eigen::VectorXd, double, State_algebra_dispatcher<Eigen::VectorXd>::algebra_type> >

Definition at line 35 of file default_integration_policy.hpp.

using BubbleProfiler::Constant_step_size_RK4 = typedef Fixed_step_size_integrator<double, Eigen::VectorXd>

Definition at line 27 of file default_integration_policy.hpp.

using BubbleProfiler::Controlled_step_size_RK4 = typedef Controlled_step_size_integrator< double, Eigen::VectorXd, boost::numeric::odeint::runge_kutta_cash_karp54< Eigen::VectorXd, double, Eigen::VectorXd, double, State_algebra_dispatcher<Eigen::VectorXd>::algebra_type> >

Definition at line 43 of file default_integration_policy.hpp.

using BubbleProfiler::Controlled_step_size_RKD5 = typedef Controlled_step_size_integrator< double, Eigen::VectorXd, boost::numeric::odeint::runge_kutta_dopri5< Eigen::VectorXd, double, Eigen::VectorXd, double, State_algebra_dispatcher<Eigen::VectorXd>::algebra_type> >

Definition at line 51 of file default_integration_policy.hpp.

Definition at line 30 of file perturbative_profiler.hpp.

Definition at line 32 of file perturbative_profiler.hpp.

Definition at line 26 of file perturbative_profiler.hpp.

using BubbleProfiler::Profiler_data = typedef std::map<std::string, std::vector<Profiler_result> >
Examples:
tabulate.cpp.

Definition at line 77 of file tabulate.cpp.

Enumeration Type Documentation

Enumerator
Runge_kutta_4 
Euler 

Definition at line 57 of file run_cmd_line_potential.cpp.

Enumerator
GK15 
GK21 
GK31 
GK41 
GK51 
GK61 

Definition at line 31 of file integration_utils.hpp.

Enumerator
SUCCESS 
FAIL 

Definition at line 25 of file root_finder.hpp.

Function Documentation

template<typename T >
T BubbleProfiler::Abs ( x)
noexcept

Definition at line 34 of file math_wrappers.hpp.

References Eigen::abs().

Referenced by BubbleProfiler::Shooting::action(), BubbleProfiler::Shooting::calculate_bubble_profile(), BubbleProfiler::Shooting_profile_guesser::calculate_potential_parameters(), BubbleProfiler::Kink_profile_guesser::calculate_potential_parameters(), BubbleProfiler::Shooting::calculate_roll_time(), BubbleProfiler::Shooting::calculate_roll_velocity(), BubbleProfiler::Shooting_profile_guesser::compute_vacuum_distance(), BubbleProfiler::Kink_profile_guesser::compute_vacuum_distance(), BubbleProfiler::Shooting::evolve(), BubbleProfiler::Logarithmic_potential::first_deriv(), BubbleProfiler::Generalized_fubini_potential::first_deriv(), BubbleProfiler::Thin_wall_potential::get_global_minimum_location(), BubbleProfiler::Shooting_profile_guesser::get_large_distance_solution(), BubbleProfiler::Logarithmic_potential::get_local_maximum_location(), BubbleProfiler::Thin_wall_potential::get_local_maximum_location(), BubbleProfiler::Thin_wall_potential::get_local_minimum_location(), BubbleProfiler::Instream_profile_guesser::get_profile_guess(), BubbleProfiler::Kink_profile_guesser::guess_domain_start(), BubbleProfiler::Shooting::interior(), BubbleProfiler::Logarithmic_potential::operator()(), BubbleProfiler::Generalized_fubini_potential::operator()(), BubbleProfiler::Shooting::quadratic_potential_solution(), BubbleProfiler::Relative_convergence_tester::relative_difference(), BubbleProfiler::Generalized_fubini_potential::second_deriv(), BubbleProfiler::Logarithmic_potential::second_deriv(), BubbleProfiler::Shooting::shoot(), and BubbleProfiler::Thin_wall_potential::Thin_wall_potential().

template<typename T >
T BubbleProfiler::Abs ( const std::complex< T > &  z)
noexcept

Definition at line 40 of file math_wrappers.hpp.

References Eigen::abs().

double BubbleProfiler::approx_root_neg_4 ( double  a)
Returns
Root of \( J_1(x) / x + a - 0.5 = 0 \)

See e.g., in Mathematica,

Numerical[a_] := x /. FindRoot[BesselJ[1, x]/x - 0.5 == -a, {x, 2}]
Taylor[a_] := 4. Sqrt[a]
Plot[{Numerical[a], Taylor[a]}, {a, 0, 1}]

Definition at line 108 of file numeric.cpp.

References Sqrt().

Referenced by BubbleProfiler::Shooting::calculate_roll_time().

double BubbleProfiler::approx_root_pos_4 ( double  a)
Returns
Root of \(I_1(x) / x - a - 0.5 = 0\)

This is used in dimension \(d=4\) case. I use Taylor expansions and a large \(x\) approximation of \(I_1(x)\).

See e.g., in Mathematica,

Numerical[a_] := x /. FindRoot[BesselI[1, x] / x - 0.5 == a, {x, 10.}]
Taylor[a_] := 4. Sqrt[a]
Lambert[a_] := -(3/2) LambertW[-1, (2 (2/\[Pi])^(1/3))/(3 (-(1 + 2 a)^2)^(1/3))] // N // Re
Plot[{Numerical[a], Taylor[a], Lambert[a]}, {a, 0, 100}]

Definition at line 96 of file numeric.cpp.

References Sqrt(), and Wm1().

Referenced by BubbleProfiler::Shooting::calculate_roll_time().

template<typename T , typename = typename std::enable_if< std::is_floating_point<T>::value,T>::type>
T BubbleProfiler::ArcCos ( x)
noexcept

Definition at line 210 of file math_wrappers.hpp.

References sm_plus_singlet::type.

template<typename T , typename = typename std::enable_if< std::is_integral<T>::value,T>::type>
double BubbleProfiler::ArcCos ( n)
noexcept

Definition at line 218 of file math_wrappers.hpp.

template<typename T >
std::complex<T> BubbleProfiler::ArcCos ( const std::complex< T > &  z)
noexcept

Definition at line 224 of file math_wrappers.hpp.

References sm_plus_singlet::type.

template<typename T , typename = typename std::enable_if< std::is_floating_point<T>::value,T>::type>
T BubbleProfiler::ArcCosh ( x)
noexcept

Definition at line 360 of file math_wrappers.hpp.

References sm_plus_singlet::type.

template<typename T , typename = typename std::enable_if< std::is_integral<T>::value,T>::type>
double BubbleProfiler::ArcCosh ( n)
noexcept

Definition at line 368 of file math_wrappers.hpp.

template<typename T >
std::complex<T> BubbleProfiler::ArcCosh ( const std::complex< T > &  z)
noexcept

Definition at line 374 of file math_wrappers.hpp.

References sm_plus_singlet::type.

template<typename T , typename = typename std::enable_if< std::is_floating_point<T>::value,T>::type>
T BubbleProfiler::ArcSin ( x)
noexcept
template<typename T , typename = typename std::enable_if< std::is_integral<T>::value,T>::type>
double BubbleProfiler::ArcSin ( n)
noexcept

Definition at line 240 of file math_wrappers.hpp.

template<typename T >
std::complex<T> BubbleProfiler::ArcSin ( const std::complex< T > &  z)
noexcept

Definition at line 246 of file math_wrappers.hpp.

References sm_plus_singlet::type.

template<typename T , typename = typename std::enable_if< std::is_floating_point<T>::value,T>::type>
T BubbleProfiler::ArcSinh ( x)
noexcept

Definition at line 382 of file math_wrappers.hpp.

References sm_plus_singlet::type.

template<typename T , typename = typename std::enable_if< std::is_integral<T>::value,T>::type>
double BubbleProfiler::ArcSinh ( n)
noexcept

Definition at line 390 of file math_wrappers.hpp.

template<typename T >
std::complex<T> BubbleProfiler::ArcSinh ( const std::complex< T > &  z)
noexcept

Definition at line 396 of file math_wrappers.hpp.

References sm_plus_singlet::type.

template<typename T , typename = typename std::enable_if< std::is_floating_point<T>::value,T>::type>
T BubbleProfiler::ArcTan ( x)
noexcept

Definition at line 254 of file math_wrappers.hpp.

References sm_plus_singlet::type.

template<typename T , typename = typename std::enable_if< std::is_integral<T>::value,T>::type>
double BubbleProfiler::ArcTan ( n)
noexcept

Definition at line 262 of file math_wrappers.hpp.

template<typename T >
std::complex<T> BubbleProfiler::ArcTan ( const std::complex< T > &  z)
noexcept

Definition at line 268 of file math_wrappers.hpp.

References sm_plus_singlet::type.

template<typename T , typename = typename std::enable_if< std::is_floating_point<T>::value,T>::type>
T BubbleProfiler::ArcTanh ( x)
noexcept
template<typename T , typename = typename std::enable_if< std::is_integral<T>::value,T>::type>
double BubbleProfiler::ArcTanh ( n)
noexcept

Definition at line 412 of file math_wrappers.hpp.

template<typename T >
std::complex<T> BubbleProfiler::ArcTanh ( const std::complex< T > &  z)
noexcept

Definition at line 418 of file math_wrappers.hpp.

References sm_plus_singlet::type.

double BubbleProfiler::area_n_sphere ( int  n)
Returns
Area of \(n\)-sphere from wiki
Parameters
n\(n\)

Definition at line 29 of file numeric.cpp.

Referenced by BubbleProfiler::Shooting::action(), calculate_full_action(), calculate_kinetic_action(), and calculate_potential_action().

double BubbleProfiler::asinc ( double  a)
Returns
Approximate inverse of the sinc function, \(\sin x / x = a \)

See e.g., in Mathematica

Asinch[a_] := InverseFunction[Sinc][a] // N // Re
f[x_] := 2*x + 3*x^3/10 + 321*x^5/2800 + 3197*x^7/56000 +
445617*x^9/13798400 + 1766784699*x^11/89689600000 +
317184685563*x^13/25113088000000 +
14328608561991*x^15/1707689984000000 +
6670995251837391*x^17/1165411287040000000 +
910588298588385889*x^19/228420612259840000000
Approx[a_] := Sqrt[3./2.] * f[Sqrt[1. - a]]
Plot[{Asinch[a], Approx[a]}, {a, 0, 1}]
Plot[{Sinc[Approx[a]]}, {a, 0, 1}]

Definition at line 92 of file numeric.cpp.

References series(), and Sqrt().

Referenced by BubbleProfiler::Shooting::calculate_roll_time().

double BubbleProfiler::asinch ( double  a)
Returns
Approximate inverse of the hyperbolic sinch function, \(\sinh x / x = a\)

I use Taylor expansions and \(\sinh x \approx e^x / 2\).

See e.g., in Mathematica,

Asinch[a_] := x /. FindRoot[Sinh[x] /x == a, {x, 5}]
Taylor[a_] := Sqrt[6. (a - 1)]
Lambert[a_] := -LambertW[-1, -0.5 / a]
Plot[{Asinch[a], Taylor[a], Lambert[a]}, {a, 1, 10}]

The Taylor expansion is more precise until about \(a \approx 1.5\).

Definition at line 63 of file numeric.cpp.

References Sqrt(), and Wm1().

Referenced by BubbleProfiler::Shooting::calculate_roll_time().

template<typename Order , typename Argument >
auto BubbleProfiler::BesselI ( Order  v,
Argument  z 
) -> decltype(boost::math::cyl_bessel_i(v, z))
template<typename Order , typename Argument >
auto BubbleProfiler::BesselJ ( Order  v,
Argument  z 
) -> decltype(boost::math::cyl_bessel_j(v, z))
template<typename Order , typename Argument >
auto BubbleProfiler::BesselK ( Order  v,
Argument  z 
) -> decltype(boost::math::cyl_bessel_k(v, z))
template<typename Order , typename Argument >
auto BubbleProfiler::BesselY ( Order  v,
Argument  z 
) -> decltype(boost::math::cyl_neumann(v, z))

Definition at line 67 of file math_wrappers.hpp.

References sm_plus_singlet::type.

double BubbleProfiler::calculate_action ( const Potential potential,
const Field_profiles profiles,
std::size_t  max_intervals = 1000,
double  rel_tol = 1.e-4,
double  abs_tol = 1.e-4,
Integration_rule  rule = Integration_rule::GK31,
bool  use_kinetic = true 
)

Calculates the action using either the kinetic or potential terms.

Parameters
potentialthe potential associated with the bounce solution
profilesthe field profiles representing the bounce solution
max_intervalsthe maximum allowed number of subintervals for adaptive quadrature routine
rel_tolthe relative error goal
abs_tolthe absolute error goal
rulethe quadrature rule to use
use_kineticflag selecting whether the kinetic or potential terms should be used to compute the action
Returns
numerical estimate for the Euclidean action

The method used to evaluate the action is chosen based on the value of use_kinetic, with a value of true indicating that the action should be evaluated using only the contribution from the field kinetic terms. If use_kinetic is false, the action is evaluated using the contribution from the potential.

Internally, this simply function calls BubbleProfiler::calculate_kinetic_action or BubbleProfiler::calculate_potential_action, as appropriate.

See also
calculate_kinetic_action, calculate_potential_action

Definition at line 257 of file euclidean_action.cpp.

References calculate_kinetic_action(), and calculate_potential_action().

Referenced by BubbleProfiler::Perturbative_profiler< Integration_policy >::calculate_bubble_profile(), BubbleProfiler::Relative_convergence_tester::is_converged(), and BubbleProfiler::Profile_observer::write_action_to_file().

double BubbleProfiler::calculate_full_action ( const Potential potential,
const Field_profiles profiles,
std::size_t  max_intervals = 1000,
double  rel_tol = 1.e-4,
double  abs_tol = 1.e-4,
Integration_rule  rule = Integration_rule::GK31 
)

Calculates the action using all terms in the action integrand.

Parameters
potentialthe potential associated with the bounce solution
profilesthe field profiles representing the bounce solution
max_intervalsthe maximum allowed number of subintervals for adaptive quadrature routine
rel_tolthe relative error goal
abs_tolthe absolute error goal
rulethe quadrature rule to use
Returns
numerical estimate for the Euclidean action

Given the bounce solution for \(n\) scalar fields in \(d\) dimensions defined on a finite domain, an approximate value for the full Euclidean action

\[ S_E \equiv S_{d-1} \int_0^\infty d\rho \rho^{d-1} \left ( \sum_{i=1}^n \frac{1}{2} \dot{\phi}_i^2 + V(\phi_i) \right ) , \]

where \(S_{d-1}\) is the surface area of the \(d\)-sphere and the potential is assumed to vanish at the false vacuum, is calculated using numerical quadrature, with the quadrature rule chosen according to the value of rule. The adaptive quadrature routine attempts to compute an estimate \(\tilde{I}\) for the integral \(I\) appearing above that satisfies

\[ | \tilde{I} - I | \leq \max ( \epsilon_{\text{abs}}, \epsilon_{\text{rel}} |I|) , \]

with \(\epsilon_{\text{abs}}\) and \(\epsilon_{\text{rel}}\) given by abs_tol and rel_tol, respectively. The maximum allowed number of subintervals used in this integration is set by the value of max_intervals.

See also
integrate_gsl_qag

Definition at line 198 of file euclidean_action.cpp.

References area_n_sphere(), BubbleProfiler::Field_profiles::derivative_at(), BubbleProfiler::Field_profiles::evaluate_at(), BubbleProfiler::Field_profiles::get_domain_end(), BubbleProfiler::Field_profiles::get_domain_start(), BubbleProfiler::Field_profiles::get_number_of_dimensions(), BubbleProfiler::Field_profiles::get_number_of_fields(), integrate_gsl_qag(), and gaussian_alpha_tests::n_fields.

double BubbleProfiler::calculate_kinetic_action ( const Potential potential,
const Field_profiles profiles,
std::size_t  max_intervals = 1000,
double  rel_tol = 1.e-4,
double  abs_tol = 1.e-4,
Integration_rule  rule = Integration_rule::GK31 
)

Calculates the action using only the kinetic term contributions.

Parameters
potentialthe potential associated with the bounce solution
profilesthe field profiles representing the bounce solution
max_intervalsthe maximum allowed number of subintervals for adaptive quadrature routine
rel_tolthe relative error goal
abs_tolthe absolute error goal
rulethe quadrature rule to use
Returns
numerical estimate for the Euclidean action

Given the bounce solution for \(n\) scalar fields in \(d\) dimensions defined on a finite domain, the contribution to the action due to the field kinetic terms,

\[ S_T \equiv S_{d-1} \int_0^\infty d\rho \rho^{d-1} \sum_{i=1}^n \frac{1}{2} \dot{\phi}_i^2 , \]

where \(S_{d-1}\) is the surface area of the \(d\)-sphere, is first calculated using numerical quadrature, with the quadrature rule chosen according to the value of rule. The adaptive quadrature routine attempts to compute an estimate \(\tilde{I}\) for the integral \(I\) appearing above that satisfies

\[ | \tilde{I} - I | \leq \max ( \epsilon_{\text{abs}}, \epsilon_{\text{rel}} |I|) , \]

with \(\epsilon_{\text{abs}}\) and \(\epsilon_{\text{rel}}\) given by abs_tol and rel_tol, respectively. The maximum allowed number of subintervals used in this integration is set by the value of max_intervals. The full value of the action to be returned is then obtained from the value of \(S_T\) using the relation

\[ S_E = \frac{2}{d} S_T . \]

See also
integrate_gsl_qag

Definition at line 63 of file euclidean_action.cpp.

References area_n_sphere(), BubbleProfiler::Field_profiles::derivative_at(), BubbleProfiler::Field_profiles::get_domain_end(), BubbleProfiler::Field_profiles::get_domain_start(), BubbleProfiler::Field_profiles::get_number_of_dimensions(), BubbleProfiler::Field_profiles::get_number_of_fields(), integrate_gsl_qag(), and gaussian_alpha_tests::n_fields.

Referenced by calculate_action().

double BubbleProfiler::calculate_potential_action ( const Potential potential,
const Field_profiles profiles,
std::size_t  max_intervals = 1000,
double  rel_tol = 1.e-4,
double  abs_tol = 1.e-4,
Integration_rule  rule = Integration_rule::GK31 
)

Calculates the action using only the potential term contributions.

Parameters
potentialthe potential associated with the bounce solution
profilesthe field profiles representing the bounce solution
max_intervalsthe maximum allowed number of subintervals for adaptive quadrature routine
rel_tolthe relative error goal
abs_tolthe absolute error goal
rulethe quadrature rule to use
Returns
numerical estimate for the Euclidean action

Given the bounce solution for \(n\) scalar fields in \(d\) dimensions defined on a finite domain, the contribution to the action due to the field potential terms,

\[ S_V \equiv S_{d-1} \int_0^\infty d\rho \rho^{d-1} V(\phi_i) , \]

where \(S_{d-1}\) is the surface area of the \(d\)-sphere and the potential is assumed to vanish at the false vacuum, is first calculated using numerical quadrature, with the quadrature rule chosen according to the value of rule. The adaptive quadrature routine attempts to compute an estimate \(\tilde{I}\) for the integral \(I\) appearing above that satisfies

\[ | \tilde{I} - I | \leq \max ( \epsilon_{\text{abs}}, \epsilon_{\text{rel}} |I|) , \]

with \(\epsilon_{\text{abs}}\) and \(\epsilon_{\text{rel}}\) given by abs_tol and rel_tol, respectively. The maximum allowed number of subintervals used in this integration is set by the value of max_intervals. The full value of the action to be returned is then obtained from the value of \(S_V\) using the relation

\[ S_E = \frac{2}{2 - d} S_V . \]

See also
integrate_gsl_qag

Definition at line 138 of file euclidean_action.cpp.

References area_n_sphere(), BubbleProfiler::Field_profiles::evaluate_at(), BubbleProfiler::Field_profiles::get_domain_end(), BubbleProfiler::Field_profiles::get_domain_start(), BubbleProfiler::Field_profiles::get_number_of_dimensions(), and integrate_gsl_qag().

Referenced by calculate_action().

Eigen::MatrixXd BubbleProfiler::calculate_rotation_to_target ( const Eigen::VectorXd &  target)

Calculate the rotation matrix for a rotation to the target vector.

Given a target vector, this method computes the rotation matrix to transform the current coordinate system to one in which the first basis vector is aligned with the target vector.

Notes:

  1. this change of coordinates is not uniquely specified, and may involve a reflection.
  2. Right multiplication will take coordinates in the rotated system to coordinates in the unrotated system. Take the transpose for the opposite effect.
Parameters
targetthe target vector to align with
Returns
the rotation matrix

Definition at line 36 of file rotation.cpp.

Referenced by BubbleProfiler::Shooting_profile_guesser::calculate_potential_parameters(), and BubbleProfiler::Kink_profile_guesser::calculate_potential_parameters().

void BubbleProfiler::configure_logging ( const Bubble_profiler_inputs input)
template<typename T , typename = typename std::enable_if< std::is_floating_point<T>::value,T>::type>
T BubbleProfiler::Cos ( x)
noexcept
template<typename T , typename = typename std::enable_if< std::is_integral<T>::value,T>::type>
double BubbleProfiler::Cos ( n)
noexcept

Definition at line 134 of file math_wrappers.hpp.

template<typename T >
std::complex<T> BubbleProfiler::Cos ( const std::complex< T > &  z)
noexcept

Definition at line 140 of file math_wrappers.hpp.

References sm_plus_singlet::type.

template<typename T , typename = typename std::enable_if< std::is_floating_point<T>::value,T>::type>
T BubbleProfiler::Cosh ( x)
noexcept
template<typename T , typename = typename std::enable_if< std::is_integral<T>::value,T>::type>
double BubbleProfiler::Cosh ( n)
noexcept

Definition at line 284 of file math_wrappers.hpp.

template<typename T >
std::complex<T> BubbleProfiler::Cosh ( const std::complex< T > &  z)
noexcept

Definition at line 290 of file math_wrappers.hpp.

References sm_plus_singlet::type.

template<typename T >
auto BubbleProfiler::Cot ( x) -> decltype(Tan(x))
noexcept

Definition at line 190 of file math_wrappers.hpp.

References Tan().

template<typename T >
auto BubbleProfiler::Coth ( x) -> decltype(Tanh(x))
noexcept

Definition at line 340 of file math_wrappers.hpp.

References Tanh().

template<typename T >
auto BubbleProfiler::Csc ( x) -> decltype(Sin(x))
noexcept

Definition at line 196 of file math_wrappers.hpp.

References Sin().

template<typename T >
auto BubbleProfiler::Csch ( x) -> decltype(Sinh(x))
noexcept

Definition at line 346 of file math_wrappers.hpp.

References Sinh().

template<typename T , typename = typename std::enable_if< std::is_floating_point<T>::value,T>::type>
T BubbleProfiler::Exp ( x)
noexcept
template<typename T , typename = typename std::enable_if< std::is_integral<T>::value,T>::type>
double BubbleProfiler::Exp ( n)
noexcept

Definition at line 84 of file math_wrappers.hpp.

template<typename T >
std::complex<T> BubbleProfiler::Exp ( const std::complex< T > &  z)
noexcept

Definition at line 90 of file math_wrappers.hpp.

References sm_plus_singlet::type.

Eigen::VectorXd BubbleProfiler::find_global_min ( const Potential potential,
double  opt_timeout 
)
Eigen::VectorXd BubbleProfiler::find_one_dimensional_barrier ( const Potential potential,
const Eigen::VectorXd &  true_vacuum_loc,
const Eigen::VectorXd &  false_vacuum_loc,
double  opt_timeout 
)

Find barrier directly between two minima of a 1D potential.

In the one-dimensional case, this function attempts to find the location of the barrier between two given minima, taken to be the location of the maximum value of the potential between the two minima.

Parameters
potentialthe potential to locate the barrier for
true_vacuum_locthe location of the global minimum
false_vacuum_locthe location of the local minimum
opt_timeoutthe allowed maximum amount of time to run for
Returns
the obtained location of the maximum

Definition at line 377 of file run_cmd_line_potential.cpp.

References BubbleProfiler::NLopt_optimizer::get_extremum_location(), BubbleProfiler::Potential::get_number_of_fields(), BubbleProfiler::NLopt_optimizer::MAX, gaussian_alpha_tests::n_fields, optimization_succeeded(), BubbleProfiler::NLopt_optimizer::optimize(), BubbleProfiler::NLopt_optimizer::set_extremum_type(), BubbleProfiler::NLopt_optimizer::set_ftol_rel(), BubbleProfiler::NLopt_optimizer::set_lower_bounds(), BubbleProfiler::NLopt_optimizer::set_max_time(), BubbleProfiler::NLopt_optimizer::set_upper_bounds(), and BubbleProfiler::NLopt_optimizer::set_xtol_rel().

Referenced by initialize_extrema().

int BubbleProfiler::get_integration_rule_key ( Integration_rule  rule)

Definition at line 22 of file integration_utils.cpp.

References GK15, GK21, GK31, GK41, GK51, and GK61.

Referenced by integrate_gsl_qag().

std::string BubbleProfiler::get_option_value ( const std::string &  option,
const std::string &  sep = "=" 
)
Examples:
general_fubini.cpp, logarithmic.cpp, and thin_wall.cpp.

Definition at line 88 of file general_fubini.cpp.

Referenced by parse_global_min(), and parse_positional_args().

template<typename T , typename = typename std::enable_if< std::is_arithmetic<T>::value,T>::type>
T BubbleProfiler::Im ( )
noexcept

Definition at line 426 of file math_wrappers.hpp.

template<typename T >
T BubbleProfiler::Im ( const std::complex< T > &  z)
noexcept

Definition at line 432 of file math_wrappers.hpp.

References sm_plus_singlet::type.

void BubbleProfiler::initialize_extrema ( const Potential potential,
const Bubble_profiler_inputs input,
Eigen::VectorXd &  true_vacuum_loc,
Eigen::VectorXd &  false_vacuum_loc 
)

Locate the requested critical points if not already provided.

If not empty, the values of input.global_min and input.local_min are used to set the values of the global minimum and local minimum, respectively.

If input.global_min is empty, the function attempts to locate the true vacuum by numerically optimizing the potential. By default, the function will not attempt to find the local minimum if input.local_min is empty, and will exit with an error. However, if input.local_min is empty and input.assume_origin_false_vacuum is set to true, then the location of the false vacuum will be set to the origin in field space.

Parameters
potentialthe potential to find the critical points of
inputsettings to be used to configure the potential and solvers
true_vacuum_locthe vector in which to store the obtained global minimum
false_vacuum_locthe vector in which to store the obtained local minimum

Definition at line 440 of file run_cmd_line_potential.cpp.

References BubbleProfiler::Bubble_profiler_inputs::assume_origin_false_vacuum, find_global_min(), BubbleProfiler::Potential::get_number_of_fields(), BubbleProfiler::Bubble_profiler_inputs::global_min, BubbleProfiler::Bubble_profiler_inputs::local_min, gaussian_alpha_tests::n_fields, and BubbleProfiler::Bubble_profiler_inputs::opt_timeout.

Referenced by initialize_extrema(), run_perturbative_profiler(), and run_shooting_profiler().

void BubbleProfiler::initialize_extrema ( const Potential potential,
const Bubble_profiler_inputs input,
Eigen::VectorXd &  true_vacuum_loc,
Eigen::VectorXd &  false_vacuum_loc,
Eigen::VectorXd &  barrier_loc 
)

Locate the necessary critical points if not already provided.

If not empty, the values of input.global_min , input.local_min , and input.barrier are used to set the values of the global minimum, local minimum, and barrier, respectively.

If input.global_min is empty, the function attempts to locate the true vacuum by numerically optimizing the potential. By default, the function will not attempt to find the local minimum if input.local_min is empty, and will exit with an error. However, if input.local_min is empty and input.assume_origin_false_vacuum is set to true, then the location of the false vacuum will be set to the origin in field space.

Locating the position of the barrier between two minima is only supported for one-field potentials. In these cases, if input.barrier is empty, the local maximum of the potential is obtained numerically; otherwise, the given value is used.

Parameters
potentialthe potential to find the critical points of
inputsettings to be used to configure the potential and solvers
true_vacuum_locthe vector in which to store the obtained global minimum
false_vacuum_locthe vector in which to store the obtained local minimum
barrier_locthe vector in which to store the obtained barrier location

Definition at line 505 of file run_cmd_line_potential.cpp.

References BubbleProfiler::Bubble_profiler_inputs::barrier, find_one_dimensional_barrier(), BubbleProfiler::Potential::get_number_of_fields(), initialize_extrema(), gaussian_alpha_tests::n_fields, and BubbleProfiler::Bubble_profiler_inputs::opt_timeout.

template<typename F >
std::tuple<int,double,double,int> BubbleProfiler::integrate_gsl_cquad ( f,
double  a,
double  b,
double  abs_tol,
double  rel_tol,
std::size_t  n_intervals 
)

Definition at line 110 of file integration_utils.hpp.

References make_raii_guard().

template<typename F >
std::tuple<int,double,double> BubbleProfiler::integrate_gsl_qag ( f,
double  a,
double  b,
double  abs_tol,
double  rel_tol,
std::size_t  max_intervals,
Integration_rule  rule 
)
template<typename F >
std::tuple<int,double,double> BubbleProfiler::integrate_gsl_qags ( f,
double  a,
double  b,
double  abs_tol,
double  rel_tol,
std::size_t  max_intervals 
)

Definition at line 75 of file integration_utils.hpp.

References make_raii_guard().

Eigen::VectorXd BubbleProfiler::interpolate_f_at ( const Eigen::VectorXd &  x_old,
const Eigen::VectorXd &  y_old,
const Eigen::VectorXd &  x_new 
)

Computes interpolated function values at a new set of grid points.

Parameters
[in]x_oldthe initial set of x values at which the function is known
[in]y_oldthe initial set of y values
[in]x_newthe x values at which to interpolate the function
Returns
a vector of the interpolated values at each point in x_new

The cubic spline used to obtain the interpolated function values is computed using all of the given (x_old, y_old) values, of which there must be at least two. The vectors x_old and y_old must have the same lengths. Default parabolic boundary conditions, with unknown derivatives, are used at the left and right ends of the domain. The obtained spline interpolant is evaluated at all points in x_new .

Definition at line 93 of file gsl_interpolator.cpp.

Referenced by BubbleProfiler::Perturbative_profiler< Integration_policy >::update_field_profiles().

Eigen::VectorXd BubbleProfiler::interpolate_f_at ( const Eigen::VectorXd &  x_old,
const Eigen::VectorXd &  y_old,
int  n_old,
const Eigen::VectorXd &  x_new,
int  n_new 
)

Computes interpolated function values at a new set of grid points.

Parameters
[in]x_oldthe initial set of x values at which the function is known
[in]y_oldthe initial set of y values
[in]n_oldthe number of values to use from x_old and y_old
[in]x_newthe x values at which to interpolate the function
[in]n_newthe number of values in x_new to evaluate at
Returns
a vector of the interpolated values at the values in x_new

The cubic spline used to obtain the interpolated function values is computed using the first n_old values from x_old and y_old, with n_old between 2 and the minimum of the number of values in x_old and y_old . The obtained spline interpolant is evaluated at the first n_new points in x_new.

Definition at line 114 of file gsl_interpolator.cpp.

References make_raii_guard().

template<typename Scalar , typename Data >
Scalar BubbleProfiler::linear_lagrange_interpolation_at_point ( Scalar  x,
const Data &  xdata,
const Data &  ydata 
)

Assumes data sorted in order of ascending x-values

Definition at line 29 of file univariate_interpolation.hpp.

template<typename T , typename = typename std::enable_if< std::is_floating_point<T>::value,T>::type>
T BubbleProfiler::Log ( x)
noexcept
template<typename T , typename = typename std::enable_if< std::is_integral<T>::value,T>::type>
double BubbleProfiler::Log ( n)
noexcept

Definition at line 106 of file math_wrappers.hpp.

template<typename T >
std::complex<T> BubbleProfiler::Log ( const std::complex< T > &  z)
noexcept

Definition at line 112 of file math_wrappers.hpp.

template<typename F >
constexpr RAII_guard<F> BubbleProfiler::make_raii_guard ( f)
bool BubbleProfiler::optimization_succeeded ( nlopt::result  r)

Helper function to check NLopt status code

Definition at line 30 of file nlopt_optimizer.cpp.

References SUCCESS.

Referenced by find_global_min(), and find_one_dimensional_barrier().

Thin_wall_options BubbleProfiler::parse_cmd_line_args ( int  argc,
const char *  argv[] 
)
Examples:
action.cpp, general_fubini.cpp, logarithmic.cpp, scale.cpp, tabulate.cpp, and thin_wall.cpp.

Definition at line 160 of file run_cmd_line_potential.cpp.

References BubbleProfiler::Bubble_profiler_inputs::action_arrived_rel, BubbleProfiler::Bubble_profiler_inputs::action_ode_abs, BubbleProfiler::Bubble_profiler_inputs::action_ode_rel, BubbleProfiler::Bubble_profiler_inputs::algorithm, BubbleProfiler::Bubble_profiler_inputs::ansatz_file, BubbleProfiler::Bubble_profiler_inputs::assume_origin_false_vacuum, BubbleProfiler::Bubble_profiler_inputs::barrier, BubbleProfiler::Bubble_profiler_inputs::bisect_lambda_max, default_max_iterations, BubbleProfiler::Bubble_profiler_inputs::domain_end, BubbleProfiler::Bubble_profiler_inputs::domain_start, BubbleProfiler::Bubble_profiler_inputs::drho_frac, BubbleProfiler::Bubble_profiler_inputs::evolve_change_rel, BubbleProfiler::Bubble_profiler_inputs::f_y_max, BubbleProfiler::Bubble_profiler_inputs::f_y_min, BubbleProfiler::Bubble_profiler_inputs::fields, BubbleProfiler::Bubble_profiler_inputs::force_output, BubbleProfiler::Bubble_profiler_inputs::global_min, BubbleProfiler::Bubble_profiler_inputs::initial_step_size, BubbleProfiler::Bubble_profiler_inputs::interpolation_fraction, BubbleProfiler::Bubble_profiler_inputs::iter_max, BubbleProfiler::Bubble_profiler_inputs::local_min, BubbleProfiler::Bubble_profiler_inputs::max_iterations, BubbleProfiler::Bubble_profiler_inputs::n_dims, BubbleProfiler::Bubble_profiler_inputs::opt_timeout, BubbleProfiler::Bubble_profiler_inputs::output_file, BubbleProfiler::Bubble_profiler_inputs::output_path, BubbleProfiler::Bubble_profiler_inputs::periods_max, BubbleProfiler::Bubble_profiler_inputs::potential, BubbleProfiler::Bubble_profiler_inputs::rtol_action, BubbleProfiler::Bubble_profiler_inputs::rtol_fields, BubbleProfiler::Bubble_profiler_inputs::shoot_bisect_bits, BubbleProfiler::Bubble_profiler_inputs::shoot_ode_abs, BubbleProfiler::Bubble_profiler_inputs::shoot_ode_rel, BubbleProfiler::Bubble_profiler_inputs::shooting_ansatz, BubbleProfiler::Bubble_profiler_inputs::use_perturbative, BubbleProfiler::Bubble_profiler_inputs::verbose, BubbleProfiler::Bubble_profiler_inputs::write_profiles, and BubbleProfiler::Bubble_profiler_inputs::y_max.

Referenced by main(), parse_global_min(), and parse_positional_args().

double BubbleProfiler::parse_global_min ( const std::string &  value)
void BubbleProfiler::parse_positional_args ( const std::vector< std::string > &  args,
Scale_options options 
)
void BubbleProfiler::parse_positional_args ( const std::vector< std::string > &  args,
Action_options options 
)
void BubbleProfiler::parse_positional_args ( const std::vector< std::string > &  args,
Tabulate_options options 
)
void BubbleProfiler::parse_positional_args ( const std::vector< std::string > &  args,
Fubini_options options 
)
void BubbleProfiler::parse_positional_args ( const std::vector< std::string > &  args,
Logarithmic_options options 
)
void BubbleProfiler::parse_positional_args ( const std::vector< std::string > &  args,
Thin_wall_options options 
)
void BubbleProfiler::print_usage ( )
template<typename Scalar , typename Data >
Scalar BubbleProfiler::quadratic_lagrange_interpolation_at_point ( Scalar  x,
const Data &  xdata,
const Data &  ydata 
)
template<typename T , typename = typename std::enable_if< std::is_arithmetic<T>::value,T>::type>
T BubbleProfiler::Re ( x)
noexcept

Definition at line 440 of file math_wrappers.hpp.

Referenced by Wm1().

template<typename T >
T BubbleProfiler::Re ( const std::complex< T > &  z)
noexcept

Definition at line 446 of file math_wrappers.hpp.

References sm_plus_singlet::type.

Profiler_result BubbleProfiler::run_perturbative_profiler ( Restricted_quartic_potential  potential,
int  dim 
)
template<class Profiler >
std::tuple<double, Field_profiles> BubbleProfiler::run_perturbative_profiler ( const Bubble_profiler_inputs input,
Profiler &  profiler 
)
std::tuple<double, Field_profiles> BubbleProfiler::run_perturbative_profiler ( const Bubble_profiler_inputs input)

Calculates the action and bubble profile using the perturbative method.

The BVP solver used in the calculation is determined from the provided setting for the integration algorithm.

Parameters
inputsettings to be used to configure the potential and solver
bvp_solverthe BVP solver to be used
Returns
a tuple containing the calculation action and, optionally, the field profiles

Definition at line 670 of file run_cmd_line_potential.cpp.

References BubbleProfiler::Bubble_profiler_inputs::algorithm, Euler, run_perturbative_profiler(), and Runge_kutta_4.

Profiler_result BubbleProfiler::run_profiler ( const std::string &  profiler,
const Restricted_quartic_potential potential,
int  dim 
)

Definition at line 284 of file tabulate.cpp.

References run_perturbative_profiler(), and run_shooting_profiler().

void BubbleProfiler::run_profiler ( const Bubble_profiler_inputs input)
Profiler_result BubbleProfiler::run_shooting_profiler ( const Restricted_quartic_potential potential,
int  dim 
)
std::tuple<double, Field_profiles> BubbleProfiler::run_shooting_profiler ( const Bubble_profiler_inputs input)

Calculates the action and bubble profile using the shooting method.

Parameters
inputsettings to be used to configure the potential and solver
Returns
a tuple containing the calculation action and, optionally, the field profiles
Examples:
tabulate.cpp.

Definition at line 542 of file run_cmd_line_potential.cpp.

References BubbleProfiler::Bubble_profiler_inputs::action_arrived_rel, BubbleProfiler::Bubble_profiler_inputs::action_ode_abs, BubbleProfiler::Bubble_profiler_inputs::action_ode_rel, BubbleProfiler::Bubble_profiler_inputs::bisect_lambda_max, BubbleProfiler::Bubble_profiler_inputs::drho_frac, BubbleProfiler::Bubble_profiler_inputs::evolve_change_rel, BubbleProfiler::Bubble_profiler_inputs::f_y_max, BubbleProfiler::Bubble_profiler_inputs::f_y_min, BubbleProfiler::Bubble_profiler_inputs::fields, BubbleProfiler::Shooting::get_bubble_profile(), BubbleProfiler::Shooting::get_euclidean_action(), BubbleProfiler::Algebraic_potential::get_number_of_fields(), initialize_extrema(), BubbleProfiler::Bubble_profiler_inputs::iter_max, BubbleProfiler::Bubble_profiler_inputs::n_dims, gaussian_alpha_tests::n_fields, BubbleProfiler::Bubble_profiler_inputs::periods_max, BubbleProfiler::Bubble_profiler_inputs::potential, BubbleProfiler::Shooting::set_action_abs_tol(), BubbleProfiler::Shooting::set_action_arrived_rel(), BubbleProfiler::Shooting::set_action_rel_tol(), BubbleProfiler::Shooting::set_bisect_lambda_max(), BubbleProfiler::Shooting::set_bisection_precision_bits(), BubbleProfiler::Shooting::set_drho_frac(), BubbleProfiler::Shooting::set_evolve_change_rel(), BubbleProfiler::Shooting::set_f_y_max(), BubbleProfiler::Shooting::set_f_y_min(), BubbleProfiler::Shooting::set_max_iterations(), BubbleProfiler::Shooting::set_max_periods(), BubbleProfiler::Shooting::set_shooting_abs_tol(), BubbleProfiler::Shooting::set_shooting_rel_tol(), BubbleProfiler::Shooting::set_y_max(), BubbleProfiler::Bubble_profiler_inputs::shoot_bisect_bits, BubbleProfiler::Bubble_profiler_inputs::shoot_ode_abs, BubbleProfiler::Bubble_profiler_inputs::shoot_ode_rel, BubbleProfiler::Shooting::solve(), BubbleProfiler::Bubble_profiler_inputs::write_profiles, and BubbleProfiler::Bubble_profiler_inputs::y_max.

Referenced by run_profiler().

template<typename T >
auto BubbleProfiler::Sec ( x) -> decltype(Cos(x))
noexcept

Definition at line 202 of file math_wrappers.hpp.

References Cos(), and sm_plus_singlet::type.

template<typename T >
auto BubbleProfiler::Sech ( x) -> decltype(Cosh(x))
noexcept
double BubbleProfiler::series ( double  a)

Helper function for approximate inverse of \(\sinc x\)

See in Mathematica,

Simplify[Normal[InverseSeries[1 - Series[Sin[x]/x, {x, 0, 20}]]]/Sqrt[3/2] /. x -> x^2, x > 0]

and this link for further information.

Definition at line 73 of file numeric.cpp.

Referenced by asinc().

template<typename T >
int BubbleProfiler::Sign ( const T &  x)
template<typename T , typename = typename std::enable_if< std::is_floating_point<T>::value,T>::type>
T BubbleProfiler::Sin ( x)
noexcept
template<typename T , typename = typename std::enable_if< std::is_integral<T>::value,T>::type>
double BubbleProfiler::Sin ( n)
noexcept

Definition at line 156 of file math_wrappers.hpp.

template<typename T >
std::complex<T> BubbleProfiler::Sin ( const std::complex< T > &  z)
noexcept

Definition at line 162 of file math_wrappers.hpp.

References sm_plus_singlet::type.

template<typename T , typename = typename std::enable_if< std::is_floating_point<T>::value,T>::type>
T BubbleProfiler::Sinh ( x)
noexcept
template<typename T , typename = typename std::enable_if< std::is_integral<T>::value,T>::type>
double BubbleProfiler::Sinh ( n)
noexcept

Definition at line 306 of file math_wrappers.hpp.

template<typename T >
std::complex<T> BubbleProfiler::Sinh ( const std::complex< T > &  z)
noexcept

Definition at line 312 of file math_wrappers.hpp.

References sm_plus_singlet::type.

std::tuple<double, Field_profiles> BubbleProfiler::solve_perturbative ( Thin_wall_potential potential,
int  dim,
const std::string &  observer_file 
)
std::tuple<double, Field_profiles> BubbleProfiler::solve_perturbative ( Logarithmic_potential potential,
double  true_min,
const std::string &  observer_file 
)
std::tuple<double, Field_profiles> BubbleProfiler::solve_perturbative ( Generalized_fubini_potential potential,
double  true_min,
const std::string &  observer_file 
)
std::tuple<double, Field_profiles> BubbleProfiler::solve_shooting ( const Thin_wall_potential potential,
int  dim,
bool  calculate_profile 
)
std::tuple<double, Field_profiles> BubbleProfiler::solve_shooting ( const Logarithmic_potential potential,
double  true_min,
bool  calculate_profile 
)
std::tuple<double, Field_profiles> BubbleProfiler::solve_shooting ( const Generalized_fubini_potential potential,
double  true_min,
bool  calculate_profile 
)
template<typename T , typename = typename std::enable_if< std::is_floating_point<T>::value,T>::type>
T BubbleProfiler::Sqrt ( x)
noexcept
template<typename T , typename = typename std::enable_if< std::is_integral<T>::value,T>::type>
double BubbleProfiler::Sqrt ( n)
noexcept

Definition at line 462 of file math_wrappers.hpp.

bool BubbleProfiler::starts_with ( const std::string &  option,
const std::string &  prefix 
)
template<typename T , typename = typename std::enable_if< std::is_floating_point<T>::value,T>::type>
T BubbleProfiler::Tan ( x)
noexcept

Definition at line 170 of file math_wrappers.hpp.

References sm_plus_singlet::type.

Referenced by Cot().

template<typename T , typename = typename std::enable_if< std::is_integral<T>::value,T>::type>
double BubbleProfiler::Tan ( n)
noexcept

Definition at line 178 of file math_wrappers.hpp.

template<typename T >
std::complex<T> BubbleProfiler::Tan ( const std::complex< T > &  z)
noexcept

Definition at line 184 of file math_wrappers.hpp.

template<typename T , typename = typename std::enable_if< std::is_floating_point<T>::value,T>::type>
T BubbleProfiler::Tanh ( x)
noexcept
template<typename T , typename = typename std::enable_if< std::is_integral<T>::value,T>::type>
double BubbleProfiler::Tanh ( n)
noexcept

Definition at line 328 of file math_wrappers.hpp.

template<typename T >
std::complex<T> BubbleProfiler::Tanh ( const std::complex< T > &  z)
noexcept

Definition at line 334 of file math_wrappers.hpp.

void BubbleProfiler::validate ( boost::any &  v,
const std::vector< std::string > &  values,
Integration_algorithm ,
int   
)

Definition at line 141 of file run_cmd_line_potential.cpp.

References Euler, and Runge_kutta_4.

double BubbleProfiler::Wm1 ( double  a)

Approximation to negative branch of the Lambert- \(W\) function.

See Eq. 4.19 of this paper

This is correct up to order \((\log \log a / \log a)^5\).

Returns
Negative branch of Lambert- \(W\) function
Parameters
a

Definition at line 41 of file numeric.cpp.

References Log(), and Re().

Referenced by approx_root_pos_4(), and asinch().

void BubbleProfiler::write_action ( double  action,
std::ostream &  ostr 
)

Write the value of the action to requested output stream.

Parameters
actionthe value of the action
ostrthe output stream to write to
Examples:
general_fubini.cpp, logarithmic.cpp, and thin_wall.cpp.

Definition at line 292 of file run_cmd_line_potential.cpp.

Referenced by main(), run_profiler(), solve_perturbative(), and write_profile().

void BubbleProfiler::write_action ( std::ostream &  ostr,
double  numerical_action,
double  exact_action 
)

Definition at line 383 of file general_fubini.cpp.

void BubbleProfiler::write_profile ( std::ostream &  ostr,
const Field_profiles numerical_profile 
)
void BubbleProfiler::write_profile ( std::ostream &  ostr,
const Field_profiles numerical_profile,
const Field_profiles exact_profile 
)
void BubbleProfiler::write_profiles ( const std::vector< std::string > &  fields,
const Field_profiles profiles,
std::ostream &  ostr 
)

Write the given field profiles to requested output stream.

Parameters
fieldsthe names of the fields
profilesthe values of the field profiles
ostrthe output stream to write to

Definition at line 304 of file run_cmd_line_potential.cpp.

References BubbleProfiler::Field_profiles::get_field_profiles(), BubbleProfiler::Field_profiles::get_spatial_grid(), and gaussian_alpha_tests::n_fields.

Referenced by run_profiler().

void BubbleProfiler::write_table ( int  dim,
double  E,
const std::vector< double > &  alpha_data,
const Profiler_data profiler_data 
)
Examples:
tabulate.cpp.

Definition at line 313 of file tabulate.cpp.

References write_table_header().

Referenced by main().

void BubbleProfiler::write_table_header ( const std::vector< std::string > &  profiler_labels)
Examples:
tabulate.cpp.

Definition at line 298 of file tabulate.cpp.

Referenced by write_table().

Variable Documentation

const double BubbleProfiler::default_action_arrived_rel = 1.e-3

Definition at line 81 of file run_cmd_line_potential.cpp.

const double BubbleProfiler::default_action_ode_abs = 1.e-6

Definition at line 84 of file run_cmd_line_potential.cpp.

const double BubbleProfiler::default_action_ode_rel = 1.e-6

Definition at line 85 of file run_cmd_line_potential.cpp.

const Integration_algorithm BubbleProfiler::default_algorithm
Initial value:
=
Integration_algorithm::Runge_kutta_4

Definition at line 64 of file run_cmd_line_potential.cpp.

const std::string BubbleProfiler::default_ansatz_file = ""

Definition at line 59 of file run_cmd_line_potential.cpp.

const bool BubbleProfiler::default_assume_origin_false_vacuum = false

Definition at line 73 of file run_cmd_line_potential.cpp.

const double BubbleProfiler::default_bisect_lambda_max = 5.0

Definition at line 88 of file run_cmd_line_potential.cpp.

const double BubbleProfiler::default_domain_end = -1.

Definition at line 61 of file run_cmd_line_potential.cpp.

const double BubbleProfiler::default_domain_start = -1.

Definition at line 60 of file run_cmd_line_potential.cpp.

const double BubbleProfiler::default_drho_frac = 1.e-3

Definition at line 86 of file run_cmd_line_potential.cpp.

const double BubbleProfiler::default_evolve_change_rel = 1.e-2

Definition at line 87 of file run_cmd_line_potential.cpp.

const double BubbleProfiler::default_f_y_max = 1.e6

Definition at line 91 of file run_cmd_line_potential.cpp.

const double BubbleProfiler::default_f_y_min = 1.e-3

Definition at line 92 of file run_cmd_line_potential.cpp.

const bool BubbleProfiler::default_force_output = false

Definition at line 70 of file run_cmd_line_potential.cpp.

const double BubbleProfiler::default_initial_step_size = 1.e-2

Definition at line 62 of file run_cmd_line_potential.cpp.

const int BubbleProfiler::default_interp_fraction = 1.0

Definition at line 63 of file run_cmd_line_potential.cpp.

const int BubbleProfiler::default_iter_max = 100000

Definition at line 89 of file run_cmd_line_potential.cpp.

const int BubbleProfiler::default_max_iterations = -1

Definition at line 72 of file run_cmd_line_potential.cpp.

Referenced by parse_cmd_line_args().

const double BubbleProfiler::default_opt_timeout = 10.

Definition at line 71 of file run_cmd_line_potential.cpp.

const std::string BubbleProfiler::default_output_file = "-"

Definition at line 68 of file run_cmd_line_potential.cpp.

const double BubbleProfiler::default_periods_max = 1.e2

Definition at line 90 of file run_cmd_line_potential.cpp.

const double BubbleProfiler::default_rtol_action = 1.e-3

Definition at line 66 of file run_cmd_line_potential.cpp.

const double BubbleProfiler::default_rtol_fields = 1.e-3

Definition at line 67 of file run_cmd_line_potential.cpp.

const int BubbleProfiler::default_shoot_bisect_bits = 5

Definition at line 80 of file run_cmd_line_potential.cpp.

const double BubbleProfiler::default_shoot_ode_abs = 1.e-4

Definition at line 82 of file run_cmd_line_potential.cpp.

const double BubbleProfiler::default_shoot_ode_rel = 1.e-4

Definition at line 83 of file run_cmd_line_potential.cpp.

const bool BubbleProfiler::default_shooting_ansatz = false

Definition at line 76 of file run_cmd_line_potential.cpp.

const bool BubbleProfiler::default_use_perturbative = false

Definition at line 74 of file run_cmd_line_potential.cpp.

const bool BubbleProfiler::default_verbose = false

Definition at line 77 of file run_cmd_line_potential.cpp.

const bool BubbleProfiler::default_write_perturbations = true

Definition at line 69 of file run_cmd_line_potential.cpp.

const bool BubbleProfiler::default_write_profiles = false

Definition at line 75 of file run_cmd_line_potential.cpp.

const double BubbleProfiler::default_y_max = 1.e1

Definition at line 93 of file run_cmd_line_potential.cpp.