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
gw.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
15#include <BSMPT/bounce_solution/bounce_solution.h> // BounceSolution
18#include <gsl/gsl_math.h>
19
20namespace BSMPT
21{
22
27{
28 bool swON = true; // enable sound wave contribution by default
29 bool turbON = true; // enable turbulence contribution by default
30 double transitionTemp = false; // transition temperature
31 double PTStrength = false; // strength of EW phase transition
32 double InvTimeScale = false; // inverse time scale beta/H
33 double vb = false; // bubble wall velocity
34
35 double kappa_sw = false; // efficiency factor for sound waves
36 double K_sw = false; // kinetic energy fraction in fluid of total bubble
37 // energy for sound waves
38 double HR = false; // time scale times max. velocity for sound waves
39 double kappa_turb = false; // efficiency factor for turbulence
40 double K_turb = false; // kinetic energy fraction in fluid of total bubble
41 // energy for turbulence
42 double gstar = false; // number of eff. d.o.f.
43 double Hstar = false; // Hubble rate at percolation temperature
44
45 double fPeakSoundWave = false; // peak frequency for sound wave
46 double h2OmegaPeakSoundWave = false; // peak amplitude for sound wave
47 double fPeakTurbulence = false; // peak frequency for turbulence
48 double h2OmegaPeakTurbulence = false; // peak amplitude for turbulence
49
50 StatusGW status = StatusGW::NotSet; // gw calculation status
51};
52
54{
55private:
56public:
59 const int &which_transition_temp = 3);
61
69 double CalcEpsTurb(double epsturb_in);
70
74 const double AbsErr = 0;
75
79 const double RelErr = 1e-6;
80
85 double h = 0.674;
86
91
96
101
106
116 double CalcGWAmplitude(double f, bool swON, bool turbON);
117
125 double GetSNR(const double fmin, const double fmax, const double T = 3);
126
130 friend double snr_integrand(double freq, void *params);
131};
132
139double SIfunc(const double f);
140
147double Rfunc(const double f);
148
155double powspec_density(const double f);
156
163double h2OmSens(const double f);
164
170struct resultErrorPair
171Nintegrate_SNR(GravitationalWave &obj, const double fmin, const double fmax);
172
181double
182Getkappa_sw(const double &alpha, const double &vwall, const double &Csound);
183
192double GetK_sw(const double &alpha, const double &vwall, const double &Csound);
193
202double
203GetHR(const double &invTimeScale, const double &vwall, const double &Csound);
204
212double GetK_turb(const double &alpha, const double &kappa);
213
221bool IsFluidTurnoverApproxOne(const double &HR, const double &K);
222
223} // namespace BSMPT
BounceSolution class that handles the calculation of the bounce solution as well as the calculation o...
Definition bounce_solution.h:41
Definition gw.h:54
double h
reduced Hubble constant
Definition gw.h:85
const double RelErr
RelErr relative error for numerical integration.
Definition gw.h:79
void CalcPeakFrequencySoundWave()
Calculate peak frequency of GW signal for sound waves.
Definition gw.cpp:67
const double AbsErr
AbsErr absolute error for numerical integration.
Definition gw.h:74
void CalcPeakFrequencyTurbulence()
Calculate peak frequency of GW signal from turbulence.
Definition gw.cpp:94
double CalcEpsTurb(double epsturb_in)
CalcEpsTurb calculate epsilon for turbulence contribution.
Definition gw.cpp:54
double GetSNR(const double fmin, const double fmax, const double T=3)
GetSNR.
Definition gw.cpp:149
void CalcPeakAmplitudeSoundWave()
Calculate peak amplitude of GW signal for sound waves.
Definition gw.cpp:75
friend double snr_integrand(double freq, void *params)
snr_integrand friend to define inner integrand of SNR integral
Definition gw.cpp:183
void CalcPeakAmplitudeTurbulence()
Calculate peak amplitude of GW signal from turbulence.
Definition gw.cpp:102
double CalcGWAmplitude(double f, bool swON, bool turbON)
Amplitude of GW signal as a function of.
Definition gw.cpp:109
This classes calculates the Bounce action of the potential with a set temperature.
Definition CalculateEtaInterface.h:24
double SIfunc(const double f)
SIfunc.
Definition gw.cpp:157
double powspec_density(const double f)
powspec_density
Definition gw.cpp:169
double h2OmSens(const double f)
return the value of LISA mission nominal sensitivity
Definition gw.cpp:176
struct resultErrorPair Nintegrate_SNR(GravitationalWave &obj, const double fmin, const double fmax)
Nintegrate_SNR Numerical integration of SNR integral.
Definition gw.cpp:194
double Getkappa_sw(const double &alpha, const double &vwall, const double &Csound)
Get efficiency factor kappa_sw for sound waves.
Definition gw.cpp:224
StatusGW
Possible results for the GW and bounce_sol class.
Definition minimum_tracer.h:172
double Rfunc(const double f)
Rfunc.
Definition gw.cpp:163
double GetK_turb(const double &alpha, const double &kappa)
Get K for turbulence.
Definition gw.cpp:269
double GetHR(const double &invTimeScale, const double &vwall, const double &Csound)
Get HR for sound waves.
Definition gw.cpp:263
bool IsFluidTurnoverApproxOne(const double &HR, const double &K)
Determine fluid turnover time regime.
Definition gw.cpp:274
double GetK_sw(const double &alpha, const double &vwall, const double &Csound)
Get K for sound waves.
Definition gw.cpp:256
struct to store all calculated GW data
Definition gw.h:27
Definition bounce_solution.h:31