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

Example thin wall potential as given by Coleman. More...

#include <thin_wall_potential.hpp>

Inheritance diagram for BubbleProfiler::Thin_wall_potential:
Inheritance graph
[legend]
Collaboration diagram for BubbleProfiler::Thin_wall_potential:
Collaboration graph
[legend]

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_potentialoperator= (const Thin_wall_potential &)=default
 
Thin_wall_potentialoperator= (Thin_wall_potential &&)=default
 
virtual Thin_wall_potentialclone () 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...
 
- Public Member Functions inherited from BubbleProfiler::Potential
virtual ~Potential ()=default
 

Private Attributes

double lambda {1.}
 
double a {1.}
 
double epsilon {0.01}
 
double origin {0.}
 
double scale {1.}
 
double offset {0.}
 

Detailed Description

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.

Examples:
thin_wall.cpp.

Definition at line 40 of file thin_wall_potential.hpp.

Constructor & Destructor Documentation

BubbleProfiler::Thin_wall_potential::Thin_wall_potential ( )
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().

virtual BubbleProfiler::Thin_wall_potential::~Thin_wall_potential ( )
virtualdefault
BubbleProfiler::Thin_wall_potential::Thin_wall_potential ( const Thin_wall_potential )
default
BubbleProfiler::Thin_wall_potential::Thin_wall_potential ( Thin_wall_potential &&  )
default

Member Function Documentation

virtual void BubbleProfiler::Thin_wall_potential::add_constant_term ( double  offset)
inlineoverridevirtual

Add a constant offset to the potential.

Parameters
offsetconstant 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().

virtual void BubbleProfiler::Thin_wall_potential::apply_basis_change ( const Eigen::MatrixXd &  cob_matrix)
inlineoverridevirtual

Apply a change of basis matrix.

Parameters
cob_matrix

Implements BubbleProfiler::Potential.

Definition at line 66 of file thin_wall_potential.hpp.

References scale.

virtual Thin_wall_potential* BubbleProfiler::Thin_wall_potential::clone ( ) const
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

Parameters
coordsthe value of the field to evaluate the derivative at
Returns
the value of \(\partial V / \partial \phi\)

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

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

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

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().

virtual std::size_t BubbleProfiler::Thin_wall_potential::get_number_of_fields ( ) const
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

Returns
value of the Euclidean 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} . \]

Examples:
thin_wall.cpp.

Definition at line 190 of file thin_wall_potential.cpp.

References a, epsilon, and lambda.

Referenced by add_constant_term(), and main().

double BubbleProfiler::Thin_wall_potential::operator() ( const Eigen::VectorXd &  coords) const
overridevirtual

Evaluate potential at point.

Parameters
coordsCoordinates at which to evaluate
Returns
Value of potential at coordinates

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

evaluates the potential for the given field value

Parameters
coordsthe value of the field to evaluate the potential at
Returns
the value of the potential

Definition at line 93 of file thin_wall_potential.cpp.

References a, epsilon, lambda, origin, and scale.

Thin_wall_potential& BubbleProfiler::Thin_wall_potential::operator= ( const Thin_wall_potential )
default
Thin_wall_potential& BubbleProfiler::Thin_wall_potential::operator= ( Thin_wall_potential &&  )
default
double BubbleProfiler::Thin_wall_potential::partial ( const Eigen::VectorXd &  coords,
int  i 
) const
overridevirtual

Partial derivative WRT coordinate i at a point.

Parameters
coordsCoordinates at which to evaluate
iIndex of coordinate to be differentiated
Returns
Value of specified partial at point

Implements BubbleProfiler::Potential.

Definition at line 53 of file thin_wall_potential.cpp.

References first_deriv().

Referenced by clone().

double BubbleProfiler::Thin_wall_potential::partial ( const Eigen::VectorXd &  coords,
int  i,
int  j 
) const
overridevirtual

Partial derivative WRT coordinates i, j at a a point.

Parameters
coordsCoordinates at which to evaluate
iIndex of first coordinate to be differentiated
jIndex of second coordinate to be differentiated
Returns
Value of specified partial at point

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

Parameters
coordsthe value of the field to evaluate the derivative at
Returns
the value of \(\partial^ V / \partial \phi^2\)

Definition at line 109 of file thin_wall_potential.cpp.

References a, lambda, origin, and scale.

Referenced by add_constant_term(), and partial().

virtual void BubbleProfiler::Thin_wall_potential::translate_origin ( const Eigen::VectorXd &  translation)
inlineoverridevirtual

Shift the location of the origin by a specified vector.

Parameters
translationshift of origin

Implements BubbleProfiler::Potential.

Definition at line 62 of file thin_wall_potential.hpp.

References origin.

Member Data Documentation

double BubbleProfiler::Thin_wall_potential::a {1.}
private
double BubbleProfiler::Thin_wall_potential::epsilon {0.01}
private
double BubbleProfiler::Thin_wall_potential::lambda {1.}
private
double BubbleProfiler::Thin_wall_potential::offset {0.}
private

Definition at line 125 of file thin_wall_potential.hpp.

Referenced by add_constant_term().

double BubbleProfiler::Thin_wall_potential::origin {0.}
private
double BubbleProfiler::Thin_wall_potential::scale {1.}
private

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