BSMPT 3.1.3
BSMPT - Beyond the Standard Model Phase Transitions : A C++ package for the computation of the EWPT in BSM models
Loading...
Searching...
No Matches
bounce_solution.h
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
16#include <BSMPT/bounce_solution/gstar.h>
17#include <BSMPT/minimum_tracer/minimum_tracer.h> // MinimumTracer
19#include <BSMPT/utility/spline/spline.h>
20#include <Eigen/Dense>
21#include <algorithm> // std::max
22#include <gsl/gsl_deriv.h> // numerical derivative
23#include <gsl/gsl_integration.h> // numerical integration
24namespace BSMPT
25{
26
31{
32 double result;
33 double error;
34};
35
41{
42public:
46 std::shared_ptr<Class_Potential_Origin> modelPointer;
47
51 std::shared_ptr<MinimumTracer> MinTracer;
52
58
62 double epsturb = 0.1;
63
67 double vwall = 0.95;
68
72 double vCJ = -1;
73
77 double Csound_false = -1;
78
82 double Csound_true = -1;
83
89
93 bool nucleation_temp_set = false;
94
99
104
108 TransitionTemperature which_transition_temp = TransitionTemperature::NotSet;
109
113 double Tc = -1;
114
118 double Tm = -1;
119
123 double Tnucl = -1;
124
128 double Tnucl_approx = -1;
129
133 double Tperc = -1;
134
138 double Tcompl = -1;
139
143 double Treh = -1;
144
148 double Tstar = -1;
149
154
158 double alpha = -1;
159
164 double betaH = -1;
165
172 double Rstar = -1;
173
177 double gstar;
178
184
191
197
203
209
215
220 std::vector<BounceActionInt> SolutionList;
221
226 std::vector<Eigen::MatrixXd> GroupElements;
227
232 Eigen::MatrixXd OptimalDiscreteSymmetry;
233
240 std::vector<double>
241 TransformIntoOptimalDiscreteSymmetry(const std::vector<double> &vev);
242
249
255 double TunnelingRate(const double &Temp);
260 double HubbleRate(const double &Temp);
261
266 friend double inner_integrand(double var, void *params);
267
272 friend double outer_integrand(double var, void *params);
273
277 double GetBounceSol(const double &Temp) const;
278
283 friend double action_ratio(double var, void *params);
284
285public:
289 const double AbsErr = 0;
290
294 const double RelErr = 1e-6;
295
302
308
314
319
324 StatusGW status_bounce_sol = StatusGW::NotSet;
325
331 const TransitionTemperature &which_transition_temp_in);
332
338 BSMPT::StatusTemperature::NotSet;
339
344 BSMPT::StatusTemperature status_nucl = BSMPT::StatusTemperature::NotSet;
345
350 BSMPT::StatusTemperature status_perc = BSMPT::StatusTemperature::NotSet;
351
356 BSMPT::StatusTemperature status_compl = BSMPT::StatusTemperature::NotSet;
357
366 double UserDefined_vwall = 0.95;
367
373
378 BounceSolution(const std::shared_ptr<Class_Potential_Origin> &pointer_in);
379
400 BounceSolution(const std::shared_ptr<Class_Potential_Origin> &pointer_in,
401 const std::shared_ptr<MinimumTracer> &MinTracer_in,
402 const CoexPhases &phase_pair_in,
403 const double &UserDefined_vwall_in,
404 const double &UserDefined_epsturb_in,
405 const int &MaxPathIntegrations_in,
406 const size_t &NumberOfInitialScanTemperatures_in,
407 const int &UserDefined_PNLO_scaling_in = 1);
408
429 BounceSolution(const std::shared_ptr<Class_Potential_Origin> &pointer_in,
430 const std::shared_ptr<MinimumTracer> &MinTracer_in,
431 const CoexPhases &phase_pair_in,
432 const double &UserDefined_vwall_in,
433 const double &UserDefined_epsturb_in,
434 const int &MaxPathIntegrations_in,
435 const size_t &NumberOfInitialScanTemperatures_in,
436 const std::vector<Eigen::MatrixXd> &GroupElements_in,
437 const int &UserDefined_PNLO_scaling_in = 1);
438
443 void GWInitialScan();
444
453 void CalculateActionAt(double T, bool smart = true);
454
461 void GWSecondaryScan();
462
468
474
479 void SetBounceSol();
480
485 double GetWallVelocity() const;
486
491 double GetChapmanJougetVelocity() const;
492
496 double GetSoundSpeedFalse() const;
497
501 double GetSoundSpeedTrue() const;
502
507 double GetEpsTurb() const;
508
512 void SetGstar(const double &gstar_in);
513
519
526 static void ConstructSplineVofT(Phase &phase, tk::spline &spline);
527
533 void InitializedVSpline();
534
541 double GetGstar(const double &T) const;
542
548 double GetGstar();
549
553 void SetCriticalTemp(const double &T_in);
554
558 double GetCriticalTemp() const;
559
563 void SetStoredTemp(const double &T_in);
567 double GetStoredTemp() const;
568
572 double GetNucleationTemp() const;
573
578 double GetNucleationTempApprox() const;
579
583 double GetPercolationTemp() const;
584
588 double GetCompletionTemp() const;
589
593 double GetTransitionTemp() const;
594
598 double GetReheatingTemp() const;
599
603 void CalcTransitionTemp();
604
612 double CalculateRhoGamma(const double &T) const;
613
617 double GetPTStrength() const;
618
623 double CalcGstarPureRad();
624
629
634
644 double FalseVacFractionExponent_I(const double &T);
645
652 double CalcTempAtFalseVacFraction(const double &false_vac_frac);
653
660 double CalcFalseVacFraction(const double &temp);
661
668 void CalculatePercolationTemp(const double &false_vac_frac = 0.71);
669
676 void CalculateCompletionTemp(const double &false_vac_frac = 0.01);
677
682
686 void CalculatePTStrength();
687
693
699 void CalculateWallVelocity(const Minimum &false_min, const Minimum &true_min);
700
707 double CalculateSoundSpeed(Phase &phase);
708
715
720
724 double GetInvTimeScale();
725
732 void CalculateRstar();
733
738 double GetRstar();
739};
740
749 const double &Tprime);
750
758
766
767} // namespace BSMPT
BounceSolution class that handles the calculation of the bounce solution as well as the calculation o...
Definition bounce_solution.h:41
double Csound_false
sound speed in false phase
Definition bounce_solution.h:77
void CalculateWallVelocity(const Minimum &false_min, const Minimum &true_min)
Calculate wall velocity.
Definition bounce_solution.cpp:1130
friend double action_ratio(double var, void *params)
action_ratio friend to define input of numerical derivative in calculation of inverse time scale
Definition bounce_solution.cpp:1216
double betaH
Inverse time scale .
Definition bounce_solution.h:164
size_t NumberOfInitialScanTemperatures
number of temperature steps in the initial scan of the bounce solver
Definition bounce_solution.h:88
double Tnucl_approx
approximate nucleation temperature
Definition bounce_solution.h:128
void CalculateOptimalDiscreteSymmetry()
Calculates which is the optimal symmetry from the group of symmetries.
Definition bounce_solution.cpp:77
double GetBounceSol(const double &Temp) const
Calculate euclidian action at temperature T.
Definition bounce_solution.cpp:1211
void SetAndCalculateGWParameters(const TransitionTemperature &which_transition_temp_in)
Set the Transition Temp object.
Definition bounce_solution.cpp:110
void CalcChapmanJougetVelocity()
Derive the Chapman-Jouget velocity from PT strength and false phase sound velocity using Eq....
Definition bounce_solution.cpp:1122
double UserDefined_vwall
defined by the user as an input parameter. If v_{\text{wall}}\f = -2$ then we use the upper bound fr...
Definition bounce_solution.h:366
double GetPTStrength() const
GetPTStrength Get PT strength alpha.
Definition bounce_solution.cpp:666
std::vector< BounceActionInt > SolutionList
Set of BounceActionInt objects with valid solutions.
Definition bounce_solution.h:220
const double RelativeErrorInCalcTempAtFalseVacFraction
Maximum relative error on the fraction of vacuum tunneled to be accepted.
Definition bounce_solution.h:307
void InitializedVSpline()
Initialize two splines for the potential across the tunneling profile. Used to improve the Hubble rat...
Definition bounce_solution.cpp:535
friend double inner_integrand(double var, void *params)
inner_integrand friend to define inner integrand of percolation temperature integral
Definition bounce_solution.cpp:990
void SetCriticalTemp(const double &T_in)
SetCriticalTemp Set critical temperature.
Definition bounce_solution.cpp:565
tk::spline GstarProfileHighT
Gstar spline, T > T_QCD (214.0 MeV)
Definition bounce_solution.h:202
double GetRstar()
Returns .
Definition bounce_solution.cpp:1251
std::vector< Eigen::MatrixXd > GroupElements
List of group elements allowed by the potential.
Definition bounce_solution.h:226
BSMPT::StatusTemperature status_compl
status of completion temperature calculation
Definition bounce_solution.h:356
tk::spline GstarProfileLowT
Gstar spline, T < T_QCD (214.0 MeV)
Definition bounce_solution.h:196
bool completion_temp_set
set to true if completion temperature is set
Definition bounce_solution.h:103
void CalcTransitionTemp()
CalcTransitionTemp Get transition temperature from int.
Definition bounce_solution.cpp:615
double TunnelingRate(const double &Temp)
Storage of the tunneling rate per volume of the transition from false to true vacuum.
Definition bounce_solution.cpp:671
void GWInitialScan()
Initially we have no idea where the transition can occur, therefore we scan the complete temperature ...
Definition bounce_solution.cpp:135
double GetInvTimeScale()
Get inverse time scale of phase transition.
Definition bounce_solution.cpp:1245
double CalcTempAtFalseVacFraction(const double &false_vac_frac)
CalcTempAtFalseVacFraction calculates the temperature at which the false vacuum fraction drops below ...
Definition bounce_solution.cpp:839
bool percolation_temp_set
set to true if percolation temperature is set
Definition bounce_solution.h:98
void GWScanTowardsLowAction()
Do linear extrapolations to calculate action at lower temperatures.
Definition bounce_solution.cpp:394
double GetChapmanJougetVelocity() const
Get Chapman-Jouget velocity.
Definition bounce_solution.cpp:485
double Tstar
transition temperature
Definition bounce_solution.h:148
double GetSoundSpeedFalse() const
Get the sound speed in the false phase.
Definition bounce_solution.cpp:490
double Csound_true
sound speed in true phase
Definition bounce_solution.h:82
void CalculateRstar()
Definition bounce_solution.cpp:1257
const double RelativeTemperatureInCalcTempAtFalseVacFraction
Maximum relative difference in temperature on the fraction of false vacuum to be accepted.
Definition bounce_solution.h:301
double HubbleRate(const double &Temp)
Storage of the temperature-dependent Hubble rate.
Definition bounce_solution.cpp:690
double CalcGstarPureRad()
CalcGstarPureRad Calculate the number of effective degrees of freedom assuming a purely radiative uni...
Definition bounce_solution.cpp:700
int MaxPathIntegrations
Number of integration of the bounce.
Definition bounce_solution.h:372
double Rstar
Definition bounce_solution.h:172
double GetPercolationTemp() const
GetPercolationTemp Get percolation temperature.
Definition bounce_solution.cpp:595
void CalculateActionAt(double T, bool smart=true)
Calculate the euclidian action of the transition from false to true phase of phase pair.
Definition bounce_solution.cpp:204
void GWSecondaryScan()
If solution were found by the GWInitialScan() then we scan temperature range in the vicinity such tha...
Definition bounce_solution.cpp:281
BSMPT::StatusTemperature status_perc
status of percolation temperature calculation
Definition bounce_solution.h:350
double Tperc
percolation temperature
Definition bounce_solution.h:133
double CalculateSoundSpeed(Phase &phase)
Calculate sound speeds at Tstar in phase.
Definition bounce_solution.cpp:1173
void SetGstar(const double &gstar_in)
SetGstar Set gstar.
Definition bounce_solution.cpp:505
friend double outer_integrand(double var, void *params)
outer_integrand friend to define outer integrand of percolation temperature integral
Definition bounce_solution.cpp:997
Eigen::MatrixXd OptimalDiscreteSymmetry
Store symmetry that produces the best tunneling rate.
Definition bounce_solution.h:232
tk::spline FalsePhaseVSpline
False V spline to interpolate.
Definition bounce_solution.h:208
std::shared_ptr< MinimumTracer > MinTracer
MinTracer object.
Definition bounce_solution.h:51
int pnlo_scaling
pressure scaling with of 1 -> N processes at NLO
Definition bounce_solution.h:57
void CalculateCompletionTemp(const double &false_vac_frac=0.01)
CalculateCompletionTemp calculation of the temperature when the false vacuum fraction drops below 1 %...
Definition bounce_solution.cpp:941
CoexPhases phase_pair
pair of coexisiting phases
Definition bounce_solution.h:318
double GetCriticalTemp() const
GetCriticalTemp Get critical temperature.
Definition bounce_solution.cpp:570
double GetNucleationTempApprox() const
GetNucleationTempApprox Get nucleation temperature via approximate method.
Definition bounce_solution.cpp:590
void CalculateNucleationTemp()
Calculation of nucleation temperature.
Definition bounce_solution.cpp:711
void InitializeGstarProfile()
Generate the spline used to interpolate the gstar SM profile.
Definition bounce_solution.cpp:510
double GetSoundSpeedTrue() const
Get the sound speed in the true phase.
Definition bounce_solution.cpp:495
double Tm
lowest temperature when a transition can occur
Definition bounce_solution.h:118
void SetBounceSol()
Set the Bounce Sol object.
Definition bounce_solution.cpp:432
const double AbsErr
AbsErr absolute error for numerical integration.
Definition bounce_solution.h:289
double Tc
critical temperature/highest temperature when transition can occur
Definition bounce_solution.h:113
int indexTrueCandidatePhase
index of the true vacuum phase candidate in the coex list
Definition bounce_solution.h:183
double store_Temp
stored temperature
Definition bounce_solution.h:153
void CalculatePTStrength()
Calculate phase transition strength alpha.
Definition bounce_solution.cpp:1071
std::vector< double > TransformIntoOptimalDiscreteSymmetry(const std::vector< double > &vev)
Transforms the vev to the optimal vev.
Definition bounce_solution.cpp:124
double GetNucleationTemp() const
GetNucleationTemp Get nucleation temperature via exact method.
Definition bounce_solution.cpp:585
double FalseVacFractionExponent_I(const double &T)
Calculate the false vacuum fraction .
Definition bounce_solution.cpp:827
static void ConstructSplineVofT(Phase &phase, tk::spline &spline)
Using the phase, constructs a spline of of that phase.
Definition bounce_solution.cpp:521
double Tnucl
nucleation temperature
Definition bounce_solution.h:123
double GetTransitionTemp() const
GetTransitionTemp Get transition temperature.
Definition bounce_solution.cpp:605
void CalculateSoundSpeeds()
Calculate sound speeds at Tstar in false and true phase.
Definition bounce_solution.cpp:1203
double GetReheatingTemp() const
GetReheatingTemp Get reheating temperature.
Definition bounce_solution.cpp:610
double Treh
reheating temperature
Definition bounce_solution.h:143
double vCJ
Chapman-Jouget velocity.
Definition bounce_solution.h:72
void GWScanTowardsHighAction()
Do linear extrapolations to calculate action at higher temperatures.
Definition bounce_solution.cpp:355
void SetStoredTemp(const double &T_in)
SetStoredTemp Set stored temperature.
Definition bounce_solution.cpp:575
double GetGstar()
Get Gstar for radiation-dominated epoch.
Definition bounce_solution.cpp:560
double Tcompl
completion temperature
Definition bounce_solution.h:138
void CalculatePercolationTemp(const double &false_vac_frac=0.71)
CalculatePercolationTemp calculation of the temperature when the false vacuum fraction drops below 71...
Definition bounce_solution.cpp:910
double vwall
wall velocity
Definition bounce_solution.h:67
void CalculateNucleationTempApprox()
Approximate calculation of nucleation temperature.
Definition bounce_solution.cpp:772
double alpha
PT strength.
Definition bounce_solution.h:158
std::shared_ptr< Class_Potential_Origin > modelPointer
modelPointer for the used parameter point
Definition bounce_solution.h:46
double GetStoredTemp() const
GetStoredTemp Get stored temperature.
Definition bounce_solution.cpp:580
BSMPT::StatusTemperature status_nucl
status of nucleation temperature calculation
Definition bounce_solution.h:344
double epsturb
epsilon of turbulence efficiency factor
Definition bounce_solution.h:62
double GetCompletionTemp() const
GetCompletionTemp Get percolation temperature.
Definition bounce_solution.cpp:600
bool nucleation_temp_set
set to true if nucleation temperature is set
Definition bounce_solution.h:93
void CalculateInvTimeScale()
Calculate inverse time scale of phase transition.
Definition bounce_solution.cpp:1239
BSMPT::StatusTemperature status_nucl_approx
status of approximate nucleation temperature calculation
Definition bounce_solution.h:337
double GetWallVelocity() const
Get the bubble wall velocity.
Definition bounce_solution.cpp:480
tk::spline S3ofT_spline
spline used to interpolate the action as a function of the temperature
Definition bounce_solution.h:190
double gstar
number of effective degrees of freedom
Definition bounce_solution.h:177
void CalculateReheatingTemp()
CalculateReheatingTemp calculation of the reheating temperature.
Definition bounce_solution.cpp:969
const double RelErr
RelErr relative error for numerical integration.
Definition bounce_solution.h:294
TransitionTemperature which_transition_temp
Temperature at which to calculate parameters.
Definition bounce_solution.h:108
double GetEpsTurb() const
Get epsturb.
Definition bounce_solution.cpp:500
tk::spline TruePhaseVSpline
True V spline to interpolate.
Definition bounce_solution.h:214
double CalculateRhoGamma(const double &T) const
Calculate .
Definition bounce_solution.cpp:1066
double CalcFalseVacFraction(const double &temp)
CalcFalseVacFraction calculates false vacuum fraction as function of temperature.
Definition bounce_solution.cpp:834
const double MarginOfCalcTempAtFalseVacFractionBeforeFailure
Additional margin of error in the while loop without admitting failure.
Definition bounce_solution.h:313
StatusGW status_bounce_sol
status of bounce solver
Definition bounce_solution.h:324
Definition spline.h:40
This classes calculates the Bounce action of the potential with a set temperature.
Definition CalculateEtaInterface.h:24
StatusGW
Possible results for the GW and bounce_sol class.
Definition minimum_tracer.h:172
TransitionTemperature
Possible transitions temperatures.
Definition minimum_tracer.h:191
struct resultErrorPair Nderive_BounceRatio(BounceSolution &obj)
Nderive_BounceRatio Numerical derivative for the inverse time scale calculation.
Definition bounce_solution.cpp:1223
StatusTemperature
Possible status for the approximated nucleation, exact nucleation, percolation and completion tempera...
Definition minimum_tracer.h:152
struct resultErrorPair Nintegrate_Inner(BounceSolution &obj, const double &Tprime)
Nintegrate_Inner Numerical integration of inner integral over inverse Hubble rate for the percolation...
Definition bounce_solution.cpp:1007
struct resultErrorPair Nintegrate_Outer(BounceSolution &obj)
Nintegrate_Outer Numerical integration of outer integral for the percolation temperature calculation.
Definition bounce_solution.cpp:1037
CoexPhases struct to save pair of coexisting phases (false and true phase)
Definition minimum_tracer.h:748
struct to store minimum and temperature
Definition minimum_tracer.h:266
Phase object.
Definition minimum_tracer.h:619
Definition bounce_solution.h:31