flexiblesusy is hosted by Hepforge, IPPP Durham
FlexibleSUSY
lowe.h
Go to the documentation of this file.
1// ====================================================================
2// This file is part of FlexibleSUSY.
3//
4// FlexibleSUSY is free software: you can redistribute it and/or modify
5// it under the terms of the GNU General Public License as published
6// by the Free Software Foundation, either version 3 of the License,
7// or (at your option) any later version.
8//
9// FlexibleSUSY is distributed in the hope that it will be useful, but
10// WITHOUT ANY WARRANTY; without even the implied warranty of
11// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
12// General Public License for more details.
13//
14// You should have received a copy of the GNU General Public License
15// along with FlexibleSUSY. If not, see
16// <http://www.gnu.org/licenses/>.
17// ====================================================================
18
28#ifndef LOWE_H
29#define LOWE_H
30
31#include "betafunction.hpp"
32#include "ckm.hpp"
33#include "pmns.hpp"
34#include <array>
35#include <iosfwd>
36#include <Eigen/Core>
37
38namespace softsusy {
39
45
58};
59
60extern const std::array<std::string, NUMBER_OF_LOW_ENERGY_INPUT_PARAMETERS> QedQcd_input_parmeter_names;
61
64{
65public:
66 using Input_t = Eigen::Array<double,NUMBER_OF_LOW_ENERGY_INPUT_PARAMETERS,1>;
67
68private:
69 Eigen::Array<double,2,1> a{Eigen::Array<double,2,1>::Zero()};
70 Eigen::Array<double,9,1> mf{Eigen::Array<double,9,1>::Zero()};
71 Input_t input{Input_t::Zero()};
72 double mbPole;
73
74 double qedBeta() const;
75 double qcdBeta() const;
76 Eigen::Array<double,9,1> massBeta() const;
77 void runto_safe(double, double eps = -1.0);
78
79 int flavours(double) const;
80
82 double extractPoleMb(double asMb);
83
84public:
85 QedQcd();
86 QedQcd(const QedQcd&) = default;
87 QedQcd(QedQcd&&) = default;
88 QedQcd& operator=(const QedQcd&) = default;
89 QedQcd& operator=(QedQcd&&) = default;
90 virtual ~QedQcd() = default;
91
92 // Beta_function interface
93 virtual Eigen::ArrayXd get() const override;
94 virtual void set(const Eigen::ArrayXd&) override;
95 virtual Eigen::ArrayXd beta() const override;
96
97 void setPoleMt(double mt) { input(Mt_pole) = mt; }
98 void setPoleMb(double mb) { mbPole = mb; }
99 void setPoleMtau(double mtau) { input(Mtau_pole) = mtau; }
100 void setPoleMmuon(double m) { input(Mm_pole) = m; }
101 void setPoleMel(double m) { input(Me_pole) = m; }
102 void setMbMb(double mb) { input(mb_mb) = mb; }
103 void setMcMc(double mc) { input(mc_mc) = mc; }
104 void setMu2GeV(double mu) { input(mu_2GeV) = mu; }
105 void setMd2GeV(double md) { input(md_2GeV) = md; }
106 void setMs2GeV(double ms) { input(ms_2GeV) = ms; }
107 void setPoleMW(double mw) { input(MW_pole) = mw; }
108 void setPoleMZ(double mz) { input(MZ_pole) = mz; }
110 void setMasses(const Eigen::Array<double,9,1>& m) { mf = m; }
112 void setMass(EMass mno, double m) { mf(mno - 1) = m; }
114 void setNeutrinoPoleMass(int i, double m) { input(Mv1_pole + i - 1) = m; }
116 void setAlphas(const Eigen::Array<double,2,1>& o) { a = o; }
118 void setAlpha(EGauge ai, double ap) { a(ai - 1) = ap; }
128 void setFermiConstant(double gf) { input(GFermi) = gf; }
130 void set_input(const Input_t& i) { input = i; }
131
133 Input_t displayInput() const { return input; }
135 double displayPoleMt() const { return input(Mt_pole); }
137 double displayPoleMtau() const { return input(Mtau_pole); }
139 double displayPoleMmuon() const { return input(Mm_pole); }
141 double displayPoleMel() const { return input(Me_pole); }
143 double displayPoleMb() const { return mbPole; }
145 double displayPoleMW() const { return input(MW_pole); }
147 double displayPoleMZ() const { return input(MZ_pole); }
149 auto displayMass() const -> decltype(mf) { return mf; }
151 double displayMass(EMass mno) const { return mf(mno - 1); }
153 double displayUpQuarkRunningMass(int) const;
155 double displayDownQuarkRunningMass(int) const;
157 double displayNeutrinoPoleMass(int i) const { return input(Mv1_pole + i - 1); }
159 double displayLeptonPoleMass(int) const;
161 double displayLeptonRunningMass(int) const;
163 double displayAlpha(EGauge ai) const { return a(ai - 1); }
165 auto displayAlphas() const -> decltype(a) { return a; }
169 double displayAlphaSInput() const { return input(alpha_s_MSbar_at_MZ); }
171 double displayFermiConstant() const { return input(GFermi); }
173 Input_t display_input() const { return input; }
175 static std::array<std::string, NUMBER_OF_LOW_ENERGY_INPUT_PARAMETERS> display_input_parameter_names();
177 double displayMbMb() const { return input(mb_mb); }
179 double displayMcMc() const { return input(mc_mc); }
181 double displayMu2GeV() const { return input(mu_2GeV); }
183 double displayMd2GeV() const { return input(md_2GeV); }
185 double displayMs2GeV() const { return input(ms_2GeV); }
189 Eigen::Matrix<double,3,3> get_real_ckm() const { return displayCKM().get_real_ckm(); }
191 Eigen::Matrix<std::complex<double>,3,3> get_complex_ckm() const { return displayCKM().get_complex_ckm(); }
195 Eigen::Matrix<double,3,3> get_real_pmns() const { return displayPMNS().get_real_pmns(); }
197 Eigen::Matrix<std::complex<double>,3,3> get_complex_pmns() const { return displayPMNS().get_complex_pmns(); }
198
200 void toMz();
202 void to(double scale, double precision_goal = 1e-5, int max_iterations = 20);
204 Eigen::Array<double,3,1> guess_alpha_SM5(double scale) const;
205};
206
208std::ostream & operator<<(std::ostream &, const QedQcd &);
209
210bool operator ==(const QedQcd&, const QedQcd&);
211
212} // namespace softsusy
213
214#endif
contains class Beta_function
beta function interface
double scale
current renormalization scale
Quark and lepton masses and gauge couplings in QEDxQCD effective theory.
Definition: lowe.h:64
static std::array< std::string, NUMBER_OF_LOW_ENERGY_INPUT_PARAMETERS > display_input_parameter_names()
returns vector of all parameter names
Definition: lowe.cpp:514
Input_t displayInput() const
Displays input parameters.
Definition: lowe.h:133
void setPMNS(const flexiblesusy::PMNS_parameters &)
sets PMNS parameters (in the MS-bar scheme at MZ)
Definition: lowe.cpp:181
double displayAlpha(EGauge ai) const
Returns a single gauge structure constant.
Definition: lowe.h:163
void runto_safe(double, double eps=-1.0)
throws if non-perturbative error occurs
Definition: lowe.cpp:149
double displayAlphaSInput() const
Returns input value alpha_s(MZ)
Definition: lowe.h:169
void setPoleMmuon(double m)
set pole muon mass
Definition: lowe.h:100
virtual Eigen::ArrayXd beta() const override
Definition: lowe.cpp:137
void setPoleMZ(double mz)
Definition: lowe.h:108
void setMs2GeV(double ms)
set ms(2 GeV)
Definition: lowe.h:106
Eigen::Matrix< std::complex< double >, 3, 3 > get_complex_ckm() const
Returns complex CKM matrix.
Definition: lowe.h:191
void setMcMc(double mc)
set mc(mc)
Definition: lowe.h:103
virtual Eigen::ArrayXd get() const override
Definition: lowe.cpp:117
flexiblesusy::CKM_parameters displayCKM() const
returns CKM parameters
Definition: lowe.cpp:191
void setMbMb(double mb)
set mb(mb)
Definition: lowe.h:102
QedQcd & operator=(const QedQcd &)=default
void setPoleMel(double m)
set pole electron mass
Definition: lowe.h:101
Eigen::Matrix< double, 3, 3 > get_real_pmns() const
Returns real PMNS matrix.
Definition: lowe.h:195
double displayPoleMt() const
Display pole top mass.
Definition: lowe.h:135
double displayFermiConstant() const
Returns Fermi constant.
Definition: lowe.h:171
double qedBeta() const
QED beta function.
Definition: lowe.cpp:281
flexiblesusy::PMNS_parameters displayPMNS() const
returns PMNS parameters
Definition: lowe.cpp:201
void toMz()
Evolves object to MZ.
Definition: lowe.cpp:408
double mbPole
pole masses of third family quarks
Definition: lowe.h:72
Input_t input
SLHA input parmeters.
Definition: lowe.h:71
auto displayMass() const -> decltype(mf)
Returns a vector of running fermion masses.
Definition: lowe.h:149
Eigen::Array< double, NUMBER_OF_LOW_ENERGY_INPUT_PARAMETERS, 1 > Input_t
Definition: lowe.h:66
double displayLeptonRunningMass(int) const
Returns a single charged lepton running mass, given a zero-based generation index i.
Definition: lowe.cpp:243
double displayMu2GeV() const
Returns mu(2 GeV)
Definition: lowe.h:181
void setNeutrinoPoleMass(int i, double m)
sets a neutrino pole mass
Definition: lowe.h:114
void setMd2GeV(double md)
set md(2 GeV)
Definition: lowe.h:105
double displayMbMb() const
Returns mb(mb) MSbar.
Definition: lowe.h:177
void setPoleMW(double mw)
set W boson pole mass
Definition: lowe.h:107
Eigen::Array< double, 9, 1 > mf
fermion running masses
Definition: lowe.h:70
double displayMs2GeV() const
Returns ms(2 GeV)
Definition: lowe.h:185
Eigen::Array< double, 9, 1 > massBeta() const
beta functions of masses
Definition: lowe.cpp:326
QedQcd & operator=(QedQcd &&)=default
QedQcd(const QedQcd &)=default
void setCKM(const flexiblesusy::CKM_parameters &)
sets CKM parameters (in the MS-bar scheme at MZ)
Definition: lowe.cpp:173
Eigen::Matrix< std::complex< double >, 3, 3 > get_complex_pmns() const
Returns complex PMNS matrix.
Definition: lowe.h:197
void setAlphaEmInput(double a)
set input value of alpha_em(MZ)
Definition: lowe.h:120
void setPoleMb(double mb)
set pole bottom mass
Definition: lowe.h:98
double displayNeutrinoPoleMass(int i) const
Returns a single neutrino pole mass.
Definition: lowe.h:157
double displayMass(EMass mno) const
Returns a single running mass.
Definition: lowe.h:151
double qcdBeta() const
QCD beta function.
Definition: lowe.cpp:296
void setMasses(const Eigen::Array< double, 9, 1 > &m)
sets running quark masses
Definition: lowe.h:110
void setMass(EMass mno, double m)
sets a running quark mass
Definition: lowe.h:112
void setAlphas(const Eigen::Array< double, 2, 1 > &o)
sets QED and QCD structure constants
Definition: lowe.h:116
Eigen::Array< double, 2, 1 > a
gauge couplings
Definition: lowe.h:69
virtual ~QedQcd()=default
double displayPoleMel() const
Display pole electron mass.
Definition: lowe.h:141
auto displayAlphas() const -> decltype(a)
Returns gauge structure constants.
Definition: lowe.h:165
double displayPoleMmuon() const
Display pole muon mass.
Definition: lowe.h:139
void setAlphaSInput(double a)
set input value of alpha_s(MZ)
Definition: lowe.h:122
int flavours(double) const
Definition: lowe.cpp:162
Eigen::Matrix< double, 3, 3 > get_real_ckm() const
Returns real CKM matrix.
Definition: lowe.h:189
double displayPoleMb() const
Returns bottom "pole" mass.
Definition: lowe.h:143
double displayDownQuarkRunningMass(int) const
Returns a single down-quark running MS-bar mass, given a zero-based generation index i.
Definition: lowe.cpp:223
double displayPoleMtau() const
Display pole tau mass.
Definition: lowe.h:137
void to(double scale, double precision_goal=1e-5, int max_iterations=20)
Evolves object to given scale.
Definition: lowe.cpp:422
double extractPoleMb(double asMb)
returns number of active flavours
Definition: lowe.cpp:379
void setFermiConstant(double gf)
sets Fermi constant
Definition: lowe.h:128
void setMu2GeV(double mu)
set mu(2 GeV)
Definition: lowe.h:104
Eigen::Array< double, 3, 1 > guess_alpha_SM5(double scale) const
guess coupling constants {alpha_1, alpha_2, alpha_3} in SM(5)
Definition: lowe.cpp:495
virtual void set(const Eigen::ArrayXd &) override
Definition: lowe.cpp:128
double displayMcMc() const
Returns mc(mc) MSbar.
Definition: lowe.h:179
double displayLeptonPoleMass(int) const
Returns a single charged lepton pole mass, given a zero-based generation index i.
Definition: lowe.cpp:233
void set_input(const Input_t &i)
sets all input parameters
Definition: lowe.h:130
double displayPoleMZ() const
Returns Z boson pole mass.
Definition: lowe.h:147
double displayUpQuarkRunningMass(int) const
Returns a single up-quark running MS-bar mass, given a zero-based generation index i.
Definition: lowe.cpp:213
void setPoleMt(double mt)
set pole top mass
Definition: lowe.h:97
QedQcd(QedQcd &&)=default
void setPoleMtau(double mtau)
set pole tau mass
Definition: lowe.h:99
double displayPoleMW() const
Returns W boson pole mass.
Definition: lowe.h:145
Input_t display_input() const
returns vector of all input parameters
Definition: lowe.h:173
double displayMd2GeV() const
Returns md(2 GeV)
Definition: lowe.h:183
double displayAlphaEmInput() const
Returns input value alpha_em(MZ)
Definition: lowe.h:167
void setAlpha(EGauge ai, double ap)
sets QED or QCD structure constant
Definition: lowe.h:118
const double mb
Definition: consts.hpp:107
const double md
Definition: consts.hpp:102
const double mc
Definition: consts.hpp:104
const double mtau
Definition: consts.hpp:61
const double ms
Definition: consts.hpp:106
Comment if you want default softsusy behaviour.
Definition: lowe.cpp:37
EGauge
order of gauge couplings stored in QedQcd
Definition: lowe.h:44
@ ALPHA
Definition: lowe.h:44
@ ALPHAS
Definition: lowe.h:44
bool operator==(const QedQcd &a, const QedQcd &b)
Definition: lowe.cpp:519
std::ostream & operator<<(std::ostream &left, const QedQcd &m)
Formatted output from QedQcd object.
Definition: lowe.cpp:253
const std::array< std::string, NUMBER_OF_LOW_ENERGY_INPUT_PARAMETERS > QedQcd_input_parmeter_names
Definition: lowe.cpp:68
EMass
used to give order of quark masses stored
Definition: lowe.h:41
@ mDown
Definition: lowe.h:41
@ mTau
Definition: lowe.h:42
@ mStrange
Definition: lowe.h:41
@ mTop
Definition: lowe.h:41
@ mUp
Definition: lowe.h:41
@ mBottom
Definition: lowe.h:41
@ mElectron
Definition: lowe.h:41
@ mMuon
Definition: lowe.h:42
@ mCharm
Definition: lowe.h:41
QedQcd_input_parmeters
Definition: lowe.h:46
@ Mv1_pole
Definition: lowe.h:51
@ PMNS_theta_23
Definition: lowe.h:56
@ PMNS_theta_13
Definition: lowe.h:56
@ Me_pole
Definition: lowe.h:52
@ PMNS_delta
Definition: lowe.h:56
@ Mm_pole
Definition: lowe.h:52
@ mc_mc
Definition: lowe.h:54
@ CKM_delta
Definition: lowe.h:55
@ ms_2GeV
Definition: lowe.h:53
@ GFermi
Definition: lowe.h:49
@ CKM_theta_12
Definition: lowe.h:55
@ Mv3_pole
Definition: lowe.h:51
@ PMNS_alpha_1
Definition: lowe.h:56
@ PMNS_theta_12
Definition: lowe.h:56
@ Mtau_pole
Definition: lowe.h:52
@ alpha_s_MSbar_at_MZ
Definition: lowe.h:48
@ Mv2_pole
Definition: lowe.h:51
@ CKM_theta_13
Definition: lowe.h:55
@ mb_mb
Definition: lowe.h:54
@ MZ_pole
Definition: lowe.h:50
@ md_2GeV
Definition: lowe.h:54
@ PMNS_alpha_2
Definition: lowe.h:56
@ CKM_theta_23
Definition: lowe.h:55
@ NUMBER_OF_LOW_ENERGY_INPUT_PARAMETERS
Definition: lowe.h:57
@ Mt_pole
Definition: lowe.h:53
@ mu_2GeV
Definition: lowe.h:53
@ alpha_em_MSbar_at_MZ
Definition: lowe.h:47
@ MW_pole
Definition: lowe.h:50
Eigen::Matrix< double, 3, 3 > get_real_ckm() const
Definition: ckm.cpp:118
Eigen::Matrix< std::complex< double >, 3, 3 > get_complex_ckm() const
Definition: ckm.cpp:146
Eigen::Matrix< std::complex< double >, 3, 3 > get_complex_pmns() const
Definition: pmns.cpp:88
Eigen::Matrix< double, 3, 3 > get_real_pmns() const
Definition: pmns.cpp:60