forked from Expander/FlexibleSUSY
/
weinberg_angle.cpp
438 lines (387 loc) · 14.3 KB
/
weinberg_angle.cpp
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
// ====================================================================
// This file is part of FlexibleSUSY.
//
// FlexibleSUSY is free software: you can redistribute it and/or modify
// it under the terms of the GNU General Public License as published
// by the Free Software Foundation, either version 3 of the License,
// or (at your option) any later version.
//
// FlexibleSUSY is distributed in the hope that it will be useful, but
// WITHOUT ANY WARRANTY; without even the implied warranty of
// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
// General Public License for more details.
//
// You should have received a copy of the GNU General Public License
// along with FlexibleSUSY. If not, see
// <http://www.gnu.org/licenses/>.
// ====================================================================
#include "weinberg_angle.hpp"
#include "ew_input.hpp"
#include "wrappers.hpp"
#include "logger.hpp"
#include "numerics.hpp"
#include "config.h"
#include "numerics.h"
#define CHECK_FOR_NON_ZERO(p,fun) \
if (is_zero(p)) \
WARNING(#fun ": " #p " is zero!");
namespace flexiblesusy {
namespace weinberg_angle {
using namespace softsusy;
Weinberg_angle::Data::Data()
: scale(0.)
, alpha_em_drbar(0.)
, fermi_contant(0.)
, self_energy_z_at_mz(0.)
, self_energy_w_at_0(0.)
, self_energy_w_at_mw(0.)
, mw_pole(0.)
, mz_pole(0.)
, mt_pole(0.)
, mh_drbar(0.)
, hmix_12(0.)
, mse_L(0.)
, msmu_L(0.)
, msnu_e(0.)
, msnu_mu(0.)
, mneut()
, mch()
, zn()
, um()
, up()
, gY(0.)
, g2(0.)
, g3(0.)
, tan_beta(0.)
, ymu(0.)
{
}
Weinberg_angle::Weinberg_angle()
: number_of_iterations(20)
, precision_goal(1.0e-8)
, rho_hat(0.)
, data()
{
}
Weinberg_angle::~Weinberg_angle()
{
}
void Weinberg_angle::set_number_of_iterations(unsigned n)
{
number_of_iterations = n;
}
void Weinberg_angle::set_data(const Data& data_)
{
data = data_;
}
void Weinberg_angle::set_precision_goal(double p)
{
precision_goal = p;
}
double Weinberg_angle::get_rho_hat() const
{
return rho_hat;
}
double Weinberg_angle::calculate() const
{
return 0.;
}
double Weinberg_angle::calculate_delta_rho(
double rho,
double sinThetaW,
const Data& data
) const
{
const double mz = data.mz_pole;
const double mw = data.mw_pole;
const double mt = data.mt_pole;
const double mh = data.mh_drbar;
const double sinb = sin(atan(data.tan_beta));
const double xt = 3.0 * data.fermi_contant * Sqr(mt) / (8.0 * Sqr(Pi) * root2);
const double alphaDRbar = data.alpha_em_drbar;
const double g3 = data.g3;
const double pizztMZ = data.self_energy_z_at_mz;
const double piwwtMW = data.self_energy_w_at_mw;
const double hmix12 = data.hmix_12;
#if defined(ENABLE_VERBOSE) || defined(ENABLE_DEBUG)
CHECK_FOR_NON_ZERO(rho, calculate_delta_rho)
CHECK_FOR_NON_ZERO(sinThetaW, calculate_delta_rho)
CHECK_FOR_NON_ZERO(mz, calculate_delta_rho)
CHECK_FOR_NON_ZERO(mw, calculate_delta_rho)
CHECK_FOR_NON_ZERO(mt, calculate_delta_rho)
CHECK_FOR_NON_ZERO(mh, calculate_delta_rho)
CHECK_FOR_NON_ZERO(sinb, calculate_delta_rho)
CHECK_FOR_NON_ZERO(xt, calculate_delta_rho)
CHECK_FOR_NON_ZERO(alphaDRbar, calculate_delta_rho)
CHECK_FOR_NON_ZERO(pizztMZ, calculate_delta_rho)
CHECK_FOR_NON_ZERO(piwwtMW, calculate_delta_rho)
CHECK_FOR_NON_ZERO(hmix12, calculate_delta_rho)
#endif
const double deltaRho2LoopSm = alphaDRbar * Sqr(g3) /
(16.0 * Pi * Sqr(Pi) * Sqr(sinThetaW)) *
(-2.145 * Sqr(mt) / Sqr(mw) + 1.262 * log(mt / mz) - 2.24
- 0.85 * Sqr(mz)
/ Sqr(mt)) + Sqr(xt) * Sqr(hmix12) / Sqr(sinb) *
rho_2(mh / mt) / 3.0;
const double deltaRhoOneLoop = pizztMZ / (rho * Sqr(mz))
- piwwtMW / Sqr(mw);
const double deltaRho = deltaRhoOneLoop + deltaRho2LoopSm;
return deltaRho;
}
double Weinberg_angle::calculate_delta_r(
double rho,
double sinThetaW,
const Data& data
) const
{
const double outcos = Cos(ArcSin(sinThetaW));
const double sinb = sin(atan(data.tan_beta));
const double mz = data.mz_pole;
const double mw = data.mw_pole;
const double mt = data.mt_pole;
const double mh = data.mh_drbar;
const double xt = 3.0 * data.fermi_contant * Sqr(mt) / (8.0 * Sqr(Pi) * root2);
const double alphaDRbar = data.alpha_em_drbar;
const double g3 = data.g3;
const double pizztMZ = data.self_energy_z_at_mz;
const double piwwt0 = data.self_energy_w_at_0;
const double hmix12 = data.hmix_12;
#if defined(ENABLE_VERBOSE) || defined(ENABLE_DEBUG)
CHECK_FOR_NON_ZERO(rho, calculate_delta_rho)
CHECK_FOR_NON_ZERO(sinThetaW, calculate_delta_rho)
CHECK_FOR_NON_ZERO(mz, calculate_delta_rho)
CHECK_FOR_NON_ZERO(mw, calculate_delta_rho)
CHECK_FOR_NON_ZERO(mt, calculate_delta_rho)
CHECK_FOR_NON_ZERO(mh, calculate_delta_rho)
CHECK_FOR_NON_ZERO(sinb, calculate_delta_rho)
CHECK_FOR_NON_ZERO(xt, calculate_delta_rho)
CHECK_FOR_NON_ZERO(alphaDRbar, calculate_delta_rho)
CHECK_FOR_NON_ZERO(g3, calculate_delta_rho)
CHECK_FOR_NON_ZERO(pizztMZ, calculate_delta_rho)
CHECK_FOR_NON_ZERO(piwwt0, calculate_delta_rho)
CHECK_FOR_NON_ZERO(hmix12, calculate_delta_rho)
#endif
const double dvb = calculate_delta_vb(rho, sinThetaW, data);
const double deltaR = rho * piwwt0 / Sqr(mw) -
pizztMZ / Sqr(mz) + dvb;
const double deltaR2LoopSm = alphaDRbar * Sqr(g3) /
(16.0 * Sqr(Pi) * Pi * Sqr(sinThetaW) * Sqr(outcos)) *
(2.145 * Sqr(mt) / Sqr(mz) + 0.575 * log(mt / mz) - 0.224
- 0.144 * Sqr(mz) / Sqr(mt)) -
Sqr(xt) * Sqr(hmix12) / Sqr(sinb) *
rho_2(mh / mt) * (1.0 - deltaR) * rho / 3.0;
const double deltaR_full = deltaR + deltaR2LoopSm;
return deltaR_full;
}
double Weinberg_angle::calculate_delta_vb(
double rho,
double sinThetaW,
const Data& data
) const
{
const double g = data.g2;
const double gp = data.gY;
const double mz = data.mz_pole;
const double mw = data.mw_pole;
const double costh = mw / mz;
const double cw2 = Sqr(costh);
const double sw2 = 1.0 - cw2;
const double sinThetaW2 = Sqr(sinThetaW);
const double outcos = Sqrt(1.0 - sinThetaW2);
const double q = data.scale;
const double alphaDRbar = data.alpha_em_drbar;
const double hmu = data.ymu;
const double mselL = data.mse_L;
const double msmuL = data.msmu_L;
const double msnue = data.msnu_e;
const double msnumu = data.msnu_mu;
const Eigen::ArrayXd& mneut(data.mneut);
const Eigen::ArrayXd& mch(data.mch);
const Eigen::MatrixXcd& n(data.zn);
const Eigen::MatrixXcd& u(data.um);
const Eigen::MatrixXcd& v(data.up);
#if defined(ENABLE_VERBOSE) || defined(ENABLE_DEBUG)
CHECK_FOR_NON_ZERO(rho, calculate_delta_rho)
CHECK_FOR_NON_ZERO(sinThetaW, calculate_delta_rho)
CHECK_FOR_NON_ZERO(g, calculate_delta_rho)
CHECK_FOR_NON_ZERO(gp, calculate_delta_rho)
CHECK_FOR_NON_ZERO(mz, calculate_delta_rho)
CHECK_FOR_NON_ZERO(mw, calculate_delta_rho)
CHECK_FOR_NON_ZERO(q, calculate_delta_rho)
CHECK_FOR_NON_ZERO(alphaDRbar, calculate_delta_rho)
CHECK_FOR_NON_ZERO(hmu, calculate_delta_rho)
CHECK_FOR_NON_ZERO(mselL, calculate_delta_rho)
CHECK_FOR_NON_ZERO(msmuL, calculate_delta_rho)
CHECK_FOR_NON_ZERO(msnue, calculate_delta_rho)
CHECK_FOR_NON_ZERO(msnumu, calculate_delta_rho)
#endif
const int dimN = mneut.rows();
const double deltaVbSm =
rho * alphaDRbar / (4.0 * Pi * sinThetaW2) *
(6.0 + log(cw2) / sw2 *
(3.5 - 2.5 * sw2 - sinThetaW2 * (5.0 - 1.5 * cw2 / Sqr(outcos))));
Eigen::VectorXd bPsi0NuNul(Eigen::VectorXd::Zero(dimN)),
bPsicNuSell(Eigen::VectorXd::Zero(2));
Eigen::VectorXd bPsi0ESell(Eigen::VectorXd::Zero(dimN)),
aPsicESnul(Eigen::VectorXd::Zero(2));
Eigen::VectorXcd bChi0NuNul(Eigen::VectorXcd::Zero(dimN)),
bChicNuSell(Eigen::VectorXcd::Zero(2));
Eigen::VectorXcd bChi0ESell(Eigen::VectorXcd::Zero(dimN)),
aChicESnul(Eigen::VectorXcd::Zero(2));
bPsicNuSell(0) = g;
bPsi0NuNul(1) = root2 * g * 0.5;
bPsi0NuNul(0) = -gp / root2;
aPsicESnul(0) = g;
bPsi0ESell(0) = -gp / root2;
bPsi0ESell(1) = -g * root2 * 0.5;
bChicNuSell = u * bPsicNuSell;
bChi0ESell = n * bPsi0ESell;
bChi0NuNul = n * bPsi0NuNul;
aChicESnul = v.conjugate() * aPsicESnul;
double deltaZnue = 0.0, deltaZe = 0.0;
for (int i = 0; i < dimN; i++) {
if (i < 2) {
deltaZnue +=
- Sqr(Abs(bChicNuSell(i))) * b1(0.0, mch(i), mselL, q);
deltaZe +=
- Sqr(Abs(aChicESnul(i))) * b1(0.0, mch(i), msnue, q);
}
deltaZnue +=
- Sqr(Abs(bChi0NuNul(i))) * b1(0.0, mneut(i), msnue, q);
deltaZe +=
- Sqr(Abs(bChi0ESell(i))) * b1(0.0, mneut(i), mselL, q);
}
Eigen::VectorXd bPsicNuSmul(Eigen::VectorXd::Zero(2));
Eigen::VectorXd bPsi0MuSmul(Eigen::VectorXd::Zero(dimN)),
aPsicMuSnul(Eigen::VectorXd::Zero(2));
Eigen::VectorXcd bChicNuSmul(Eigen::VectorXcd::Zero(2));
Eigen::VectorXcd bChi0MuSmul(Eigen::VectorXcd::Zero(dimN)),
aChicMuSnul(Eigen::VectorXcd::Zero(2));
bPsicNuSmul(0) = g;
bPsicNuSmul(1) = -hmu;
aPsicMuSnul(0) = g;
aPsicMuSnul(1) = -hmu;
bPsi0MuSmul(0) = -gp / root2;
bPsi0MuSmul(1) = -g * root2 * 0.5;
bChicNuSmul = u * bPsicNuSmul;
bChi0MuSmul = n * bPsi0MuSmul;
bChi0NuNul = n * bPsi0NuNul;
aChicMuSnul = v.conjugate() * aPsicMuSnul;
double deltaZnumu = 0.0, deltaZmu = 0.0;
for(int i = 0; i < dimN; i++) {
if (i < 2) {
deltaZnumu +=
- Sqr(Abs(bChicNuSmul(i))) * b1(0.0, mch(i), msmuL, q);
deltaZmu +=
- Sqr(Abs(aChicMuSnul(i))) * b1(0.0, mch(i), msnumu, q);
}
deltaZnumu +=
- Sqr(Abs(bChi0NuNul(i))) * b1(0.0, mneut(i), msnumu, q);
deltaZmu +=
- Sqr(Abs(bChi0MuSmul(i))) * b1(0.0, mneut(i), msmuL, q);
}
Eigen::MatrixXd aPsi0PsicW(Eigen::MatrixXd::Zero(dimN,2)),
bPsi0PsicW(Eigen::MatrixXd::Zero(dimN,2)),
fW(Eigen::MatrixXd::Zero(dimN,2)),
gW(Eigen::MatrixXd::Zero(dimN,2));
Eigen::MatrixXcd aChi0ChicW(Eigen::MatrixXcd::Zero(dimN,2)),
bChi0ChicW(Eigen::MatrixXcd::Zero(dimN,2));
aPsi0PsicW(1, 0) = - g;
bPsi0PsicW(1, 0) = - g;
aPsi0PsicW(3, 1) = g / root2;
bPsi0PsicW(2, 1) = -g / root2;
/// These ought to be in physpars
aChi0ChicW = n.conjugate() * aPsi0PsicW * v.transpose();
bChi0ChicW = n * bPsi0PsicW * u.adjoint();
std::complex<double> deltaVE;
for(int i = 0; i < 2; i++) {
for(int j = 0; j < dimN; j++) {
deltaVE += bChicNuSell(i) * Conj(bChi0ESell(j)) *
(- root2 / g * aChi0ChicW(j, i) * mch(i) * mneut(j) *
c0(mselL, mch(i), mneut(j)) + 1.0 / (root2 * g) *
bChi0ChicW(j, i) *
(b0(0.0, mch(i), mneut(j), q) + Sqr(mselL) *
c0(mselL, mch(i), mneut(j)) - 0.5));
deltaVE += - aChicESnul(i) * bChi0NuNul(j) *
(- root2 / g * bChi0ChicW(j, i) * mch(i) * mneut(j) *
c0(msnue, mch(i), mneut(j)) + 1.0 / (root2 * g) *
aChi0ChicW(j, i) *
(b0(0.0, mch(i), mneut(j), q) + Sqr(msnue) *
c0(msnue, mch(i), mneut(j)) - 0.5));
if (i == 0) {
deltaVE +=
0.5 * Conj(bChi0ESell(j)) * bChi0NuNul(j) *
(b0(0.0, mselL, msnue, q) + Sqr(mneut(j)) *
c0(mneut(j), mselL, msnue) + 0.5);
}
}
}
std::complex<double> deltaVMu;
for(int i = 0; i < 2; i++) {
for(int j = 0; j < dimN; j++) {
deltaVMu += bChicNuSmul(i) * Conj(bChi0MuSmul(j)) *
(- root2 / g * aChi0ChicW(j, i) * mch(i) * mneut(j) *
c0(msmuL, mch(i), mneut(j)) + 1.0 / (root2 * g) *
bChi0ChicW(j, i) *
(b0(0.0, mch(i), mneut(j), q) + Sqr(msmuL) *
c0(msmuL, mch(i), mneut(j)) - 0.5));
deltaVMu += - aChicMuSnul(i) * bChi0NuNul(j) *
(- root2 / g * bChi0ChicW(j, i) * mch(i) * mneut(j) *
c0(msnumu, mch(i), mneut(j)) + 1.0 / (root2 * g) *
aChi0ChicW(j, i) *
(b0(0.0, mch(i), mneut(j), q) + Sqr(msnumu) *
c0(msnumu, mch(i), mneut(j)) - 0.5));
if (i == 0) {
deltaVMu +=
0.5 * Conj(bChi0MuSmul(j)) * bChi0NuNul(j) *
(b0(0.0, msmuL, msnumu, q) + Sqr(mneut(j)) * c0(mneut(j), msmuL,
msnumu) + 0.5);
}
}
}
std::complex<double> a1;
for(int i = 0; i < 2; i++) {
for(int j = 0; j < dimN; j++) {
a1 += 0.5 * aChicMuSnul(i) * Conj(bChicNuSell(i)) *
bChi0NuNul(j) * bChi0ESell(j) * mch(i) * mneut(j) *
d0(mselL, msnumu, mch(i), mneut(j));
a1 += 0.5 * Conj(aChicESnul(i)) * bChicNuSmul(i) *
Conj(bChi0NuNul(j)) * Conj(bChi0MuSmul(j)) * mch(i) * mneut(j) *
d0(msmuL, msnue, mch(i), mneut(j));
a1 += bChicNuSmul(i) * Conj(bChicNuSell(i)) *
Conj(bChi0MuSmul(j)) * bChi0ESell(j) *
d27(msmuL, mselL, mch(i), mneut(j));
a1 += Conj(aChicMuSnul(i)) * aChicESnul(i) *
bChi0NuNul(j) * Conj(bChi0NuNul(j)) *
d27(msnumu, msnue, mch(i), mneut(j));
}
}
const double deltaVbSusy =
(-sinThetaW2 * Sqr(outcos) / (2.0 * Pi * alphaDRbar) * Sqr(mz)
* a1.real() + deltaVE.real() + deltaVMu.real() +
0.5 * (deltaZe + deltaZnue + deltaZmu + deltaZnumu) ) * oneOver16PiSqr;
const double deltaVb = deltaVbSusy + deltaVbSm;
return deltaVb;
}
double Weinberg_angle::rho_2(double r) const
{
if (r <= 1.9) {
return 19.0 - 16.5 * r + 43.0 * Sqr(r) / 12.0 + 7.0 / 120.0 * Sqr(r) * r -
Pi * Sqrt(r) * (4.0 - 1.5 * r + 3.0 / 32.0 * Sqr(r) + Sqr(r) * r /
256.0) - Sqr(Pi) * (2.0 - 2.0 * r + 0.5 * Sqr(r)) -
Log(r) * (3.0 * r - 0.5 * Sqr(r));
} else {
const double rm1 = 1.0 / r, rm2 = Sqr(rm1), rm3 = rm2 * rm1,
rm4 = rm3 * rm1, rm5 = rm4 * rm1;
return Sqr(Log(r)) * (1.5 - 9.0 * rm1 - 15.0 * rm2 - 48.0 * rm3 - 168.0
* rm4 - 612.0 * rm5) -
Log(r) * (13.5 + 4.0 * rm1 - 125.0 / 4.0 * rm2 - 558.0 / 5.0 * rm3 -
8307.0 / 20.0 * rm4 - 109321.0 / 70.0 * rm5)
+ Sqr(Pi) * (1.0 - 4.0 * rm1 - 5.0 * rm2 - 16.0 * rm3 -
56.0 * rm4 - 204.0 * rm5)
+ 49.0 / 4.0 + 2.0 / 3.0 * rm1 + 1613.0 / 48.0 * rm2 + 87.57 * rm3 +
341959.0 / 1200.0 * rm4 + 9737663.0 / 9800.0 * rm5;
}
}
} // namespace weinberg_angle
} // namespace flexiblesusy