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 double Delta_C7_U;
359 double Delta_C9_U;
360
361 /*LATTICE fit parameters*/
362 double b_0_fplus;
363 double b_1_fplus;
364 double b_2_fplus;
365 double m_fit2_fplus_lat;
366 double b_0_fT;
367 double b_1_fT;
368 double b_2_fT;
369 double m_fit2_fT_lat;
370 double b_0_f0;
371 double b_1_f0;
372 double b_2_f0;
373 double m_fit2_f0_lat;
374 /*LCSR fit parameters*/
375 double r_1_fplus;
376 double r_2_fplus;
377 double m_fit2_fplus;
378 double r_1_fT;
379 double r_2_fT;
380 double m_fit2_fT;
381 double r_2_f0;
382 double m_fit2_f0;
385 gslpp::vector<gslpp::complex> ** allcoeff;
386 gslpp::vector<gslpp::complex> ** allcoeffh;
387 gslpp::vector<gslpp::complex> ** allcoeffprime;
389 gslpp::vector<gslpp::complex> ** allcoeff_noSM;
391 gslpp::vector<gslpp::complex> ** allcoeff_nu;
393 gslpp::vector<gslpp::complex> ** allcoeff_noSM_nu;
395 gslpp::complex C_1;
396 gslpp::complex C_1L_bar;
397 gslpp::complex C_1Lh_bar;
398 gslpp::complex C_2;
399 gslpp::complex C_2L_bar;
400 gslpp::complex C_2Lh_bar;
401 gslpp::complex C_3;
402 gslpp::complex C_4;
403 gslpp::complex C_5;
404 gslpp::complex C_6;
405 gslpp::complex C_7;
406 gslpp::complex C_8;
407 gslpp::complex C_8L;
408 gslpp::complex C_8Lh;
409 gslpp::complex C_9;
410 gslpp::complex C_10;
411 gslpp::complex C_S;
412 gslpp::complex C_P;
414 gslpp::complex C_7p;
415 gslpp::complex C_9p;
416 gslpp::complex C_10p;
417 gslpp::complex C_Sp;
418 gslpp::complex C_Pp;
420 gslpp::complex C_L_nunu_e;
421 gslpp::complex C_R_nunu_e;
422 gslpp::complex C_L_nunu_mu;
423 gslpp::complex C_R_nunu_mu;
424 gslpp::complex C_L_nunu_tau;
425 gslpp::complex C_R_nunu_tau;
427 gsl_interp_accel *acc_Re_deltaC7_QCDF;
428 gsl_interp_accel *acc_Im_deltaC7_QCDF;
429 gsl_interp_accel *acc_Re_deltaC9_QCDF;
430 gsl_interp_accel *acc_Im_deltaC9_QCDF;
431
432 gsl_spline *spline_Re_deltaC7_QCDF;
433 gsl_spline *spline_Im_deltaC7_QCDF;
434 gsl_spline *spline_Re_deltaC9_QCDF;
435 gsl_spline *spline_Im_deltaC9_QCDF;
436
437
438 std::vector<double> ReDeltaC9;
439 std::vector<double> ImDeltaC9;
440 std::vector<double> myq2;
441 TFitResultPtr refres;
442 TFitResultPtr imfres;
444 TGraph gr1;
445 TGraph gr2;
447 TF1 reffit;
448 TF1 imffit;
450 double tmpq2;
452 gslpp::complex H_0_pre;
453 gslpp::complex H_0_WC;
454 gslpp::complex H_c_WC;
455 gslpp::complex H_b_WC;
457 gslpp::complex ihalfMPI;
458 double fournineth;
459 double half;
460 double twothird;
461 double Mc2;
462 double Mb2;
463 double logMc;
464 double logMb;
465 double mu_b2;
466 double fourMc2;
467 double fourMb2;
468 double Mlep2;
469 double NN;
470 double MM2;
471 double MM4;
472 double MP2;
473 double MP4;
474 double MM2mMP2;
475 double twoMP2;
476 double twoMM;
477 double twoMM2;
478 double twoMM2_MMpMP;
479 double twoMM_MbpMs;
480 double S_L_pre;
481 double fourMM2;
482 double twoMboMM;
483 double sixteenM_PI2;
484 double ninetysixM_PI3MM3;
485 double MboMW;
486 double MboMM;
487 double MsoMb;
488 double twoMlepMb;
489 double threeGegen0;
490 double threeGegen1otwo;
491 double M_PI2osix;
492 double twoMc2;
493 double sixMMoMb;
494 double CF;
495 double deltaT_0;
496 double deltaT_1par;
498 gslpp::complex ubar;
499 gslpp::complex arg1;
500 gslpp::complex B01;
501 gslpp::complex B00;
502 gslpp::complex xp;
503 gslpp::complex xm;
504 gslpp::complex yp;
505 gslpp::complex ym;
506 gslpp::complex L1xp;
507 gslpp::complex L1xm;
508 gslpp::complex L1yp;
509 gslpp::complex L1ym;
510 gslpp::complex F87_0;
511 gslpp::complex F87_1;
512 gslpp::complex F87_2;
513 gslpp::complex F87_3;
514 gslpp::complex F29_0;
515 gslpp::complex F29_L1;
516 gslpp::complex F29_1;
517 gslpp::complex F29_2;
518 gslpp::complex F29_3;
519 gslpp::complex F29_L1_1;
520 gslpp::complex F29_L1_2;
521 gslpp::complex F29_L1_3;
522 gslpp::complex F27_0;
523 gslpp::complex F27_1;
524 gslpp::complex F27_2;
525 gslpp::complex F27_3;
526 gslpp::complex F27_L1_1;
527 gslpp::complex F27_L1_2;
528 gslpp::complex F27_L1_3;
529 double F89_0;
530 double F89_1;
531 double F89_2;
532 double F89_3;
533 double Ee;
535 //additional variables for B to K nu nu
536 const double M_PI2 = M_PI * M_PI;
537 double GF4;
538 double MM3;
539 double fM2;
540 double fP2;
541
542 double mtau;
543 double mtau2;
544 double Gammatau;
545 double VusVub_abs2;
546
547 unsigned int fplus_lat_updated;
548 gslpp::vector<double> fplus_lat_cache;
550 unsigned int fT_lat_updated;
551 gslpp::vector<double> fT_lat_cache;
553 unsigned int f0_lat_updated;
554 gslpp::vector<double> f0_lat_cache;
556 unsigned int fplus_updated;
557 gslpp::vector<double> fplus_cache;
559 unsigned int fT_updated;
560 gslpp::vector<double> fT_cache;
562 unsigned int f0_updated;
563 double f0_cache;
565 unsigned int k2_updated;
566 gslpp::vector<double> k2_cache;
568 unsigned int beta_updated;
569 double beta_cache;
571 unsigned int lambda_updated;
573 unsigned int F_updated;
575 unsigned int VL_updated;
577 unsigned int TL_updated;
579 unsigned int SL_updated;
580 gslpp::vector<double> SL_cache;
582 unsigned int N_updated;
583 gslpp::vector<double> N_cache;
584 gslpp::complex Nc_cache;
586 unsigned int C_1_updated;
587 gslpp::complex C_1_cache;
589 unsigned int C_2_updated;
590 gslpp::complex C_2_cache;
592 unsigned int C_3_updated;
593 gslpp::complex C_3_cache;
595 unsigned int C_4_updated;
596 gslpp::complex C_4_cache;
598 unsigned int C_5_updated;
599 gslpp::complex C_5_cache;
601 unsigned int C_6_updated;
602 gslpp::complex C_6_cache;
604 unsigned int C_7_updated;
605 gslpp::complex C_7_cache;
607 unsigned int C_9_updated;
608 gslpp::complex C_9_cache;
610 unsigned int C_10_updated;
611 gslpp::complex C_10_cache;
613 unsigned int C_7p_updated;
614 gslpp::complex C_7p_cache;
616 unsigned int C_9p_updated;
617 gslpp::complex C_9p_cache;
619 unsigned int C_10p_updated;
620 gslpp::complex C_10p_cache;
622 unsigned int C_S_updated;
623 gslpp::complex C_S_cache;
625 unsigned int C_P_updated;
626 gslpp::complex C_P_cache;
628 unsigned int C_Sp_updated;
629 gslpp::complex C_Sp_cache;
631 unsigned int C_Pp_updated;
632 gslpp::complex C_Pp_cache;
634 unsigned int C_2Lh_updated;
635 gslpp::complex C_2Lh_cache;
637 unsigned int C_8Lh_updated;
638 gslpp::complex C_8Lh_cache;
640 unsigned int C_L_nunu_e_updated;
641 unsigned int C_L_nunu_mu_updated;
642 unsigned int C_L_nunu_tau_updated;
643 gslpp::complex C_L_nunu_e_cache;
644 gslpp::complex C_L_nunu_mu_cache;
645 gslpp::complex C_L_nunu_tau_cache;
647 unsigned int C_R_nunu_e_updated;
648 unsigned int C_R_nunu_mu_updated;
649 unsigned int C_R_nunu_tau_updated;
650 gslpp::complex C_R_nunu_e_cache;
651 gslpp::complex C_R_nunu_mu_cache;
652 gslpp::complex C_R_nunu_tau_cache;
654 unsigned int Yupdated;
655 gslpp::vector<double> Ycache;
657 unsigned int H_V0updated;
658 gslpp::vector<double> H_V0cache;
659 gslpp::complex H_V0Ccache[3];
660 gslpp::complex H_V0Ccache_dispersion[4];
662 unsigned int H_A0updated;
664 unsigned int H_Supdated;
665 gslpp::vector<double> H_Scache;
667 unsigned int H_P_updated;
668 gslpp::vector<double> H_P_cache;
670 unsigned int I0_updated;
671 unsigned int I2_updated;
672 unsigned int I8_updated;
674 unsigned int Itree_updated;
675 gslpp::vector<double> Itree_cache;
677 std::map<std::pair<double, double>, unsigned int > I1Cached;
679 std::map<std::pair<double, double>, unsigned int > sigma0Cached;
680 std::map<std::pair<double, double>, unsigned int > sigma2Cached;
682 std::map<std::pair<double, double>, unsigned int > delta0Cached;
683 std::map<std::pair<double, double>, unsigned int > delta2Cached;
685 std::map<std::pair<double, double>, unsigned int > sigmaTreeCached;
687 double avaSigma;
688 double avaDelta;
689 double avaSigmaTree;
690 double avaDTPPR;
692 double errSigma;
693 double errDelta;
694 double errSigmaTree;
695 double errDTPPR;
697 gsl_function FS;
698 gsl_function FD;
699 gsl_function DTPPR;
701 gsl_integration_cquad_workspace * w_sigma;
702 gsl_integration_cquad_workspace * w_delta;
703 gsl_integration_cquad_workspace * w_sigmaTree;
704 gsl_integration_cquad_workspace * w_DTPPR;
706 gsl_error_handler_t * old_handler;
708 std::map<std::pair<double, double>, gslpp::complex > cacheI1;
710 std::map<std::pair<double, double>, double > cacheSigma0;
711 std::map<std::pair<double, double>, double > cacheSigma2;
713 std::map<std::pair<double, double>, double > cacheDelta0;
714 std::map<std::pair<double, double>, double > cacheDelta2;
716 std::map<std::pair<double, double>, double > cacheSigmaTree;
718 std::map<double, unsigned int> deltaTparpCached;
719 std::map<double, unsigned int> deltaTparmCached;
721 std::map<double, gslpp::complex> cacheDeltaTparp;
722 std::map<double, gslpp::complex> cacheDeltaTparm;
724 unsigned int deltaTparpupdated;
725 unsigned int deltaTparmupdated;
727 unsigned int T_updated;
728 gslpp::vector<double> T_cache;
735 void updateParameters();
736
740 void checkCache();
741
750 double LCSR_fit1(double q2, double r_1, double r_2, double m_fit2);
751
752
760 double LCSR_fit2(double q2, double r_2, double m_fit2);
761
762
772 double LCSR_fit3(double q2, double b_0, double b_1, double b_2, double m_fit2);
773
774
780 double zeta(double q2);
781
782
788 double zeta_DM(double q2);
789
796 double phiplus_DM(double q2, double m_fit2);
797
804 double phi0_DM(double q2);
805
812 double phiT_DM(double q2, double m_fit2);
813
823 double fplus_DM(double q2, double b_0, double b_1, double b_2, double m_fit2);
824
834 double f0_DM(double q2, double b_0, double b_1, double b_2);
835
845 double fT_DM(double q2, double b_0, double b_1, double b_2, double m_fit2);
846
847
857 double LATTICE_fit1(double q2, double b_0, double b_1, double b_2, double m_fit2);
858
859
869 double LATTICE_fit2(double q2, double b_0, double b_1, double b_2, double m_fit2);
870
871
877 double f_plus(double q2);
878
879
885 double f_T(double q2);
886
887
893 double f_0(double q2);
894
895
901 gslpp::complex H_0(double q2);
902
909 gslpp::complex H_c(double q2, double mu2);
910
917 gslpp::complex H_b(double q2, double mu2);
918
919
925 gslpp::complex Y(double q2);
926
927
933 double k2 (double q2);
934
935
941 double beta (double q2);
942
948 double beta2 (double q2);
949
950
956 double lambda(double q2);
957
958
964 double F(double q2);
965
971 double I_1c(double q2);
972
978 double I_2c(double q2);
979
985 double I_6c(double q2);
986
987
994 double Delta(int i, double q2);
995
1001 double SigmaTree(double q2);
1002
1007 double getintegratedSigmaTree();
1008
1014 double getSigma1c(double q2)
1015 {
1016 return I_1c(q2);
1017 };
1018
1024 double getSigma2c(double q2)
1025 {
1026 return I_2c(q2);
1027 };
1028
1034 double getSigma6c(double q2)
1035 {
1036 return I_6c(q2);
1037 };
1038
1044 double getDelta1c(double q2)
1045 {
1046 return Delta(0, q2);
1047 };
1048
1054 double getDelta2c(double q2)
1055 {
1056 return Delta(2, q2);
1057 };
1058
1065 gslpp::complex I1(double u, double q2);
1066
1073 gslpp::complex Tparplus(double u, double q2);
1074
1081 gslpp::complex Tparminus(double u, double q2);
1082
1088 double Integrand_ReTparplus(double up);
1089
1095 double Integrand_ImTparplus(double up);
1096
1102 double Integrand_ReTparminus(double up);
1103
1109 double Integrand_ImTparminus(double up);
1110
1116 double Integrand_ReTpar_pm(double up);
1117
1123 double Integrand_ImTpar_pm(double up);
1124
1130 gslpp::complex F19(double q2);
1131
1137 gslpp::complex F27(double q2);
1138
1144 gslpp::complex F29(double q2);
1145
1151 gslpp::complex F87(double q2);
1152
1158 double F89(double q2);
1159
1165 gslpp::complex Cpar(double q2);
1166
1172 gslpp::complex deltaTpar(double q2);
1173
1180 double reDC9fit(double* x, double* p);
1181
1188 double imDC9fit(double* x, double* p);
1189
1193 void fit_DeltaC9_mumu();
1194
1200 gslpp::complex fDeltaC9(double q2);
1201
1207 gslpp::complex DeltaC9(double q2);
1208
1214 gslpp::complex deltaC7_QCDF(double q2, bool spline = false);
1215
1221 gslpp::complex deltaC9_QCDF(double q2, bool spline = false);
1222
1223 void spline_QCDF_func();
1224
1225};
1226
1227#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:1057
double getSigma1c(double q2)
The CP average .
Definition: MPll.h:1014
bool FixedWCbtos
Definition: MPll.h:323
double SigmaTree(double q2)
Definition: MPll.cpp:1549
bool MPll_GRvDV_flag
Definition: MPll.h:325
gslpp::complex H_A(double q2)
The helicity amplitude .
Definition: MPll.cpp:1344
double Mb
Definition: MPll.h:335
gslpp::complex h_lambda(double q2)
The non-pertubative ccbar contributions to the helicity amplitudes.
Definition: MPll.cpp:1320
double F89(double q2)
The correction from .
Definition: MPll.cpp:1101
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:1222
gslpp::complex DeltaC9_KD(double q2)
Definition: MPll.cpp:1314
gslpp::complex H_P(double q2)
The helicity amplitude .
Definition: MPll.cpp:1354
QCD::meson pseudoscalar
Definition: MPll.h:318
gslpp::complex H_S(double q2)
The helicity amplitude .
Definition: MPll.cpp:1349
gslpp::complex F87(double q2)
The correction from .
Definition: MPll.cpp:1094
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:1041
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:972
double Integrand_ReTparminus(double up)
The real part of the integral involving at fixed , according to .
Definition: MPll.cpp:1030
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:1151
gslpp::complex Tparminus(double u, double q2)
The function from .
Definition: MPll.cpp:1007
double imDC9fit(double *x, double *p)
The fit function for the imaginary part of the QCDF correction .
Definition: MPll.cpp:1144
double getDelta2c(double q2)
The CP asymmetry .
Definition: MPll.h:1054
QCD::meson meson
Definition: MPll.h:317
double width
Definition: MPll.h:343
gslpp::complex F19(double q2)
The correction from .
Definition: MPll.cpp:1062
double integrateSigma(int i, double q_min, double q_max)
The integral of from to .
Definition: MPll.cpp:1424
gslpp::complex Cpar(double q2)
The correction from .
Definition: MPll.cpp:1108
double reDC9fit(double *x, double *p)
The fit function for the real part of the QCDF correction .
Definition: MPll.cpp:1139
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:953
gslpp::complex H_nunu(double q2, QCD::lepton lep)
The helicity amplitude for the invisible decay.
Definition: MPll.cpp:1359
double getSigma(int i, double q_2)
The value of from to .
Definition: MPll.cpp:1461
QCD::lepton lep
Definition: MPll.h:316
virtual ~MPll()
Destructor.
Definition: MPll.cpp:76
gslpp::complex funct_g(double q2)
Definition: MPll.cpp:1306
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:1052
double Integrand_ReTparplus(double up)
The real part of the integral involving at fixed , according to .
Definition: MPll.cpp:1014
gslpp::complex T_L(double q2)
The helicity form factor .
Definition: MPll.cpp:958
double alpha_s_mub
Definition: MPll.h:344
double getSigma2c(double q2)
The CP average .
Definition: MPll.h:1024
gslpp::complex Tparplus(double u, double q2)
The function from .
Definition: MPll.cpp:994
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:1519
double getintegratedSigmaTree()
The integral of from 0 to .
Definition: MPll.cpp:1554
double Ms
Definition: MPll.h:341
double getDelta1c(double q2)
The CP asymmetry .
Definition: MPll.h:1044
double MP
Definition: MPll.h:334
gslpp::complex H_V(double q2)
The helicity amplitude .
Definition: MPll.cpp:1337
gslpp::complex DeltaC9(double q2)
The total QCDF correction computed integrating over .
Definition: MPll.cpp:1192
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:1197
bool dispersion
Definition: MPll.h:322
double getSigma6c(double q2)
The CP average .
Definition: MPll.h:1034
double mb_pole
Definition: MPll.h:339
gslpp::complex F27(double q2)
The correction from .
Definition: MPll.cpp:1076
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:1183
double S_L(double q2)
The helicity form factor .
Definition: MPll.cpp:963
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:1114
double integrateDelta(int i, double q_min, double q_max)
The integral of from to .
Definition: MPll.cpp:1482
gslpp::complex F29(double q2)
The correction from .
Definition: MPll.cpp:1085
void spline_QCDF_func()
Definition: MPll.cpp:1247
double Integrand_ImTparplus(double up)
The imaginary part of the integral involving at fixed , according to .
Definition: MPll.cpp:1022
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 .