BSMPT 3.0.7
BSMPT - Beyond the Standard Model Phase Transitions : A C++ package for the computation of the EWPT in BSM models
Loading...
Searching...
No Matches
action_calculation.h
Go to the documentation of this file.
1// Copyright (C) 2024 Lisa Biermann, Margarete Mühlleitner, Rui Santos, João
2// Viana
3//
4// SPDX-FileCopyrightText: 2024 Lisa Biermann, Margarete Mühlleitner, Rui
5// Santos, João Viana
6//
7// SPDX-License-Identifier: GPL-3.0-or-later
8
9#pragma once
10
15#include <BSMPT/utility/const_velocity_spline.h>
17#include <Eigen/Dense>
18#include <gsl/gsl_math.h>
19#include <numeric>
20#include <optional> // std::optional
21#include <sys/stat.h>
22#include <sys/types.h>
23
24using Eigen::MatrixXd;
25using Eigen::VectorXd;
26
32namespace BSMPT
33{
34
36{
37private:
38 // These two variables check if the method undershot and overshot at least
39 // once. If only one type of convergence was achieved it is an indication that
40 // no solution has been found.
44 bool UndershotOnce = false;
48 bool OvershotOnce = false;
49
55
61
67 double Vfalse;
68
74
81
87 std::optional<double> ExactSolutionThreshold;
88
95
96public:
101 int dim = -1;
102
107 enum class ActionStatus
108 {
109 Success,
110 NotCalculated,
111 Integration1DFailed,
112 PathDeformationNotConverged,
113 PathDeformationCrashed,
114 FalseVacuumNotMinimum,
115 BackwardsPropagationFailed,
116 NeverUndershootOvershoot,
117 UndershootOvershootNegativeGrad,
118 NotEnoughPointsForSpline
119 };
120
126 {
127 Converged,
128 Undershoot,
129 Overshoot
130 };
131
137 {
138 Converged,
139 NotConverged,
140 };
141
147 {
148 Converged,
149 NotConverged,
150 };
151
156 double Action = -1;
157
162 ActionStatus StateOfBounceActionInt = ActionStatus::NotCalculated;
163
168 Integration1DStatus StateOf1DIntegration = Integration1DStatus::NotConverged;
169
175 PathDeformationStatus::NotConverged;
176
185 double Alpha = 2; // alpha = 2 if T > 0 | alpha = 3 if T = 0
186
191 double T = -1;
192
198
204
208 std::vector<double> rho_sol;
209
214 std::vector<double> l_sol;
215
220 std::vector<double> dldrho_sol;
221
226 double eps = 0.01;
227
233
238 double NumberPathKnots = 50;
243 std::vector<double> TrueVacuum;
248 std::vector<double> FalseVacuum;
254 std::function<double(std::vector<double>)> V; // Potential
260 std::function<std::vector<double>(std::vector<double>)>
261 dV; // Potential gradient
267 std::function<std::vector<std::vector<double>>(std::vector<double>)> Hessian;
272 std::vector<std::vector<double>> InitPath; // Initial path given to class
277 std::vector<std::vector<double>> Path; // Initial path given to class
283 cvspline Spline; // Constant velocity spline object
288 tk::spline RasterizeddVdl; // RasterizeddVdl
289
295
306 std::vector<std::vector<double>> InitPath_In,
307 std::vector<double> TrueVacuum_In,
308 std::vector<double> FalseVacuum_In,
309 std::function<double(std::vector<double>)> &V_In,
310 std::function<std::vector<double>(std::vector<double>)> &dV_In,
311 double T_In,
312 int MaxPathIntegrations_in);
313
322 BounceActionInt(std::vector<std::vector<double>> InitPath_In,
323 std::vector<double> TrueVacuum_In,
324 std::vector<double> FalseVacuum_In,
325 std::function<double(std::vector<double>)> &V_In,
326 double T_In,
327 int MaxPathIntegrations_in);
328
334 void SetPath(std::vector<std::vector<double>> InitPath_In);
335
343 void RasterizedVdl(double l_start = 0);
344
350 void PrintVector(std::vector<double> vec);
351
364 double Calc_dVdl(double l);
365
376 double Calc_d2Vdl2(double l);
377
386 std::vector<double> NormalForce(const double &l,
387 const double &dldrho,
388 const std::vector<double> &gradient);
389
398 void AuxFunctionDev(const double &rho,
399 const std::vector<double> &dvs,
400 std::vector<double> &aks);
401
429 void RK5_step(const std::vector<double> &y,
430 const std::vector<double> &dydx,
431 int n,
432 float rho,
433 float h,
434 std::vector<double> &yout,
435 std::vector<double> &yerr);
436
452 double BesselI(double alpha, double x, int terms = 100);
453
469 double BesselJ(double x, int terms = 100);
470
494 std::vector<double>
495 ExactSolutionLin(double l0, double l, double dVdl, double d2Vdl2);
496
505 std::vector<double> ExactSolutionFromMinimum(double l);
506
517 double LogisticFunction(const double &x);
518
527 void CalculateExactSolutionThreshold(double MinError = 1e100);
528
537 std::vector<double> ExactSolution(double l0);
538
547
563 void IntegrateBounce(double l0,
565 std::vector<double> &rho,
566 std::vector<double> &l,
567 std::vector<double> &dl_drho,
568 std::vector<double> &d2l_drho2,
569 int maxiter,
570 double error,
571 double eps_abs,
572 double max_step = 0);
587 void Solve1DBounce(std::vector<double> &rho,
588 std::vector<double> &l,
589 std::vector<double> &dl_drho,
590 std::vector<double> &d2l_drho2,
591 double error = 1e-7,
592 int maxiter = 100);
593
603 double ReductorCalculator(const double &MaximumGradient);
604
614 bool PathDeformationCheck(std::vector<double> &l, tk::spline &rho_l_spl);
615
635 void SinglePathDeformation(double &stepsize,
636 double &reductor,
637 std::vector<double> &l,
638 tk::spline &rho_l_spl,
639 std::vector<double> &l_fornextpath,
640 std::vector<std::vector<double>> &best_path,
641 std::vector<std::vector<double>> &next_path,
642 double &MaximumGradient,
643 double &MaximumForce,
644 double &MaximumRelativeError,
645 double &Maximum_dldrho,
646 double &PerpendicularGradient,
647 MatrixXd &inverseK,
648 std::vector<std::vector<double>> &forces);
649
657 void PathDeformation(std::vector<double> &l, tk::spline &rho_l_spl);
667 unsigned nChoosek(unsigned n, unsigned k);
668
680 double Bernstein(int n, int nu, double x);
681
697 std::vector<double> NormalForceBernstein(const double &dldrho,
698 const std::vector<double> &gradient,
699 const std::vector<double> &dphidl,
700 const std::vector<double> &d2phidl2);
701
713 double d2ldrho2(double l, double rho, double dldrho);
714
722 double CalculateKineticTermAction(const std::vector<double> &rho,
723 const tk::spline &dl_drho_spl);
724
732 double CalculatePotentialTermAction(const std::vector<double> &rho,
733 const tk::spline &l_rho_spl);
734
743 void CalculateAction(double error = 1e-6);
744};
745} // namespace BSMPT
Definition action_calculation.h:36
int dim
Dimension of the VEV space.
Definition action_calculation.h:101
UndershootOvershootStatus
Possible results of the undershoot/overshoot algorithm.
Definition action_calculation.h:126
void CalculateExactSolutionThreshold(double MinError=1e100)
Calculate that splits the two branches of analytical integration. This function works recursively un...
Definition action_calculation.cpp:536
std::vector< double > NormalForce(const double &l, const double &dldrho, const std::vector< double > &gradient)
Calculated the normal force on a spline point.
Definition action_calculation.cpp:125
std::optional< double > ExactSolutionThreshold
Lmin - L0 that splits between the two branches of the H > 0 analytical solution. If unset,...
Definition action_calculation.h:87
std::function< std::vector< double >(std::vector< double >)> dV
Potential gradient of the class, can be either numerical or analytical.
Definition action_calculation.h:261
std::vector< std::vector< double > > InitPath
First path given to class.
Definition action_calculation.h:272
double NumberPathKnots
// Number of knots in the new path
Definition action_calculation.h:238
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.
Definition action_calculation.cpp:1011
double BesselI(double alpha, double x, int terms=100)
Modified Bessel function of the first kind.
Definition action_calculation.cpp:232
void PathDeformation(std::vector< double > &l, tk::spline &rho_l_spl)
Deforms the path minimizing the force without solving.
Definition action_calculation.cpp:1309
bool UndershotOnce
Records if overshoot/undershoot method undershots at least once.
Definition action_calculation.h:44
std::vector< double > rho_sol
list of of the solution
Definition action_calculation.h:208
Integration1DStatus
Possible status of the 1D bounce solver.
Definition action_calculation.h:137
PathDeformationStatus StateOfPathDeformation
Status of the path deformation algorithm.
Definition action_calculation.h:174
std::vector< double > dldrho_sol
list of of the solution
Definition action_calculation.h:220
std::function< double(std::vector< double >)> V
Potential of the class.
Definition action_calculation.h:254
double TrueVacuumHessian
Value of d2Vdl^2 near the true vacuum.
Definition action_calculation.h:54
double BesselJ(double x, int terms=100)
Modified Bessel function i of the first kind.
Definition action_calculation.cpp:249
double Action
either returns a -1 (if failed) or the value of the action
Definition action_calculation.h:156
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.
Definition action_calculation.cpp:1092
double LogisticFunction(const double &x)
Logistic function with patched edges to account for numerical instability/nans.
Definition action_calculation.cpp:529
std::vector< double > ExactSolution(double l0)
Calculates the 1D solution by comparing the ExactSolutionCons and ExactSolutionLin so that the analyt...
Definition action_calculation.cpp:572
unsigned nChoosek(unsigned n, unsigned k)
Number of combinations of chosing k in n.
Definition action_calculation.cpp:1458
void PrintVector(std::vector< double > vec)
Prints a vector.
Definition action_calculation.cpp:114
double l0_minus_lmin
l0 - Initial_lmin for solutions starting very near the true vacuum
Definition action_calculation.h:73
ActionStatus StateOfBounceActionInt
Status of the Action calculation.
Definition action_calculation.h:162
void BackwardsPropagation()
This is done to make sure that we can still find solution after path deformation. This propagates the...
Definition action_calculation.cpp:731
double Vfalse
Potential at the false vacuum in the unshifted potential (in the shifted Vfalse = 0)
Definition action_calculation.h:67
double ReductorCalculator(const double &MaximumGradient)
Calculates the normalization of the force . We have that , where is a small parameter.
Definition action_calculation.cpp:1004
bool OvershotOnce
Records if overshoot/undershoot method overshots at least once.
Definition action_calculation.h:48
double eps
Step for the numerical derivative.
Definition action_calculation.h:226
int BernsteinDegree
Number of basis function that are used + 1.
Definition action_calculation.h:232
void SetPath(std::vector< std::vector< double > > InitPath_In)
Used set the path of the class.
Definition action_calculation.cpp:71
double CalculatePotentialTermAction(const std::vector< double > &rho, const tk::spline &l_rho_spl)
Calculate potential term of the action.
Definition action_calculation.cpp:1512
double d2ldrho2(double l, double rho, double dldrho)
Calculates .
Definition action_calculation.cpp:153
int MaxSinglePathDeformations
Number of path deformations before integrating again.
Definition action_calculation.h:203
std::vector< std::vector< double > > Path
Current class path, can be changed by PathDeformation.
Definition action_calculation.h:277
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 b...
Definition action_calculation.cpp:781
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.
Definition action_calculation.cpp:172
double T
Temperature of the potential. Irrelevant but helpful.
Definition action_calculation.h:191
Integration1DStatus StateOf1DIntegration
Status of the 1D bounce solver.
Definition action_calculation.h:168
BounceActionInt()
Default constructor (unit tests)
Definition action_calculation.cpp:16
void CalculateAction(double error=1e-6)
Calculates the action of the bounce equation by deforming the given path and minimizing the normal fo...
Definition action_calculation.cpp:1546
double FractionOfThePathExact
First guess of l(rho) - l0 used to to integrate in the analytical solution. Can be decreased if the e...
Definition action_calculation.h:94
double Bernstein(int n, int nu, double x)
Calculates the Bernstein polynomial of degree at .
Definition action_calculation.cpp:987
std::vector< double > l_sol
list of of the solution
Definition action_calculation.h:214
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.
Definition action_calculation.cpp:628
int MaxPathIntegrations
Number of integration of the bounce.
Definition action_calculation.h:197
cvspline Spline
We describe the tunneling path using a cubid spline. The parameterization is the length along the spl...
Definition action_calculation.h:283
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.
Definition action_calculation.cpp:162
std::vector< double > ExactSolutionLin(double l0, double l, double dVdl, double d2Vdl2)
Integrates the 1D profile assuming is a linear in l, i.e. . The solution is for is.
Definition action_calculation.cpp:337
double CalculateKineticTermAction(const std::vector< double > &rho, const tk::spline &dl_drho_spl)
Calculate kinect term of the action.
Definition action_calculation.cpp:1475
PathDeformationStatus
Possible status of the path deformation algorithm.
Definition action_calculation.h:147
double Calc_d2Vdl2(double l)
Calculates using the spline and potential derivatives.
Definition action_calculation.cpp:522
std::vector< double > FalseVacuum
False vacuum. Should be the same as the last path knot.
Definition action_calculation.h:248
std::vector< double > ExactSolutionFromMinimum(double l)
Integrates the 1D profile assuming is a linear in l, i.e. . This correspondes to a purely qudratic p...
Definition action_calculation.cpp:266
std::function< std::vector< std::vector< double > >(std::vector< double >)> Hessian
Potential hessian of the class, completly numerical.
Definition action_calculation.h:267
ActionStatus
Possible status of the Action calculation.
Definition action_calculation.h:108
tk::spline RasterizeddVdl
Spline used to save .
Definition action_calculation.h:288
double Initial_lmin
Store the value of backwards propagation.
Definition action_calculation.h:60
bool PathDeformationConvergedWithout1D
True if path deformation reached the desired results without solving the 1D equation one more time.
Definition action_calculation.h:80
void RasterizedVdl(double l_start=0)
Precalculates dVdl and creates a spline with the result. This is done to increase the runtime in larg...
Definition action_calculation.cpp:97
double Calc_dVdl(double l)
Calculates using the spline and potential derivatives.
Definition action_calculation.cpp:516
double Alpha
Factor produced by the spherical symmetry of the potential.
Definition action_calculation.h:185
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 of multiple path knots at the same time. Use in the path deformation alg...
Definition action_calculation.cpp:141
std::vector< double > TrueVacuum
True vacuum candidate.
Definition action_calculation.h:243
Constructs a spline with constant velocity, i.e. . acts as the length alongisde the spline....
Definition const_velocity_spline.h:34
Definition spline.h:40
This classes calculates the Bounce action of the potential with a set temperature.
Definition CalculateEtaInterface.h:24