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

Potential composed of two multivariate Gaussian functions. More...

#include <gaussian_potential.hpp>

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

Public Member Functions

 Gaussian_potential (double gamma_, double lambda_, std::size_t n_fields_)
 instantiate a new Gaussian potential More...
 
virtual ~Gaussian_potential ()=default
 
Gaussian_potentialclone () const
 return a copy of this potential More...
 
virtual double operator() (const Eigen::VectorXd &coords) const override
 Evaluate potential at point. More...
 
virtual double partial (const Eigen::VectorXd &coords, int i) const override
 Partial derivative WRT coordinate i at a point. More...
 
virtual double partial (const Eigen::VectorXd &coords, int i, int j) 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 &) override
 Shift the location of the origin by a specified vector. More...
 
virtual void apply_basis_change (const Eigen::MatrixXd &) override
 Apply a change of basis matrix. More...
 
virtual void add_constant_term (double) override
 Add a constant offset to the potential. More...
 
- Public Member Functions inherited from BubbleProfiler::Potential
virtual ~Potential ()=default
 

Private Member Functions

double gaussian (const Eigen::VectorXd &, const Eigen::VectorXd &) const
 
double gaussian_deriv (const Eigen::VectorXd &, const Eigen::VectorXd &, int i) const
 
double gaussian_deriv2 (const Eigen::VectorXd &, const Eigen::VectorXd &, int i, int j) const
 

Private Attributes

double gamma {2.}
 
double lambda {1.}
 
std::size_t n_fields {0}
 
double constant_term {0.}
 
Eigen::VectorXd origin {}
 
Eigen::VectorXd origin_translation {}
 
Eigen::VectorXd mu {}
 
Eigen::MatrixXd basis_transform {}
 

Detailed Description

Potential composed of two multivariate Gaussian functions.

This class implements a potential made of two (multivariate) Gaussian functions. The functional form of the potential is

\[ V(\phi) = -(N(\phi, 0) + \gamma N(\phi, \mu)), \]

where \(N(\phi,\mu) = \frac{1}{(2\pi)^{n/2}} \exp(-\frac{1}{2} | \phi - \mu|^2)\) is a unit Gaussian, \(\gamma\) controls the relative depth of the two minima, and \(\mu = \frac{1}{\sqrt{n_f}} (\lambda, \ldots, \lambda)\), where \(n_f\) is the number of fields and the parameter \(\lambda\) controls the geometric distance between the minima.

Note that we enforce \(\gamma > 1\), so that the Gaussian at the origin is a false vacuum.

Definition at line 53 of file gaussian_potential.hpp.

Constructor & Destructor Documentation

BubbleProfiler::Gaussian_potential::Gaussian_potential ( double  gamma_,
double  lambda_,
std::size_t  n_fields_ 
)

instantiate a new Gaussian potential

Parameters
gamma_relative depth of the two minima
lambda_geometric distance between minima ( \(\mu = \frac{1}{\sqrt{n_fields}}(\lambda, ..., \lambda)\))
n_fields_number of fields

Definition at line 31 of file gaussian_potential.cpp.

References basis_transform, gamma, lambda, mu, n_fields, origin, and origin_translation.

Referenced by clone().

virtual BubbleProfiler::Gaussian_potential::~Gaussian_potential ( )
virtualdefault

Member Function Documentation

void BubbleProfiler::Gaussian_potential::add_constant_term ( double  offset)
overridevirtual

Add a constant offset to the potential.

Parameters
offsetconstant term to add to potential

Implements BubbleProfiler::Potential.

Definition at line 137 of file gaussian_potential.cpp.

References constant_term.

Referenced by clone().

void BubbleProfiler::Gaussian_potential::apply_basis_change ( const Eigen::MatrixXd &  cob_matrix)
overridevirtual

Apply a change of basis matrix.

Parameters
cob_matrix

Implements BubbleProfiler::Potential.

Definition at line 131 of file gaussian_potential.cpp.

References basis_transform.

Referenced by clone().

Gaussian_potential* BubbleProfiler::Gaussian_potential::clone ( ) const
inlinevirtual

return a copy of this potential

Returns
a copy of this object

Implements BubbleProfiler::Potential.

Definition at line 70 of file gaussian_potential.hpp.

References add_constant_term(), apply_basis_change(), Gaussian_potential(), get_number_of_fields(), operator()(), partial(), and translate_origin().

double BubbleProfiler::Gaussian_potential::gaussian ( const Eigen::VectorXd &  coords,
const Eigen::VectorXd &  mu 
) const
private

Definition at line 53 of file gaussian_potential.cpp.

References n_fields.

Referenced by operator()().

double BubbleProfiler::Gaussian_potential::gaussian_deriv ( const Eigen::VectorXd &  coords,
const Eigen::VectorXd &  mu,
int  i 
) const
private

Definition at line 60 of file gaussian_potential.cpp.

References mu, and n_fields.

Referenced by partial().

double BubbleProfiler::Gaussian_potential::gaussian_deriv2 ( const Eigen::VectorXd &  coords,
const Eigen::VectorXd &  mu,
int  i,
int  j 
) const
private

Definition at line 68 of file gaussian_potential.cpp.

References mu, and n_fields.

Referenced by partial().

std::size_t BubbleProfiler::Gaussian_potential::get_number_of_fields ( ) const
overridevirtual

Implements BubbleProfiler::Potential.

Definition at line 142 of file gaussian_potential.cpp.

References n_fields.

Referenced by clone().

double BubbleProfiler::Gaussian_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 78 of file gaussian_potential.cpp.

References basis_transform, constant_term, gamma, gaussian(), mu, origin, and origin_translation.

Referenced by clone().

double BubbleProfiler::Gaussian_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 86 of file gaussian_potential.cpp.

References basis_transform, gamma, gaussian_deriv(), mu, n_fields, origin, and origin_translation.

Referenced by clone().

double BubbleProfiler::Gaussian_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 104 of file gaussian_potential.cpp.

References basis_transform, gamma, gaussian_deriv2(), mu, n_fields, origin, and origin_translation.

void BubbleProfiler::Gaussian_potential::translate_origin ( const Eigen::VectorXd &  translation)
overridevirtual

Shift the location of the origin by a specified vector.

Parameters
translationshift of origin

Implements BubbleProfiler::Potential.

Definition at line 126 of file gaussian_potential.cpp.

References origin_translation.

Referenced by clone().

Member Data Documentation

Eigen::MatrixXd BubbleProfiler::Gaussian_potential::basis_transform {}
private

Definition at line 92 of file gaussian_potential.hpp.

Referenced by apply_basis_change(), Gaussian_potential(), operator()(), and partial().

double BubbleProfiler::Gaussian_potential::constant_term {0.}
private

Definition at line 88 of file gaussian_potential.hpp.

Referenced by add_constant_term(), and operator()().

double BubbleProfiler::Gaussian_potential::gamma {2.}
private

Definition at line 85 of file gaussian_potential.hpp.

Referenced by Gaussian_potential(), operator()(), and partial().

double BubbleProfiler::Gaussian_potential::lambda {1.}
private

Definition at line 86 of file gaussian_potential.hpp.

Referenced by Gaussian_potential().

Eigen::VectorXd BubbleProfiler::Gaussian_potential::mu {}
private
std::size_t BubbleProfiler::Gaussian_potential::n_fields {0}
private
Eigen::VectorXd BubbleProfiler::Gaussian_potential::origin {}
private

Definition at line 89 of file gaussian_potential.hpp.

Referenced by Gaussian_potential(), operator()(), and partial().

Eigen::VectorXd BubbleProfiler::Gaussian_potential::origin_translation {}
private

Definition at line 90 of file gaussian_potential.hpp.

Referenced by Gaussian_potential(), operator()(), partial(), and translate_origin().


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