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
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/minimum_tracer/minimum_tracer.h> // MinimumTracer
18#include <BSMPT/utility/spline/spline.h>
19#include <Eigen/Dense>
20#include <algorithm> // std::max
21#include <gsl/gsl_deriv.h> // numerical derivative
22#include <gsl/gsl_integration.h> // numerical integration
23
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
56 double epsturb = 0.1;
57
61 double vwall = 0.95;
62
68
72 bool nucleation_temp_set = false;
73
78
82 bool completion_temp_set = false;
83
87 double Tc = -1;
88
92 double Tm = -1;
93
97 double Tnucl = -1;
98
102 double Tnucl_approx = -1;
103
107 double Tperc = -1;
108
112 double Tcompl = -1;
113
118
122 double alpha = -1;
123
128 double betaH = -1;
129
133 double gstar;
134
140
147
152 std::vector<BounceActionInt> SolutionList;
153
158 std::vector<Eigen::MatrixXd> GroupElements;
159
164 Eigen::MatrixXd OptimalDiscreteSymmetry;
165
172 std::vector<double>
173 TransformIntoOptimalDiscreteSymmetry(const std::vector<double> &vev);
174
181
187 double TunnelingRate(const double &Temp);
192 double HubbleRate(const double &Temp);
193
198 friend double inner_integrand(double var, void *params);
199
204 friend double outer_integrand(double var, void *params);
205
209 double GetBounceSol(const double &Temp) const;
210
215 friend double action_ratio(double var, void *params);
216
217public:
221 const double AbsErr = 0;
222
226 const double RelErr = 1e-6;
227
234
240
246
251
256 StatusGW status_bounce_sol = StatusGW::NotSet;
257
263 BSMPT::StatusTemperature::NotSet;
264
269 BSMPT::StatusTemperature status_nucl = BSMPT::StatusTemperature::NotSet;
270
275 BSMPT::StatusTemperature status_perc = BSMPT::StatusTemperature::NotSet;
276
281 BSMPT::StatusTemperature status_compl = BSMPT::StatusTemperature::NotSet;
282
291 double UserDefined_vwall = 0.95;
292
298
303 BounceSolution(const std::shared_ptr<Class_Potential_Origin> &pointer_in);
304
323 BounceSolution(const std::shared_ptr<Class_Potential_Origin> &pointer_in,
324 const std::shared_ptr<MinimumTracer> &MinTracer_in,
325 const CoexPhases &phase_pair_in,
326 const double &UserDefined_vwall_in,
327 const double &UserDefined_epsturb_in,
328 const int &MaxPathIntegrations_in,
329 const size_t &NumberOfInitialScanTemperatures_in);
330
349 BounceSolution(const std::shared_ptr<Class_Potential_Origin> &pointer_in,
350 const std::shared_ptr<MinimumTracer> &MinTracer_in,
351 const CoexPhases &phase_pair_in,
352 const double &UserDefined_vwall_in,
353 const double &UserDefined_epsturb_in,
354 const int &MaxPathIntegrations_in,
355 const size_t &NumberOfInitialScanTemperatures_in,
356 std::vector<Eigen::MatrixXd> GroupElements_in);
357
362 void GWInitialScan();
363
372 void CalculateActionAt(double T, bool smart = true);
373
380 void GWSecondaryScan();
381
387
393
398 void SetBounceSol();
399
404 double GetWallVelocity() const;
405
410 double GetEpsTurb() const;
411
415 void SetGstar(const double &gstar_in);
416
420 double GetGstar() const;
421
425 void SetCriticalTemp(const double &T_in);
426
430 double GetCriticalTemp() const;
431
435 void SetStoredTemp(const double &T_in);
439 double GetStoredTemp() const;
440
444 double GetNucleationTemp() const;
445
450 double GetNucleationTempApprox() const;
451
455 double GetPercolationTemp() const;
456
460 double GetCompletionTemp() const;
461
465 double CalcTransitionTemp(const int &which_transition_temp);
466
470 double GetPTStrength() const;
471
476 void CalcGstarPureRad();
477
482
487
494 double CalcTempAtFalseVacFraction(const double &false_vac_frac);
495
502 double CalcFalseVacFraction(const double &temp);
503
510 void CalculatePercolationTemp(const double &false_vac_frac = 0.71);
511
518 void CalculateCompletionTemp(const double &false_vac_frac = 0.01);
519
523 void CalculatePTStrength();
524
530 void CalculateWallVelocity(const Minimum &false_min, const Minimum &true_min);
531
536
540 double GetInvTimeScale();
541};
542
551 const double &Tprime);
552
560
568
569} // namespace BSMPT
BounceSolution class that handles the calculation of the bounce solution as well as the calculation o...
Definition bounce_solution.h:41
void CalculateWallVelocity(const Minimum &false_min, const Minimum &true_min)
Calculate wall velocity.
Definition bounce_solution.cpp:992
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:1060
double betaH
Inverse time scale .
Definition bounce_solution.h:128
size_t NumberOfInitialScanTemperatures
number of temperature steps in the initial scan of the bounce solver
Definition bounce_solution.h:67
double Tnucl_approx
approximate nucleation temperature
Definition bounce_solution.h:102
void CalculateOptimalDiscreteSymmetry()
Calculates which is the optimal symmetry from the group of symmetries.
Definition bounce_solution.cpp:74
double GetBounceSol(const double &Temp) const
Calculate euclidian action at temperature T.
Definition bounce_solution.cpp:1055
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:291
double GetPTStrength() const
GetPTStrength Get PT strength alpha.
Definition bounce_solution.cpp:559
std::vector< BounceActionInt > SolutionList
Set of BounceActionInt objects with valid solutions.
Definition bounce_solution.h:152
const double RelativeErrorInCalcTempAtFalseVacFraction
Maximum relative error on the fraction of vacuum tunneled to be accepted.
Definition bounce_solution.h:239
friend double inner_integrand(double var, void *params)
inner_integrand friend to define inner integrand of percolation temperature integral
Definition bounce_solution.cpp:858
void SetCriticalTemp(const double &T_in)
SetCriticalTemp Set critical temperature.
Definition bounce_solution.cpp:480
double CalcTransitionTemp(const int &which_transition_temp)
CalcTransitionTemp Get transition temperature from int.
Definition bounce_solution.cpp:520
std::vector< Eigen::MatrixXd > GroupElements
List of group elements allowed by the potential.
Definition bounce_solution.h:158
BSMPT::StatusTemperature status_compl
status of completion temperature calculation
Definition bounce_solution.h:281
bool completion_temp_set
set to true if completion temperature is set
Definition bounce_solution.h:82
double TunnelingRate(const double &Temp)
Storage of the tunneling rate per volume of the transition from false to true vacuum.
Definition bounce_solution.cpp:564
void GWInitialScan()
Initially we have no idea where the transition can occur, therefore we scan the complete temperature ...
Definition bounce_solution.cpp:118
void CalcGstarPureRad()
CalcGstarPureRad Calculate the number of effective degrees of freedom assuming a purely radiative uni...
Definition bounce_solution.cpp:590
double GetInvTimeScale()
Get inverse time scale of phase transition.
Definition bounce_solution.cpp:1095
double CalcTempAtFalseVacFraction(const double &false_vac_frac)
CalcTempAtFalseVacFraction calculates the temperature at which the false vacuum fraction drops below ...
Definition bounce_solution.cpp:724
bool percolation_temp_set
set to true if percolation temperature is set
Definition bounce_solution.h:77
void GWScanTowardsLowAction()
Do linear extrapolations to calculate action at lower temperatures.
Definition bounce_solution.cpp:377
const double RelativeTemperatureInCalcTempAtFalseVacFraction
Maximum relative difference in temperature on the fraction of false vacuum to be accepted.
Definition bounce_solution.h:233
double HubbleRate(const double &Temp)
Storage of the temperature-dependent Hubble rate.
Definition bounce_solution.cpp:583
int MaxPathIntegrations
Number of integration of the bounce.
Definition bounce_solution.h:297
double GetPercolationTemp() const
GetPercolationTemp Get percolation temperature.
Definition bounce_solution.cpp:510
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:187
void GWSecondaryScan()
If solution were found by the GWInitialScan() then we scan temperature range in the vicinity such tha...
Definition bounce_solution.cpp:264
BSMPT::StatusTemperature status_perc
status of percolation temperature calculation
Definition bounce_solution.h:275
double Tperc
percolation temperature
Definition bounce_solution.h:107
void SetGstar(const double &gstar_in)
SetGstar Set gstar.
Definition bounce_solution.cpp:470
friend double outer_integrand(double var, void *params)
outer_integrand friend to define outer integrand of percolation temperature integral
Definition bounce_solution.cpp:865
Eigen::MatrixXd OptimalDiscreteSymmetry
Store symmetry that produces the best tunneling rate.
Definition bounce_solution.h:164
std::shared_ptr< MinimumTracer > MinTracer
MinTracer object.
Definition bounce_solution.h:51
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:830
CoexPhases phase_pair
pair of coexisiting phases
Definition bounce_solution.h:250
double GetCriticalTemp() const
GetCriticalTemp Get critical temperature.
Definition bounce_solution.cpp:485
double GetGstar() const
GetGstar Get gstar.
Definition bounce_solution.cpp:475
double GetNucleationTempApprox() const
GetNucleationTempApprox Get nucleation temperature via approximate method.
Definition bounce_solution.cpp:505
void CalculateNucleationTemp()
Calculation of nucleation temperature.
Definition bounce_solution.cpp:601
double Tm
lowest temperature when a transition can occur
Definition bounce_solution.h:92
void SetBounceSol()
Set the Bounce Sol object.
Definition bounce_solution.cpp:415
const double AbsErr
AbsErr absolute error for numerical integration.
Definition bounce_solution.h:221
double Tc
critical temperature/highest temperature when transition can occur
Definition bounce_solution.h:87
int indexTrueCandidatePhase
index of the true vacuum phase candidate in the coex list
Definition bounce_solution.h:139
double store_Temp
stored temperature
Definition bounce_solution.h:117
void CalculatePTStrength()
Calculate phase transition strength alpha at percolation temperature.
Definition bounce_solution.cpp:934
std::vector< double > TransformIntoOptimalDiscreteSymmetry(const std::vector< double > &vev)
Transforms the vev to the optimal vev.
Definition bounce_solution.cpp:107
double GetNucleationTemp() const
GetNucleationTemp Get nucleation temperature via exact method.
Definition bounce_solution.cpp:500
double Tnucl
nucleation temperature
Definition bounce_solution.h:97
void GWScanTowardsHighAction()
Do linear extrapolations to calculate action at higher temperatures.
Definition bounce_solution.cpp:338
void SetStoredTemp(const double &T_in)
SetStoredTemp Set stored temperature.
Definition bounce_solution.cpp:490
double Tcompl
completion temperature
Definition bounce_solution.h:112
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:799
double vwall
wall velocity
Definition bounce_solution.h:61
void CalculateNucleationTempApprox()
Approximate calculation of nucleation temperature.
Definition bounce_solution.cpp:662
double alpha
PT strength.
Definition bounce_solution.h:122
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:495
BSMPT::StatusTemperature status_nucl
status of nucleation temperature calculation
Definition bounce_solution.h:269
double epsturb
epsilon of turbulence efficiency factor
Definition bounce_solution.h:56
double GetCompletionTemp() const
GetCompletionTemp Get percolation temperature.
Definition bounce_solution.cpp:515
bool nucleation_temp_set
set to true if nucleation temperature is set
Definition bounce_solution.h:72
void CalculateInvTimeScale()
Calculate inverse time scale of phase transition.
Definition bounce_solution.cpp:1083
BSMPT::StatusTemperature status_nucl_approx
status of approximate nucleation temperature calculation
Definition bounce_solution.h:262
double GetWallVelocity() const
Get the bubble wall velocity.
Definition bounce_solution.cpp:460
tk::spline S3ofT_spline
spline used to interpolate the action as a function of the temperature
Definition bounce_solution.h:146
double gstar
number of effective degrees of freedom
Definition bounce_solution.h:133
const double RelErr
RelErr relative error for numerical integration.
Definition bounce_solution.h:226
double GetEpsTurb() const
Get epsturb.
Definition bounce_solution.cpp:465
double CalcFalseVacFraction(const double &temp)
CalcFalseVacFraction calculates false vacuum fraction as function of temperature.
Definition bounce_solution.cpp:717
const double MarginOfCalcTempAtFalseVacFractionBeforeFailure
Additional margin of error in the while loop without admitting failure.
Definition bounce_solution.h:245
StatusGW status_bounce_sol
status of bounce solver
Definition bounce_solution.h:256
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
struct resultErrorPair Nderive_BounceRatio(BounceSolution &obj)
Nderive_BounceRatio Numerical derivative for the inverse time scale calculation.
Definition bounce_solution.cpp:1067
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:875
struct resultErrorPair Nintegrate_Outer(BounceSolution &obj)
Nintegrate_Outer Numerical integration of outer integral for the percolation temperature calculation.
Definition bounce_solution.cpp:905
CoexPhases struct to save pair of coexisting phases (false and true phase)
Definition minimum_tracer.h:735
struct to store minimum and temperature
Definition minimum_tracer.h:253
Definition bounce_solution.h:31