BSMPT 3.2.1
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 TransitionTemperature which_transition_temp = TransitionTemperature::NotSet;
94
98 double Tc = -1;
99
103 double Tm = -1;
104
108 double Tnucl = -1;
109
113 double Tnucl_approx = -1;
114
118 double Tperc = -1;
119
123 double Tcompl = -1;
124
128 double Treh = -1;
129
133 double Tstar = -1;
134
139
143 double alpha = -1;
144
149 double betaH = -1;
150
157 double Rstar = -1;
158
162 double gstar;
163
169
176
182
188
194
200
205 std::vector<BounceActionInt> SolutionList;
206
211 std::vector<Eigen::MatrixXd> GroupElements;
212
217 Eigen::MatrixXd OptimalDiscreteSymmetry;
218
225 std::vector<double>
226 TransformIntoOptimalDiscreteSymmetry(const std::vector<double> &vev);
227
234
240 double TunnelingRate(const double &Temp);
245 double HubbleRate(const double &Temp);
246
251 friend double inner_integrand(double var, void *params);
252
257 friend double outer_integrand(double var, void *params);
258
262 double GetBounceSol(const double &Temp) const;
263
268 friend double action_ratio(double var, void *params);
269
270public:
274 const double AbsErr = 0;
275
279 const double RelErr = 1e-6;
280
287
293
298
303 StatusGW status_bounce_sol = StatusGW::NotSet;
304
310 const TransitionTemperature &which_transition_temp_in);
311
317 BSMPT::StatusTemperature::NotSet;
318
323 BSMPT::StatusTemperature status_nucl = BSMPT::StatusTemperature::NotSet;
324
329 BSMPT::StatusTemperature status_perc = BSMPT::StatusTemperature::NotSet;
330
335 BSMPT::StatusTemperature status_compl = BSMPT::StatusTemperature::NotSet;
336
345 double UserDefined_vwall = 0.95;
346
352
357 BounceSolution(const std::shared_ptr<Class_Potential_Origin> &pointer_in);
358
379 BounceSolution(const std::shared_ptr<Class_Potential_Origin> &pointer_in,
380 const std::shared_ptr<MinimumTracer> &MinTracer_in,
381 const CoexPhases &phase_pair_in,
382 const double &UserDefined_vwall_in,
383 const double &UserDefined_epsturb_in,
384 const int &MaxPathIntegrations_in,
385 const size_t &NumberOfInitialScanTemperatures_in,
386 const int &UserDefined_PNLO_scaling_in = 1);
387
408 BounceSolution(const std::shared_ptr<Class_Potential_Origin> &pointer_in,
409 const std::shared_ptr<MinimumTracer> &MinTracer_in,
410 const CoexPhases &phase_pair_in,
411 const double &UserDefined_vwall_in,
412 const double &UserDefined_epsturb_in,
413 const int &MaxPathIntegrations_in,
414 const size_t &NumberOfInitialScanTemperatures_in,
415 const std::vector<Eigen::MatrixXd> &GroupElements_in,
416 const int &UserDefined_PNLO_scaling_in = 1);
417
422 void GWInitialScan();
423
432 void CalculateActionAt(double T, bool smart = true);
433
440 void GWSecondaryScan();
441
447
453
458 void SetBounceSol();
459
464 double GetWallVelocity() const;
465
470 double GetChapmanJougetVelocity() const;
471
475 double GetSoundSpeedFalse() const;
476
480 double GetSoundSpeedTrue() const;
481
486 double GetEpsTurb() const;
487
491 void SetGstar(const double &gstar_in);
492
498
505 static void ConstructSplineVofT(Phase &phase, tk::spline &spline);
506
512 void InitializedVSpline();
513
520 double GetGstar(const double &T) const;
521
527 double GetGstar();
528
532 void SetCriticalTemp(const double &T_in);
533
537 double GetCriticalTemp() const;
538
542 void SetStoredTemp(const double &T_in);
546 double GetStoredTemp() const;
547
551 double GetNucleationTemp() const;
552
557 double GetNucleationTempApprox() const;
558
562 double GetPercolationTemp() const;
563
567 double GetCompletionTemp() const;
568
572 double GetTransitionTemp() const;
573
577 double GetReheatingTemp() const;
578
582 void CalcTransitionTemp();
583
591 double CalculateRhoGamma(const double &T) const;
592
596 double GetPTStrength() const;
597
602 double CalcGstarPureRad();
603
608
613
623 double FalseVacFractionExponent_I(const double &T);
624
631 double CalcTempAtFalseVacFraction(const double &false_vac_frac);
632
639 double CalcFalseVacFraction(const double &temp);
640
647 void CalculatePercolationTemp(const double &false_vac_frac = 0.71);
648
655 void CalculateCompletionTemp(const double &false_vac_frac = 0.01);
656
661
665 void CalculatePTStrength();
666
672
678 void CalculateWallVelocity(const Minimum &false_min, const Minimum &true_min);
679
686 double CalculateSoundSpeed(Phase &phase);
687
694
699
703 double GetInvTimeScale();
704
711 void CalculateRstar();
712
717 double GetRstar();
718};
719
728 const double &Tprime);
729
737
745
746} // 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:1122
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:1214
double betaH
Inverse time scale .
Definition bounce_solution.h:149
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:113
void CalculateOptimalDiscreteSymmetry()
Calculates which is the optimal symmetry from the group of symmetries.
Definition bounce_solution.cpp:78
double GetBounceSol(const double &Temp) const
Calculate euclidian action at temperature T.
Definition bounce_solution.cpp:1209
void SetAndCalculateGWParameters(const TransitionTemperature &which_transition_temp_in)
Set the Transition Temp object.
Definition bounce_solution.cpp:111
void CalcChapmanJougetVelocity()
Derive the Chapman-Jouget velocity from PT strength and false phase sound velocity using Eq....
Definition bounce_solution.cpp:1114
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:345
double GetPTStrength() const
GetPTStrength Get PT strength alpha.
Definition bounce_solution.cpp:667
std::vector< BounceActionInt > SolutionList
Set of BounceActionInt objects with valid solutions.
Definition bounce_solution.h:205
const double RelativeErrorInCalcTempAtFalseVacFraction
Maximum relative error on the fraction of vacuum tunneled to be accepted.
Definition bounce_solution.h:292
void InitializedVSpline()
Initialize two splines for the potential across the tunneling profile. Used to improve the Hubble rat...
Definition bounce_solution.cpp:536
friend double inner_integrand(double var, void *params)
inner_integrand friend to define inner integrand of percolation temperature integral
Definition bounce_solution.cpp:978
void SetCriticalTemp(const double &T_in)
SetCriticalTemp Set critical temperature.
Definition bounce_solution.cpp:566
tk::spline GstarProfileHighT
Gstar spline, T > T_QCD (214.0 MeV)
Definition bounce_solution.h:187
double GetRstar()
Returns .
Definition bounce_solution.cpp:1249
std::vector< Eigen::MatrixXd > GroupElements
List of group elements allowed by the potential.
Definition bounce_solution.h:211
BSMPT::StatusTemperature status_compl
status of completion temperature calculation
Definition bounce_solution.h:335
tk::spline GstarProfileLowT
Gstar spline, T < T_QCD (214.0 MeV)
Definition bounce_solution.h:181
void CalcTransitionTemp()
CalcTransitionTemp Get transition temperature from int.
Definition bounce_solution.cpp:616
double TunnelingRate(const double &Temp)
Storage of the tunneling rate per volume of the transition from false to true vacuum.
Definition bounce_solution.cpp:672
void GWInitialScan()
Initially we have no idea where the transition can occur, therefore we scan the complete temperature ...
Definition bounce_solution.cpp:136
double GetInvTimeScale()
Get inverse time scale of phase transition.
Definition bounce_solution.cpp:1243
double CalcTempAtFalseVacFraction(const double &false_vac_frac)
CalcTempAtFalseVacFraction calculates the temperature at which the false vacuum fraction drops below ...
Definition bounce_solution.cpp:839
void GWScanTowardsLowAction()
Do linear extrapolations to calculate action at lower temperatures.
Definition bounce_solution.cpp:395
double GetChapmanJougetVelocity() const
Get Chapman-Jouget velocity.
Definition bounce_solution.cpp:486
double Tstar
transition temperature
Definition bounce_solution.h:133
double GetSoundSpeedFalse() const
Get the sound speed in the false phase.
Definition bounce_solution.cpp:491
double Csound_true
sound speed in true phase
Definition bounce_solution.h:82
void CalculateRstar()
Definition bounce_solution.cpp:1255
const double RelativeTemperatureInCalcTempAtFalseVacFraction
Maximum relative difference in temperature on the fraction of false vacuum to be accepted.
Definition bounce_solution.h:286
double HubbleRate(const double &Temp)
Storage of the temperature-dependent Hubble rate.
Definition bounce_solution.cpp:691
double CalcGstarPureRad()
CalcGstarPureRad Calculate the number of effective degrees of freedom assuming a purely radiative uni...
Definition bounce_solution.cpp:701
int MaxPathIntegrations
Number of integration of the bounce.
Definition bounce_solution.h:351
double Rstar
Definition bounce_solution.h:157
double GetPercolationTemp() const
GetPercolationTemp Get percolation temperature.
Definition bounce_solution.cpp:596
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:205
void GWSecondaryScan()
If solution were found by the GWInitialScan() then we scan temperature range in the vicinity such tha...
Definition bounce_solution.cpp:282
BSMPT::StatusTemperature status_perc
status of percolation temperature calculation
Definition bounce_solution.h:329
double Tperc
percolation temperature
Definition bounce_solution.h:118
double CalculateSoundSpeed(Phase &phase)
Calculate sound speeds at Tstar in phase.
Definition bounce_solution.cpp:1165
void SetGstar(const double &gstar_in)
SetGstar Set gstar.
Definition bounce_solution.cpp:506
friend double outer_integrand(double var, void *params)
outer_integrand friend to define outer integrand of percolation temperature integral
Definition bounce_solution.cpp:985
Eigen::MatrixXd OptimalDiscreteSymmetry
Store symmetry that produces the best tunneling rate.
Definition bounce_solution.h:217
tk::spline FalsePhaseVSpline
False V spline to interpolate.
Definition bounce_solution.h:193
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:926
CoexPhases phase_pair
pair of coexisiting phases
Definition bounce_solution.h:297
double GetCriticalTemp() const
GetCriticalTemp Get critical temperature.
Definition bounce_solution.cpp:571
double GetNucleationTempApprox() const
GetNucleationTempApprox Get nucleation temperature via approximate method.
Definition bounce_solution.cpp:591
void CalculateNucleationTemp()
Calculation of nucleation temperature.
Definition bounce_solution.cpp:712
void InitializeGstarProfile()
Generate the spline used to interpolate the gstar SM profile.
Definition bounce_solution.cpp:511
double GetSoundSpeedTrue() const
Get the sound speed in the true phase.
Definition bounce_solution.cpp:496
double Tm
lowest temperature when a transition can occur
Definition bounce_solution.h:103
void SetBounceSol()
Set the Bounce Sol object.
Definition bounce_solution.cpp:433
const double AbsErr
AbsErr absolute error for numerical integration.
Definition bounce_solution.h:274
double Tc
critical temperature/highest temperature when transition can occur
Definition bounce_solution.h:98
int indexTrueCandidatePhase
index of the true vacuum phase candidate in the coex list
Definition bounce_solution.h:168
double store_Temp
stored temperature
Definition bounce_solution.h:138
void CalculatePTStrength()
Calculate phase transition strength alpha.
Definition bounce_solution.cpp:1059
std::vector< double > TransformIntoOptimalDiscreteSymmetry(const std::vector< double > &vev)
Transforms the vev to the optimal vev.
Definition bounce_solution.cpp:125
double GetNucleationTemp() const
GetNucleationTemp Get nucleation temperature via exact method.
Definition bounce_solution.cpp:586
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:522
double Tnucl
nucleation temperature
Definition bounce_solution.h:108
double GetTransitionTemp() const
GetTransitionTemp Get transition temperature.
Definition bounce_solution.cpp:606
void CalculateSoundSpeeds()
Calculate sound speeds at Tstar in false and true phase.
Definition bounce_solution.cpp:1201
double GetReheatingTemp() const
GetReheatingTemp Get reheating temperature.
Definition bounce_solution.cpp:611
double Treh
reheating temperature
Definition bounce_solution.h:128
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:356
void SetStoredTemp(const double &T_in)
SetStoredTemp Set stored temperature.
Definition bounce_solution.cpp:576
double GetGstar()
Get Gstar for radiation-dominated epoch.
Definition bounce_solution.cpp:561
double Tcompl
completion temperature
Definition bounce_solution.h:123
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:895
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:143
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:581
BSMPT::StatusTemperature status_nucl
status of nucleation temperature calculation
Definition bounce_solution.h:323
double epsturb
epsilon of turbulence efficiency factor
Definition bounce_solution.h:62
double GetCompletionTemp() const
GetCompletionTemp Get percolation temperature.
Definition bounce_solution.cpp:601
void CalculateInvTimeScale()
Calculate inverse time scale of phase transition.
Definition bounce_solution.cpp:1237
BSMPT::StatusTemperature status_nucl_approx
status of approximate nucleation temperature calculation
Definition bounce_solution.h:316
double GetWallVelocity() const
Get the bubble wall velocity.
Definition bounce_solution.cpp:481
tk::spline S3ofT_spline
spline used to interpolate the action as a function of the temperature
Definition bounce_solution.h:175
double gstar
number of effective degrees of freedom
Definition bounce_solution.h:162
void CalculateReheatingTemp()
CalculateReheatingTemp calculation of the reheating temperature.
Definition bounce_solution.cpp:957
const double RelErr
RelErr relative error for numerical integration.
Definition bounce_solution.h:279
TransitionTemperature which_transition_temp
Temperature at which to calculate parameters.
Definition bounce_solution.h:93
double GetEpsTurb() const
Get epsturb.
Definition bounce_solution.cpp:501
tk::spline TruePhaseVSpline
True V spline to interpolate.
Definition bounce_solution.h:199
double CalculateRhoGamma(const double &T) const
Calculate .
Definition bounce_solution.cpp:1054
double CalcFalseVacFraction(const double &temp)
CalcFalseVacFraction calculates false vacuum fraction as function of temperature.
Definition bounce_solution.cpp:834
StatusGW status_bounce_sol
status of bounce solver
Definition bounce_solution.h:303
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:162
TransitionTemperature
Possible transitions temperatures.
Definition minimum_tracer.h:181
struct resultErrorPair Nderive_BounceRatio(BounceSolution &obj)
Nderive_BounceRatio Numerical derivative for the inverse time scale calculation.
Definition bounce_solution.cpp:1221
StatusTemperature
Possible status for the approximated nucleation, exact nucleation, percolation and completion tempera...
Definition minimum_tracer.h:142
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:995
struct resultErrorPair Nintegrate_Outer(BounceSolution &obj)
Nintegrate_Outer Numerical integration of outer integral for the percolation temperature calculation.
Definition bounce_solution.cpp:1025
CoexPhases struct to save pair of coexisting phases (false and true phase)
Definition minimum_tracer.h:745
struct to store minimum and temperature
Definition minimum_tracer.h:278
Phase object.
Definition minimum_tracer.h:616
Definition bounce_solution.h:31