softsusy
is hosted by
Hepforge
,
IPPP Durham
SOFTSUSY
4.1
src
lowe.h
Go to the documentation of this file.
1
12
#ifndef LOWE_H
13
#define LOWE_H
14
15
#include <iostream>
16
#include <fstream>
17
#include <sstream>
18
#include <cstring>
19
#include "
utils.h
"
20
#include "
linalg.h
"
21
#include "
rge.h
"
22
24
namespace
softsusy
{
25
const
double
MUP
= 2.2e-3;
26
const
double
MDOWN
= 4.7e-3;
27
const
double
MSTRANGE
= 0.096;
28
const
double
MCHARM
= 1.27;
29
const
double
MBOTTOM
= 4.18;
30
const
double
MTOP
= 165.0;
32
const
double
MELECTRON
= 5.109989461e-4;
33
const
double
MMUON
= 1.0565835745e-1;
34
const
double
MTAU
= 1.77686;
35
const
double
ALPHASMZ
= 0.1181;
36
const
double
ALPHAMZ
= 1.0 / 127.950;
37
38
const
double
PMTOP
= 173.21;
39
const
double
PMBOTTOM
= 4.9;
41
const
double
THETA12CKM
= 0.229206;
42
const
double
THETA13CKM
= 0.003960;
43
const
double
THETA23CKM
= 0.042223;
44
46
typedef
enum
{mUp=1, mCharm, mTop, mDown, mStrange, mBottom, mElectron,
47
mMuon, mTau}
mass
;
49
typedef
enum
{ALPHA=1, ALPHAS}
leGauge
;
50
52
DoubleVector
gaugeDerivs
(
double
,
const
DoubleVector
&);
53
55
class
QedQcd
:
public
RGE
,
public
Approx
{
56
private
:
57
DoubleVector
a;
58
DoubleVector
mf;
59
double
mtPole, mbPole;
60
double
mbMb;
61
double
mtauPole;
62
DoubleVector
aMz;
63
64
public
:
65
QedQcd
();
66
QedQcd
(
const
QedQcd
&);
67
const
QedQcd
&
operator=
(
const
QedQcd
& m);
68
virtual
~
QedQcd
() {};
69
70
void
setPoleMt(
double
mt) { mtPole = mt; };
71
void
setPoleMb
(
double
mb) { mbPole = mb; };
72
void
setPoleMtau
(
double
mtau) { mtauPole = mtau; };
73
void
setMbMb
(
double
mb) { mbMb = mb; };
75
void
setMass
(
mass
mno,
double
m) { mf(mno) = m; };
77
void
setAlpha
(
leGauge
ai,
double
ap) { a(ai) = ap; };
78
void
setAlphaMz(
leGauge
ai,
double
ap) { aMz(ai) = ap; };
80
void
set
(
const
DoubleVector
&);
81
83
double
displayPoleMt
()
const
{
return
mtPole; };
85
double
displayPoleMtau
()
const
{
return
mtauPole; };
87
double
displayPoleMb
()
const
{
return
mbPole; };
89
const
DoubleVector
&
displayMass
()
const
{
return
mf; };
91
double
displayMass
(
mass
mno)
const
{
return
mf.
display
(mno); };
93
double
displayAlpha
(
leGauge
ai)
const
{
return
a.
display
(ai); };
94
double
displayAlphaMz(
leGauge
ai)
const
{
return
aMz.
display
(ai); };
96
const
DoubleVector
display
()
const
;
98
double
displayMbMb
()
const
{
return
mbMb; }
99
100
int
flavours(
double
)
const
;
101
102
double
qedBeta
()
const
;
103
double
qcdBeta
()
const
;
104
void
massBeta
(
DoubleVector
&)
const
;
106
DoubleVector
beta
()
const
;
107
109
void
runGauge
(
double
start,
double
end);
111
double
extractPoleMb
(
double
asMb);
113
double
extractRunningMb
(
double
asMb);
115
void
calcRunningMb
();
118
void
calcPoleMb
();
119
121
void
toMt
();
123
void
toMz
();
128
// alpha1 is in the GUT normalisation. sinth = sin^2 thetaW(Q) in MSbar
129
// scheme
130
DoubleVector
getGaugeMu
(
const
double
m2,
const
131
double
sinth)
const
;
132
};
133
135
ostream &
operator <<
(ostream &,
const
QedQcd
&);
137
istream &
operator >>
(istream &left,
QedQcd
&m);
138
141
double
getRunMt
(
double
poleMt,
double
asmt);
144
double
getAsmt
(
double
mtop,
double
alphasMz);
146
double
getRunMtFromMz
(
double
poleMt,
double
asMZ);
147
148
inline
QedQcd::QedQcd
(
const
QedQcd
&m)
149
:
RGE
(),
Approx
(m.displayApprox()), a(m.a), mf(m.mf), mtPole(m.mtPole),
150
mbPole(m.mbPole), mbMb(m.mbMb),
151
mtauPole(m.mtauPole), aMz(m.aMz) {
152
setPars
(11);
153
setMu
(m.
displayMu
());
154
}
155
157
void
massFermions
(
const
QedQcd
& r,
DoubleMatrix
& mDon,
158
DoubleMatrix
& mUpq,
DoubleMatrix
& mEle);
159
160
}
// namespace softsusy
161
162
#endif
163
DoubleMatrix
Matrix from 1..rows, 1..cols of double values.
Definition:
linalg.h:214
DoubleVector
DoubleVector is of variable length, and contains double precision.
Definition:
linalg.h:35
DoubleVector::display
double display(int i) const
display one element
Definition:
linalg.h:90
softsusy::Approx
Definition:
rge.h:23
softsusy::QedQcd
Quark and lepton masses and gauge couplings in QEDxQCD effective theory.
Definition:
lowe.h:55
softsusy::QedQcd::setMass
void setMass(mass mno, double m)
sets a running quark mass
Definition:
lowe.h:75
softsusy::QedQcd::beta
DoubleVector beta() const
Beta functions of both beta-functions and all MSbar masses.
Definition:
lowe.cpp:240
softsusy::QedQcd::toMt
void toMt()
Evolves object to running top mass.
Definition:
lowe.cpp:360
softsusy::QedQcd::calcRunningMb
void calcRunningMb()
calculates running bottom mass given alpha_s(Mb)^{MSbar} from pole m_b
Definition:
lowe.cpp:317
softsusy::QedQcd::setMbMb
void setMbMb(double mb)
set pole tau mass
Definition:
lowe.h:73
softsusy::QedQcd::displayAlpha
double displayAlpha(leGauge ai) const
Returns a single gauge structure constant.
Definition:
lowe.h:93
softsusy::QedQcd::displayPoleMt
double displayPoleMt() const
Display pole top mass.
Definition:
lowe.h:83
softsusy::QedQcd::qedBeta
double qedBeta() const
returns number of active flavours
Definition:
lowe.cpp:161
softsusy::QedQcd::toMz
void toMz()
Evolves object to MZ.
Definition:
lowe.cpp:382
softsusy::QedQcd::setAlpha
void setAlpha(leGauge ai, double ap)
sets QED or QCD structure constant
Definition:
lowe.h:77
softsusy::QedQcd::set
void set(const DoubleVector &)
For exporting beta functions to Runge-Kutta.
Definition:
lowe.cpp:49
softsusy::QedQcd::massBeta
void massBeta(DoubleVector &) const
Definition:
lowe.cpp:207
softsusy::QedQcd::runGauge
void runGauge(double start, double end)
Does not run the masses, just gauge couplings from start to end.
Definition:
lowe.cpp:251
softsusy::QedQcd::displayMbMb
double displayMbMb() const
Returns mb(mb) MSbar.
Definition:
lowe.h:98
softsusy::QedQcd::setPoleMb
void setPoleMb(double mb)
set pole top mass
Definition:
lowe.h:71
softsusy::QedQcd::qcdBeta
double qcdBeta() const
QCD beta function.
Definition:
lowe.cpp:178
softsusy::QedQcd::operator=
const QedQcd & operator=(const QedQcd &m)
Sets two objects equal.
Definition:
lowe.cpp:33
softsusy::QedQcd::calcPoleMb
void calcPoleMb()
Definition:
lowe.cpp:343
softsusy::QedQcd::displayPoleMb
double displayPoleMb() const
Returns bottom "pole" mass.
Definition:
lowe.h:87
softsusy::QedQcd::extractRunningMb
double extractRunningMb(double asMb)
Done at pole mb: extracts running mb(polemb)
Definition:
lowe.cpp:266
softsusy::QedQcd::displayPoleMtau
double displayPoleMtau() const
Display pole tau mass.
Definition:
lowe.h:85
softsusy::QedQcd::extractPoleMb
double extractPoleMb(double asMb)
calculates pole bottom mass given alpha_s(Mb)^{MSbar} from running b mass
Definition:
lowe.cpp:293
softsusy::QedQcd::displayMass
const DoubleVector & displayMass() const
Returns a vector of running fermion masses.
Definition:
lowe.h:89
softsusy::QedQcd::displayMass
double displayMass(mass mno) const
Returns a single running mass.
Definition:
lowe.h:91
softsusy::QedQcd::QedQcd
QedQcd()
Initialises with default values defined in lowe.h.
Definition:
lowe.cpp:17
softsusy::QedQcd::display
const DoubleVector display() const
Obgligatory: returns vector of all running parameters.
Definition:
lowe.cpp:56
softsusy::QedQcd::getGaugeMu
DoubleVector getGaugeMu(const double m2, const double sinth) const
Definition:
lowe.cpp:407
softsusy::QedQcd::setPoleMtau
void setPoleMtau(double mtau)
set pole bottom mass
Definition:
lowe.h:72
softsusy::RGE
Describes a set of parameters and RGEs in general.
Definition:
rge.h:49
softsusy::RGE::setPars
void setPars(int i)
Set number of parameters in RGE object.
Definition:
rge.h:64
softsusy::RGE::displayMu
double displayMu() const
Return renomalisation scale.
Definition:
rge.h:67
softsusy::RGE::setMu
void setMu(double e)
Sets renormalisation scale to e.
Definition:
rge.h:59
linalg.h
DoubleVector and DoubleMatrix classes of doubles and operations between them, complexified copies als...
softsusy
global variable declaration
Definition:
def.cpp:13
softsusy::MTAU
const double MTAU
default pole lepton mass from PDG
Definition:
lowe.h:34
softsusy::gaugeDerivs
DoubleVector gaugeDerivs(double x, const DoubleVector &y)
Returns beta functions of alpha, alpha_s only.
Definition:
lowe.cpp:440
softsusy::mass
mass
used to give order of quark masses stored
Definition:
lowe.h:46
softsusy::getAsmt
double getAsmt(double mtop, double alphasMz)
Definition:
lowe.cpp:464
softsusy::MSTRANGE
const double MSTRANGE
default running quark mass from PDG at 2 GeV
Definition:
lowe.h:27
softsusy::getRunMt
double getRunMt(double poleMt, double asmt)
Definition:
lowe.cpp:458
softsusy::MMUON
const double MMUON
default pole lepton mass from PDG
Definition:
lowe.h:33
softsusy::ALPHASMZ
const double ALPHASMZ
default running value from PDG
Definition:
lowe.h:35
softsusy::MELECTRON
const double MELECTRON
default pole lepton mass from PDG
Definition:
lowe.h:32
softsusy::ALPHAMZ
const double ALPHAMZ
default running alpha(MZ) from PDG
Definition:
lowe.h:36
softsusy::getRunMtFromMz
double getRunMtFromMz(double poleMt, double asMZ)
Given pole mass and alphaS(MZ), returns running top mass – one loop qcd.
Definition:
lowe.cpp:452
softsusy::THETA13CKM
const double THETA13CKM
From Vub in global CKM fit, PDG.
Definition:
lowe.h:42
softsusy::MDOWN
const double MDOWN
default running quark mass from PDG at 2 GeV
Definition:
lowe.h:26
softsusy::MUP
const double MUP
default running quark mass from PDG at 2 GeV
Definition:
lowe.h:25
softsusy::THETA23CKM
const double THETA23CKM
From Vcb/Vtb in global CKM fit, PDG.
Definition:
lowe.h:43
softsusy::PMBOTTOM
const double PMBOTTOM
Definition:
lowe.h:39
softsusy::MBOTTOM
const double MBOTTOM
default running quark mass from PDG
Definition:
lowe.h:29
softsusy::MCHARM
const double MCHARM
default running quark mass from PDG at 2 GeV
Definition:
lowe.h:28
softsusy::operator<<
ostream & operator<<(ostream &left, const FlavourMssmSoftsusy &m)
Formatted output.
Definition:
flavoursoft.cpp:217
softsusy::PMTOP
const double PMTOP
PDG 2014.
Definition:
lowe.h:38
softsusy::massFermions
void massFermions(const QedQcd &r, DoubleMatrix &mDon, DoubleMatrix &mUpq, DoubleMatrix &mEle)
Returns diagonal fermion mass matrices given input object r.
Definition:
lowe.cpp:470
softsusy::leGauge
leGauge
order of gauge couplings stored in QedQcd
Definition:
lowe.h:49
softsusy::operator>>
istream & operator>>(istream &left, flavourPhysical &s)
Formatted input of physical parameters.
Definition:
flavoursoft.cpp:1220
softsusy::MTOP
const double MTOP
Definition:
lowe.h:30
softsusy::THETA12CKM
const double THETA12CKM
default central values of CKM matrix elements from PDG 2006 in radians
Definition:
lowe.h:41
rge.h
RGE objects consisting of energy scale and parameters and loops (order in perturbation theory) and th...
utils.h
A few handy bits and pieces - little mathematical functions and the like.
Generated by
1.9.1