Loading [MathJax]/extensions/tex2jax.js
BSMPT 3.0.7
BSMPT - Beyond the Standard Model Phase Transitions : A C++ package for the computation of the EWPT in BSM models
All Classes Namespaces Files Functions Variables Typedefs Enumerations Friends Pages
BSMPT::BounceActionInt Class Reference

Public Types

enum class  ActionStatus {
  Success , NotCalculated , Integration1DFailed , PathDeformationNotConverged ,
  PathDeformationCrashed , FalseVacuumNotMinimum , BackwardsPropagationFailed , NeverUndershootOvershoot ,
  UndershootOvershootNegativeGrad , NotEnoughPointsForSpline
}
 Possible status of the Action calculation.
 
enum class  UndershootOvershootStatus { Converged , Undershoot , Overshoot }
 Possible results of the undershoot/overshoot algorithm.
 
enum class  Integration1DStatus { Converged , NotConverged }
 Possible status of the 1D bounce solver.
 
enum class  PathDeformationStatus { Converged , NotConverged }
 Possible status of the path deformation algorithm.
 

Public Member Functions

 BounceActionInt ()
 Default constructor (unit tests)
 
 BounceActionInt (std::vector< std::vector< double > > InitPath_In, std::vector< double > TrueVacuum_In, std::vector< double > FalseVacuum_In, std::function< double(std::vector< double >)> &V_In, std::function< std::vector< double >(std::vector< double >)> &dV_In, double T_In, int MaxPathIntegrations_in)
 Construct a new Bounce Action Int object.
 
 BounceActionInt (std::vector< std::vector< double > > InitPath_In, std::vector< double > TrueVacuum_In, std::vector< double > FalseVacuum_In, std::function< double(std::vector< double >)> &V_In, double T_In, int MaxPathIntegrations_in)
 Construct a new Bounce Action Int object.
 
void SetPath (std::vector< std::vector< double > > InitPath_In)
 Used set the path of the class.
 
void RasterizedVdl (double l_start=0)
 Precalculates dVdl and creates a spline with the result. This is done to increase the runtime in large dimensional models.
 
void PrintVector (std::vector< double > vec)
 Prints a vector.
 
double Calc_dVdl (double l)
 Calculates \( \frac{dV}{dl} \) using the spline and potential derivatives.
 
double Calc_d2Vdl2 (double l)
 Calculates \( \frac{d^2V}{dl^2} \) using the spline and potential derivatives.
 
std::vector< double > NormalForce (const double &l, const double &dldrho, const std::vector< double > &gradient)
 Calculated the normal force \( \vec{N} \) on a spline point.
 
void AuxFunctionDev (const double &rho, const std::vector< double > &dvs, std::vector< double > &aks)
 Auxiliary function used in the Runge-Kutta 5th order adaptative step RK5_step.
 
void RK5_step (const std::vector< double > &y, const std::vector< double > &dydx, int n, float rho, float h, std::vector< double > &yout, std::vector< double > &yerr)
 Runge-Kutta 5th order step.
 
double BesselI (double alpha, double x, int terms=100)
 Modified Bessel function \(I_\alpha (x) \) of the first kind.
 
double BesselJ (double x, int terms=100)
 Modified Bessel function i \(J_1 (i x) \) of the first kind.
 
std::vector< double > ExactSolutionLin (double l0, double l, double dVdl, double d2Vdl2)
 Integrates the 1D profile assuming \( \frac{dV}{dl} \) is a linear in l, i.e. \( \frac{dV}{dl} \approx dV + H (l - l_0) \). The solution is for \( \alpha = 2 \) is.
 
std::vector< double > ExactSolutionFromMinimum (double l)
 Integrates the 1D profile assuming \( \frac{dV}{dl} \) is a linear in l, i.e. \( \frac{dV}{dl} \approx H (l - l_{min}) \). This correspondes to a purely qudratic potential.
 
double LogisticFunction (const double &x)
 Logistic function with patched edges to account for numerical instability/nans.
 
void CalculateExactSolutionThreshold (double MinError=1e100)
 Calculate \(l_\text{threshold}\) that splits the two branches of analytical integration. This function works recursively until the error MinError is small enough or the integration step is too small to be numerically stable.
 
std::vector< double > ExactSolution (double l0)
 Calculates the 1D solution by comparing the ExactSolutionCons and ExactSolutionLin so that the analytical step is appropriate.
 
void BackwardsPropagation ()
 This is done to make sure that we can still find solution after path deformation. This propagates the spline into negative values (by extrapolating).
 
void IntegrateBounce (double l0, UndershootOvershootStatus &conv, std::vector< double > &rho, std::vector< double > &l, std::vector< double > &dl_drho, std::vector< double > &d2l_drho2, int maxiter, double error, double eps_abs, double max_step=0)
 Integrates 1D bounce equation once.
 
void Solve1DBounce (std::vector< double > &rho, std::vector< double > &l, std::vector< double > &dl_drho, std::vector< double > &d2l_drho2, double error=1e-7, int maxiter=100)
 Performs a binary search using the overshooting/undershooting method to find the solution to the 1D bounce equation.
 
double ReductorCalculator (const double &MaximumGradient)
 Calculates the normalization of the force \( \vec{\phi} \rightarrow \vec{\phi} + \vec{N}/reductor \). We have that \( reductor = \varepsilon \max{\nabla V}/L \), where \( 10^{-4} \le \varepsilon \le 10^{-1} \) is a small parameter.
 
bool PathDeformationCheck (std::vector< double > &l, tk::spline &rho_l_spl)
 Check if the force in each point is sufficient small compared to the gradient in each point.
 
void SinglePathDeformation (double &stepsize, double &reductor, std::vector< double > &l, tk::spline &rho_l_spl, std::vector< double > &l_fornextpath, std::vector< std::vector< double > > &best_path, std::vector< std::vector< double > > &next_path, double &MaximumGradient, double &MaximumForce, double &MaximumRelativeError, double &Maximum_dldrho, double &PerpendicularGradient, MatrixXd &inverseK, std::vector< std::vector< double > > &forces)
 Takes a single path deformation step.
 
void PathDeformation (std::vector< double > &l, tk::spline &rho_l_spl)
 Deforms the path minimizing the force \( \vec{N} \) without solving.
 
unsigned nChoosek (unsigned n, unsigned k)
 Number of combinations of chosing k in n.
 
double Bernstein (int n, int nu, double x)
 Calculates the \( k^{th}\) Bernstein polynomial of degree \( n \) at \( x \).
 
std::vector< double > NormalForceBernstein (const double &dldrho, const std::vector< double > &gradient, const std::vector< double > &dphidl, const std::vector< double > &d2phidl2)
 Calculates the force vector \( \vec{N} \) of multiple path knots at the same time. Use in the path deformation algorithm.
 
double d2ldrho2 (double l, double rho, double dldrho)
 Calculates \( \frac{d^2l}{d\rho^2} \).
 
double CalculateKineticTermAction (const std::vector< double > &rho, const tk::spline &dl_drho_spl)
 Calculate kinect term of the action.
 
double CalculatePotentialTermAction (const std::vector< double > &rho, const tk::spline &l_rho_spl)
 Calculate potential term of the action.
 
void CalculateAction (double error=1e-6)
 Calculates the action of the bounce equation by deforming the given path and minimizing the normal force \( \vec{N} \) until it gets sufficiently small.
 

Public Attributes

int dim = -1
 Dimension of the VEV space.
 
double Action = -1
 either returns a -1 (if failed) or the value of the action
 
ActionStatus StateOfBounceActionInt = ActionStatus::NotCalculated
 Status of the Action calculation.
 
Integration1DStatus StateOf1DIntegration = Integration1DStatus::NotConverged
 Status of the 1D bounce solver.
 
PathDeformationStatus StateOfPathDeformation
 Status of the path deformation algorithm.
 
double Alpha = 2
 Factor produced by the spherical symmetry of the potential.
 
double T = -1
 Temperature of the potential. Irrelevant but helpful.
 
int MaxPathIntegrations
 Number of integration of the bounce.
 
int MaxSinglePathDeformations = 200
 Number of path deformations before integrating again.
 
std::vector< double > rho_sol
 list of \( \rho \) of the solution
 
std::vector< double > l_sol
 list of \( l(\rho) \) of the solution
 
std::vector< double > dldrho_sol
 list of \( \frac{l}{\rho} \) of the solution
 
double eps = 0.01
 Step for the numerical derivative.
 
int BernsteinDegree = 10
 Number of basis function that are used + 1.
 
double NumberPathKnots = 50
 // Number of knots in the new path
 
std::vector< double > TrueVacuum
 True vacuum candidate.
 
std::vector< double > FalseVacuum
 False vacuum. Should be the same as the last path knot.
 
std::function< double(std::vector< double >)> V
 Potential of the class.
 
std::function< std::vector< double >(std::vector< double >)> dV
 Potential gradient of the class, can be either numerical or analytical.
 
std::function< std::vector< std::vector< double > >(std::vector< double >)> Hessian
 Potential hessian of the class, completly numerical.
 
std::vector< std::vector< double > > InitPath
 First path given to class.
 
std::vector< std::vector< double > > Path
 Current class path, can be changed by PathDeformation.
 
cvspline Spline
 We describe the tunneling path using a cubid spline. The parameterization is the length along the spline.
 
tk::spline RasterizeddVdl
 Spline used to save \( \frac{dV}{dl} \).
 

Private Attributes

bool UndershotOnce = false
 Records if overshoot/undershoot method undershots at least once.
 
bool OvershotOnce = false
 Records if overshoot/undershoot method overshots at least once.
 
double TrueVacuumHessian
 Value of d2Vdl^2 near the true vacuum.
 
double Initial_lmin
 Store the value of backwards propagation.
 
double Vfalse
 Potential at the false vacuum in the unshifted potential (in the shifted Vfalse = 0)
 
double l0_minus_lmin
 l0 - Initial_lmin for solutions starting very near the true vacuum
 
bool PathDeformationConvergedWithout1D = false
 True if path deformation reached the desired results without solving the 1D equation one more time.
 
std::optional< double > ExactSolutionThreshold
 Lmin - L0 that splits between the two branches of the H > 0 analytical solution. If unset, we only use the linear approximation.
 
double FractionOfThePathExact = 1e-4
 First guess of l(rho) - l0 used to to integrate in the analytical solution. Can be decreased if the error is too high.
 

Constructor & Destructor Documentation

◆ BounceActionInt() [1/2]

BSMPT::BounceActionInt::BounceActionInt ( std::vector< std::vector< double > >  InitPath_In,
std::vector< double >  TrueVacuum_In,
std::vector< double >  FalseVacuum_In,
std::function< double(std::vector< double >)> &  V_In,
std::function< std::vector< double >(std::vector< double >)> &  dV_In,
double  T_In,
int  MaxPathIntegrations_in 
)

Construct a new Bounce Action Int object.

Parameters
init_pathis the initial path guess
TrueVacuumInis the true vacuum candidate of the potential
FalseVacuumInis the false vacuum
Vis the class potential
dVis the class gradient

◆ BounceActionInt() [2/2]

BSMPT::BounceActionInt::BounceActionInt ( std::vector< std::vector< double > >  InitPath_In,
std::vector< double >  TrueVacuum_In,
std::vector< double >  FalseVacuum_In,
std::function< double(std::vector< double >)> &  V_In,
double  T_In,
int  MaxPathIntegrations_in 
)

Construct a new Bounce Action Int object.

Parameters
init_pathis the initial path guess
TrueVacuumInis the true vacuum candidate of the potential
FalseVacuumInis the false vacuum
Vis the class potential

Member Function Documentation

◆ AuxFunctionDev()

void BSMPT::BounceActionInt::AuxFunctionDev ( const double &  rho,
const std::vector< double > &  dvs,
std::vector< double > &  aks 
)

Auxiliary function used in the Runge-Kutta 5th order adaptative step RK5_step.

Parameters
rhois the integration variable \( \rho \).
dvsare the functions values.
aksare used to return the function derivatives.

◆ BackwardsPropagation()

void BSMPT::BounceActionInt::BackwardsPropagation ( )

This is done to make sure that we can still find solution after path deformation. This propagates the spline into negative values (by extrapolating).

Returns
double negative or zero value.

◆ Bernstein()

double BSMPT::BounceActionInt::Bernstein ( int  n,
int  nu,
double  x 
)

Calculates the \( k^{th}\) Bernstein polynomial of degree \( n \) at \( x \).

\( B_{\nu,n}(x) = {n \choose \nu} x^\nu (1-x)^{n-\nu} \)

Parameters
nBernstein degree.
nuBernstein basis function, \( 0 \leq \nu \leq n \).
xis the independent parameter.
Returns
double \( B_{\nu,n}(x) \).

◆ BesselI()

double BSMPT::BounceActionInt::BesselI ( double  alpha,
double  x,
int  terms = 100 
)

Modified Bessel function \(I_\alpha (x) \) of the first kind.

\(I_\alpha (x) = \sum_{m=0}^\infty \frac{1}{m! \Gamma(m + \alpha + 1)} \left(\frac{x}{2}\right)^{2m + \alpha} \)

The convergence should be quite fast. We use as many terms to get the last term to impact the result in \( 10^{-15}\) (relatively). The default maximum number of terms is 50.

Parameters
alphais an complex number.
xwhere to calculate the Bessel function.
termsare the maximum number of terms allowed in the calculation.
Returns
double returns \(I_\alpha (x) \)

◆ BesselJ()

double BSMPT::BounceActionInt::BesselJ ( double  x,
int  terms = 100 
)

Modified Bessel function i \(J_1 (i x) \) of the first kind.

\(i J_1 (i x) = i \sum_{m=0}^\infty (-1)^{m + 1} \frac{1}{m! \Gamma(m + \alpha + 1)} \left(\frac{x}{2}\right)^{2m + \alpha} \)

The convergence should be quite fast. We use as many terms to get the last term to impact the result in \( 10^{-15}\) (relatively). The default maximum number of terms is 50.

Parameters
alphais an complex number.
xwhere to calculate the Bessel function.
termsare the maximum number of terms allowed in the calculation.
Returns
double returns \(I_\alpha (x) \)

◆ Calc_d2Vdl2()

double BSMPT::BounceActionInt::Calc_d2Vdl2 ( double  l)

Calculates \( \frac{d^2V}{dl^2} \) using the spline and potential derivatives.

\( \frac{d^2V}{dl^2} = \nabla V(\vec{\phi}) \cdot \frac{d^2 \vec{\phi}}{dl^2} + \left(\frac{d\vec{\phi}}{dl}\right)^T H(\vec{\phi}) \frac{d\vec{\phi}}{dl} \)

Returns
double Returns \( \frac{d^2V}{dl^2} \)

◆ Calc_dVdl()

double BSMPT::BounceActionInt::Calc_dVdl ( double  l)

Calculates \( \frac{dV}{dl} \) using the spline and potential derivatives.

\( \frac{dV}{dl} = \frac{\partial V}{\partial \vec{\phi}} \cdot \frac{d \vec{\phi}}{dl} = \nabla V(\vec{\phi}) \cdot \frac{d \vec{\phi}}{dl} \)

Parameters
lSpline parameterization point where \( \frac{dV}{dl} \) is calculated
Returns
double Returns \( \frac{dV}{dl} \)

◆ CalculateAction()

void BSMPT::BounceActionInt::CalculateAction ( double  error = 1e-6)

Calculates the action of the bounce equation by deforming the given path and minimizing the normal force \( \vec{N} \) until it gets sufficiently small.

Parameters
erroris the acceptance of undershoot/overshoot method.

◆ CalculateExactSolutionThreshold()

void BSMPT::BounceActionInt::CalculateExactSolutionThreshold ( double  MinError = 1e100)

Calculate \(l_\text{threshold}\) that splits the two branches of analytical integration. This function works recursively until the error MinError is small enough or the integration step is too small to be numerically stable.

Parameters
MinErrorLowest error found

◆ CalculateKineticTermAction()

double BSMPT::BounceActionInt::CalculateKineticTermAction ( const std::vector< double > &  rho,
const tk::spline dl_drho_spl 
)

Calculate kinect term of the action.

Parameters
rholist of rho coming from the integration
dl_drho_splSpline of \(\frac{dl(\rho)}{d\rho}\)
Returns
double kinetic part of the action

◆ CalculatePotentialTermAction()

double BSMPT::BounceActionInt::CalculatePotentialTermAction ( const std::vector< double > &  rho,
const tk::spline l_rho_spl 
)

Calculate potential term of the action.

Parameters
rholist of rho coming from the integration
l_rho_splSpline of \(l(\rho)\)
Returns
double potential part of the action

◆ d2ldrho2()

double BSMPT::BounceActionInt::d2ldrho2 ( double  l,
double  rho,
double  dldrho 
)

Calculates \( \frac{d^2l}{d\rho^2} \).

\( \frac{d^2l}{d\rho^2} = \frac{dV}{dl} - \frac{\alpha}{\rho}\frac{dl}{d\rho}\)

Parameters
l\( = l(\rho) \)
rhowhere we want to calculate.
dldrho\( = \frac{dl(\rho)}{d\rho} \)
Returns
double \( \frac{d^2l}{d\rho^2} \) at \( \rho \)

◆ ExactSolution()

std::vector< double > BSMPT::BounceActionInt::ExactSolution ( double  l0)

Calculates the 1D solution by comparing the ExactSolutionCons and ExactSolutionLin so that the analytical step is appropriate.

Parameters
l0starting point
Returns
std::vector<double>

◆ ExactSolutionFromMinimum()

std::vector< double > BSMPT::BounceActionInt::ExactSolutionFromMinimum ( double  l)

Integrates the 1D profile assuming \( \frac{dV}{dl} \) is a linear in l, i.e. \( \frac{dV}{dl} \approx H (l - l_{min}) \). This correspondes to a purely qudratic potential.

Parameters
lis the integration final point.
Returns
std::vector<double>

◆ ExactSolutionLin()

std::vector< double > BSMPT::BounceActionInt::ExactSolutionLin ( double  l0,
double  l,
double  dVdl,
double  d2Vdl2 
)

Integrates the 1D profile assuming \( \frac{dV}{dl} \) is a linear in l, i.e. \( \frac{dV}{dl} \approx dV + H (l - l_0) \). The solution is for \( \alpha = 2 \) is.

\( l(\rho) = l_0 - \frac{dV}{H} + \frac{dV}{H^{3/2}}\frac{\sinh(\rho\sqrt{H})}{\rho} \)

and the solution for \( \alpha = 3 \) is

\( l(\rho) = l_0 - \frac{dV}{H} + \frac{2 dV}{H^{3/2}} \frac{I_\alpha(\rho \sqrt{H})}{\rho}\)

where \( I_\alpha(\rho) \) is the modified Bessel function of the first kind.

Parameters
l0is the integration starting point.
lis the integration final point.
dVdl\( \frac{dV}{dl} \) at \( l_0 \).
d2Vdl2\( \frac{d^2V}{dl^2} \) at \( l_0 \).
maxiterMaximum iteration when calculating \( \rho(l) \).
Returns
std::vector<double>

◆ IntegrateBounce()

void BSMPT::BounceActionInt::IntegrateBounce ( double  l0,
UndershootOvershootStatus conv,
std::vector< double > &  rho,
std::vector< double > &  l,
std::vector< double > &  dl_drho,
std::vector< double > &  d2l_drho2,
int  maxiter,
double  error,
double  eps_abs,
double  max_step = 0 
)

Integrates 1D bounce equation once.

Parameters
l0is the starting value.
convchecks type of convergence. Converged. Undershoot. Overshoot.
rhovector of integration variable \( \rho \) steps.
lvector of variable \( l \) steps.
dl_drhovector of variable \( \frac{dl}{d\rho} \) steps.
d2l_drho2vector of variable \( \frac{d^2l}{d\rho^2} \) steps.
maxiteris the maximum integration steps.
erroris the acceptance of undershoot/overshoot.
eps_absis used to control the step size error (RK4 vs RK5).
max_stepin the case you want to set a maximum step size in \( \rho \).

◆ LogisticFunction()

double BSMPT::BounceActionInt::LogisticFunction ( const double &  x)

Logistic function with patched edges to account for numerical instability/nans.

\( \text{LogisticFunction}(x)=\begin{cases}1,\quad x\ge10\\0,\quad x\le -10\\\frac{1}{1+e^{-x}},\quad\text{otherwise}\end{cases} \)

Parameters
xindependent variable
Returns
double logistic function at at x

◆ nChoosek()

unsigned BSMPT::BounceActionInt::nChoosek ( unsigned  n,
unsigned  k 
)

Number of combinations of chosing k in n.

\( {n \choose k} = \frac{n!}{k!(n-k)!}\)

Parameters
n
k
Returns
unsigned

◆ NormalForce()

std::vector< double > BSMPT::BounceActionInt::NormalForce ( const double &  l,
const double &  dldrho,
const std::vector< double > &  gradient 
)

Calculated the normal force \( \vec{N} \) on a spline point.

Parameters
lis the spline parameter where the force is calculated.
dldrhois \( \frac{dl}{d\rho}\).
gradientis the gradient evaluated at spline parameter l.
Returns
std::vector<double> is the \( \vec{N} \) at spline parameter l.

◆ NormalForceBernstein()

std::vector< double > BSMPT::BounceActionInt::NormalForceBernstein ( const double &  dldrho,
const std::vector< double > &  gradient,
const std::vector< double > &  dphidl,
const std::vector< double > &  d2phidl2 
)

Calculates the force vector \( \vec{N} \) of multiple path knots at the same time. Use in the path deformation algorithm.

Parameters
dldrhois the vector of \( \frac{dl}{d\rho} \) values from the solution.
gradientis the vector of \( \nabla(\vec{\phi}(l)) \) values from the solution.
dphidlis the vector of \( \frac{d\vec{\phi}}{dl} \) values from the solution.
d2phidl2is the vector of \( \frac{d^2\vec{\phi}}{dl^2} \) values from the solution.
Returns
std::vector<double> list of force vectors \( \vec{N} \) applied to each path knot.

◆ PathDeformation()

void BSMPT::BounceActionInt::PathDeformation ( std::vector< double > &  l,
tk::spline rho_l_spl 
)

Deforms the path minimizing the force \( \vec{N} \) without solving.

Parameters
lis the vector of \( \rho \) values from the solution
rho_l_spltk::spline of \( \rho(l) \)

◆ PathDeformationCheck()

bool BSMPT::BounceActionInt::PathDeformationCheck ( std::vector< double > &  l,
tk::spline rho_l_spl 
)

Check if the force in each point is sufficient small compared to the gradient in each point.

Parameters
llist of spline parameter \( l \) of the knots
rho_l_spllist of \( \frac{dl}{d\rho}\) of the knots
Returns
true if converged
false if not converged

◆ PrintVector()

void BSMPT::BounceActionInt::PrintVector ( std::vector< double >  vec)

Prints a vector.

Parameters
vecVector to be printed

◆ RasterizedVdl()

void BSMPT::BounceActionInt::RasterizedVdl ( double  l_start = 0)

Precalculates dVdl and creates a spline with the result. This is done to increase the runtime in large dimensional models.

  • Parameters
    l_startis the starting position produced in the backwardspropagation part

◆ ReductorCalculator()

double BSMPT::BounceActionInt::ReductorCalculator ( const double &  MaximumGradient)

Calculates the normalization of the force \( \vec{\phi} \rightarrow \vec{\phi} + \vec{N}/reductor \). We have that \( reductor = \varepsilon \max{\nabla V}/L \), where \( 10^{-4} \le \varepsilon \le 10^{-1} \) is a small parameter.

Parameters
MaximumGradient
Returns
double

◆ RK5_step()

void BSMPT::BounceActionInt::RK5_step ( const std::vector< double > &  y,
const std::vector< double > &  dydx,
int  n,
float  rho,
float  h,
std::vector< double > &  yout,
std::vector< double > &  yerr 
)

Runge-Kutta 5th order step.

Although the Runge-Kutta methods are valid for \( y' = f(t,y) \) integration one can generalize the method for higher order ODE using auxiliary functions. One can write \( y'' = f(t,y, y') \) as

  • \( y' = m \)
  • \( m' = f(t, y, m) \)

linearizing the ODE system allowing the application of Runge-Kutta method to each one.

This method takes a 5th order Runge-Kutta step that can be embedded into a 4th order Runge-Kutta step. By comparing the result of the 4th and 5th order we can control the error within our calculation by adjusting the step size accordingly.

Parameters
yfunction values \( \{y(\rho_i), m(\rho_i) \} \).
dydxfunction derivatives \( \{y'(\rho_i), m'(\rho_i) \} \).
nis the number of linearizations, i.e. 2.
rhois the integration variable \( \rho \).
his the step size.
youtis the 5th order Runge-Kutta integration result.
yerris the difference between the 4th order and 5th order Runge-Kutta result.

◆ SetPath()

void BSMPT::BounceActionInt::SetPath ( std::vector< std::vector< double > >  InitPath_In)

Used set the path of the class.

Parameters
init_pathKnots that are used to describe the path

◆ SinglePathDeformation()

void BSMPT::BounceActionInt::SinglePathDeformation ( double &  stepsize,
double &  reductor,
std::vector< double > &  l,
tk::spline rho_l_spl,
std::vector< double > &  l_fornextpath,
std::vector< std::vector< double > > &  best_path,
std::vector< std::vector< double > > &  next_path,
double &  MaximumGradient,
double &  MaximumForce,
double &  MaximumRelativeError,
double &  Maximum_dldrho,
double &  PerpendicularGradient,
MatrixXd &  inverseK,
std::vector< std::vector< double > > &  forces 
)

Takes a single path deformation step.

Parameters
stepsize\( \varepsilon \)
reductoris the reductor
llist of \( l \) at the knots of the old solution
rho_l_spllist of \( \frac{dl}{d\rho} \) at the knots of the old solution
l_fornextpathlist of new \( l \) at the new path iteration
best_pathsave the best path
savesthe current iteration on the fly
MaximumGradientmaximum \( \nabla V \)
MaximumForcemaximum \( \vec{N} \)
MaximumRelativeErrormaximum \( \frac{|\vec{N}|}{|\nabla V|} \)
Maximum_dldrhomaximum \( \frac{dl}{d\rho} \)
PerpendicularGradientmaximum \( \nabla_\perp V \)
inverseKinverse of the kernel matrix
forceslist of forces to check if path is converging or not

◆ Solve1DBounce()

void BSMPT::BounceActionInt::Solve1DBounce ( std::vector< double > &  rho,
std::vector< double > &  l,
std::vector< double > &  dl_drho,
std::vector< double > &  d2l_drho2,
double  error = 1e-7,
int  maxiter = 100 
)

Performs a binary search using the overshooting/undershooting method to find the solution to the 1D bounce equation.

Parameters
rhois the vector of \( \rho \) values from the solution
lis the vector of \( l \) values from the solution
dl_drhois the vector of \( \frac{dl}{d\rho} \) values from the solution
d2l_drho2isthe vector of \( \frac{d^2l}{d\rho^2} \) values from the solution
erroris the IntegrateBounce error (acceptance of undershoot/overshoot)
maxiteris the maximum number of binary searches

Member Data Documentation

◆ Alpha

double BSMPT::BounceActionInt::Alpha = 2

Factor produced by the spherical symmetry of the potential.

  • = 2 if \( T > 0\) ( \(O(3)\) symmetry).
  • = 3 if \( T = 0\) ( \(O(4)\) symmetry).

Default value is \( 2 \) since most of our calculation are done at finite temperature.

◆ Hessian

std::function<std::vector<std::vector<double> >(std::vector<double>)> BSMPT::BounceActionInt::Hessian

Potential hessian of the class, completly numerical.

TODO: Calculate hessian from analytical gradient

◆ StateOfPathDeformation

PathDeformationStatus BSMPT::BounceActionInt::StateOfPathDeformation
Initial value:
=
PathDeformationStatus::NotConverged

Status of the path deformation algorithm.


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