a Code for the Combination of Indirect and Direct Constraints on High Energy Physics Models Logo
MPll.h
Go to the documentation of this file.
1/*
2 * Copyright (C) 2014 HEPfit Collaboration
3 *
4 *
5 * For the licensing terms see doc/COPYING.
6 */
7
8#ifndef MPLL_H
9#define MPLL_H
10
11class StandardModel;
12class F_1;
13class F_2;
14#include <gsl/gsl_integration.h>
15#include <TF1.h>
16#include <TGraph.h>
17#include <TFitResultPtr.h>
18#include <gsl/gsl_spline.h>
19#include <memory>
20
21#define MPllSWITCH 8.2
22#define LATTICE true
23#define BSZ false
24#define GSL_INTERP_DIM_DC 10
25#define SPLINE true
26
27#define NFPOLARBASIS_MPLL false
28
173class MPll{
174public:
182 MPll(const StandardModel& SM_i, QCD::meson meson_i, QCD::meson pseudoscalar_i, QCD::lepton lep_i);
183
187 virtual ~MPll();
188
189 gslpp::complex funct_g(double q2);
190
191 gslpp::complex DeltaC9_KD(double q2);
192
199 gslpp::complex h_lambda(double q2);
200
206 gslpp::complex V_L(double q2);
207
208
214 gslpp::complex T_L(double q2);
215
216
222 double S_L(double q2);
223
229 gslpp::complex H_V(double q2);
230
231
237 gslpp::complex H_A(double q2);
238
239
245 gslpp::complex H_S(double q2);
246
247
253 gslpp::complex H_P(double q2);
254
261 gslpp::complex H_nunu(double q2, QCD::lepton lep);
262
263
271 double integrateSigma(int i, double q_min, double q_max);
272
279 double getSigma(int i, double q_2);
280
288 double integrateDelta(int i, double q_min, double q_max);
289
296 double integrateSigmaTree(double q_min, double q_max);
297
302 double getwidth(){
303 updateParameters();
304 return width;
305 }
306
311 std::vector<std::string> initializeMPllParameters();
312
313
314private:
319 std::vector<std::string> mpllParameters;
320 std::unique_ptr<F_1> myF_1;
321 std::unique_ptr<F_2> myF_2;
322 bool dispersion = false;
323 bool FixedWCbtos = false;
328 double mJ2;
329
330 double GF;
331 double ale;
332 double Mlep;
333 double MM;
334 double MP;
335 double Mb;
336 double mu_b;
337 double mu_h;
338 double Mc;
339 double mb_pole;
340 double mc_pole;
341 double Ms;
343 double width;
344 double alpha_s_mub;
345 double angmomP;
346 int etaP;
348 double MW;
349 gslpp::complex lambda_t;
350 gslpp::complex h_0;
351 gslpp::complex h_1;
352 gslpp::complex h_2;
353 gslpp::complex r_1;
354 gslpp::complex r_2;
355 gslpp::complex Delta_C9;
356 gslpp::complex exp_Phase;
357
358 /*LATTICE fit parameters*/
359 double b_0_fplus;
360 double b_1_fplus;
361 double b_2_fplus;
362 double m_fit2_fplus_lat;
363 double b_0_fT;
364 double b_1_fT;
365 double b_2_fT;
366 double m_fit2_fT_lat;
367 double b_0_f0;
368 double b_1_f0;
369 double b_2_f0;
370 double m_fit2_f0_lat;
371 /*LCSR fit parameters*/
372 double r_1_fplus;
373 double r_2_fplus;
374 double m_fit2_fplus;
375 double r_1_fT;
376 double r_2_fT;
377 double m_fit2_fT;
378 double r_2_f0;
379 double m_fit2_f0;
382 gslpp::vector<gslpp::complex> ** allcoeff;
383 gslpp::vector<gslpp::complex> ** allcoeffh;
384 gslpp::vector<gslpp::complex> ** allcoeffprime;
386 gslpp::vector<gslpp::complex> ** allcoeff_noSM;
388 gslpp::vector<gslpp::complex> ** allcoeff_nu;
390 gslpp::vector<gslpp::complex> ** allcoeff_noSM_nu;
392 gslpp::complex C_1;
393 gslpp::complex C_1L_bar;
394 gslpp::complex C_1Lh_bar;
395 gslpp::complex C_2;
396 gslpp::complex C_2L_bar;
397 gslpp::complex C_2Lh_bar;
398 gslpp::complex C_3;
399 gslpp::complex C_4;
400 gslpp::complex C_5;
401 gslpp::complex C_6;
402 gslpp::complex C_7;
403 gslpp::complex C_8;
404 gslpp::complex C_8L;
405 gslpp::complex C_8Lh;
406 gslpp::complex C_9;
407 gslpp::complex C_10;
408 gslpp::complex C_S;
409 gslpp::complex C_P;
411 gslpp::complex C_7p;
412 gslpp::complex C_9p;
413 gslpp::complex C_10p;
414 gslpp::complex C_Sp;
415 gslpp::complex C_Pp;
417 gslpp::complex C_L_nunu_e;
418 gslpp::complex C_R_nunu_e;
419 gslpp::complex C_L_nunu_mu;
420 gslpp::complex C_R_nunu_mu;
421 gslpp::complex C_L_nunu_tau;
422 gslpp::complex C_R_nunu_tau;
424 gsl_interp_accel *acc_Re_deltaC7_QCDF;
425 gsl_interp_accel *acc_Im_deltaC7_QCDF;
426 gsl_interp_accel *acc_Re_deltaC9_QCDF;
427 gsl_interp_accel *acc_Im_deltaC9_QCDF;
428
429 gsl_spline *spline_Re_deltaC7_QCDF;
430 gsl_spline *spline_Im_deltaC7_QCDF;
431 gsl_spline *spline_Re_deltaC9_QCDF;
432 gsl_spline *spline_Im_deltaC9_QCDF;
433
434
435 std::vector<double> ReDeltaC9;
436 std::vector<double> ImDeltaC9;
437 std::vector<double> myq2;
438 TFitResultPtr refres;
439 TFitResultPtr imfres;
441 TGraph gr1;
442 TGraph gr2;
444 TF1 reffit;
445 TF1 imffit;
447 double tmpq2;
449 gslpp::complex H_0_pre;
450 gslpp::complex H_0_WC;
451 gslpp::complex H_c_WC;
452 gslpp::complex H_b_WC;
454 gslpp::complex ihalfMPI;
455 double fournineth;
456 double half;
457 double twothird;
458 double Mc2;
459 double Mb2;
460 double logMc;
461 double logMb;
462 double mu_b2;
463 double fourMc2;
464 double fourMb2;
465 double Mlep2;
466 double NN;
467 double MM2;
468 double MM4;
469 double MP2;
470 double MP4;
471 double MM2mMP2;
472 double twoMP2;
473 double twoMM;
474 double twoMM2;
475 double twoMM2_MMpMP;
476 double twoMM_MbpMs;
477 double S_L_pre;
478 double fourMM2;
479 double twoMboMM;
480 double sixteenM_PI2;
481 double ninetysixM_PI3MM3;
482 double MboMW;
483 double MboMM;
484 double MsoMb;
485 double twoMlepMb;
486 double threeGegen0;
487 double threeGegen1otwo;
488 double M_PI2osix;
489 double twoMc2;
490 double sixMMoMb;
491 double CF;
492 double deltaT_0;
493 double deltaT_1par;
495 gslpp::complex ubar;
496 gslpp::complex arg1;
497 gslpp::complex B01;
498 gslpp::complex B00;
499 gslpp::complex xp;
500 gslpp::complex xm;
501 gslpp::complex yp;
502 gslpp::complex ym;
503 gslpp::complex L1xp;
504 gslpp::complex L1xm;
505 gslpp::complex L1yp;
506 gslpp::complex L1ym;
507 gslpp::complex F87_0;
508 gslpp::complex F87_1;
509 gslpp::complex F87_2;
510 gslpp::complex F87_3;
511 gslpp::complex F29_0;
512 gslpp::complex F29_L1;
513 gslpp::complex F29_1;
514 gslpp::complex F29_2;
515 gslpp::complex F29_3;
516 gslpp::complex F29_L1_1;
517 gslpp::complex F29_L1_2;
518 gslpp::complex F29_L1_3;
519 gslpp::complex F27_0;
520 gslpp::complex F27_1;
521 gslpp::complex F27_2;
522 gslpp::complex F27_3;
523 gslpp::complex F27_L1_1;
524 gslpp::complex F27_L1_2;
525 gslpp::complex F27_L1_3;
526 double F89_0;
527 double F89_1;
528 double F89_2;
529 double F89_3;
530 double Ee;
532 //additional variables for B to K nu nu
533 const double M_PI2 = M_PI * M_PI;
534 double GF4;
535 double MM3;
536 double fM2;
537 double fP2;
538
539 double mtau;
540 double mtau2;
541 double Gammatau;
542 double VusVub_abs2;
543
544 unsigned int fplus_lat_updated;
545 gslpp::vector<double> fplus_lat_cache;
547 unsigned int fT_lat_updated;
548 gslpp::vector<double> fT_lat_cache;
550 unsigned int f0_lat_updated;
551 gslpp::vector<double> f0_lat_cache;
553 unsigned int fplus_updated;
554 gslpp::vector<double> fplus_cache;
556 unsigned int fT_updated;
557 gslpp::vector<double> fT_cache;
559 unsigned int f0_updated;
560 double f0_cache;
562 unsigned int k2_updated;
563 gslpp::vector<double> k2_cache;
565 unsigned int beta_updated;
566 double beta_cache;
568 unsigned int lambda_updated;
570 unsigned int F_updated;
572 unsigned int VL_updated;
574 unsigned int TL_updated;
576 unsigned int SL_updated;
577 gslpp::vector<double> SL_cache;
579 unsigned int N_updated;
580 gslpp::vector<double> N_cache;
581 gslpp::complex Nc_cache;
583 unsigned int C_1_updated;
584 gslpp::complex C_1_cache;
586 unsigned int C_2_updated;
587 gslpp::complex C_2_cache;
589 unsigned int C_3_updated;
590 gslpp::complex C_3_cache;
592 unsigned int C_4_updated;
593 gslpp::complex C_4_cache;
595 unsigned int C_5_updated;
596 gslpp::complex C_5_cache;
598 unsigned int C_6_updated;
599 gslpp::complex C_6_cache;
601 unsigned int C_7_updated;
602 gslpp::complex C_7_cache;
604 unsigned int C_9_updated;
605 gslpp::complex C_9_cache;
607 unsigned int C_10_updated;
608 gslpp::complex C_10_cache;
610 unsigned int C_7p_updated;
611 gslpp::complex C_7p_cache;
613 unsigned int C_9p_updated;
614 gslpp::complex C_9p_cache;
616 unsigned int C_10p_updated;
617 gslpp::complex C_10p_cache;
619 unsigned int C_S_updated;
620 gslpp::complex C_S_cache;
622 unsigned int C_P_updated;
623 gslpp::complex C_P_cache;
625 unsigned int C_Sp_updated;
626 gslpp::complex C_Sp_cache;
628 unsigned int C_Pp_updated;
629 gslpp::complex C_Pp_cache;
631 unsigned int C_2Lh_updated;
632 gslpp::complex C_2Lh_cache;
634 unsigned int C_8Lh_updated;
635 gslpp::complex C_8Lh_cache;
637 unsigned int C_L_nunu_e_updated;
638 unsigned int C_L_nunu_mu_updated;
639 unsigned int C_L_nunu_tau_updated;
640 gslpp::complex C_L_nunu_e_cache;
641 gslpp::complex C_L_nunu_mu_cache;
642 gslpp::complex C_L_nunu_tau_cache;
644 unsigned int C_R_nunu_e_updated;
645 unsigned int C_R_nunu_mu_updated;
646 unsigned int C_R_nunu_tau_updated;
647 gslpp::complex C_R_nunu_e_cache;
648 gslpp::complex C_R_nunu_mu_cache;
649 gslpp::complex C_R_nunu_tau_cache;
651 unsigned int Yupdated;
652 gslpp::vector<double> Ycache;
654 unsigned int H_V0updated;
655 gslpp::vector<double> H_V0cache;
656 gslpp::complex H_V0Ccache[3];
657 gslpp::complex H_V0Ccache_dispersion[4];
659 unsigned int H_A0updated;
661 unsigned int H_Supdated;
662 gslpp::vector<double> H_Scache;
664 unsigned int H_P_updated;
665 gslpp::vector<double> H_P_cache;
667 unsigned int I0_updated;
668 unsigned int I2_updated;
669 unsigned int I8_updated;
671 unsigned int Itree_updated;
672 gslpp::vector<double> Itree_cache;
674 std::map<std::pair<double, double>, unsigned int > I1Cached;
676 std::map<std::pair<double, double>, unsigned int > sigma0Cached;
677 std::map<std::pair<double, double>, unsigned int > sigma2Cached;
679 std::map<std::pair<double, double>, unsigned int > delta0Cached;
680 std::map<std::pair<double, double>, unsigned int > delta2Cached;
682 std::map<std::pair<double, double>, unsigned int > sigmaTreeCached;
684 double avaSigma;
685 double avaDelta;
686 double avaSigmaTree;
687 double avaDTPPR;
689 double errSigma;
690 double errDelta;
691 double errSigmaTree;
692 double errDTPPR;
694 gsl_function FS;
695 gsl_function FD;
696 gsl_function DTPPR;
698 gsl_integration_cquad_workspace * w_sigma;
699 gsl_integration_cquad_workspace * w_delta;
700 gsl_integration_cquad_workspace * w_sigmaTree;
701 gsl_integration_cquad_workspace * w_DTPPR;
703 gsl_error_handler_t * old_handler;
705 std::map<std::pair<double, double>, gslpp::complex > cacheI1;
707 std::map<std::pair<double, double>, double > cacheSigma0;
708 std::map<std::pair<double, double>, double > cacheSigma2;
710 std::map<std::pair<double, double>, double > cacheDelta0;
711 std::map<std::pair<double, double>, double > cacheDelta2;
713 std::map<std::pair<double, double>, double > cacheSigmaTree;
715 std::map<double, unsigned int> deltaTparpCached;
716 std::map<double, unsigned int> deltaTparmCached;
718 std::map<double, gslpp::complex> cacheDeltaTparp;
719 std::map<double, gslpp::complex> cacheDeltaTparm;
721 unsigned int deltaTparpupdated;
722 unsigned int deltaTparmupdated;
724 unsigned int T_updated;
725 gslpp::vector<double> T_cache;
732 void updateParameters();
733
737 void checkCache();
738
747 double LCSR_fit1(double q2, double r_1, double r_2, double m_fit2);
748
749
757 double LCSR_fit2(double q2, double r_2, double m_fit2);
758
759
769 double LCSR_fit3(double q2, double b_0, double b_1, double b_2, double m_fit2);
770
771
777 double zeta(double q2);
778
779
785 double zeta_DM(double q2);
786
793 double phiplus_DM(double q2, double m_fit2);
794
801 double phi0_DM(double q2);
802
809 double phiT_DM(double q2, double m_fit2);
810
820 double fplus_DM(double q2, double b_0, double b_1, double b_2, double m_fit2);
821
831 double f0_DM(double q2, double b_0, double b_1, double b_2);
832
842 double fT_DM(double q2, double b_0, double b_1, double b_2, double m_fit2);
843
844
854 double LATTICE_fit1(double q2, double b_0, double b_1, double b_2, double m_fit2);
855
856
866 double LATTICE_fit2(double q2, double b_0, double b_1, double b_2, double m_fit2);
867
868
874 double f_plus(double q2);
875
876
882 double f_T(double q2);
883
884
890 double f_0(double q2);
891
892
898 gslpp::complex H_0(double q2);
899
906 gslpp::complex H_c(double q2, double mu2);
907
914 gslpp::complex H_b(double q2, double mu2);
915
916
922 gslpp::complex Y(double q2);
923
924
930 double k2 (double q2);
931
932
938 double beta (double q2);
939
945 double beta2 (double q2);
946
947
953 double lambda(double q2);
954
955
961 double F(double q2);
962
968 double I_1c(double q2);
969
975 double I_2c(double q2);
976
982 double I_6c(double q2);
983
984
991 double Delta(int i, double q2);
992
998 double SigmaTree(double q2);
999
1004 double getintegratedSigmaTree();
1005
1011 double getSigma1c(double q2)
1012 {
1013 return I_1c(q2);
1014 };
1015
1021 double getSigma2c(double q2)
1022 {
1023 return I_2c(q2);
1024 };
1025
1031 double getSigma6c(double q2)
1032 {
1033 return I_6c(q2);
1034 };
1035
1041 double getDelta1c(double q2)
1042 {
1043 return Delta(0, q2);
1044 };
1045
1051 double getDelta2c(double q2)
1052 {
1053 return Delta(2, q2);
1054 };
1055
1062 gslpp::complex I1(double u, double q2);
1063
1070 gslpp::complex Tparplus(double u, double q2);
1071
1078 gslpp::complex Tparminus(double u, double q2);
1079
1085 double Integrand_ReTparplus(double up);
1086
1092 double Integrand_ImTparplus(double up);
1093
1099 double Integrand_ReTparminus(double up);
1100
1106 double Integrand_ImTparminus(double up);
1107
1113 double Integrand_ReTpar_pm(double up);
1114
1120 double Integrand_ImTpar_pm(double up);
1121
1127 gslpp::complex F19(double q2);
1128
1134 gslpp::complex F27(double q2);
1135
1141 gslpp::complex F29(double q2);
1142
1148 gslpp::complex F87(double q2);
1149
1155 double F89(double q2);
1156
1162 gslpp::complex Cpar(double q2);
1163
1169 gslpp::complex deltaTpar(double q2);
1170
1177 double reDC9fit(double* x, double* p);
1178
1185 double imDC9fit(double* x, double* p);
1186
1190 void fit_DeltaC9_mumu();
1191
1197 gslpp::complex fDeltaC9(double q2);
1198
1204 gslpp::complex DeltaC9(double q2);
1205
1211 gslpp::complex deltaC7_QCDF(double q2, bool spline = false);
1212
1218 gslpp::complex deltaC9_QCDF(double q2, bool spline = false);
1219
1220 void spline_QCDF_func();
1221
1222};
1223
1224#endif /* MPLL_H */
Definition: F_1.h:15
Definition: F_2.h:15
A class for the decay.
Definition: MPll.h:173
std::vector< std::string > initializeMPllParameters()
A method for initializing the parameters necessary for MPll.
Definition: MPll.cpp:80
double Integrand_ReTpar_pm(double up)
The sum of Integrand_ReTparplus() and Integrand_ReTparminus().
Definition: MPll.cpp:1043
double getSigma1c(double q2)
The CP average .
Definition: MPll.h:1011
bool FixedWCbtos
Definition: MPll.h:323
double SigmaTree(double q2)
Definition: MPll.cpp:1532
bool MPll_GRvDV_flag
Definition: MPll.h:325
gslpp::complex H_A(double q2)
The helicity amplitude .
Definition: MPll.cpp:1328
double Mb
Definition: MPll.h:335
gslpp::complex h_lambda(double q2)
The non-pertubative ccbar contributions to the helicity amplitudes.
Definition: MPll.cpp:1306
double F89(double q2)
The correction from .
Definition: MPll.cpp:1087
double mc_pole
Definition: MPll.h:340
gslpp::complex deltaC9_QCDF(double q2, bool spline=false)
QCDF Correction from various BFS papers (hep-ph/0403185, hep-ph/0412400) and Greub et....
Definition: MPll.cpp:1208
gslpp::complex DeltaC9_KD(double q2)
Definition: MPll.cpp:1300
gslpp::complex H_P(double q2)
The helicity amplitude .
Definition: MPll.cpp:1338
QCD::meson pseudoscalar
Definition: MPll.h:318
gslpp::complex H_S(double q2)
The helicity amplitude .
Definition: MPll.cpp:1333
gslpp::complex F87(double q2)
The correction from .
Definition: MPll.cpp:1080
std::vector< std::string > mpllParameters
Definition: MPll.h:319
double Integrand_ImTparminus(double up)
The imaginary part of the integral involving at fixed , according to .
Definition: MPll.cpp:1027
std::unique_ptr< F_2 > myF_2
Definition: MPll.h:321
bool MPll_Lattice_flag
Definition: MPll.h:324
double Mlep
Definition: MPll.h:332
double mu_h
Definition: MPll.h:337
gslpp::complex I1(double u, double q2)
The function from .
Definition: MPll.cpp:958
double Integrand_ReTparminus(double up)
The real part of the integral involving at fixed , according to .
Definition: MPll.cpp:1016
bool MPll_DM_flag
Definition: MPll.h:326
void fit_DeltaC9_mumu()
The fitting routine for the QCDF correction in the muon channel.
Definition: MPll.cpp:1137
gslpp::complex Tparminus(double u, double q2)
The function from .
Definition: MPll.cpp:993
double imDC9fit(double *x, double *p)
The fit function for the imaginary part of the QCDF correction .
Definition: MPll.cpp:1130
double getDelta2c(double q2)
The CP asymmetry .
Definition: MPll.h:1051
QCD::meson meson
Definition: MPll.h:317
double width
Definition: MPll.h:343
gslpp::complex F19(double q2)
The correction from .
Definition: MPll.cpp:1048
double integrateSigma(int i, double q_min, double q_max)
The integral of from to .
Definition: MPll.cpp:1407
gslpp::complex Cpar(double q2)
The correction from .
Definition: MPll.cpp:1094
double reDC9fit(double *x, double *p)
The fit function for the real part of the QCDF correction .
Definition: MPll.cpp:1125
double Mc
Definition: MPll.h:338
const StandardModel & mySM
Definition: MPll.h:315
gslpp::complex V_L(double q2)
The helicity form factor .
Definition: MPll.cpp:939
gslpp::complex H_nunu(double q2, QCD::lepton lep)
The helicity amplitude for the invisible decay.
Definition: MPll.cpp:1343
double getSigma(int i, double q_2)
The value of from to .
Definition: MPll.cpp:1444
QCD::lepton lep
Definition: MPll.h:316
virtual ~MPll()
Destructor.
Definition: MPll.cpp:76
gslpp::complex funct_g(double q2)
Definition: MPll.cpp:1292
double mJ2
Definition: MPll.h:328
double GF
Definition: MPll.h:330
double Integrand_ImTpar_pm(double up)
The sum of Integrand_ImTparplus() and Integrand_ImTparminus().
Definition: MPll.cpp:1038
double Integrand_ReTparplus(double up)
The real part of the integral involving at fixed , according to .
Definition: MPll.cpp:1000
gslpp::complex T_L(double q2)
The helicity form factor .
Definition: MPll.cpp:944
double alpha_s_mub
Definition: MPll.h:344
double getSigma2c(double q2)
The CP average .
Definition: MPll.h:1021
gslpp::complex Tparplus(double u, double q2)
The function from .
Definition: MPll.cpp:980
double MM
Definition: MPll.h:333
std::unique_ptr< F_1 > myF_1
Definition: MPll.h:320
bool NeutrinoTree_flag
Definition: MPll.h:327
double integrateSigmaTree(double q_min, double q_max)
The integral of from to (arxiv/2301.06990)
Definition: MPll.cpp:1502
double getintegratedSigmaTree()
The integral of from 0 to .
Definition: MPll.cpp:1537
double Ms
Definition: MPll.h:341
double getDelta1c(double q2)
The CP asymmetry .
Definition: MPll.h:1041
double MP
Definition: MPll.h:334
gslpp::complex H_V(double q2)
The helicity amplitude .
Definition: MPll.cpp:1321
gslpp::complex DeltaC9(double q2)
The total QCDF correction computed integrating over .
Definition: MPll.cpp:1178
gslpp::complex deltaC7_QCDF(double q2, bool spline=false)
QCDF Correction from various BFS papers (hep-ph/0403185, hep-ph/0412400) and Greub et....
Definition: MPll.cpp:1183
bool dispersion
Definition: MPll.h:322
double getSigma6c(double q2)
The CP average .
Definition: MPll.h:1031
double mb_pole
Definition: MPll.h:339
gslpp::complex F27(double q2)
The correction from .
Definition: MPll.cpp:1062
double ale
Definition: MPll.h:331
double getwidth()
The width of the meson M.
Definition: MPll.h:302
double mu_b
Definition: MPll.h:336
gslpp::complex fDeltaC9(double q2)
The total QCDF correction computed fitting over .
Definition: MPll.cpp:1169
double S_L(double q2)
The helicity form factor .
Definition: MPll.cpp:949
double spectator_charge
Definition: MPll.h:342
MPll(const StandardModel &SM_i, QCD::meson meson_i, QCD::meson pseudoscalar_i, QCD::lepton lep_i)
Constructor.
Definition: MPll.cpp:20
gslpp::complex deltaTpar(double q2)
The total correction from .
Definition: MPll.cpp:1100
double integrateDelta(int i, double q_min, double q_max)
The integral of from to .
Definition: MPll.cpp:1465
gslpp::complex F29(double q2)
The correction from .
Definition: MPll.cpp:1071
void spline_QCDF_func()
Definition: MPll.cpp:1233
double Integrand_ImTparplus(double up)
The imaginary part of the integral involving at fixed , according to .
Definition: MPll.cpp:1008
meson
An enum type for mesons.
Definition: QCD.h:336
lepton
An enum type for leptons.
Definition: QCD.h:310
A model class for the Standard Model.
A class for the correction in .