a Code for the Combination of Indirect and Direct Constraints on High Energy Physics Models Logo
MVll.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 MVLL_H
9#define MVLL_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 SWITCH 8.2
22#define NFPOLARBASIS_MVLL false
23#define COMPUTECP false
24#define GSL_INTERP_DIM 10
25#define GSL_INTERP_DIM_DC 10
26#define SPLINE true
27#define FULLNLOQCDF_MVll false
28
308class MVll {
309public:
310
318 MVll(const StandardModel& SM_i, QCD::meson meson_i, QCD::meson vector_i, QCD::lepton lep_i);
319
323 virtual ~MVll();
324
332 double integrateSigma(int i, double q_min, double q_max);
333
340 double getSigma(int i, double q_2);
341
349 double integrateDelta(int i, double q_min, double q_max);
350
357 double integrateSigmaTree(double q_min, double q_max);
358
363 double getwidth(){
364 updateParameters();
365 return width;
366 }
367
372 double getMlep(){
373 updateParameters();
374 return Mlep;
375 }
376
382 double beta (double q2);
383
389 double getV0(double q2)
390 {
391 updateParameters();
392 return (2. * MM * sqrt(q2))/sqrt(lambda(q2)) * V_0t(q2);
393 };
394
400 double getVp(double q2)
401 {
402 updateParameters();
403 return V_p(q2);
404 };
405
411 double getVm(double q2)
412 {
413 return V_m(q2);
414 };
415
421 double getT0(double q2)
422 {
423 updateParameters();
424 return twoMM3/sqrt(q2 * lambda(q2)) * T_0t(q2);
425 };
426
432 double getTp(double q2)
433 {
434 updateParameters();
435 return T_p(q2);
436 };
437
443 double getTm(double q2)
444 {
445 updateParameters();
446 return T_m(q2);
447 };
448
454 double getS(double q2)
455 {
456 updateParameters();
457 return S_L_pre/sqrt(lambda(q2)) * S_L(q2);
458 };
459
465 {
466 updateParameters();
468 };
469
475 {
476 updateParameters();
477 return unitarity_bound_g;
478 };
479
485 {
486 updateParameters();
487 return unitarity_bound_F2;
488 };
489
495 {
496 updateParameters();
497 return unitarity_bound_T1;
498 };
499
505 {
506 updateParameters();
508 };
509
516 gslpp::complex h_lambda(int hel, double q2);
517
523 double Delta_C9_zExp(int hel);
524
530 {
531 updateParameters();
532 return Delta_C9_zExp(0);
533 };
534
540 {
541 updateParameters();
542 return Delta_C9_zExp(1);
543 };
544
550 {
551 updateParameters();
552 return Delta_C9_zExp(2);
553 };
554
560 gslpp::complex H_V_0(double q2, bool bar);
561
567 gslpp::complex H_V_p(double q2, bool bar);
568
574 gslpp::complex H_V_m(double q2, bool bar);
575
581 gslpp::complex H_A_0(double q2, bool bar);
582
588 gslpp::complex H_A_p(double q2, bool bar);
589
595 gslpp::complex H_A_m(double q2, bool bar);
596
602 gslpp::complex H_S(double q2, bool bar);
603
609 gslpp::complex H_P(double q2, bool bar);
610
617 gslpp::complex AmpMVpsi_zExpansion(double mpsi, int tran);
618
619 gslpp::complex getQCDf_1(double q2)
620 {
621 updateParameters();
622// return (gtilde_1_pre/(sqrt(lambda(q2)) * V(q2)) * 1./(16. * M_PI * M_PI * MM*MM) * (fDeltaC9_m(q2) * V_m(q2) - fDeltaC9_p(q2) * V_p(q2)));
623 return 0.;
624 }
625
626 gslpp::complex getQCDf_2(double q2)
627 {
628 updateParameters();
629// return (gtilde_2_pre/A_1(q2) * 1./(16. * M_PI * M_PI * MM*MM) * (fDeltaC9_m(q2) * V_m(q2) + fDeltaC9_p(q2) * V_p(q2)));
630 return 0.;
631 }
632
633 gslpp::complex getQCDf_3(double q2)
634 {
635 updateParameters();
636// return (gtilde_3_pre/(lambda(q2) * A_2(q2)) * (sqrt(q2) * 1./(16. * M_PI * M_PI * MM*MM) * fDeltaC9_0(q2) * V_0t(q2) - (MM2mMV2 - q2)/(4.*MV) * 1./(16. * M_PI * M_PI * MM*MM) * (fDeltaC9_m(q2) * V_m(q2) + fDeltaC9_p(q2) * V_p(q2))));
637 return 0.;
638 }
639
640 double getQCDfC9_1(double q2, double cutoff)
641 {
642 updateParameters();
643 return (getQCDf_1(q2) - getQCDf_1(cutoff) * cutoff/q2).abs();
644 }
645
646 double getQCDfC9_2(double q2, double cutoff)
647 {
648 updateParameters();
649 return (getQCDf_2(q2) - getQCDf_2(cutoff) * cutoff/q2).abs();
650 }
651
652 double getQCDfC9_3(double q2, double cutoff)
653 {
654 updateParameters();
655 return (getQCDf_3(q2) - getQCDf_3(cutoff) * cutoff/q2).abs();
656 }
657
658 double getQCDfC9p_1(double cutoff)
659 {
660 updateParameters();
661 return (getQCDf_1(cutoff) * cutoff).abs();
662 }
663
664 double getQCDfC9p_2(double cutoff)
665 {
666 updateParameters();
667 return (getQCDf_2(cutoff) * cutoff).abs();
668 }
669
670 double getQCDfC9p_3(double cutoff)
671 {
672 updateParameters();
673 return (getQCDf_3(cutoff) * cutoff).abs();
674 }
675
681 double getgtilde_1_re(double q2)
682 {
683 updateParameters();
684 return C2_inv * (gtilde_1_pre/(sqrt(lambda(q2)) * V(q2)) * (h_lambda(2,q2)-h_lambda(1,q2))).real()/q2;
685 }
686
692 double getgtilde_1_im(double q2)
693 {
694 updateParameters();
695 return C2_inv * (gtilde_1_pre/(sqrt(lambda(q2)) * V(q2)) * (h_lambda(2,q2)-h_lambda(1,q2))).imag()/q2;
696 }
697
703 double getgtilde_2_re(double q2)
704 {
705 updateParameters();
706 return C2_inv * (gtilde_2_pre/A_1(q2) * (h_lambda(1,q2)+h_lambda(2,q2))).real()/q2;
707 }
708
714 double getgtilde_2_im(double q2)
715 {
716 updateParameters();
717 return C2_inv * (gtilde_2_pre/A_1(q2) * (h_lambda(1,q2)+h_lambda(2,q2))).imag()/q2;
718 }
719
725 double getgtilde_3_re(double q2)
726 {
727 updateParameters();
728 return C2_inv * (gtilde_3_pre/(lambda(q2) * A_2(q2)) * (sqrt(q2)*h_lambda(0,q2)/q2-
729 (MM2mMV2 - q2)/(4.*MV) * (h_lambda(1,q2)+h_lambda(2,q2))/q2)).real();
730 }
731
737 double getgtilde_3_im(double q2)
738 {
739 updateParameters();
740 return C2_inv * (gtilde_3_pre/(lambda(q2) * A_2(q2)) * (sqrt(q2)*h_lambda(0,q2)/q2-
741 (MM2mMV2 - q2)/(4.*MV) * (h_lambda(1,q2)+h_lambda(2,q2))/q2)).imag();
742 }
743
749 double geth_0_re(double q2)
750 {
751 return (sixteenM_PI2MM2 * h_lambda(0,q2)/q2).real();
752 }
753
759 double geth_0_im(double q2)
760 {
761 return (sixteenM_PI2MM2 * h_lambda(0,q2)/q2).imag();
762 }
763
768 gslpp::complex geth_p_0()
769 {
770 return h_lambda(1,0.);
771 }
772
778 double geth_p_re(double q2)
779 {
780 return (sixteenM_PI2MM2 * h_lambda(1,q2)/q2).real();
781 }
782
788 double geth_p_im(double q2)
789 {
790 return (sixteenM_PI2MM2 * h_lambda(1,q2)/q2).imag();
791 }
792
797 gslpp::complex geth_m_0()
798 {
799 return h_lambda(2,0.);
800 }
801
807 double geth_m_re(double q2)
808 {
809 return (sixteenM_PI2MM2 * h_lambda(2,q2)/q2).real();
810 }
811
817 double geth_m_im(double q2)
818 {
819 return (sixteenM_PI2MM2 * h_lambda(2,q2)/q2).imag();
820 }
821
826 std::vector<std::string> initializeMVllParameters();
827
828private:
833 std::vector<std::string> mvllParameters;
834 std::unique_ptr<F_1> myF_1;
835 std::unique_ptr<F_2> myF_2;
841 double mJpsi, mJ2;
843 double mD2;
844 gslpp::complex exp_Phase[3];
845
846 double GF;
847 double ale;
848 double Mlep;
849 double MM;
850 double MV;
851 double Mb;
852 double mu_b;
853 double mu_h;
854 double Mc;
855 double mb_pole;
856 double mc_pole;
857 double Ms;
859 double width;
860 double ys;
861 double xs;
862 double angmomV;
863 int etaV;
864 double alpha_s_mub;
865 double fB;
866 double fpara;
867 double fperp;
868
869 double MW;
870 gslpp::complex lambda_t;
871 gslpp::complex lambda_u;
872 double b;
873 gslpp::complex h_0[3];
874 gslpp::complex h_1[3];
875 gslpp::complex h_2[3];
876 gslpp::complex SU3_breaking;
877
878 gslpp::complex beta_0[7];
879 gslpp::complex beta_1[7];
880 gslpp::complex beta_2[7];
882 double t_p;
883 double t_m;
884 double t_0;
885 double z_0;
886 double s_p;
887 double s_0;
888 double Q2;
889 double chiOPE;
890 double twoalphaBtoKst;
891 double rho_0;
892 double rho_1;
893 double rho_2;
894 double rho_3;
895 double rho_4;
896 double rho_5;
897 double onemrho_0_2;
898 double onemrho_1_2;
899 double onemrho_2_2;
900 double onemrho_3_2;
901 double onemrho_4_2;
902 double onemrho_5_2;
903 double DeltaC9;
904 double DeltaC10;
905 double MMpMV;
906 double MMpMV2;
907 double MMmMV;
908 double MMmMV2;
909 double rV;
910 double Chi1minus;
911 double Chi1plus;
912 double Chi0plus;
913 double Chi0minus;
914 double ChiTT;
915 double ChiBB;
916 double MM2;
917 double MM4;
918 double MV2;
919 double MV4;
920 double MMMV;
921 double MM2mMV2;
922 double MM2pMV2;
923 double fourMV;
924 double onepMMoMV;
925 double MM_MMpMV;
926 double twoMM2;
927 double twoMV2;
928 double twoMM_mbpms;
929 double fourMM2;
930 double Mlep2;
931 double twoMlepMb;
932 double MboMW;
933 double MsoMb;
934 double M_PI2osix;
935 double N_QCDF;
936 double twoMM;
937 double ninetysixM_PI3MM3;
938 double M_PI2;
939 double sixteenM_PI2;
940 double sixteenM_PI2MM2;
941 double twoMboMM;
942 gslpp::complex H_0_pre;
943 double mu_b2;
944 double Mc2;
945 double Mb2;
946 double fourMc2;
947 double fourMb2;
948 double logMc;
949 double logMb;
950 gslpp::complex H_0_WC;
951 gslpp::complex H_c_WC;
952 gslpp::complex H_b_WC;
953 double fournineth;
954 double half;
955 double twothird;
956 double sqrt3;
957 gslpp::complex ihalfMPI;
958 double twoMM3;
959 double gtilde_1_pre;
960 double gtilde_2_pre;
961 double gtilde_3_pre;
962 double C2_inv;
963 double S_L_pre;
964 gslpp::complex NN;
965 gslpp::complex NN_conjugate;
966 double CF;
967 double deltaT_0;
968 double deltaT_1par;
969 double deltaT_1perp;
971 bool h_pole;
972
973 gslpp::complex ubar;
974 gslpp::complex arg1;
975 gslpp::complex B01;
976 gslpp::complex B00;
977 gslpp::complex xp;
978 gslpp::complex xm;
979 gslpp::complex yp;
980 gslpp::complex ym;
981 gslpp::complex L1xp;
982 gslpp::complex L1xm;
983 gslpp::complex L1yp;
984 gslpp::complex L1ym;
985 gslpp::complex F87_0;
986 gslpp::complex F87_1;
987 gslpp::complex F87_2;
988 gslpp::complex F87_3;
989 double F89_0;
990 double F89_1;
991 double F89_2;
992 double F89_3;
994 double a_0V;
995 double a_1V;
996 double a_2V;
997 double MRV_2;
998 double a_0A0;
999 double a_1A0;
1000 double a_2A0;
1001 double MRA0_2;
1002 double a_0A1;
1003 double a_1A1;
1004 double a_2A1;
1005 double MRA1_2;
1006 double a_0A12;
1007 double a_1A12;
1008 double a_2A12;
1009 double MRA12_2;
1010 double a_0T1;
1011 double a_1T1;
1012 double a_2T1;
1013 double MRT1_2;
1014 double a_0T2;
1015 double a_1T2;
1016 double a_2T2;
1017 double MRT2_2;
1018 double a_0T23;
1019 double a_1T23;
1020 double a_2T23;
1021 double MRT23_2;
1023 double a_0f;
1024 double a_1f;
1025 double a_2f;
1026 double MRf_2;
1027 double a_0g;
1028 double a_1g;
1029 double a_2g;
1030 double MRg_2;
1031 double a_0F1;
1032 double a_1F1;
1033 double a_2F1;
1034 double MRF1_2;
1035 double a_0F2;
1036 double a_1F2;
1037 double a_2F2;
1038 double MRF2_2;
1039 double a_0T0;
1040 double a_1T0;
1041 double a_2T0;
1042 double MRT0_2;
1044 double unitarity_bound_f_F1;
1045 double unitarity_bound_g;
1046 double unitarity_bound_F2;
1047 double unitarity_bound_T1;
1048 double unitarity_bound_T2_T0;
1050 //additional variables for B to K nu nu
1051 double GF4;
1052 double MM3;
1053 double fM2;
1054 double fV2;
1055
1056 double mtau;
1057 double mtau2;
1058 double Gammatau;
1059 double VusVub_abs2;
1060
1061 gslpp::vector<gslpp::complex> ** allcoeff;
1062 gslpp::vector<gslpp::complex> ** allcoeffh;
1063 gslpp::vector<gslpp::complex> ** allcoeffprime;
1065 gslpp::vector<gslpp::complex> ** allcoeff_noSM;
1067 gslpp::vector<gslpp::complex> ** allcoeff_nu;
1069 gslpp::vector<gslpp::complex> ** allcoeff_noSM_nu;
1071 gslpp::complex C_1;
1072 gslpp::complex C_1L_bar;
1073 gslpp::complex C_1Lh_bar;
1074 gslpp::complex C_2;
1075 gslpp::complex C_2L_bar;
1076 gslpp::complex C_2Lh_bar;
1077 gslpp::complex C_3;
1078 gslpp::complex C_4;
1079 gslpp::complex C_5;
1080 gslpp::complex C_6;
1081 gslpp::complex C_7;
1082 gslpp::complex C_8;
1083 gslpp::complex C_8L;
1084 gslpp::complex C_8Lh;
1085 gslpp::complex C_9;
1086 gslpp::complex C_10;
1087 gslpp::complex C_S;
1088 gslpp::complex C_P;
1090 gslpp::complex C_7p;
1091 gslpp::complex C_9p;
1092 gslpp::complex C_10p;
1093 gslpp::complex C_Sp;
1094 gslpp::complex C_Pp;
1096 gslpp::complex C_L_nunu;
1097 gslpp::complex C_R_nunu;
1099 gslpp::complex C_L_nunu_e;
1100 gslpp::complex C_R_nunu_e;
1101 gslpp::complex C_L_nunu_mu;
1102 gslpp::complex C_R_nunu_mu;
1103 gslpp::complex C_L_nunu_tau;
1104 gslpp::complex C_R_nunu_tau;
1106 std::vector<double> Re_T_perp;
1107 std::vector<double> Im_T_perp;
1108 std::vector<double> Re_T_para;
1109 std::vector<double> Im_T_para;
1111 gsl_interp_accel *acc_Re_T_perp;
1112 gsl_interp_accel *acc_Im_T_perp;
1113 gsl_interp_accel *acc_Re_T_para;
1114 gsl_interp_accel *acc_Im_T_para;
1115
1116 gsl_spline *spline_Re_T_perp;
1117 gsl_spline *spline_Im_T_perp;
1118 gsl_spline *spline_Re_T_para;
1119 gsl_spline *spline_Im_T_para;
1120
1121 gsl_interp_accel *acc_Re_deltaC7_QCDF;
1122 gsl_interp_accel *acc_Im_deltaC7_QCDF;
1123 gsl_interp_accel *acc_Re_deltaC9_QCDF;
1124 gsl_interp_accel *acc_Im_deltaC9_QCDF;
1125
1126 gsl_spline *spline_Re_deltaC7_QCDF;
1127 gsl_spline *spline_Im_deltaC7_QCDF;
1128 gsl_spline *spline_Re_deltaC9_QCDF;
1129 gsl_spline *spline_Im_deltaC9_QCDF;
1130
1131#if COMPUTECP
1132 std::vector<double> Re_T_perp_conj;
1133 std::vector<double> Im_T_perp_conj;
1134 std::vector<double> Re_T_para_conj;
1135 std::vector<double> Im_T_para_conj;
1137 gsl_interp_accel *acc_Re_T_perp_conj;
1138 gsl_interp_accel *acc_Im_T_perp_conj;
1139 gsl_interp_accel *acc_Re_T_para_conj;
1140 gsl_interp_accel *acc_Im_T_para_conj;
1141
1142 gsl_interp_accel *acc_Re_deltaC7_QCDF_conj;
1143 gsl_interp_accel *acc_Im_deltaC7_QCDF_conj;
1144 gsl_interp_accel *acc_Re_deltaC9_QCDF_conj;
1145 gsl_interp_accel *acc_Im_deltaC9_QCDF_conj;
1146
1147 gsl_spline *spline_Re_T_perp_conj;
1148 gsl_spline *spline_Im_T_perp_conj;
1149 gsl_spline *spline_Re_T_para_conj;
1150 gsl_spline *spline_Im_T_para_conj;
1151
1152 gsl_spline *spline_Re_deltaC7_QCDF_conj;
1153 gsl_spline *spline_Im_deltaC7_QCDF_conj;
1154 gsl_spline *spline_Re_deltaC9_QCDF_conj;
1155 gsl_spline *spline_Im_deltaC9_QCDF_conj;
1156#endif
1157
1158 std::vector<double> myq2;
1160 TFitResultPtr Re_T_perp_res;
1161 TFitResultPtr Im_T_perp_res;
1162 TFitResultPtr Re_T_para_res;
1163 TFitResultPtr Im_T_para_res;
1165 TFitResultPtr Re_T_perp_res_conj;
1166 TFitResultPtr Im_T_perp_res_conj;
1167 TFitResultPtr Re_T_para_res_conj;
1168 TFitResultPtr Im_T_para_res_conj;
1170 TGraph gr1;
1171 TGraph gr2;
1173 TF1 QCDFfit;
1175 TF1 reffit;
1176 TF1 imffit;
1178 double avaSigma;
1179 double avaDelta;
1180 double avaSigmaTree;
1182 double errSigma;
1183 double errDelta;
1184 double errSigmaTree;
1186 gsl_function FS;
1187 gsl_function FD;
1189 gsl_integration_cquad_workspace * w_sigma;
1190 gsl_integration_cquad_workspace * w_delta;
1191 gsl_integration_cquad_workspace * w_sigmaTree;
1193 gsl_error_handler_t * old_handler;
1195 std::map<std::pair<double, double>, gslpp::complex > cacheI1;
1197 std::map<std::pair<double, double>, double > cacheSigma0;
1198 std::map<std::pair<double, double>, double > cacheSigma1;
1199 std::map<std::pair<double, double>, double > cacheSigma2;
1200 std::map<std::pair<double, double>, double > cacheSigma3;
1201 std::map<std::pair<double, double>, double > cacheSigma4;
1202 std::map<std::pair<double, double>, double > cacheSigma5;
1203 std::map<std::pair<double, double>, double > cacheSigma6;
1204 std::map<std::pair<double, double>, double > cacheSigma7;
1205 std::map<std::pair<double, double>, double > cacheSigma8;
1206 std::map<std::pair<double, double>, double > cacheSigma9;
1207 std::map<std::pair<double, double>, double > cacheSigma10;
1208 std::map<std::pair<double, double>, double > cacheSigma11;
1210 std::map<std::pair<double, double>, double > cacheDelta0;
1211 std::map<std::pair<double, double>, double > cacheDelta1;
1212 std::map<std::pair<double, double>, double > cacheDelta2;
1213 std::map<std::pair<double, double>, double > cacheDelta3;
1214 std::map<std::pair<double, double>, double > cacheDelta6;
1215 std::map<std::pair<double, double>, double > cacheDelta7;
1216 std::map<std::pair<double, double>, double > cacheDelta8;
1217 std::map<std::pair<double, double>, double > cacheDelta10;
1218 std::map<std::pair<double, double>, double > cacheDelta11;
1220 std::map<std::pair<double, double>, double > cacheSigmaTree;
1222 unsigned int N_updated;
1223 gslpp::vector<double> N_cache;
1224 gslpp::complex Nc_cache;
1226 unsigned int V_updated;
1227 gslpp::vector<double> V_cache;
1229 unsigned int A0_updated;
1230 gslpp::vector<double> A0_cache;
1232 unsigned int A1_updated;
1233 gslpp::vector<double> A1_cache;
1235 unsigned int T1_updated;
1236 gslpp::vector<double> T1_cache;
1238 unsigned int T2_updated;
1239 gslpp::vector<double> T2_cache;
1241 unsigned int k2_updated;
1242 gslpp::vector<double> k2_cache;
1244 unsigned int z_updated;
1246 unsigned int lambda_updated;
1248 unsigned int beta_updated;
1249 double beta_cache;
1251 unsigned int F_updated;
1253 unsigned int VL1_updated;
1254 unsigned int VL2_updated;
1256 unsigned int TL1_updated;
1257 unsigned int TL2_updated;
1259 unsigned int VR1_updated;
1260 unsigned int VR2_updated;
1262 unsigned int TR1_updated;
1263 unsigned int TR2_updated;
1265 unsigned int VL0_updated;
1266 gslpp::vector<double> VL0_cache;
1268 unsigned int TL0_updated;
1269 gslpp::vector<double> TL0_cache;
1271 unsigned int VR0_updated;
1273 unsigned int TR0_updated;
1275 unsigned int Mb_Ms_updated;
1277 unsigned int SL_updated;
1278 gslpp::vector<double> SL_cache;
1280 unsigned int SR_updated;
1282 unsigned int C_1_updated;
1283 gslpp::complex C_1_cache;
1285 unsigned int C_2_updated;
1286 gslpp::complex C_2_cache;
1288 unsigned int C_3_updated;
1289 gslpp::complex C_3_cache;
1291 unsigned int C_4_updated;
1292 gslpp::complex C_4_cache;
1294 unsigned int C_5_updated;
1295 gslpp::complex C_5_cache;
1297 unsigned int C_6_updated;
1298 gslpp::complex C_6_cache;
1300 unsigned int C_7_updated;
1301 gslpp::complex C_7_cache;
1303 unsigned int C_9_updated;
1304 gslpp::complex C_9_cache;
1306 unsigned int C_10_updated;
1307 gslpp::complex C_10_cache;
1309 unsigned int C_7p_updated;
1310 gslpp::complex C_7p_cache;
1312 unsigned int C_9p_updated;
1313 gslpp::complex C_9p_cache;
1315 unsigned int C_10p_updated;
1316 gslpp::complex C_10p_cache;
1318 unsigned int C_S_updated;
1319 gslpp::complex C_S_cache;
1321 unsigned int C_P_updated;
1322 gslpp::complex C_P_cache;
1324 unsigned int C_Sp_updated;
1325 gslpp::complex C_Sp_cache;
1327 unsigned int C_Pp_updated;
1328 gslpp::complex C_Pp_cache;
1330 unsigned int C_2Lh_updated;
1331 gslpp::complex C_2Lh_cache;
1333 unsigned int C_8Lh_updated;
1334 gslpp::complex C_8Lh_cache;
1336 unsigned int C_L_nunu_updated;
1337 gslpp::complex C_L_nunu_cache;
1339 unsigned int C_R_nunu_updated;
1340 gslpp::complex C_R_nunu_cache;
1342 unsigned int Yupdated;
1343 gslpp::vector<double> Ycache;
1345 gslpp::complex h0Ccache[4];
1346 gslpp::complex h1Ccache[4];
1347 gslpp::complex h2Ccache[4];
1349 gslpp::complex beta0Ccache[8];
1350 gslpp::complex beta1Ccache[8];
1351 gslpp::complex beta2Ccache[8];
1353 unsigned int h0_updated;
1354 unsigned int h1_updated;
1355 unsigned int h2_updated;
1357 unsigned int H_V0updated;
1358 gslpp::vector<double> H_V0cache;
1360 unsigned int H_V1updated;
1361 gslpp::vector<double> H_V1cache;
1363 unsigned int H_V2updated;
1364 gslpp::vector<double> H_V2cache;
1366 unsigned int H_A0updated;
1367 unsigned int H_A1updated;
1368 unsigned int H_A2updated;
1370 unsigned int H_Supdated;
1371 gslpp::vector<double> H_Scache;
1373 unsigned int H_Pupdated;
1374 gslpp::vector<double> H_Pcache;
1376 unsigned int I0_updated;
1377 unsigned int I1_updated;
1378 unsigned int I2_updated;
1379 unsigned int I3_updated;
1380 unsigned int I4_updated;
1381 unsigned int I5_updated;
1382 unsigned int I6_updated;
1383 unsigned int I7_updated;
1384 unsigned int I8_updated;
1385 unsigned int I9_updated;
1386 unsigned int I10_updated;
1387 unsigned int I11_updated;
1389 unsigned int Itree_updated;
1390 gslpp::vector<double> Itree_cache;
1392 std::map<std::pair<double, double>, unsigned int > I1Cached;
1394 std::map<std::pair<double, double>, unsigned int > sigma0Cached;
1395 std::map<std::pair<double, double>, unsigned int > sigma1Cached;
1396 std::map<std::pair<double, double>, unsigned int > sigma2Cached;
1397 std::map<std::pair<double, double>, unsigned int > sigma3Cached;
1398 std::map<std::pair<double, double>, unsigned int > sigma4Cached;
1399 std::map<std::pair<double, double>, unsigned int > sigma5Cached;
1400 std::map<std::pair<double, double>, unsigned int > sigma6Cached;
1401 std::map<std::pair<double, double>, unsigned int > sigma7Cached;
1402 std::map<std::pair<double, double>, unsigned int > sigma8Cached;
1403 std::map<std::pair<double, double>, unsigned int > sigma9Cached;
1404 std::map<std::pair<double, double>, unsigned int > sigma10Cached;
1405 std::map<std::pair<double, double>, unsigned int > sigma11Cached;
1407 std::map<std::pair<double, double>, unsigned int > delta0Cached;
1408 std::map<std::pair<double, double>, unsigned int > delta1Cached;
1409 std::map<std::pair<double, double>, unsigned int > delta2Cached;
1410 std::map<std::pair<double, double>, unsigned int > delta3Cached;
1411 std::map<std::pair<double, double>, unsigned int > delta6Cached;
1412 std::map<std::pair<double, double>, unsigned int > delta7Cached;
1413 std::map<std::pair<double, double>, unsigned int > delta8Cached;
1414 std::map<std::pair<double, double>, unsigned int > delta10Cached;
1415 std::map<std::pair<double, double>, unsigned int > delta11Cached;
1417 std::map<std::pair<double, double>, unsigned int > sigmaTreeCached;
1419 std::map<double, unsigned int> deltaTparpCached;
1420 std::map<double, unsigned int> deltaTparmCached;
1421 std::map<double, unsigned int> deltaTperpCached;
1423 std::map<double, gslpp::complex> cacheDeltaTparp;
1424 std::map<double, gslpp::complex> cacheDeltaTparm;
1425 std::map<double, gslpp::complex> cacheDeltaTperp;
1427 unsigned int deltaTparpupdated;
1428 unsigned int deltaTparmupdated;
1429 unsigned int deltaTperpupdated;
1431 unsigned int T_updated;
1432 gslpp::vector<double> T_cache;
1437 void updateParameters();
1438
1442 void checkCache();
1443
1449 double z(double q2);
1450
1456 double z_DM(double q2);
1457
1464 double phi_f(double q2, double MRf_2);
1465
1472 double phi_g(double q2, double MRg_2);
1473
1480 double phi_F1(double q2, double MRF1_2);
1481
1488 double phi_F2(double q2, double MRF2_2);
1489
1496 double phi_T0(double q2, double MRT0_2);
1497
1504 double phi_T1(double q2, double MRT1_2);
1505
1512 double phi_T2(double q2, double MRT2_2);
1513
1523 double f_DM(double q2, double a_0f, double a_1f, double a_2f, double MRf_2);
1524
1534 double g_DM(double q2, double a_0g, double a_1g, double a_2g, double MRg_2);
1535
1545 double F1_DM(double q2, double a_0F1, double a_1F1, double a_2F1, double MRF1_2);
1546
1556 double F2_DM(double q2, double a_0F2, double a_1F2, double a_2F2, double MRF2_2);
1557
1567 double T0_DM(double q2, double a_0T0, double a_1T0, double a_2T0, double MRT0_2);
1568
1578 double T1_DM(double q2, double a_0T1, double a_1T1, double a_2T1, double MRT1_2);
1579
1589 double T2_DM(double q2, double a_0T2, double a_1T2, double a_2T2, double MRT2_2);
1590
1596 double V(double q2);
1597
1598
1604 double A_0(double q2);
1605
1606
1612 double A_1(double q2);
1613
1619 double A_2(double q2);
1620
1626 double T_1(double q2);
1627
1628
1634 double T_2(double q2);
1635
1641 double V_0t(double q2);
1642
1648 double V_p(double q2);
1649
1655 double V_m(double q2);
1656
1662 double T_0t(double q2);
1663
1669 double T_p(double q2);
1670
1676 double T_m(double q2);
1677
1683 double S_L(double q2);
1684
1690 gslpp::complex H_0(double q2);
1691
1699 gslpp::complex H(double q2, double m2, double mu2);
1700
1706 gslpp::complex Y(double q2);
1707
1708 gslpp::complex funct_g(double q2);
1709
1710 gslpp::complex DeltaC9_KD(double q2, int com);
1711
1717 gslpp::complex zh(double q2);
1718
1724 gslpp::complex P(double q2);
1725
1731 gslpp::complex Phi_1(double q2);
1732
1738 gslpp::complex Phi_1_st(double q2);
1739
1745 gslpp::complex Phi_2(double q2);
1746
1752 gslpp::complex Phi_2_st(double q2);
1753
1759 gslpp::complex Phi_3(double q2);
1760
1766 gslpp::complex Phi_3_st(double q2);
1767
1773 gslpp::complex Phi_4(double q2);
1774
1780 gslpp::complex Phi_4_st(double q2);
1781
1787 gslpp::complex Phi_5(double q2);
1788
1794 gslpp::complex Phi_5_st(double q2);
1795
1801 gslpp::complex Phi_6(double q2);
1802
1808 gslpp::complex Phi_6_st(double q2);
1809
1814 gslpp::complex p0();
1815
1821 gslpp::complex p1(double q2);
1822
1828 gslpp::complex p2(double q2);
1829
1835 gslpp::complex p3(double q2);
1836
1842 gslpp::complex p4(double q2);
1843
1849 gslpp::complex p5(double q2);
1850
1856 gslpp::complex p6(double q2);
1857
1863 gslpp::complex phi_1(double q2);
1864
1870 gslpp::complex phi_2(double q2);
1871
1877 gslpp::complex phi_3(double q2);
1878
1884 gslpp::complex phi_4(double q2);
1885
1892 gslpp::complex DeltaC9_zExpansion(double q2, int tran);
1893
1899 double k2 (double q2);
1900
1906 double beta2 (double q2);
1907
1913 double lambda(double q2);
1914
1921 double F(double q2, double b_i);
1922
1929 double I_1c(double q2, bool bar);
1930
1937 double I_1s(double q2, bool bar);
1938
1945 double I_2c(double q2, bool bar);
1946
1953 double I_2s(double q2, bool bar);
1954
1961 double I_3(double q2, bool bar);
1962
1969 double I_4(double q2, bool bar);
1970
1977 double I_5(double q2, bool bar);
1978
1985 double I_6c(double q2, bool bar);
1986
1993 double I_6s(double q2, bool bar);
1994
2001 double I_7(double q2, bool bar);
2002
2009 double I_8(double q2, bool bar);
2010
2017 double I_9(double q2, bool bar);
2018
2025 double h_1s(double q2, bool bar);
2026
2033 double h_1c(double q2, bool bar);
2034
2041 double h_2s(double q2, bool bar);
2042
2049 double h_2c(double q2, bool bar);
2050
2057 double h_3(double q2, bool bar);
2058
2065 double h_4(double q2, bool bar);
2066
2073 double h_7(double q2, bool bar);
2074
2081 double s_5(double q2, bool bar);
2082
2089 double s_6s(double q2, bool bar);
2090
2097 double s_6c(double q2, bool bar);
2098
2105 double s_8(double q2, bool bar);
2106
2113 double s_9(double q2, bool bar);
2114
2120 double getSigma1c(double q2)
2121 {
2122 switch(vectorM){
2123 case QCD::K_star:
2124 return (I_1c(q2, 0) + I_1c(q2, 1))/2.;
2125 break;
2126 case QCD::K_star_P:
2127 return (I_1c(q2, 0) + I_1c(q2, 1))/2.;
2128 break;
2129 case QCD::PHI:
2130 return (I_1c(q2, 0) + I_1c(q2, 1) - ys * h_1c(q2, 0) )/2.;
2131 break;
2132 default:
2133 std::stringstream out;
2134 out << vectorM;
2135 throw std::runtime_error("MVll::getSigma1c : vector " + out.str() + " not implemented");
2136 }
2137 };
2138
2144 double getSigma1s(double q2)
2145 {
2146 switch(vectorM){
2147 case QCD::K_star:
2148 return (I_1s(q2, 0) + I_1s(q2, 1))/2.;
2149 break;
2150 case QCD::K_star_P:
2151 return (I_1s(q2, 0) + I_1s(q2, 1))/2.;
2152 break;
2153 case QCD::PHI:
2154 return (I_1s(q2, 0) + I_1s(q2, 1) - ys * h_1s(q2, 0))/2.;
2155 break;
2156 default:
2157 std::stringstream out;
2158 out << vectorM;
2159 throw std::runtime_error("MVll::getSigma1s : vector " + out.str() + " not implemented");
2160 }
2161 };
2162
2168 double getSigma2c(double q2)
2169 {
2170 switch(vectorM){
2171 case QCD::K_star:
2172 return (I_2c(q2, 0) + I_2c(q2, 1))/2.;
2173 break;
2174 case QCD::K_star_P:
2175 return (I_2c(q2, 0) + I_2c(q2, 1))/2.;
2176 break;
2177 case QCD::PHI:
2178 return (I_2c(q2, 0) + I_2c(q2, 1) - ys * h_2c(q2, 0))/2.;
2179 break;
2180 default:
2181 std::stringstream out;
2182 out << vectorM;
2183 throw std::runtime_error("MVll::getSigma2c : vector " + out.str() + " not implemented");
2184 }
2185 };
2186
2192 double getSigma2s(double q2)
2193 {
2194 switch(vectorM){
2195 case QCD::K_star:
2196 return (I_2s(q2, 0) + I_2s(q2, 1))/2.;
2197 break;
2198 case QCD::K_star_P:
2199 return (I_2s(q2, 0) + I_2s(q2, 1))/2.;
2200 break;
2201 case QCD::PHI:
2202 return (I_2s(q2, 0) + I_2s(q2, 1) - ys * h_2s(q2, 0))/2.;
2203 break;
2204 default:
2205 std::stringstream out;
2206 out << vectorM;
2207 throw std::runtime_error("MVll::getSigma2s : vector " + out.str() + " not implemented");
2208 }
2209 };
2210
2216 double getSigma3(double q2)
2217 {
2218 switch(vectorM){
2219 case QCD::K_star:
2220 return (I_3(q2, 0) + I_3(q2, 1))/2.;
2221 break;
2222 case QCD::K_star_P:
2223 return (I_3(q2, 0) + I_3(q2, 1))/2.;
2224 break;
2225 case QCD::PHI:
2226 return (I_3(q2, 0) + I_3(q2, 1) - ys * h_3(q2, 0))/2.;
2227 break;
2228 default:
2229 std::stringstream out;
2230 out << vectorM;
2231 throw std::runtime_error("MVll::getSigma3 : vector " + out.str() + " not implemented");
2232 }
2233 };
2234
2240 double getSigma4(double q2)
2241 {
2242 switch(vectorM){
2243 case QCD::K_star:
2244 return (I_4(q2, 0) + I_4(q2, 1))/2.;
2245 break;
2246 case QCD::K_star_P:
2247 return (I_4(q2, 0) + I_4(q2, 1))/2.;
2248 break;
2249 case QCD::PHI:
2250 return (I_4(q2, 0) + I_4(q2, 1) - ys * h_4(q2, 0))/2.;
2251 break;
2252 default:
2253 std::stringstream out;
2254 out << vectorM;
2255 throw std::runtime_error("MVll::getSigma4 : vector " + out.str() + " not implemented");
2256 }
2257 };
2258
2264 double getSigma5(double q2)
2265 {
2266 return (I_5(q2, 0) + I_5(q2, 1))/2.;
2267 };
2268
2274 double getSigma6s(double q2)
2275 {
2276 return (I_6s(q2, 0) + I_6s(q2, 1))/2.;
2277 };
2278
2284 double getSigma6c(double q2)
2285 {
2286 return (I_6c(q2, 0) + I_6c(q2, 1))/2.;
2287 };
2288
2294 double getSigma7(double q2)
2295 {
2296 switch(vectorM){
2297 case QCD::K_star:
2298 return (I_7(q2, 0) + I_7(q2, 1))/2.;
2299 break;
2300 case QCD::K_star_P:
2301 return (I_7(q2, 0) + I_7(q2, 1))/2.;
2302 break;
2303 case QCD::PHI:
2304 return (I_7(q2, 0) + I_7(q2, 1) - ys * h_7(q2, 0))/2.;
2305 break;
2306 default:
2307 std::stringstream out;
2308 out << vectorM;
2309 throw std::runtime_error("MVll::getSigma7 : vector " + out.str() + " not implemented");
2310 }
2311 };
2312
2318 double getSigma8(double q2)
2319 {
2320 return (I_8(q2, 0) + I_8(q2, 1))/2.;
2321 };
2322
2328 double getSigma9(double q2)
2329 {
2330 return (I_9(q2, 0) + I_9(q2, 1))/2.;
2331 };
2332
2338 double getDelta1c(double q2)
2339 {
2340 return (I_1c(q2, 0) - I_1c(q2, 1)) / 2.;
2341 };
2342
2348 double getDelta1s(double q2)
2349 {
2350 return (I_1s(q2, 0) - I_1s(q2, 1)) / 2.;
2351 };
2352
2358 double getDelta2c(double q2)
2359 {
2360 return (I_2c(q2, 0) - I_2c(q2, 1)) / 2.;
2361 };
2362
2368 double getDelta2s(double q2)
2369 {
2370 return (I_2s(q2, 0) - I_2s(q2, 1))/2.;
2371 };
2372
2378 double getDelta3(double q2)
2379 {
2380 return (I_3(q2, 0) - I_3(q2, 1))/2.;
2381 };
2382
2388 double getDelta4(double q2)
2389 {
2390 return (I_4(q2, 0) - I_4(q2, 1))/2.;
2391 };
2392
2398 double getDelta5(double q2)
2399 {
2400 switch(vectorM){
2401 case QCD::K_star:
2402 return (I_5(q2, 0) - I_5(q2, 1))/2.;
2403 break;
2404 case QCD::K_star_P:
2405 return (I_5(q2, 0) - I_5(q2, 1))/2.;
2406 break;
2407 case QCD::PHI:
2408 return (1. - ys*ys)/(1. + xs*xs) * (I_5(q2, 0) - I_5(q2, 1) - xs * s_5(q2, 0))/2.;
2409 break;
2410 default:
2411 std::stringstream out;
2412 out << vectorM;
2413 throw std::runtime_error("MVll::getDelta5 : vector " + out.str() + " not implemented");
2414 }
2415 };
2416
2422 double getDelta6s(double q2)
2423 {
2424 switch(vectorM){
2425 case QCD::K_star:
2426 return (I_6s(q2, 0) - I_6s(q2, 1))/2.;
2427 break;
2428 case QCD::K_star_P:
2429 return (I_6s(q2, 0) - I_6s(q2, 1))/2.;
2430 break;
2431 case QCD::PHI:
2432 return (1. - ys*ys)/(1. + xs*xs) * (I_6s(q2, 0) - I_6s(q2, 1) - xs * s_6s(q2, 0))/2.;
2433 break;
2434 default:
2435 std::stringstream out;
2436 out << vectorM;
2437 throw std::runtime_error("MVll::getDelta6s : vector " + out.str() + " not implemented");
2438 }
2439 };
2440
2446 double getDelta6c(double q2)
2447 {
2448 switch(vectorM){
2449 case QCD::K_star:
2450 return (I_6c(q2, 0) - I_6c(q2, 1))/2.;
2451 break;
2452 case QCD::K_star_P:
2453 return (I_6c(q2, 0) - I_6c(q2, 1))/2.;
2454 break;
2455 case QCD::PHI:
2456 return (1. - ys*ys)/(1. + xs*xs) * (I_6c(q2, 0) - I_6c(q2, 1) - xs * s_6c(q2, 0))/2.;
2457 break;
2458 default:
2459 std::stringstream out;
2460 out << vectorM;
2461 throw std::runtime_error("MVll::getDelta6c : vector " + out.str() + " not implemented");
2462 }
2463 };
2464
2470 double getDelta7(double q2)
2471 {
2472 return (I_7(q2, 0) - I_7(q2, 1))/2.;
2473 };
2474
2480 double getDelta8(double q2)
2481 {
2482 switch(vectorM){
2483 case QCD::K_star:
2484 return (I_8(q2, 0) - I_8(q2, 1))/2.;
2485 break;
2486 case QCD::K_star_P:
2487 return (I_8(q2, 0) - I_8(q2, 1))/2.;
2488 break;
2489 case QCD::PHI:
2490 return (1. - ys*ys)/(1. + xs*xs) * (I_8(q2, 0) - I_8(q2, 1) - xs * s_8(q2, 0))/2.;
2491 break;
2492 default:
2493 std::stringstream out;
2494 out << vectorM;
2495 throw std::runtime_error("MVll::getDelta8 : vector " + out.str() + " not implemented");
2496 }
2497 };
2498
2504 double getDelta9(double q2)
2505 {
2506 switch(vectorM){
2507 case QCD::K_star:
2508 return (I_9(q2, 0) - I_9(q2, 1))/2.;
2509 break;
2510 case QCD::K_star_P:
2511 return (I_9(q2, 0) - I_9(q2, 1))/2.;
2512 break;
2513 case QCD::PHI:
2514 return (1. - ys*ys)/(1. + xs*xs) * (I_9(q2, 0) - I_9(q2, 1) - xs * s_9(q2, 0))/2.;
2515 break;
2516 default:
2517 std::stringstream out;
2518 out << vectorM;
2519 throw std::runtime_error("MVll::getDelta9 : vector " + out.str() + " not implemented");
2520 }
2521 };
2522
2528 double SigmaTree(double q2);
2529
2534 double getintegratedSigmaTree();
2535
2536 gslpp::complex A_Seidel(double q2, double mb2);
2537
2538 gslpp::complex B_Seidel(double q2, double mb2);
2539
2540 gslpp::complex C_Seidel(double q2);
2541
2547 gslpp::complex deltaC7_QCDF(double q2, bool conjugate, bool spline = false);
2548
2554 gslpp::complex deltaC9_QCDF(double q2, bool conjugate, bool spline = false);
2555
2561 gslpp::complex Cq34(bool conjugate);
2562
2567 gslpp::complex T_para_minus_WA(bool conjugate);
2568
2573 gslpp::complex T_perp_WA_1();
2574
2580 gslpp::complex T_perp_WA_2(bool conjugate);
2581
2587 gslpp::complex T_perp_plus_O8(double q2, double u);
2588
2594 gslpp::complex T_para_minus_O8(double q2, double u);
2595
2602 gslpp::complex t_perp(double q2, double u, double m2);
2603
2610 gslpp::complex t_para(double q2, double u, double m2);
2611
2612 gslpp::complex I1(double q2, double u, double m2);
2613
2614 gslpp::complex B0diff(double q2, double u, double m2);
2615
2616 gslpp::complex B0(double s, double m2);
2617
2618 gslpp::complex h_func(double s, double m2);
2619
2626 gslpp::complex T_perp_plus_QSS(double q2, double u, bool conjugate);
2627
2634 gslpp::complex T_para_plus_QSS(double q2, double u, bool conjugate);
2635
2642 gslpp::complex T_para_minus_QSS(double q2, double u, bool conjugate);
2643
2649 double phi_V(double u);
2650
2651 gslpp::complex lambda_B_minus(double q2);
2652
2660 double T_perp_real(double q2, double u, bool conjugate);
2661
2669 double T_perp_imag(double q2, double u, bool conjugate);
2670
2678 double T_para_real(double q2, double u, bool conjugate);
2679
2687 double T_para_imag(double q2, double u, bool conjugate);
2688
2695 double T_perp_real(double q2, bool conjugate);
2696
2703 double T_perp_imag(double q2, bool conjugate);
2704
2711 double T_para_real(double q2, bool conjugate);
2712
2719 double T_para_imag(double q2, bool conjugate);
2720
2721 double QCDF_fit_func(double* x, double* p);
2722
2723 void fit_QCDF_func();
2724
2725 void spline_QCDF_func();
2726
2727 gslpp::complex T_minus(double q2, bool conjugate);
2728
2729 gslpp::complex T_0(double q2, bool conjugate);
2730
2740 double FF_fit(double q2, double a_0, double a_1, double a_2, double MR2);
2741
2742};
2743
2744#endif /* MVLL_H */
2745
Definition: F_1.h:15
Definition: F_2.h:15
A class for the decay.
Definition: MVll.h:308
gslpp::complex T_para_minus_WA(bool conjugate)
QCDF Correction from various BFS paper (hep-ph/0412400). Weak Annihilation.
Definition: MVll.cpp:1837
gslpp::complex getQCDf_1(double q2)
Definition: MVll.h:619
gslpp::complex deltaC7_QCDF(double q2, bool conjugate, bool spline=false)
QCDF Correction from various BFS papers (hep-ph/0403185, hep-ph/0412400) and Greub et....
Definition: MVll.cpp:1730
bool FixedWCbtos
Definition: MVll.h:838
std::vector< std::string > mvllParameters
Definition: MVll.h:833
const StandardModel & mySM
Definition: MVll.h:829
double xs
Definition: MVll.h:861
double mu_h
Definition: MVll.h:853
double getDelta_C9_zExp_p()
The non-pertubative ccbar contributions to the helicity amplitudes.
Definition: MVll.h:539
bool zExpansion
Definition: MVll.h:837
double getgtilde_2_re(double q2)
The real part of .
Definition: MVll.h:703
double getS(double q2)
The form factor .
Definition: MVll.h:454
double phi_V(double u)
QCDF Correction from various BFS paper (hep-ph/0106067).Vector meson distribution amplitude.
Definition: MVll.cpp:1994
double getQCDfC9_2(double q2, double cutoff)
Definition: MVll.h:646
void spline_QCDF_func()
Definition: MVll.cpp:2156
double getQCDfC9p_1(double cutoff)
Definition: MVll.h:658
double getDelta_C9_zExp_m()
The non-pertubative ccbar contributions to the helicity amplitudes.
Definition: MVll.h:549
gslpp::complex t_para(double q2, double u, double m2)
QCDF Correction from various BFS paper (hep-ph/0106067). Part of 4 quark operator contribution.
Definition: MVll.cpp:1876
gslpp::complex B_Seidel(double q2, double mb2)
Definition: MVll.cpp:1692
double geth_0_re(double q2)
The real part of .
Definition: MVll.h:749
bool MVll_DM_flag
Definition: MVll.h:840
double geth_p_im(double q2)
The imaginary part of .
Definition: MVll.h:788
gslpp::complex H_A_p(double q2, bool bar)
The helicity amplitude .
Definition: MVll.cpp:2561
gslpp::complex T_perp_plus_QSS(double q2, double u, bool conjugate)
QCDF Correction from various BFS paper (hep-ph/0106067). 4 quark operator contribution.
Definition: MVll.cpp:1928
double ale
Definition: MVll.h:847
double T_para_real(double q2, double u, bool conjugate)
QCDF Correction from various BFS papers (hep-ph/0106067, hep-ph/0412400). Total.
Definition: MVll.cpp:2031
gslpp::complex T_perp_WA_1()
QCDF Correction from various BFS paper (hep-ph/0412400). Weak Annihilation.
Definition: MVll.cpp:1842
double getMlep()
The mass of the lepton l.
Definition: MVll.h:372
gslpp::complex deltaC9_QCDF(double q2, bool conjugate, bool spline=false)
QCDF Correction from various BFS papers (hep-ph/0403185, hep-ph/0412400) and Greub et....
Definition: MVll.cpp:1778
double Mb
Definition: MVll.h:851
std::unique_ptr< F_2 > myF_2
Definition: MVll.h:835
gslpp::complex Cq34(bool conjugate)
QCDF Correction from various BFS paper (hep-ph/0412400). Part of Weak Annihilation.
Definition: MVll.cpp:1827
double QCDF_fit_func(double *x, double *p)
Definition: MVll.cpp:2087
double mPsi2S2
Definition: MVll.h:842
double MM
Definition: MVll.h:849
double T_perp_real(double q2, double u, bool conjugate)
QCDF Correction from various BFS papers (hep-ph/0106067, hep-ph/0412400). Total.
Definition: MVll.cpp:2005
double getgtilde_2_im(double q2)
The immaginary part of .
Definition: MVll.h:714
double geth_m_im(double q2)
The imaginary part of .
Definition: MVll.h:817
gslpp::complex T_para_plus_QSS(double q2, double u, bool conjugate)
QCDF Correction from various BFS paper (hep-ph/0106067). 4 quark operator contribution.
Definition: MVll.cpp:1950
gslpp::complex getQCDf_2(double q2)
Definition: MVll.h:626
double getVm(double q2)
The form factor .
Definition: MVll.h:411
double get_unitarity_bound_g()
The unitarity constraints on form factors .
Definition: MVll.h:474
gslpp::complex T_para_minus_O8(double q2, double u)
QCDF Correction from various BFS paper (hep-ph/0106067). Chromomagnetic dipole contribution contribut...
Definition: MVll.cpp:1860
gslpp::complex C_Seidel(double q2)
Definition: MVll.cpp:1724
double getVp(double q2)
The form factor .
Definition: MVll.h:400
gslpp::complex H_S(double q2, bool bar)
The helicity amplitude .
Definition: MVll.cpp:2581
double mD2
Definition: MVll.h:843
std::vector< std::string > initializeMVllParameters()
A method for initializing the parameters necessary for MVll.
Definition: MVll.cpp:160
std::unique_ptr< F_1 > myF_1
Definition: MVll.h:834
double integrateDelta(int i, double q_min, double q_max)
The integral of from to .
Definition: MVll.cpp:2961
double width
Definition: MVll.h:859
double getgtilde_3_im(double q2)
The imaginary part of .
Definition: MVll.h:737
double alpha_s_mub
Definition: MVll.h:864
double getTm(double q2)
The form factor .
Definition: MVll.h:443
gslpp::complex H_V_0(double q2, bool bar)
The helicity amplitude .
Definition: MVll.cpp:2520
QCD::meson meson
Definition: MVll.h:831
double T_para_imag(double q2, double u, bool conjugate)
QCDF Correction from various BFS papers (hep-ph/0106067, hep-ph/0412400). Total.
Definition: MVll.cpp:2043
double getT0(double q2)
The form factor .
Definition: MVll.h:421
virtual ~MVll()
Destructor.
Definition: MVll.cpp:156
void fit_QCDF_func()
Definition: MVll.cpp:2092
double T_perp_imag(double q2, double u, bool conjugate)
QCDF Correction from various BFS papers (hep-ph/0106067, hep-ph/0412400). Total.
Definition: MVll.cpp:2018
bool dispersion
Definition: MVll.h:836
gslpp::complex h_func(double s, double m2)
Definition: MVll.cpp:1914
double GF
Definition: MVll.h:846
gslpp::complex T_minus(double q2, bool conjugate)
Definition: MVll.cpp:2225
double getSigma(int i, double q_2)
The value of from to .
Definition: MVll.cpp:2913
gslpp::complex geth_p_0()
.
Definition: MVll.h:768
int etaV
Definition: MVll.h:863
gslpp::complex H_V_p(double q2, bool bar)
The helicity amplitude .
Definition: MVll.cpp:2531
double getgtilde_3_re(double q2)
The real part of .
Definition: MVll.h:725
gslpp::complex getQCDf_3(double q2)
Definition: MVll.h:633
gslpp::complex lambda_B_minus(double q2)
Definition: MVll.cpp:1999
double getwidth()
The width of the meson M.
Definition: MVll.h:363
gslpp::complex T_0(double q2, bool conjugate)
Definition: MVll.cpp:2247
double get_unitarity_bound_T1()
The unitarity constraints on form factors .
Definition: MVll.h:494
double getV0(double q2)
The form factor .
Definition: MVll.h:389
double Ms
Definition: MVll.h:857
double get_unitarity_bound_T2_T0()
The unitarity constraints on form factors and .
Definition: MVll.h:504
double mPsi2S
Definition: MVll.h:842
gslpp::complex h_lambda(int hel, double q2)
The non-pertubative ccbar contributions to the helicity amplitudes.
Definition: MVll.cpp:2463
gslpp::complex exp_Phase[3]
Definition: MVll.h:844
double getDelta_C9_zExp_0()
The non-pertubative ccbar contributions to the helicity amplitudes.
Definition: MVll.h:529
double mJpsi
Definition: MVll.h:841
double MV
Definition: MVll.h:850
gslpp::complex geth_m_0()
.
Definition: MVll.h:797
double geth_p_re(double q2)
The real part of .
Definition: MVll.h:778
double getintegratedSigmaTree()
The integral of from 0 to .
Definition: MVll.cpp:3095
double integrateSigmaTree(double q_min, double q_max)
The integral of from to (arxiv/2301.06990)
Definition: MVll.cpp:3060
gslpp::complex T_para_minus_QSS(double q2, double u, bool conjugate)
QCDF Correction from various BFS paper (hep-ph/0106067). 4 quark operator contribution.
Definition: MVll.cpp:1972
double geth_0_im(double q2)
The imaginary part of .
Definition: MVll.h:759
double mc_pole
Definition: MVll.h:856
double getTp(double q2)
The form factor .
Definition: MVll.h:432
double angmomV
Definition: MVll.h:862
gslpp::complex T_perp_WA_2(bool conjugate)
QCDF Correction from various BFS paper (hep-ph/0412400). Weak Annihilation.
Definition: MVll.cpp:1847
double FF_fit(double q2, double a_0, double a_1, double a_2, double MR2)
The fit function from , .
Definition: MVll.cpp:1455
gslpp::complex t_perp(double q2, double u, double m2)
QCDF Correction from various BFS paper (hep-ph/0106067). Part of 4 quark operator contribution.
Definition: MVll.cpp:1867
MVll(const StandardModel &SM_i, QCD::meson meson_i, QCD::meson vector_i, QCD::lepton lep_i)
Constructor.
Definition: MVll.cpp:22
QCD::meson vectorM
Definition: MVll.h:832
gslpp::complex H_V_m(double q2, bool bar)
The helicity amplitude .
Definition: MVll.cpp:2541
gslpp::complex A_Seidel(double q2, double mb2)
Definition: MVll.cpp:1678
gslpp::complex AmpMVpsi_zExpansion(double mpsi, int tran)
Polarization amplitudes for M to V psi, Eq. B.16 of arXiv:2206.03797.
Definition: MVll.cpp:2597
double spectator_charge
Definition: MVll.h:858
double Mlep
Definition: MVll.h:848
double get_unitarity_bound_F2()
The unitarity constraints on form factors .
Definition: MVll.h:484
gslpp::complex B0diff(double q2, double u, double m2)
Definition: MVll.cpp:1900
double Delta_C9_zExp(int hel)
The non-pertubative ccbar contributions to the helicity amplitudes.
Definition: MVll.cpp:2493
double geth_m_re(double q2)
The real part of .
Definition: MVll.h:807
gslpp::complex H_A_0(double q2, bool bar)
The helicity amplitude .
Definition: MVll.cpp:2551
double getQCDfC9p_2(double cutoff)
Definition: MVll.h:664
double SigmaTree(double q2)
Definition: MVll.cpp:3090
double getQCDfC9p_3(double cutoff)
Definition: MVll.h:670
QCD::lepton lep
Definition: MVll.h:830
gslpp::complex I1(double q2, double u, double m2)
Definition: MVll.cpp:1883
double getgtilde_1_im(double q2)
The immaginary part of .
Definition: MVll.h:692
double getQCDfC9_1(double q2, double cutoff)
Definition: MVll.h:640
gslpp::complex B0(double s, double m2)
Definition: MVll.cpp:1908
gslpp::complex H_A_m(double q2, bool bar)
The helicity amplitude .
Definition: MVll.cpp:2571
gslpp::complex T_perp_plus_O8(double q2, double u)
QCDF Correction from various BFS paper (hep-ph/0106067). Chromomagnetic dipole contribution contribut...
Definition: MVll.cpp:1852
gslpp::complex H_P(double q2, bool bar)
The helicity amplitude .
Definition: MVll.cpp:2589
double Mc
Definition: MVll.h:854
double integrateSigma(int i, double q_min, double q_max)
The integral of from to .
Definition: MVll.cpp:2786
bool NeutrinoTree_flag
Definition: MVll.h:839
double mb_pole
Definition: MVll.h:855
double getQCDfC9_3(double q2, double cutoff)
Definition: MVll.h:652
double beta(double q2)
The factor used in the angular coefficients .
Definition: MVll.cpp:2629
double mu_b
Definition: MVll.h:852
double ys
Definition: MVll.h:860
double mJ2
Definition: MVll.h:841
double getgtilde_1_re(double q2)
The real part of .
Definition: MVll.h:681
double get_unitarity_bound_f_F1()
The unitarity constraints on form factors and .
Definition: MVll.h:464
meson
An enum type for mesons.
Definition: QCD.h:336
@ PHI
Definition: QCD.h:348
@ K_star
Definition: QCD.h:349
@ K_star_P
Definition: QCD.h:350
lepton
An enum type for leptons.
Definition: QCD.h:310
A model class for the Standard Model.
A class for the correction in .
Test Observable.
A class for the unitarity constraints on form factors .
A class for the unitarity constraints on form factors .
A class for the unitarity constraints on form factors and .
A class for the unitarity constraints on form factors and .
A class for the unitarity constraints on form factors .