BubbleProfiler
0.3.0
by Peter Athron, Csaba Balazs, Michael Bardsley, Andrew Fowlie, Dylan Harries & Graham White
|
Example thin wall potential as given by Coleman. More...
#include <thin_wall_potential.hpp>
Public Member Functions | |
Thin_wall_potential ()=default | |
Thin_wall_potential (double, double, double) | |
virtual | ~Thin_wall_potential ()=default |
Thin_wall_potential (const Thin_wall_potential &)=default | |
Thin_wall_potential (Thin_wall_potential &&)=default | |
Thin_wall_potential & | operator= (const Thin_wall_potential &)=default |
Thin_wall_potential & | operator= (Thin_wall_potential &&)=default |
virtual Thin_wall_potential * | clone () const override |
Subclasses must implement a clone method. More... | |
virtual double | operator() (const Eigen::VectorXd &) const override |
Evaluate potential at point. More... | |
virtual double | partial (const Eigen::VectorXd &, int) const override |
Partial derivative WRT coordinate i at a point. More... | |
virtual double | partial (const Eigen::VectorXd &, int, int) const override |
Partial derivative WRT coordinates i, j at a a point. More... | |
virtual std::size_t | get_number_of_fields () const override |
virtual void | translate_origin (const Eigen::VectorXd &translation) override |
Shift the location of the origin by a specified vector. More... | |
virtual void | apply_basis_change (const Eigen::MatrixXd &cob_matrix) override |
Apply a change of basis matrix. More... | |
virtual void | add_constant_term (double _offset) override |
Add a constant offset to the potential. More... | |
double | operator() (double) const |
evaluates the potential for the given field value More... | |
double | first_deriv (double) const |
evaluates the first derivative of the potential More... | |
double | second_deriv (double) const |
evaluates the second derivative of the potential More... | |
double | get_local_minimum_location () const |
returns the location of the local minimum More... | |
double | get_local_maximum_location () const |
returns the location of the local maximum More... | |
double | get_global_minimum_location () const |
returns the location of the global minimum More... | |
double | get_thin_wall_action () const |
returns the approximate action in the thin-wall approximation More... | |
![]() | |
virtual | ~Potential ()=default |
Private Attributes | |
double | lambda {1.} |
double | a {1.} |
double | epsilon {0.01} |
double | origin {0.} |
double | scale {1.} |
double | offset {0.} |
Example thin wall potential as given by Coleman.
This class implements the example potential
\[ V(\phi) = \frac{\lambda}{8} (\phi^2 - a^2)^2 + \frac{\epsilon}{2 a} (\phi - a) . \]
and provides analytic expressions for the tunnelling action in the thin-wall approximation.
Definition at line 40 of file thin_wall_potential.hpp.
|
default |
Referenced by clone().
BubbleProfiler::Thin_wall_potential::Thin_wall_potential | ( | double | lambda_, |
double | a_, | ||
double | epsilon_ | ||
) |
Definition at line 24 of file thin_wall_potential.cpp.
References a, BubbleProfiler::Abs(), epsilon, lambda, and BubbleProfiler::Sqrt().
|
virtualdefault |
|
default |
|
default |
|
inlineoverridevirtual |
Add a constant offset to the potential.
offset | constant term to add to potential |
Implements BubbleProfiler::Potential.
Definition at line 70 of file thin_wall_potential.hpp.
References first_deriv(), get_global_minimum_location(), get_local_maximum_location(), get_local_minimum_location(), get_thin_wall_action(), offset, operator()(), and second_deriv().
|
inlineoverridevirtual |
Apply a change of basis matrix.
cob_matrix |
Implements BubbleProfiler::Potential.
Definition at line 66 of file thin_wall_potential.hpp.
References scale.
|
inlineoverridevirtual |
Subclasses must implement a clone method.
Implements BubbleProfiler::Potential.
Definition at line 52 of file thin_wall_potential.hpp.
References operator()(), partial(), and Thin_wall_potential().
double BubbleProfiler::Thin_wall_potential::first_deriv | ( | double | coords | ) | const |
evaluates the first derivative of the potential
coords | the value of the field to evaluate the derivative at |
Definition at line 101 of file thin_wall_potential.cpp.
References a, epsilon, lambda, origin, and scale.
Referenced by add_constant_term(), and partial().
double BubbleProfiler::Thin_wall_potential::get_global_minimum_location | ( | ) | const |
returns the location of the global minimum
The location of the global minimum is computed exactly using the solution to the cubic equation that results from imposing \( \partial V / \partial \phi = 0\). For \(\epsilon / (\lambda a^2) \ll 1\), one finds that the global minimum is located approximately at \(\phi_t \approx -a\).
Definition at line 166 of file thin_wall_potential.cpp.
References a, BubbleProfiler::Abs(), BubbleProfiler::ArcSin(), epsilon, lambda, origin, scale, BubbleProfiler::Sin(), and BubbleProfiler::Sqrt().
Referenced by add_constant_term(), BubbleProfiler::solve_perturbative(), and BubbleProfiler::solve_shooting().
double BubbleProfiler::Thin_wall_potential::get_local_maximum_location | ( | ) | const |
returns the location of the local maximum
The location of the local maximum (i.e., the potential barrier) is computed exactly using the solution to the cubic equation that results from imposing \( \partial V / \partial \phi = 0\). For \(\epsilon / (\lambda a^2) \ll 1\), one finds that the barrier is approximately located at \(\phi_b \approx \epsilon / (\lambda a^3)\).
Definition at line 148 of file thin_wall_potential.cpp.
References a, BubbleProfiler::Abs(), BubbleProfiler::ArcSin(), epsilon, lambda, origin, scale, BubbleProfiler::Sin(), and BubbleProfiler::Sqrt().
Referenced by add_constant_term(), and BubbleProfiler::solve_shooting().
double BubbleProfiler::Thin_wall_potential::get_local_minimum_location | ( | ) | const |
returns the location of the local minimum
The location of the local minimum is computed exactly using the solution to the cubic equation that results from imposing \( \partial V / \partial \phi = 0\). For \(\epsilon / (\lambda a^2) \ll 1\), one finds that the local minimum is located approximately at \(\phi_f \approx a\).
Definition at line 124 of file thin_wall_potential.cpp.
References a, BubbleProfiler::Abs(), BubbleProfiler::ArcSin(), epsilon, lambda, origin, scale, BubbleProfiler::Sin(), and BubbleProfiler::Sqrt().
Referenced by add_constant_term(), BubbleProfiler::solve_perturbative(), and BubbleProfiler::solve_shooting().
|
inlineoverridevirtual |
Implements BubbleProfiler::Potential.
Definition at line 60 of file thin_wall_potential.hpp.
double BubbleProfiler::Thin_wall_potential::get_thin_wall_action | ( | ) | const |
returns the approximate action in the thin-wall approximation
An approximate value for the Euclidean action is calculated using the analytic expression for the action derived in the thin-wall approximation,
\[ S \approx \frac{8 \pi^2 \lambda^2 a^{12}}{3 \epsilon^3} . \]
Definition at line 190 of file thin_wall_potential.cpp.
References a, epsilon, and lambda.
Referenced by add_constant_term(), and main().
|
overridevirtual |
Evaluate potential at point.
coords | Coordinates at which to evaluate |
Implements BubbleProfiler::Potential.
Definition at line 42 of file thin_wall_potential.cpp.
Referenced by add_constant_term(), and clone().
double BubbleProfiler::Thin_wall_potential::operator() | ( | double | coords | ) | const |
|
default |
|
default |
|
overridevirtual |
Partial derivative WRT coordinate i at a point.
coords | Coordinates at which to evaluate |
i | Index of coordinate to be differentiated |
Implements BubbleProfiler::Potential.
Definition at line 53 of file thin_wall_potential.cpp.
References first_deriv().
Referenced by clone().
|
overridevirtual |
Partial derivative WRT coordinates i, j at a a point.
coords | Coordinates at which to evaluate |
i | Index of first coordinate to be differentiated |
j | Index of second coordinate to be differentiated |
Implements BubbleProfiler::Potential.
Definition at line 70 of file thin_wall_potential.cpp.
References second_deriv().
double BubbleProfiler::Thin_wall_potential::second_deriv | ( | double | coords | ) | const |
evaluates the second derivative of the potential
coords | the value of the field to evaluate the derivative at |
Definition at line 109 of file thin_wall_potential.cpp.
References a, lambda, origin, and scale.
Referenced by add_constant_term(), and partial().
|
inlineoverridevirtual |
Shift the location of the origin by a specified vector.
translation | shift of origin |
Implements BubbleProfiler::Potential.
Definition at line 62 of file thin_wall_potential.hpp.
References origin.
|
private |
Definition at line 121 of file thin_wall_potential.hpp.
Referenced by first_deriv(), get_global_minimum_location(), get_local_maximum_location(), get_local_minimum_location(), get_thin_wall_action(), operator()(), second_deriv(), and Thin_wall_potential().
|
private |
Definition at line 122 of file thin_wall_potential.hpp.
Referenced by first_deriv(), get_global_minimum_location(), get_local_maximum_location(), get_local_minimum_location(), get_thin_wall_action(), operator()(), and Thin_wall_potential().
|
private |
Definition at line 120 of file thin_wall_potential.hpp.
Referenced by first_deriv(), get_global_minimum_location(), get_local_maximum_location(), get_local_minimum_location(), get_thin_wall_action(), operator()(), second_deriv(), and Thin_wall_potential().
|
private |
Definition at line 125 of file thin_wall_potential.hpp.
Referenced by add_constant_term().
|
private |
Definition at line 123 of file thin_wall_potential.hpp.
Referenced by first_deriv(), get_global_minimum_location(), get_local_maximum_location(), get_local_minimum_location(), operator()(), second_deriv(), and translate_origin().
|
private |
Definition at line 124 of file thin_wall_potential.hpp.
Referenced by apply_basis_change(), first_deriv(), get_global_minimum_location(), get_local_maximum_location(), get_local_minimum_location(), operator()(), and second_deriv().