a Code for the Combination of Indirect and Direct Constraints on High Energy Physics Models Logo
BXqll.h
Go to the documentation of this file.
1/*
2 * Copyright (C) 2016 HEPfit Collaboration
3 *
4 *
5 * For the licensing terms see doc/COPYING.
6 */
7
8#ifndef BXQLL_H
9#define BXQLL_H
10
11class StandardModel;
12
13#include "QCD.h"
14#include "ThObservable.h"
15//#include "Particle.h"
16#include "gslpp.h"
17#include "Expanded.h"
18#include <gsl/gsl_integration.h>
19#include "HeffDF1.h"
20#include "F_1.h"
21#include "F_2.h"
22
35class BXqll {
36public:
37
44 BXqll(const StandardModel& SM_i, QCD::quark quark_i, QCD::lepton lep_i);
45
49 virtual ~BXqll();
50
55 std::vector<std::string> initializeBXqllParameters();
56
61 double getR_LOWQ2(double sh);
62
68 double getH(std::string obs, double sh);
69
76 double integrateH(std::string obs, double q_min, double q_max);
77
78private:
91 unsigned int QCD_max, QED_max;
92
93 std::vector< gslpp::vector<gslpp::complex> > M_7;
94 std::vector< gslpp::vector<gslpp::complex> > M_9;
95 std::vector< gslpp::vector<gslpp::complex> > M_10;
96
97 std::vector< gslpp::matrix<gslpp::complex> > Hij_T;
98 std::vector< gslpp::matrix<gslpp::complex> > Hij_L;
99 std::vector< gslpp::matrix<gslpp::complex> > Hij_A;
100
101 std::vector<std::string> BXqllParameters;
103// gslpp::vector<gslpp::complex> ** allcoeff_smm;/**<Vector that contains the Wilson coeffients */
104 Expanded<gslpp::vector<gslpp::complex> > allcoeff;
106 gslpp::matrix<gslpp::complex> WC;
108 double aveH;
109 double errH;
110 gsl_function FH;
111 gsl_integration_cquad_workspace * w_H;
112 gsl_error_handler_t * old_handler;
117 void updateParameters();
118
119 double F_17re(double muh, double z, double sh, int maxpow=20);
120 double F_17im(double muh, double z, double sh, int maxpow=20);
121 double F_19re(double muh, double z, double sh, int maxpow=20);
122 double F_19im(double muh, double z, double sh, int maxpow=20);
123 double F_27re(double muh, double z, double sh, int maxpow=20);
124 double F_27im(double muh, double z, double sh, int maxpow=20);
125 double F_29re(double muh, double z, double sh, int maxpow=20);
126 double F_29im(double muh, double z, double sh, int maxpow=20);
127
128 double DeltaF_19re(double muh, double z, double sh, int maxpow=20);
129 double DeltaF_19im(double muh, double z, double sh, int maxpow=20);
130 double DeltaF_29re(double muh, double z, double sh, int maxpow=20);
131 double DeltaF_29im(double muh, double z, double sh, int maxpow=20);
137 gslpp::complex F17(double sh);
138
144 gslpp::complex F27(double sh);
145
151 gslpp::complex F19(double sh);
152
158 gslpp::complex F29(double sh);
159
165 gslpp::complex F87(double sh);
166
172 double F89(double sh);
173
179 double S77_T(double sh, orders order);
180 double S79_T(double sh, orders order);
181 double S99_T(double sh, orders order);
182 double S1010_T(double sh, orders order);
183
189 double S77_L(double sh, orders order);
190 double S79_L(double sh, orders order);
191 double S99_L(double sh, orders order);
192 double S1010_L(double sh, orders order);
193
199 double S710_A(double sh, orders order);
200 double S910_A(double sh, orders order);
201
208 gslpp::complex cij_T(unsigned int i, unsigned int j, double sh, unsigned int ord);
209 gslpp::complex cij_L(unsigned int i, unsigned int j, double sh, unsigned int ord);
210 gslpp::complex cij_A(unsigned int i, unsigned int j, double sh);
211
217 gslpp::complex eij_T(unsigned int i, unsigned int j, double sh);
218 gslpp::complex eij_L(unsigned int i, unsigned int j, double sh);
219 gslpp::complex eij_A(unsigned int i, unsigned int j, double sh);
220
225 double omega77_T(double sh);
226 double omega79_T(double sh);
227 double omega99_T(double sh);
228
233 double omega77_L(double sh);
234 double omega79_L(double sh);
235 double omega99_L(double sh);
236
241 double omega710_A(double sh);
242 double omega910_A(double sh);
243
248 double omega77em_T(double sh);
249 double omega79em_T(double sh);
250 double omega99em_T(double sh);
251 double omega22em_T(double sh);
252 gslpp::complex omega27em_T(double sh);
253 gslpp::complex omega29em_T(double sh);
254
259 double omega77em_L(double sh);
260 double omega79em_L(double sh);
261 double omega99em_L(double sh);
262 double omega22em_L(double sh);
263 gslpp::complex omega27em_L(double sh);
264 gslpp::complex omega29em_L(double sh);
265
270 double omega710em_A(double sh);
271 double omega910em_A(double sh);
272 gslpp::complex omega210em_A(double sh);
273
280 gslpp::complex f_Huber(double sh, double gamma_9, double rho_c, double rho_b, double rho_0, double rho_num);
281
286 gslpp::complex f9pen_Huber(double sh);
287
292 gslpp::complex g_Huber(double y);
293
298 gslpp::complex KS_cc(double sh);
299
308 gslpp::complex KS_aux(double sh, double m, double Gamma, double Br_ll, double Br_had);
309
314 gslpp::complex F_BIR(double r);
315
321 void computeMi(double sh);
322
327 double H_T (double sh);
328
333 double H_L (double sh);
334
339 double H_A (double sh);
340
345 double Phi_u(orders ord);
346 double Phi_u(orders_qed ord_qed);
347
352 double Phi_u_inv(unsigned int ord_qcd, unsigned int ord_qed);
353
358 double PhiTL_brems(double sh);
359
364 double PhiA_brems(double sh);
365
371 gslpp::complex C9mod(double sh);
372
379 gslpp::complex h_z(double zed, double sh);
380
386 double tau78(double sh);
387
393 double tau89(double sh);
394
400 double tau88(double sh);
401
407 double tau22fit(double sh);
408
414 double tau27fit_Re(double sh);
415
421 double tau27fit_Im(double sh);
422
428 double tau28fit_Re(double sh);
429
435 double tau28fit_Im(double sh);
436
442 double tau29fit_Re(double sh);
443
449 double tau29fit_Im(double sh);
450
456 double t810(double sh);
457
463 double t210fit(double sh);
464
469 unsigned int int_qed(orders_qed order_qed);
470
475 double CCH_multiplication(std::vector< gslpp::matrix<gslpp::complex> >& Hij);
476
480 void Test_WC_DF1();
481
482 friend double gslpp_special_functions::dilog(double x);
483};
484#endif /* BXqLL_H */
A class for the decay.
Definition: BXqll.h:35
double tau88(double sh)
The finite bremsstrahlung correction from .
Definition: BXqll.cpp:2387
double getH(std::string obs, double sh)
Method to obtain each observable as defined in .
Definition: BXqll.cpp:277
double omega22em_T(double sh)
Definition: BXqll.cpp:1389
double phinv_00_01
Definition: BXqll.h:89
double z
Definition: BXqll.h:87
double tau28fit_Re(double sh)
The fit of the real part of finite bremsstrahlung correction from .
Definition: BXqll.cpp:2468
double Mb
Definition: BXqll.h:86
gsl_function FH
Definition: BXqll.h:110
gslpp::complex KS_cc(double sh)
Kruger-Sehgal factorizable non-perturbative charm contributions following .
Definition: BXqll.cpp:1561
double omega79_L(double sh)
Definition: BXqll.cpp:1298
double phinv11
Definition: BXqll.h:89
unsigned int QED_max
Definition: BXqll.h:91
double PhiTL_brems(double sh)
The finite bremsstrahlung corrections to dGamma/ds for from .
Definition: BXqll.cpp:2281
gslpp::complex h_z(double zed, double sh)
Auxiliary function from .
Definition: BXqll.cpp:2339
std::vector< gslpp::vector< gslpp::complex > > M_10
Definition: BXqll.h:95
double Phi_u(orders ord)
Normalization function for from eq. (4.8) of 1503.04849.
Definition: BXqll.cpp:1610
double Mc
Definition: BXqll.h:86
double H_A(double sh)
Angular observable as defined in .
Definition: BXqll.cpp:694
double phi00
Definition: BXqll.h:88
double alstilde
Definition: BXqll.h:85
double S79_T(double sh, orders order)
Definition: BXqll.cpp:889
gslpp::complex F29(double sh)
The correction from .
Definition: BXqll.cpp:153
double omega77em_T(double sh)
Auxiliary functions from .
Definition: BXqll.cpp:1365
gslpp::complex eij_L(unsigned int i, unsigned int j, double sh)
Definition: BXqll.cpp:1177
double phinv20
Definition: BXqll.h:89
gslpp::complex f_Huber(double sh, double gamma_9, double rho_c, double rho_b, double rho_0, double rho_num)
Auxiliary function from .
Definition: BXqll.cpp:1522
BXqll(const StandardModel &SM_i, QCD::quark quark_i, QCD::lepton lep_i)
Constructor.
Definition: BXqll.cpp:21
double F_27re(double muh, double z, double sh, int maxpow=20)
Definition: BXqll.cpp:215
double integrateH(std::string obs, double q_min, double q_max)
The integral of each observable as defined in .
Definition: BXqll.cpp:260
double omega79_T(double sh)
Definition: BXqll.cpp:1257
double F_19re(double muh, double z, double sh, int maxpow=20)
Definition: BXqll.cpp:205
double S910_A(double sh, orders order)
Definition: BXqll.cpp:1021
double DeltaF_29im(double muh, double z, double sh, int maxpow=20)
Definition: BXqll.cpp:250
double H_T(double sh)
Angular observable as defined in .
Definition: BXqll.cpp:296
double alsmuc
Definition: BXqll.h:85
double kappa
Definition: BXqll.h:85
gslpp::complex omega27em_T(double sh)
Definition: BXqll.cpp:1401
std::vector< gslpp::vector< gslpp::complex > > M_9
Definition: BXqll.h:94
double S77_L(double sh, orders order)
Auxiliary functions from .
Definition: BXqll.cpp:936
double S77_T(double sh, orders order)
Auxiliary functions from .
Definition: BXqll.cpp:870
gslpp::complex cij_A(unsigned int i, unsigned int j, double sh)
Definition: BXqll.cpp:1120
double Mlep
Definition: BXqll.h:86
double H_L(double sh)
Angular observable as defined in .
Definition: BXqll.cpp:495
gslpp::complex F87(double sh)
The correction from .
Definition: BXqll.cpp:174
gslpp::complex g_Huber(double y)
Auxiliary function from .
Definition: BXqll.cpp:1545
void computeMi(double sh)
Vectors of auxiliary functions from Table 6 of .
Definition: BXqll.cpp:803
double lambda_1
Definition: BXqll.h:87
double S99_L(double sh, orders order)
Definition: BXqll.cpp:974
double ale
Definition: BXqll.h:85
double phi2
Definition: BXqll.h:88
std::vector< gslpp::matrix< gslpp::complex > > Hij_L
Definition: BXqll.h:98
double tau27fit_Re(double sh)
The fit of the real part of finite bremsstrahlung correction from .
Definition: BXqll.cpp:2422
double C_ratio
Definition: BXqll.h:90
double tau27fit_Im(double sh)
The fit of the imaginary part of finite bremsstrahlung correction from .
Definition: BXqll.cpp:2445
double alsmu
Definition: BXqll.h:85
unsigned int QCD_max
Definition: BXqll.h:91
double aletilde
Definition: BXqll.h:85
double omega77_L(double sh)
Auxiliary functions from .
Definition: BXqll.cpp:1284
void Test_WC_DF1()
Temporary method to test Wilson coefficients with C10_OS1 matching and HeffDF1 evolution.
Definition: BXqll.cpp:1712
gslpp::complex F17(double sh)
The correction from .
Definition: BXqll.cpp:160
gslpp::complex eij_T(unsigned int i, unsigned int j, double sh)
Log-enhanced electromagnetic corrections as defined in .
Definition: BXqll.cpp:1137
double mu_c
Definition: BXqll.h:86
double tau28fit_Im(double sh)
The fit of the imaginary part of finite bremsstrahlung correction from .
Definition: BXqll.cpp:2491
double S99_T(double sh, orders order)
Definition: BXqll.cpp:908
double Mtau
Definition: BXqll.h:86
void updateParameters()
The update parameter method for BXqll.
Definition: BXqll.cpp:48
double phinv00
Definition: BXqll.h:89
double F_27im(double muh, double z, double sh, int maxpow=20)
Definition: BXqll.cpp:220
double t810(double sh)
The finite bremsstrahlung correction from .
Definition: BXqll.cpp:2560
gslpp::complex cij_L(unsigned int i, unsigned int j, double sh, unsigned int ord)
Definition: BXqll.cpp:1082
gslpp::complex omega29em_L(double sh)
Definition: BXqll.cpp:1476
gsl_integration_cquad_workspace * w_H
Definition: BXqll.h:111
gslpp::complex C9mod(double sh)
The modified coefficient from .
Definition: BXqll.cpp:2325
double omega79em_T(double sh)
Definition: BXqll.cpp:1373
double phinv01
Definition: BXqll.h:89
QCD::lepton lep
Definition: BXqll.h:83
double phi1
Definition: BXqll.h:88
double omega79em_L(double sh)
Definition: BXqll.cpp:1435
double omega22em_L(double sh)
Definition: BXqll.cpp:1451
gslpp::complex cij_T(unsigned int i, unsigned int j, double sh, unsigned int ord)
contributions as defined in
Definition: BXqll.cpp:1044
double Vts_over_Vcb
Definition: BXqll.h:87
double omega77em_L(double sh)
Auxiliary functions from .
Definition: BXqll.cpp:1427
gslpp::complex omega210em_A(double sh)
Definition: BXqll.cpp:1505
double CCH_multiplication(std::vector< gslpp::matrix< gslpp::complex > > &Hij)
Auxiliary function that performs the multiplication of Wilson coefficients and matrix elements.
Definition: BXqll.cpp:1668
double lambda_2
Definition: BXqll.h:87
double MW
Definition: BXqll.h:86
double omega710_A(double sh)
Auxiliary functions from .
Definition: BXqll.cpp:1325
double tau78(double sh)
The finite bremsstrahlung correction from .
Definition: BXqll.cpp:2362
gslpp::complex KS_aux(double sh, double m, double Gamma, double Br_ll, double Br_had)
Auxiliary function for the Kruger-Sehgal charm contributions.
Definition: BXqll.cpp:1582
double F_29im(double muh, double z, double sh, int maxpow=20)
Definition: BXqll.cpp:230
gslpp::complex F19(double sh)
The correction from .
Definition: BXqll.cpp:146
double phinv21
Definition: BXqll.h:89
double pre
Definition: BXqll.h:90
F_1 myF_1
Definition: BXqll.h:80
std::vector< gslpp::matrix< gslpp::complex > > Hij_T
Definition: BXqll.h:97
double DeltaF_19re(double muh, double z, double sh, int maxpow=20)
Definition: BXqll.cpp:235
double F_17re(double muh, double z, double sh, int maxpow=20)
Definition: BXqll.cpp:195
HeffDF1 myHeff
Definition: BXqll.h:82
double tau22fit(double sh)
The fit of the finite bremsstrahlung correction from .
Definition: BXqll.cpp:2399
gslpp::matrix< gslpp::complex > WC
Definition: BXqll.h:106
double tau89(double sh)
The finite bremsstrahlung correction from .
Definition: BXqll.cpp:2375
double S1010_L(double sh, orders order)
Definition: BXqll.cpp:997
double omega99em_T(double sh)
Definition: BXqll.cpp:1381
double S79_L(double sh, orders order)
Definition: BXqll.cpp:955
double phi01
Definition: BXqll.h:88
double abslambdat_over_Vcb
Definition: BXqll.h:87
double mu_b
Definition: BXqll.h:86
std::vector< std::string > BXqllParameters
Definition: BXqll.h:101
double Lbl
Definition: BXqll.h:87
std::vector< std::string > initializeBXqllParameters()
A method for initializing the parameters necessary for BXqll.
Definition: BXqll.cpp:41
unsigned int int_qed(orders_qed order_qed)
Auxiliary function that matches orders_qed to an integer.
Definition: BXqll.cpp:1662
const StandardModel & mySM
Definition: BXqll.h:79
double omega77_T(double sh)
Auxiliary functions from .
Definition: BXqll.cpp:1243
double DeltaF_19im(double muh, double z, double sh, int maxpow=20)
Definition: BXqll.cpp:240
gslpp::complex omega27em_L(double sh)
Definition: BXqll.cpp:1463
double tau29fit_Re(double sh)
The fit of the real part of finite bremsstrahlung correction from .
Definition: BXqll.cpp:2514
double phi20
Definition: BXqll.h:88
double BR_BXcenu
Definition: BXqll.h:90
double omega910em_A(double sh)
Definition: BXqll.cpp:1497
double CF
Definition: BXqll.h:85
double Mc_pole
Definition: BXqll.h:86
gslpp::complex F27(double sh)
The correction from .
Definition: BXqll.cpp:167
double F_19im(double muh, double z, double sh, int maxpow=20)
Definition: BXqll.cpp:210
double F89(double sh)
The correction from .
Definition: BXqll.cpp:185
std::vector< gslpp::vector< gslpp::complex > > M_7
Definition: BXqll.h:93
double F_17im(double muh, double z, double sh, int maxpow=20)
Definition: BXqll.cpp:200
double tau29fit_Im(double sh)
The fit of the imaginary part of finite bremsstrahlung correction from .
Definition: BXqll.cpp:2537
double S710_A(double sh, orders order)
Auxiliary functions from .
Definition: BXqll.cpp:1002
gslpp::complex F_BIR(double r)
Auxiliary function from .
Definition: BXqll.cpp:1597
friend double gslpp_special_functions::dilog(double x)
virtual ~BXqll()
Destructor.
Definition: BXqll.cpp:37
double F_29re(double muh, double z, double sh, int maxpow=20)
Definition: BXqll.cpp:225
double getR_LOWQ2(double sh)
dGamma/ds for in the low dilepton invariant mass region.
Definition: BXqll.cpp:134
gslpp::complex eij_A(unsigned int i, unsigned int j, double sh)
Definition: BXqll.cpp:1217
double aveH
Definition: BXqll.h:108
double PhiA_brems(double sh)
The finite bremsstrahlung corrections to dAFB/ds for from .
Definition: BXqll.cpp:2315
double omega910_A(double sh)
Definition: BXqll.cpp:1351
gsl_error_handler_t * old_handler
Definition: BXqll.h:112
Expanded< gslpp::vector< gslpp::complex > > allcoeff
Definition: BXqll.h:104
gslpp::complex f9pen_Huber(double sh)
Auxiliary function from .
Definition: BXqll.cpp:1534
double muh
Definition: BXqll.h:87
std::vector< gslpp::matrix< gslpp::complex > > Hij_A
Definition: BXqll.h:99
double t210fit(double sh)
The fit of finite bremsstrahlung correction from .
Definition: BXqll.cpp:2582
double omega99em_L(double sh)
Definition: BXqll.cpp:1443
double DeltaF_29re(double muh, double z, double sh, int maxpow=20)
Definition: BXqll.cpp:245
double Mb_pole
Definition: BXqll.h:86
double phi00_2
Definition: BXqll.h:88
F_2 myF_2
Definition: BXqll.h:81
QCD::quark quark
Definition: BXqll.h:84
double GF
Definition: BXqll.h:85
double omega99_L(double sh)
Definition: BXqll.cpp:1311
double omega99_T(double sh)
Definition: BXqll.cpp:1270
double phinv10
Definition: BXqll.h:89
double omega710em_A(double sh)
Auxiliary functions from .
Definition: BXqll.cpp:1489
gslpp::complex omega29em_T(double sh)
Definition: BXqll.cpp:1414
double Phi_u_inv(unsigned int ord_qcd, unsigned int ord_qed)
Inverse of the normalization function for from eq. (4.8) of 1503.04849.
Definition: BXqll.cpp:1644
double errH
Definition: BXqll.h:109
double Ms
Definition: BXqll.h:86
double S1010_T(double sh, orders order)
Definition: BXqll.cpp:931
Definition: F_1.h:15
Definition: F_2.h:15
quark
An enum type for quarks.
Definition: QCD.h:323
lepton
An enum type for leptons.
Definition: QCD.h:310
A model class for the Standard Model.
orders
An enum type for orders in QCD.
Definition: OrderScheme.h:33
orders_qed
An enum type for orders in electroweak.
Definition: OrderScheme.h:56