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
255
263 double integrateSigma(int i, double q_min, double q_max);
264
271 double getSigma(int i, double q_2);
272
280 double integrateDelta(int i, double q_min, double q_max);
281
288 double integrateSigmaTree(double q_min, double q_max);
289
294 double getwidth(){
295 updateParameters();
296 return width;
297 }
298
303 std::vector<std::string> initializeMPllParameters();
304
305
306private:
311 std::vector<std::string> mpllParameters;
312 std::unique_ptr<F_1> myF_1;
313 std::unique_ptr<F_2> myF_2;
314 bool dispersion = false;
315 bool FixedWCbtos = false;
320 double mJ2;
321
322 double GF;
323 double ale;
324 double Mlep;
325 double MM;
326 double MP;
327 double Mb;
328 double mu_b;
329 double mu_h;
330 double Mc;
331 double mb_pole;
332 double mc_pole;
333 double Ms;
335 double width;
336 double alpha_s_mub;
337 double angmomP;
338 int etaP;
340 double MW;
341 gslpp::complex lambda_t;
342 gslpp::complex h_0;
343 gslpp::complex h_1;
344 gslpp::complex h_2;
345 gslpp::complex r_1;
346 gslpp::complex r_2;
347 gslpp::complex Delta_C9;
348 gslpp::complex exp_Phase;
349
350 /*LATTICE fit parameters*/
351 double b_0_fplus;
352 double b_1_fplus;
353 double b_2_fplus;
354 double m_fit2_fplus_lat;
355 double b_0_fT;
356 double b_1_fT;
357 double b_2_fT;
358 double m_fit2_fT_lat;
359 double b_0_f0;
360 double b_1_f0;
361 double b_2_f0;
362 double m_fit2_f0_lat;
363 /*LCSR fit parameters*/
364 double r_1_fplus;
365 double r_2_fplus;
366 double m_fit2_fplus;
367 double r_1_fT;
368 double r_2_fT;
369 double m_fit2_fT;
370 double r_2_f0;
371 double m_fit2_f0;
374 gslpp::vector<gslpp::complex> ** allcoeff;
375 gslpp::vector<gslpp::complex> ** allcoeffh;
376 gslpp::vector<gslpp::complex> ** allcoeffprime;
378 gslpp::vector<gslpp::complex> ** allcoeff_noSM;
380 gslpp::vector<gslpp::complex> ** allcoeff_nu;
382 gslpp::vector<gslpp::complex> ** allcoeff_noSM_nu;
384 gslpp::complex C_1;
385 gslpp::complex C_1L_bar;
386 gslpp::complex C_1Lh_bar;
387 gslpp::complex C_2;
388 gslpp::complex C_2L_bar;
389 gslpp::complex C_2Lh_bar;
390 gslpp::complex C_3;
391 gslpp::complex C_4;
392 gslpp::complex C_5;
393 gslpp::complex C_6;
394 gslpp::complex C_7;
395 gslpp::complex C_8;
396 gslpp::complex C_8L;
397 gslpp::complex C_8Lh;
398 gslpp::complex C_9;
399 gslpp::complex C_10;
400 gslpp::complex C_S;
401 gslpp::complex C_P;
403 gslpp::complex C_7p;
404 gslpp::complex C_9p;
405 gslpp::complex C_10p;
406 gslpp::complex C_Sp;
407 gslpp::complex C_Pp;
409 gslpp::complex C_L_nunu;
410 gslpp::complex C_R_nunu;
412 gslpp::complex C_L_nunu_e;
413 gslpp::complex C_R_nunu_e;
414 gslpp::complex C_L_nunu_mu;
415 gslpp::complex C_R_nunu_mu;
416 gslpp::complex C_L_nunu_tau;
417 gslpp::complex C_R_nunu_tau;
419 gsl_interp_accel *acc_Re_deltaC7_QCDF;
420 gsl_interp_accel *acc_Im_deltaC7_QCDF;
421 gsl_interp_accel *acc_Re_deltaC9_QCDF;
422 gsl_interp_accel *acc_Im_deltaC9_QCDF;
423
424 gsl_spline *spline_Re_deltaC7_QCDF;
425 gsl_spline *spline_Im_deltaC7_QCDF;
426 gsl_spline *spline_Re_deltaC9_QCDF;
427 gsl_spline *spline_Im_deltaC9_QCDF;
428
429
430 std::vector<double> ReDeltaC9;
431 std::vector<double> ImDeltaC9;
432 std::vector<double> myq2;
433 TFitResultPtr refres;
434 TFitResultPtr imfres;
436 TGraph gr1;
437 TGraph gr2;
439 TF1 reffit;
440 TF1 imffit;
442 double tmpq2;
444 gslpp::complex H_0_pre;
445 gslpp::complex H_0_WC;
446 gslpp::complex H_c_WC;
447 gslpp::complex H_b_WC;
449 gslpp::complex ihalfMPI;
450 double fournineth;
451 double half;
452 double twothird;
453 double Mc2;
454 double Mb2;
455 double logMc;
456 double logMb;
457 double mu_b2;
458 double fourMc2;
459 double fourMb2;
460 double Mlep2;
461 double NN;
462 double MM2;
463 double MM4;
464 double MP2;
465 double MP4;
466 double MM2mMP2;
467 double twoMP2;
468 double twoMM;
469 double twoMM2;
470 double twoMM2_MMpMP;
471 double twoMM_MbpMs;
472 double S_L_pre;
473 double fourMM2;
474 double twoMboMM;
475 double sixteenM_PI2;
476 double ninetysixM_PI3MM3;
477 double MboMW;
478 double MboMM;
479 double MsoMb;
480 double twoMlepMb;
481 double threeGegen0;
482 double threeGegen1otwo;
483 double M_PI2osix;
484 double twoMc2;
485 double sixMMoMb;
486 double CF;
487 double deltaT_0;
488 double deltaT_1par;
490 gslpp::complex ubar;
491 gslpp::complex arg1;
492 gslpp::complex B01;
493 gslpp::complex B00;
494 gslpp::complex xp;
495 gslpp::complex xm;
496 gslpp::complex yp;
497 gslpp::complex ym;
498 gslpp::complex L1xp;
499 gslpp::complex L1xm;
500 gslpp::complex L1yp;
501 gslpp::complex L1ym;
502 gslpp::complex F87_0;
503 gslpp::complex F87_1;
504 gslpp::complex F87_2;
505 gslpp::complex F87_3;
506 gslpp::complex F29_0;
507 gslpp::complex F29_L1;
508 gslpp::complex F29_1;
509 gslpp::complex F29_2;
510 gslpp::complex F29_3;
511 gslpp::complex F29_L1_1;
512 gslpp::complex F29_L1_2;
513 gslpp::complex F29_L1_3;
514 gslpp::complex F27_0;
515 gslpp::complex F27_1;
516 gslpp::complex F27_2;
517 gslpp::complex F27_3;
518 gslpp::complex F27_L1_1;
519 gslpp::complex F27_L1_2;
520 gslpp::complex F27_L1_3;
521 double F89_0;
522 double F89_1;
523 double F89_2;
524 double F89_3;
525 double Ee;
527 //additional variables for B to K nu nu
528 const double M_PI2 = M_PI * M_PI;
529 double GF4;
530 double MM3;
531 double fM2;
532 double fP2;
533
534 double mtau;
535 double mtau2;
536 double Gammatau;
537 double VusVub_abs2;
538
539 unsigned int fplus_lat_updated;
540 gslpp::vector<double> fplus_lat_cache;
542 unsigned int fT_lat_updated;
543 gslpp::vector<double> fT_lat_cache;
545 unsigned int f0_lat_updated;
546 gslpp::vector<double> f0_lat_cache;
548 unsigned int fplus_updated;
549 gslpp::vector<double> fplus_cache;
551 unsigned int fT_updated;
552 gslpp::vector<double> fT_cache;
554 unsigned int f0_updated;
555 double f0_cache;
557 unsigned int k2_updated;
558 gslpp::vector<double> k2_cache;
560 unsigned int beta_updated;
561 double beta_cache;
563 unsigned int lambda_updated;
565 unsigned int F_updated;
567 unsigned int VL_updated;
569 unsigned int TL_updated;
571 unsigned int SL_updated;
572 gslpp::vector<double> SL_cache;
574 unsigned int N_updated;
575 gslpp::vector<double> N_cache;
576 gslpp::complex Nc_cache;
578 unsigned int C_1_updated;
579 gslpp::complex C_1_cache;
581 unsigned int C_2_updated;
582 gslpp::complex C_2_cache;
584 unsigned int C_3_updated;
585 gslpp::complex C_3_cache;
587 unsigned int C_4_updated;
588 gslpp::complex C_4_cache;
590 unsigned int C_5_updated;
591 gslpp::complex C_5_cache;
593 unsigned int C_6_updated;
594 gslpp::complex C_6_cache;
596 unsigned int C_7_updated;
597 gslpp::complex C_7_cache;
599 unsigned int C_9_updated;
600 gslpp::complex C_9_cache;
602 unsigned int C_10_updated;
603 gslpp::complex C_10_cache;
605 unsigned int C_7p_updated;
606 gslpp::complex C_7p_cache;
608 unsigned int C_9p_updated;
609 gslpp::complex C_9p_cache;
611 unsigned int C_10p_updated;
612 gslpp::complex C_10p_cache;
614 unsigned int C_S_updated;
615 gslpp::complex C_S_cache;
617 unsigned int C_P_updated;
618 gslpp::complex C_P_cache;
620 unsigned int C_Sp_updated;
621 gslpp::complex C_Sp_cache;
623 unsigned int C_Pp_updated;
624 gslpp::complex C_Pp_cache;
626 unsigned int C_2Lh_updated;
627 gslpp::complex C_2Lh_cache;
629 unsigned int C_8Lh_updated;
630 gslpp::complex C_8Lh_cache;
632 unsigned int C_L_nunu_updated;
633 gslpp::complex C_L_nunu_cache;
635 unsigned int C_R_nunu_updated;
636 gslpp::complex C_R_nunu_cache;
638 unsigned int Yupdated;
639 gslpp::vector<double> Ycache;
641 unsigned int H_V0updated;
642 gslpp::vector<double> H_V0cache;
643 gslpp::complex H_V0Ccache[3];
644 gslpp::complex H_V0Ccache_dispersion[4];
646 unsigned int H_A0updated;
648 unsigned int H_Supdated;
649 gslpp::vector<double> H_Scache;
651 unsigned int H_P_updated;
652 gslpp::vector<double> H_P_cache;
654 unsigned int I0_updated;
655 unsigned int I2_updated;
656 unsigned int I8_updated;
658 unsigned int Itree_updated;
659 gslpp::vector<double> Itree_cache;
661 std::map<std::pair<double, double>, unsigned int > I1Cached;
663 std::map<std::pair<double, double>, unsigned int > sigma0Cached;
664 std::map<std::pair<double, double>, unsigned int > sigma2Cached;
666 std::map<std::pair<double, double>, unsigned int > delta0Cached;
667 std::map<std::pair<double, double>, unsigned int > delta2Cached;
669 std::map<std::pair<double, double>, unsigned int > sigmaTreeCached;
671 double avaSigma;
672 double avaDelta;
673 double avaSigmaTree;
674 double avaDTPPR;
676 double errSigma;
677 double errDelta;
678 double errSigmaTree;
679 double errDTPPR;
681 gsl_function FS;
682 gsl_function FD;
683 gsl_function DTPPR;
685 gsl_integration_cquad_workspace * w_sigma;
686 gsl_integration_cquad_workspace * w_delta;
687 gsl_integration_cquad_workspace * w_sigmaTree;
688 gsl_integration_cquad_workspace * w_DTPPR;
690 gsl_error_handler_t * old_handler;
692 std::map<std::pair<double, double>, gslpp::complex > cacheI1;
694 std::map<std::pair<double, double>, double > cacheSigma0;
695 std::map<std::pair<double, double>, double > cacheSigma2;
697 std::map<std::pair<double, double>, double > cacheDelta0;
698 std::map<std::pair<double, double>, double > cacheDelta2;
700 std::map<std::pair<double, double>, double > cacheSigmaTree;
702 std::map<double, unsigned int> deltaTparpCached;
703 std::map<double, unsigned int> deltaTparmCached;
705 std::map<double, gslpp::complex> cacheDeltaTparp;
706 std::map<double, gslpp::complex> cacheDeltaTparm;
708 unsigned int deltaTparpupdated;
709 unsigned int deltaTparmupdated;
711 unsigned int T_updated;
712 gslpp::vector<double> T_cache;
719 void updateParameters();
720
724 void checkCache();
725
734 double LCSR_fit1(double q2, double r_1, double r_2, double m_fit2);
735
736
744 double LCSR_fit2(double q2, double r_2, double m_fit2);
745
746
756 double LCSR_fit3(double q2, double b_0, double b_1, double b_2, double m_fit2);
757
758
764 double zeta(double q2);
765
766
772 double zeta_DM(double q2);
773
780 double phiplus_DM(double q2, double m_fit2);
781
788 double phi0_DM(double q2);
789
796 double phiT_DM(double q2, double m_fit2);
797
807 double fplus_DM(double q2, double b_0, double b_1, double b_2, double m_fit2);
808
818 double f0_DM(double q2, double b_0, double b_1, double b_2);
819
829 double fT_DM(double q2, double b_0, double b_1, double b_2, double m_fit2);
830
831
841 double LATTICE_fit1(double q2, double b_0, double b_1, double b_2, double m_fit2);
842
843
853 double LATTICE_fit2(double q2, double b_0, double b_1, double b_2, double m_fit2);
854
855
861 double f_plus(double q2);
862
863
869 double f_T(double q2);
870
871
877 double f_0(double q2);
878
879
885 gslpp::complex H_0(double q2);
886
893 gslpp::complex H_c(double q2, double mu2);
894
901 gslpp::complex H_b(double q2, double mu2);
902
903
909 gslpp::complex Y(double q2);
910
911
917 double k2 (double q2);
918
919
925 double beta (double q2);
926
932 double beta2 (double q2);
933
934
940 double lambda(double q2);
941
942
948 double F(double q2);
949
955 double I_1c(double q2);
956
962 double I_2c(double q2);
963
969 double I_6c(double q2);
970
971
978 double Delta(int i, double q2);
979
985 double SigmaTree(double q2);
986
991 double getintegratedSigmaTree();
992
998 double getSigma1c(double q2)
999 {
1000 return I_1c(q2);
1001 };
1002
1008 double getSigma2c(double q2)
1009 {
1010 return I_2c(q2);
1011 };
1012
1018 double getSigma6c(double q2)
1019 {
1020 return I_6c(q2);
1021 };
1022
1028 double getDelta1c(double q2)
1029 {
1030 return Delta(0, q2);
1031 };
1032
1038 double getDelta2c(double q2)
1039 {
1040 return Delta(2, q2);
1041 };
1042
1049 gslpp::complex I1(double u, double q2);
1050
1057 gslpp::complex Tparplus(double u, double q2);
1058
1065 gslpp::complex Tparminus(double u, double q2);
1066
1072 double Integrand_ReTparplus(double up);
1073
1079 double Integrand_ImTparplus(double up);
1080
1086 double Integrand_ReTparminus(double up);
1087
1093 double Integrand_ImTparminus(double up);
1094
1100 double Integrand_ReTpar_pm(double up);
1101
1107 double Integrand_ImTpar_pm(double up);
1108
1114 gslpp::complex F19(double q2);
1115
1121 gslpp::complex F27(double q2);
1122
1128 gslpp::complex F29(double q2);
1129
1135 gslpp::complex F87(double q2);
1136
1142 double F89(double q2);
1143
1149 gslpp::complex Cpar(double q2);
1150
1156 gslpp::complex deltaTpar(double q2);
1157
1164 double reDC9fit(double* x, double* p);
1165
1172 double imDC9fit(double* x, double* p);
1173
1177 void fit_DeltaC9_mumu();
1178
1184 gslpp::complex fDeltaC9(double q2);
1185
1191 gslpp::complex DeltaC9(double q2);
1192
1198 gslpp::complex deltaC7_QCDF(double q2, bool spline = false);
1199
1205 gslpp::complex deltaC9_QCDF(double q2, bool spline = false);
1206
1207 void spline_QCDF_func();
1208
1209};
1210
1211#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:1018
double getSigma1c(double q2)
The CP average .
Definition: MPll.h:998
bool FixedWCbtos
Definition: MPll.h:315
double SigmaTree(double q2)
Definition: MPll.cpp:1499
bool MPll_GRvDV_flag
Definition: MPll.h:317
gslpp::complex H_A(double q2)
The helicity amplitude .
Definition: MPll.cpp:1306
double Mb
Definition: MPll.h:327
gslpp::complex h_lambda(double q2)
The non-pertubative ccbar contributions to the helicity amplitudes.
Definition: MPll.cpp:1281
double F89(double q2)
The correction from .
Definition: MPll.cpp:1062
double mc_pole
Definition: MPll.h:332
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:1183
gslpp::complex DeltaC9_KD(double q2)
Definition: MPll.cpp:1275
gslpp::complex H_P(double q2)
The helicity amplitude .
Definition: MPll.cpp:1319
QCD::meson pseudoscalar
Definition: MPll.h:310
gslpp::complex H_S(double q2)
The helicity amplitude .
Definition: MPll.cpp:1314
gslpp::complex F87(double q2)
The correction from .
Definition: MPll.cpp:1055
std::vector< std::string > mpllParameters
Definition: MPll.h:311
double Integrand_ImTparminus(double up)
The imaginary part of the integral involving at fixed , according to .
Definition: MPll.cpp:1002
std::unique_ptr< F_2 > myF_2
Definition: MPll.h:313
bool MPll_Lattice_flag
Definition: MPll.h:316
double Mlep
Definition: MPll.h:324
double mu_h
Definition: MPll.h:329
gslpp::complex I1(double u, double q2)
The function from .
Definition: MPll.cpp:933
double Integrand_ReTparminus(double up)
The real part of the integral involving at fixed , according to .
Definition: MPll.cpp:991
bool MPll_DM_flag
Definition: MPll.h:318
void fit_DeltaC9_mumu()
The fitting routine for the QCDF correction in the muon channel.
Definition: MPll.cpp:1112
gslpp::complex Tparminus(double u, double q2)
The function from .
Definition: MPll.cpp:968
double imDC9fit(double *x, double *p)
The fit function for the imaginary part of the QCDF correction .
Definition: MPll.cpp:1105
double getDelta2c(double q2)
The CP asymmetry .
Definition: MPll.h:1038
QCD::meson meson
Definition: MPll.h:309
double width
Definition: MPll.h:335
gslpp::complex F19(double q2)
The correction from .
Definition: MPll.cpp:1023
double integrateSigma(int i, double q_min, double q_max)
The integral of from to .
Definition: MPll.cpp:1374
gslpp::complex Cpar(double q2)
The correction from .
Definition: MPll.cpp:1069
double reDC9fit(double *x, double *p)
The fit function for the real part of the QCDF correction .
Definition: MPll.cpp:1100
double Mc
Definition: MPll.h:330
const StandardModel & mySM
Definition: MPll.h:307
gslpp::complex V_L(double q2)
The helicity form factor .
Definition: MPll.cpp:914
double getSigma(int i, double q_2)
The value of from to .
Definition: MPll.cpp:1411
QCD::lepton lep
Definition: MPll.h:308
virtual ~MPll()
Destructor.
Definition: MPll.cpp:76
gslpp::complex funct_g(double q2)
Definition: MPll.cpp:1267
double mJ2
Definition: MPll.h:320
double GF
Definition: MPll.h:322
double Integrand_ImTpar_pm(double up)
The sum of Integrand_ImTparplus() and Integrand_ImTparminus().
Definition: MPll.cpp:1013
double Integrand_ReTparplus(double up)
The real part of the integral involving at fixed , according to .
Definition: MPll.cpp:975
gslpp::complex T_L(double q2)
The helicity form factor .
Definition: MPll.cpp:919
double alpha_s_mub
Definition: MPll.h:336
double getSigma2c(double q2)
The CP average .
Definition: MPll.h:1008
gslpp::complex Tparplus(double u, double q2)
The function from .
Definition: MPll.cpp:955
double MM
Definition: MPll.h:325
std::unique_ptr< F_1 > myF_1
Definition: MPll.h:312
bool NeutrinoTree_flag
Definition: MPll.h:319
double integrateSigmaTree(double q_min, double q_max)
The integral of from to (arxiv/2301.06990)
Definition: MPll.cpp:1469
double getintegratedSigmaTree()
The integral of from 0 to .
Definition: MPll.cpp:1504
double Ms
Definition: MPll.h:333
double getDelta1c(double q2)
The CP asymmetry .
Definition: MPll.h:1028
double MP
Definition: MPll.h:326
gslpp::complex H_V(double q2)
The helicity amplitude .
Definition: MPll.cpp:1296
gslpp::complex DeltaC9(double q2)
The total QCDF correction computed integrating over .
Definition: MPll.cpp:1153
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:1158
bool dispersion
Definition: MPll.h:314
double getSigma6c(double q2)
The CP average .
Definition: MPll.h:1018
double mb_pole
Definition: MPll.h:331
gslpp::complex F27(double q2)
The correction from .
Definition: MPll.cpp:1037
double ale
Definition: MPll.h:323
double getwidth()
The width of the meson M.
Definition: MPll.h:294
double mu_b
Definition: MPll.h:328
gslpp::complex fDeltaC9(double q2)
The total QCDF correction computed fitting over .
Definition: MPll.cpp:1144
double S_L(double q2)
The helicity form factor .
Definition: MPll.cpp:924
double spectator_charge
Definition: MPll.h:334
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:1075
double integrateDelta(int i, double q_min, double q_max)
The integral of from to .
Definition: MPll.cpp:1432
gslpp::complex F29(double q2)
The correction from .
Definition: MPll.cpp:1046
void spline_QCDF_func()
Definition: MPll.cpp:1208
double Integrand_ImTparplus(double up)
The imaginary part of the integral involving at fixed , according to .
Definition: MPll.cpp:983
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 .