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 H_0_nunu(double q2, bool bar, QCD::lepton lep);
618
625 gslpp::complex H_p_nunu(double q2, bool bar, QCD::lepton lep);
626
633 gslpp::complex H_m_nunu(double q2, bool bar, QCD::lepton lep);
634
641 gslpp::complex AmpMVpsi_zExpansion(double mpsi, int tran);
642
643 gslpp::complex getQCDf_1(double q2)
644 {
645 updateParameters();
646// 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)));
647 return 0.;
648 }
649
650 gslpp::complex getQCDf_2(double q2)
651 {
652 updateParameters();
653// 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)));
654 return 0.;
655 }
656
657 gslpp::complex getQCDf_3(double q2)
658 {
659 updateParameters();
660// 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))));
661 return 0.;
662 }
663
664 double getQCDfC9_1(double q2, double cutoff)
665 {
666 updateParameters();
667 return (getQCDf_1(q2) - getQCDf_1(cutoff) * cutoff/q2).abs();
668 }
669
670 double getQCDfC9_2(double q2, double cutoff)
671 {
672 updateParameters();
673 return (getQCDf_2(q2) - getQCDf_2(cutoff) * cutoff/q2).abs();
674 }
675
676 double getQCDfC9_3(double q2, double cutoff)
677 {
678 updateParameters();
679 return (getQCDf_3(q2) - getQCDf_3(cutoff) * cutoff/q2).abs();
680 }
681
682 double getQCDfC9p_1(double cutoff)
683 {
684 updateParameters();
685 return (getQCDf_1(cutoff) * cutoff).abs();
686 }
687
688 double getQCDfC9p_2(double cutoff)
689 {
690 updateParameters();
691 return (getQCDf_2(cutoff) * cutoff).abs();
692 }
693
694 double getQCDfC9p_3(double cutoff)
695 {
696 updateParameters();
697 return (getQCDf_3(cutoff) * cutoff).abs();
698 }
699
705 double getgtilde_1_re(double q2)
706 {
707 updateParameters();
708 return C2_inv * (gtilde_1_pre/(sqrt(lambda(q2)) * V(q2)) * (h_lambda(2,q2)-h_lambda(1,q2))).real()/q2;
709 }
710
716 double getgtilde_1_im(double q2)
717 {
718 updateParameters();
719 return C2_inv * (gtilde_1_pre/(sqrt(lambda(q2)) * V(q2)) * (h_lambda(2,q2)-h_lambda(1,q2))).imag()/q2;
720 }
721
727 double getgtilde_2_re(double q2)
728 {
729 updateParameters();
730 return C2_inv * (gtilde_2_pre/A_1(q2) * (h_lambda(1,q2)+h_lambda(2,q2))).real()/q2;
731 }
732
738 double getgtilde_2_im(double q2)
739 {
740 updateParameters();
741 return C2_inv * (gtilde_2_pre/A_1(q2) * (h_lambda(1,q2)+h_lambda(2,q2))).imag()/q2;
742 }
743
749 double getgtilde_3_re(double q2)
750 {
751 updateParameters();
752 return C2_inv * (gtilde_3_pre/(lambda(q2) * A_2(q2)) * (sqrt(q2)*h_lambda(0,q2)/q2-
753 (MM2mMV2 - q2)/(4.*MV) * (h_lambda(1,q2)+h_lambda(2,q2))/q2)).real();
754 }
755
761 double getgtilde_3_im(double q2)
762 {
763 updateParameters();
764 return C2_inv * (gtilde_3_pre/(lambda(q2) * A_2(q2)) * (sqrt(q2)*h_lambda(0,q2)/q2-
765 (MM2mMV2 - q2)/(4.*MV) * (h_lambda(1,q2)+h_lambda(2,q2))/q2)).imag();
766 }
767
773 double geth_0_re(double q2)
774 {
775 return (sixteenM_PI2MM2 * h_lambda(0,q2)/q2).real();
776 }
777
783 double geth_0_im(double q2)
784 {
785 return (sixteenM_PI2MM2 * h_lambda(0,q2)/q2).imag();
786 }
787
792 gslpp::complex geth_p_0()
793 {
794 return h_lambda(1,0.);
795 }
796
802 double geth_p_re(double q2)
803 {
804 return (sixteenM_PI2MM2 * h_lambda(1,q2)/q2).real();
805 }
806
812 double geth_p_im(double q2)
813 {
814 return (sixteenM_PI2MM2 * h_lambda(1,q2)/q2).imag();
815 }
816
821 gslpp::complex geth_m_0()
822 {
823 return h_lambda(2,0.);
824 }
825
831 double geth_m_re(double q2)
832 {
833 return (sixteenM_PI2MM2 * h_lambda(2,q2)/q2).real();
834 }
835
841 double geth_m_im(double q2)
842 {
843 return (sixteenM_PI2MM2 * h_lambda(2,q2)/q2).imag();
844 }
845
850 std::vector<std::string> initializeMVllParameters();
851
852private:
857 std::vector<std::string> mvllParameters;
858 std::unique_ptr<F_1> myF_1;
859 std::unique_ptr<F_2> myF_2;
865 double mJpsi, mJ2;
867 double mD2;
868 gslpp::complex exp_Phase[3];
869
870 double GF;
871 double ale;
872 double Mlep;
873 double MM;
874 double MV;
875 double Mb;
876 double mu_b;
877 double mu_h;
878 double Mc;
879 double mb_pole;
880 double mc_pole;
881 double Ms;
883 double width;
884 double ys;
885 double xs;
886 double angmomV;
887 int etaV;
888 double alpha_s_mub;
889 double fB;
890 double fpara;
891 double fperp;
892
893 double MW;
894 gslpp::complex lambda_t;
895 gslpp::complex lambda_u;
896 double b;
897 gslpp::complex h_0[3];
898 gslpp::complex h_1[3];
899 gslpp::complex h_2[3];
900 gslpp::complex SU3_breaking;
901
902 gslpp::complex beta_0[7];
903 gslpp::complex beta_1[7];
904 gslpp::complex beta_2[7];
906 double Delta_C7_U;
907 double Delta_C9_U;
908
909 double t_p;
910 double t_m;
911 double t_0;
912 double z_0;
913 double s_p;
914 double s_0;
915 double Q2;
916 double chiOPE;
917 double twoalphaBtoKst;
918 double rho_0;
919 double rho_1;
920 double rho_2;
921 double rho_3;
922 double rho_4;
923 double rho_5;
924 double onemrho_0_2;
925 double onemrho_1_2;
926 double onemrho_2_2;
927 double onemrho_3_2;
928 double onemrho_4_2;
929 double onemrho_5_2;
930 double DeltaC9;
931 double DeltaC10;
932 double MMpMV;
933 double MMpMV2;
934 double MMmMV;
935 double MMmMV2;
936 double rV;
937 double Chi1minus;
938 double Chi1plus;
939 double Chi0plus;
940 double Chi0minus;
941 double ChiTT;
942 double ChiBB;
943 double MM2;
944 double MM4;
945 double MV2;
946 double MV4;
947 double MMMV;
948 double MM2mMV2;
949 double MM2pMV2;
950 double fourMV;
951 double onepMMoMV;
952 double MM_MMpMV;
953 double twoMM2;
954 double twoMV2;
955 double twoMM_mbpms;
956 double fourMM2;
957 double Mlep2;
958 double twoMlepMb;
959 double MboMW;
960 double MsoMb;
961 double M_PI2osix;
962 double N_QCDF;
963 double twoMM;
964 double ninetysixM_PI3MM3;
965 double M_PI2;
966 double sixteenM_PI2;
967 double sixteenM_PI2MM2;
968 double twoMboMM;
969 gslpp::complex H_0_pre;
970 double mu_b2;
971 double Mc2;
972 double Mb2;
973 double fourMc2;
974 double fourMb2;
975 double logMc;
976 double logMb;
977 gslpp::complex H_0_WC;
978 gslpp::complex H_c_WC;
979 gslpp::complex H_b_WC;
980 double fournineth;
981 double half;
982 double twothird;
983 double sqrt3;
984 gslpp::complex ihalfMPI;
985 double twoMM3;
986 double gtilde_1_pre;
987 double gtilde_2_pre;
988 double gtilde_3_pre;
989 double C2_inv;
990 double S_L_pre;
991 gslpp::complex NN;
992 gslpp::complex NN_conjugate;
993 double CF;
994 double deltaT_0;
995 double deltaT_1par;
996 double deltaT_1perp;
998 bool h_pole;
999
1000 gslpp::complex ubar;
1001 gslpp::complex arg1;
1002 gslpp::complex B01;
1003 gslpp::complex B00;
1004 gslpp::complex xp;
1005 gslpp::complex xm;
1006 gslpp::complex yp;
1007 gslpp::complex ym;
1008 gslpp::complex L1xp;
1009 gslpp::complex L1xm;
1010 gslpp::complex L1yp;
1011 gslpp::complex L1ym;
1012 gslpp::complex F87_0;
1013 gslpp::complex F87_1;
1014 gslpp::complex F87_2;
1015 gslpp::complex F87_3;
1016 double F89_0;
1017 double F89_1;
1018 double F89_2;
1019 double F89_3;
1021 double a_0V;
1022 double a_1V;
1023 double a_2V;
1024 double MRV_2;
1025 double a_0A0;
1026 double a_1A0;
1027 double a_2A0;
1028 double MRA0_2;
1029 double a_0A1;
1030 double a_1A1;
1031 double a_2A1;
1032 double MRA1_2;
1033 double a_0A12;
1034 double a_1A12;
1035 double a_2A12;
1036 double MRA12_2;
1037 double a_0T1;
1038 double a_1T1;
1039 double a_2T1;
1040 double MRT1_2;
1041 double a_0T2;
1042 double a_1T2;
1043 double a_2T2;
1044 double MRT2_2;
1045 double a_0T23;
1046 double a_1T23;
1047 double a_2T23;
1048 double MRT23_2;
1050 double a_0f;
1051 double a_1f;
1052 double a_2f;
1053 double MRf_2;
1054 double a_0g;
1055 double a_1g;
1056 double a_2g;
1057 double MRg_2;
1058 double a_0F1;
1059 double a_1F1;
1060 double a_2F1;
1061 double MRF1_2;
1062 double a_0F2;
1063 double a_1F2;
1064 double a_2F2;
1065 double MRF2_2;
1066 double a_0T0;
1067 double a_1T0;
1068 double a_2T0;
1069 double MRT0_2;
1071 double unitarity_bound_f_F1;
1072 double unitarity_bound_g;
1073 double unitarity_bound_F2;
1074 double unitarity_bound_T1;
1075 double unitarity_bound_T2_T0;
1077 //additional variables for B to K nu nu
1078 double GF4;
1079 double MM3;
1080 double fM2;
1081 double fV2;
1082
1083 double mtau;
1084 double mtau2;
1085 double Gammatau;
1086 double VusVub_abs2;
1087
1088 gslpp::vector<gslpp::complex> ** allcoeff;
1089 gslpp::vector<gslpp::complex> ** allcoeffh;
1090 gslpp::vector<gslpp::complex> ** allcoeffprime;
1092 gslpp::vector<gslpp::complex> ** allcoeff_noSM;
1094 gslpp::vector<gslpp::complex> ** allcoeff_nu;
1096 gslpp::vector<gslpp::complex> ** allcoeff_noSM_nu;
1098 gslpp::complex C_1;
1099 gslpp::complex C_1L_bar;
1100 gslpp::complex C_1Lh_bar;
1101 gslpp::complex C_2;
1102 gslpp::complex C_2L_bar;
1103 gslpp::complex C_2Lh_bar;
1104 gslpp::complex C_3;
1105 gslpp::complex C_4;
1106 gslpp::complex C_5;
1107 gslpp::complex C_6;
1108 gslpp::complex C_7;
1109 gslpp::complex C_8;
1110 gslpp::complex C_8L;
1111 gslpp::complex C_8Lh;
1112 gslpp::complex C_9;
1113 gslpp::complex C_10;
1114 gslpp::complex C_S;
1115 gslpp::complex C_P;
1117 gslpp::complex C_7p;
1118 gslpp::complex C_9p;
1119 gslpp::complex C_10p;
1120 gslpp::complex C_Sp;
1121 gslpp::complex C_Pp;
1123 gslpp::complex C_L_nunu_e;
1124 gslpp::complex C_R_nunu_e;
1125 gslpp::complex C_L_nunu_mu;
1126 gslpp::complex C_R_nunu_mu;
1127 gslpp::complex C_L_nunu_tau;
1128 gslpp::complex C_R_nunu_tau;
1130 std::vector<double> Re_T_perp;
1131 std::vector<double> Im_T_perp;
1132 std::vector<double> Re_T_para;
1133 std::vector<double> Im_T_para;
1135 gsl_interp_accel *acc_Re_T_perp;
1136 gsl_interp_accel *acc_Im_T_perp;
1137 gsl_interp_accel *acc_Re_T_para;
1138 gsl_interp_accel *acc_Im_T_para;
1139
1140 gsl_spline *spline_Re_T_perp;
1141 gsl_spline *spline_Im_T_perp;
1142 gsl_spline *spline_Re_T_para;
1143 gsl_spline *spline_Im_T_para;
1144
1145 gsl_interp_accel *acc_Re_deltaC7_QCDF;
1146 gsl_interp_accel *acc_Im_deltaC7_QCDF;
1147 gsl_interp_accel *acc_Re_deltaC9_QCDF;
1148 gsl_interp_accel *acc_Im_deltaC9_QCDF;
1149
1150 gsl_spline *spline_Re_deltaC7_QCDF;
1151 gsl_spline *spline_Im_deltaC7_QCDF;
1152 gsl_spline *spline_Re_deltaC9_QCDF;
1153 gsl_spline *spline_Im_deltaC9_QCDF;
1154
1155#if COMPUTECP
1156 std::vector<double> Re_T_perp_conj;
1157 std::vector<double> Im_T_perp_conj;
1158 std::vector<double> Re_T_para_conj;
1159 std::vector<double> Im_T_para_conj;
1161 gsl_interp_accel *acc_Re_T_perp_conj;
1162 gsl_interp_accel *acc_Im_T_perp_conj;
1163 gsl_interp_accel *acc_Re_T_para_conj;
1164 gsl_interp_accel *acc_Im_T_para_conj;
1165
1166 gsl_interp_accel *acc_Re_deltaC7_QCDF_conj;
1167 gsl_interp_accel *acc_Im_deltaC7_QCDF_conj;
1168 gsl_interp_accel *acc_Re_deltaC9_QCDF_conj;
1169 gsl_interp_accel *acc_Im_deltaC9_QCDF_conj;
1170
1171 gsl_spline *spline_Re_T_perp_conj;
1172 gsl_spline *spline_Im_T_perp_conj;
1173 gsl_spline *spline_Re_T_para_conj;
1174 gsl_spline *spline_Im_T_para_conj;
1175
1176 gsl_spline *spline_Re_deltaC7_QCDF_conj;
1177 gsl_spline *spline_Im_deltaC7_QCDF_conj;
1178 gsl_spline *spline_Re_deltaC9_QCDF_conj;
1179 gsl_spline *spline_Im_deltaC9_QCDF_conj;
1180#endif
1181
1182 std::vector<double> myq2;
1184 TFitResultPtr Re_T_perp_res;
1185 TFitResultPtr Im_T_perp_res;
1186 TFitResultPtr Re_T_para_res;
1187 TFitResultPtr Im_T_para_res;
1189 TFitResultPtr Re_T_perp_res_conj;
1190 TFitResultPtr Im_T_perp_res_conj;
1191 TFitResultPtr Re_T_para_res_conj;
1192 TFitResultPtr Im_T_para_res_conj;
1194 TGraph gr1;
1195 TGraph gr2;
1197 TF1 QCDFfit;
1199 TF1 reffit;
1200 TF1 imffit;
1202 double avaSigma;
1203 double avaDelta;
1204 double avaSigmaTree;
1206 double errSigma;
1207 double errDelta;
1208 double errSigmaTree;
1210 gsl_function FS;
1211 gsl_function FD;
1213 gsl_integration_cquad_workspace * w_sigma;
1214 gsl_integration_cquad_workspace * w_delta;
1215 gsl_integration_cquad_workspace * w_sigmaTree;
1217 gsl_error_handler_t * old_handler;
1219 std::map<std::pair<double, double>, gslpp::complex > cacheI1;
1221 std::map<std::pair<double, double>, double > cacheSigma0;
1222 std::map<std::pair<double, double>, double > cacheSigma1;
1223 std::map<std::pair<double, double>, double > cacheSigma2;
1224 std::map<std::pair<double, double>, double > cacheSigma3;
1225 std::map<std::pair<double, double>, double > cacheSigma4;
1226 std::map<std::pair<double, double>, double > cacheSigma5;
1227 std::map<std::pair<double, double>, double > cacheSigma6;
1228 std::map<std::pair<double, double>, double > cacheSigma7;
1229 std::map<std::pair<double, double>, double > cacheSigma8;
1230 std::map<std::pair<double, double>, double > cacheSigma9;
1231 std::map<std::pair<double, double>, double > cacheSigma10;
1232 std::map<std::pair<double, double>, double > cacheSigma11;
1234 std::map<std::pair<double, double>, double > cacheDelta0;
1235 std::map<std::pair<double, double>, double > cacheDelta1;
1236 std::map<std::pair<double, double>, double > cacheDelta2;
1237 std::map<std::pair<double, double>, double > cacheDelta3;
1238 std::map<std::pair<double, double>, double > cacheDelta6;
1239 std::map<std::pair<double, double>, double > cacheDelta7;
1240 std::map<std::pair<double, double>, double > cacheDelta8;
1241 std::map<std::pair<double, double>, double > cacheDelta10;
1242 std::map<std::pair<double, double>, double > cacheDelta11;
1244 std::map<std::pair<double, double>, double > cacheSigmaTree;
1246 unsigned int N_updated;
1247 gslpp::vector<double> N_cache;
1248 gslpp::complex Nc_cache;
1250 unsigned int V_updated;
1251 gslpp::vector<double> V_cache;
1253 unsigned int A0_updated;
1254 gslpp::vector<double> A0_cache;
1256 unsigned int A1_updated;
1257 gslpp::vector<double> A1_cache;
1259 unsigned int T1_updated;
1260 gslpp::vector<double> T1_cache;
1262 unsigned int T2_updated;
1263 gslpp::vector<double> T2_cache;
1265 unsigned int k2_updated;
1266 gslpp::vector<double> k2_cache;
1268 unsigned int z_updated;
1270 unsigned int lambda_updated;
1272 unsigned int beta_updated;
1273 double beta_cache;
1275 unsigned int F_updated;
1277 unsigned int VL1_updated;
1278 unsigned int VL2_updated;
1280 unsigned int TL1_updated;
1281 unsigned int TL2_updated;
1283 unsigned int VR1_updated;
1284 unsigned int VR2_updated;
1286 unsigned int TR1_updated;
1287 unsigned int TR2_updated;
1289 unsigned int VL0_updated;
1290 gslpp::vector<double> VL0_cache;
1292 unsigned int TL0_updated;
1293 gslpp::vector<double> TL0_cache;
1295 unsigned int VR0_updated;
1297 unsigned int TR0_updated;
1299 unsigned int Mb_Ms_updated;
1301 unsigned int SL_updated;
1302 gslpp::vector<double> SL_cache;
1304 unsigned int SR_updated;
1306 unsigned int C_1_updated;
1307 gslpp::complex C_1_cache;
1309 unsigned int C_2_updated;
1310 gslpp::complex C_2_cache;
1312 unsigned int C_3_updated;
1313 gslpp::complex C_3_cache;
1315 unsigned int C_4_updated;
1316 gslpp::complex C_4_cache;
1318 unsigned int C_5_updated;
1319 gslpp::complex C_5_cache;
1321 unsigned int C_6_updated;
1322 gslpp::complex C_6_cache;
1324 unsigned int C_7_updated;
1325 gslpp::complex C_7_cache;
1327 unsigned int C_9_updated;
1328 gslpp::complex C_9_cache;
1330 unsigned int C_10_updated;
1331 gslpp::complex C_10_cache;
1333 unsigned int C_7p_updated;
1334 gslpp::complex C_7p_cache;
1336 unsigned int C_9p_updated;
1337 gslpp::complex C_9p_cache;
1339 unsigned int C_10p_updated;
1340 gslpp::complex C_10p_cache;
1342 unsigned int C_S_updated;
1343 gslpp::complex C_S_cache;
1345 unsigned int C_P_updated;
1346 gslpp::complex C_P_cache;
1348 unsigned int C_Sp_updated;
1349 gslpp::complex C_Sp_cache;
1351 unsigned int C_Pp_updated;
1352 gslpp::complex C_Pp_cache;
1354 unsigned int C_2Lh_updated;
1355 gslpp::complex C_2Lh_cache;
1357 unsigned int C_8Lh_updated;
1358 gslpp::complex C_8Lh_cache;
1360 unsigned int C_L_nunu_e_updated;
1361 unsigned int C_L_nunu_mu_updated;
1362 unsigned int C_L_nunu_tau_updated;
1363 gslpp::complex C_L_nunu_e_cache;
1364 gslpp::complex C_L_nunu_mu_cache;
1365 gslpp::complex C_L_nunu_tau_cache;
1367 unsigned int C_R_nunu_e_updated;
1368 unsigned int C_R_nunu_mu_updated;
1369 unsigned int C_R_nunu_tau_updated;
1370 gslpp::complex C_R_nunu_e_cache;
1371 gslpp::complex C_R_nunu_mu_cache;
1372 gslpp::complex C_R_nunu_tau_cache;
1374 unsigned int Yupdated;
1375 gslpp::vector<double> Ycache;
1377 gslpp::complex h0Ccache[4];
1378 gslpp::complex h1Ccache[4];
1379 gslpp::complex h2Ccache[4];
1381 gslpp::complex beta0Ccache[8];
1382 gslpp::complex beta1Ccache[8];
1383 gslpp::complex beta2Ccache[8];
1385 unsigned int h0_updated;
1386 unsigned int h1_updated;
1387 unsigned int h2_updated;
1389 unsigned int H_V0updated;
1390 gslpp::vector<double> H_V0cache;
1392 unsigned int H_V1updated;
1393 gslpp::vector<double> H_V1cache;
1395 unsigned int H_V2updated;
1396 gslpp::vector<double> H_V2cache;
1398 unsigned int H_A0updated;
1399 unsigned int H_A1updated;
1400 unsigned int H_A2updated;
1402 unsigned int H_Supdated;
1403 gslpp::vector<double> H_Scache;
1405 unsigned int H_Pupdated;
1406 gslpp::vector<double> H_Pcache;
1408 unsigned int I0_updated;
1409 unsigned int I1_updated;
1410 unsigned int I2_updated;
1411 unsigned int I3_updated;
1412 unsigned int I4_updated;
1413 unsigned int I5_updated;
1414 unsigned int I6_updated;
1415 unsigned int I7_updated;
1416 unsigned int I8_updated;
1417 unsigned int I9_updated;
1418 unsigned int I10_updated;
1419 unsigned int I11_updated;
1421 unsigned int Itree_updated;
1422 gslpp::vector<double> Itree_cache;
1424 std::map<std::pair<double, double>, unsigned int > I1Cached;
1426 std::map<std::pair<double, double>, unsigned int > sigma0Cached;
1427 std::map<std::pair<double, double>, unsigned int > sigma1Cached;
1428 std::map<std::pair<double, double>, unsigned int > sigma2Cached;
1429 std::map<std::pair<double, double>, unsigned int > sigma3Cached;
1430 std::map<std::pair<double, double>, unsigned int > sigma4Cached;
1431 std::map<std::pair<double, double>, unsigned int > sigma5Cached;
1432 std::map<std::pair<double, double>, unsigned int > sigma6Cached;
1433 std::map<std::pair<double, double>, unsigned int > sigma7Cached;
1434 std::map<std::pair<double, double>, unsigned int > sigma8Cached;
1435 std::map<std::pair<double, double>, unsigned int > sigma9Cached;
1436 std::map<std::pair<double, double>, unsigned int > sigma10Cached;
1437 std::map<std::pair<double, double>, unsigned int > sigma11Cached;
1439 std::map<std::pair<double, double>, unsigned int > delta0Cached;
1440 std::map<std::pair<double, double>, unsigned int > delta1Cached;
1441 std::map<std::pair<double, double>, unsigned int > delta2Cached;
1442 std::map<std::pair<double, double>, unsigned int > delta3Cached;
1443 std::map<std::pair<double, double>, unsigned int > delta6Cached;
1444 std::map<std::pair<double, double>, unsigned int > delta7Cached;
1445 std::map<std::pair<double, double>, unsigned int > delta8Cached;
1446 std::map<std::pair<double, double>, unsigned int > delta10Cached;
1447 std::map<std::pair<double, double>, unsigned int > delta11Cached;
1449 std::map<std::pair<double, double>, unsigned int > sigmaTreeCached;
1451 std::map<double, unsigned int> deltaTparpCached;
1452 std::map<double, unsigned int> deltaTparmCached;
1453 std::map<double, unsigned int> deltaTperpCached;
1455 std::map<double, gslpp::complex> cacheDeltaTparp;
1456 std::map<double, gslpp::complex> cacheDeltaTparm;
1457 std::map<double, gslpp::complex> cacheDeltaTperp;
1459 unsigned int deltaTparpupdated;
1460 unsigned int deltaTparmupdated;
1461 unsigned int deltaTperpupdated;
1463 unsigned int T_updated;
1464 gslpp::vector<double> T_cache;
1469 void updateParameters();
1470
1474 void checkCache();
1475
1481 double z(double q2);
1482
1488 double z_DM(double q2);
1489
1496 double phi_f(double q2, double MRf_2);
1497
1504 double phi_g(double q2, double MRg_2);
1505
1512 double phi_F1(double q2, double MRF1_2);
1513
1520 double phi_F2(double q2, double MRF2_2);
1521
1528 double phi_T0(double q2, double MRT0_2);
1529
1536 double phi_T1(double q2, double MRT1_2);
1537
1544 double phi_T2(double q2, double MRT2_2);
1545
1555 double f_DM(double q2, double a_0f, double a_1f, double a_2f, double MRf_2);
1556
1566 double g_DM(double q2, double a_0g, double a_1g, double a_2g, double MRg_2);
1567
1577 double F1_DM(double q2, double a_0F1, double a_1F1, double a_2F1, double MRF1_2);
1578
1588 double F2_DM(double q2, double a_0F2, double a_1F2, double a_2F2, double MRF2_2);
1589
1599 double T0_DM(double q2, double a_0T0, double a_1T0, double a_2T0, double MRT0_2);
1600
1610 double T1_DM(double q2, double a_0T1, double a_1T1, double a_2T1, double MRT1_2);
1611
1621 double T2_DM(double q2, double a_0T2, double a_1T2, double a_2T2, double MRT2_2);
1622
1628 double V(double q2);
1629
1630
1636 double A_0(double q2);
1637
1638
1644 double A_1(double q2);
1645
1651 double A_2(double q2);
1652
1658 double T_1(double q2);
1659
1660
1666 double T_2(double q2);
1667
1673 double V_0t(double q2);
1674
1680 double V_p(double q2);
1681
1687 double V_m(double q2);
1688
1694 double T_0t(double q2);
1695
1701 double T_p(double q2);
1702
1708 double T_m(double q2);
1709
1715 double S_L(double q2);
1716
1722 gslpp::complex H_0(double q2);
1723
1731 gslpp::complex H(double q2, double m2, double mu2);
1732
1738 gslpp::complex Y(double q2);
1739
1740 gslpp::complex funct_g(double q2);
1741
1742 gslpp::complex DeltaC9_KD(double q2, int com);
1743
1749 gslpp::complex zh(double q2);
1750
1756 gslpp::complex P(double q2);
1757
1763 gslpp::complex Phi_1(double q2);
1764
1770 gslpp::complex Phi_1_st(double q2);
1771
1777 gslpp::complex Phi_2(double q2);
1778
1784 gslpp::complex Phi_2_st(double q2);
1785
1791 gslpp::complex Phi_3(double q2);
1792
1798 gslpp::complex Phi_3_st(double q2);
1799
1805 gslpp::complex Phi_4(double q2);
1806
1812 gslpp::complex Phi_4_st(double q2);
1813
1819 gslpp::complex Phi_5(double q2);
1820
1826 gslpp::complex Phi_5_st(double q2);
1827
1833 gslpp::complex Phi_6(double q2);
1834
1840 gslpp::complex Phi_6_st(double q2);
1841
1846 gslpp::complex p0();
1847
1853 gslpp::complex p1(double q2);
1854
1860 gslpp::complex p2(double q2);
1861
1867 gslpp::complex p3(double q2);
1868
1874 gslpp::complex p4(double q2);
1875
1881 gslpp::complex p5(double q2);
1882
1888 gslpp::complex p6(double q2);
1889
1895 gslpp::complex phi_1(double q2);
1896
1902 gslpp::complex phi_2(double q2);
1903
1909 gslpp::complex phi_3(double q2);
1910
1916 gslpp::complex phi_4(double q2);
1917
1924 gslpp::complex DeltaC9_zExpansion(double q2, int tran);
1925
1931 double k2 (double q2);
1932
1938 double beta2 (double q2);
1939
1945 double lambda(double q2);
1946
1953 double F(double q2, double b_i);
1954
1961 double I_1c(double q2, bool bar);
1962
1969 double I_1s(double q2, bool bar);
1970
1977 double I_2c(double q2, bool bar);
1978
1985 double I_2s(double q2, bool bar);
1986
1993 double I_3(double q2, bool bar);
1994
2001 double I_4(double q2, bool bar);
2002
2009 double I_5(double q2, bool bar);
2010
2017 double I_6c(double q2, bool bar);
2018
2025 double I_6s(double q2, bool bar);
2026
2033 double I_7(double q2, bool bar);
2034
2041 double I_8(double q2, bool bar);
2042
2049 double I_9(double q2, bool bar);
2050
2057 double h_1s(double q2, bool bar);
2058
2065 double h_1c(double q2, bool bar);
2066
2073 double h_2s(double q2, bool bar);
2074
2081 double h_2c(double q2, bool bar);
2082
2089 double h_3(double q2, bool bar);
2090
2097 double h_4(double q2, bool bar);
2098
2105 double h_7(double q2, bool bar);
2106
2113 double s_5(double q2, bool bar);
2114
2121 double s_6s(double q2, bool bar);
2122
2129 double s_6c(double q2, bool bar);
2130
2137 double s_8(double q2, bool bar);
2138
2145 double s_9(double q2, bool bar);
2146
2152 double getSigma1c(double q2)
2153 {
2154 switch(vectorM){
2155 case QCD::K_star:
2156 return (I_1c(q2, 0) + I_1c(q2, 1))/2.;
2157 break;
2158 case QCD::K_star_P:
2159 return (I_1c(q2, 0) + I_1c(q2, 1))/2.;
2160 break;
2161 case QCD::PHI:
2162 return (I_1c(q2, 0) + I_1c(q2, 1) - ys * h_1c(q2, 0) )/2.;
2163 break;
2164 default:
2165 std::stringstream out;
2166 out << vectorM;
2167 throw std::runtime_error("MVll::getSigma1c : vector " + out.str() + " not implemented");
2168 }
2169 };
2170
2176 double getSigma1s(double q2)
2177 {
2178 switch(vectorM){
2179 case QCD::K_star:
2180 return (I_1s(q2, 0) + I_1s(q2, 1))/2.;
2181 break;
2182 case QCD::K_star_P:
2183 return (I_1s(q2, 0) + I_1s(q2, 1))/2.;
2184 break;
2185 case QCD::PHI:
2186 return (I_1s(q2, 0) + I_1s(q2, 1) - ys * h_1s(q2, 0))/2.;
2187 break;
2188 default:
2189 std::stringstream out;
2190 out << vectorM;
2191 throw std::runtime_error("MVll::getSigma1s : vector " + out.str() + " not implemented");
2192 }
2193 };
2194
2200 double getSigma2c(double q2)
2201 {
2202 switch(vectorM){
2203 case QCD::K_star:
2204 return (I_2c(q2, 0) + I_2c(q2, 1))/2.;
2205 break;
2206 case QCD::K_star_P:
2207 return (I_2c(q2, 0) + I_2c(q2, 1))/2.;
2208 break;
2209 case QCD::PHI:
2210 return (I_2c(q2, 0) + I_2c(q2, 1) - ys * h_2c(q2, 0))/2.;
2211 break;
2212 default:
2213 std::stringstream out;
2214 out << vectorM;
2215 throw std::runtime_error("MVll::getSigma2c : vector " + out.str() + " not implemented");
2216 }
2217 };
2218
2224 double getSigma2s(double q2)
2225 {
2226 switch(vectorM){
2227 case QCD::K_star:
2228 return (I_2s(q2, 0) + I_2s(q2, 1))/2.;
2229 break;
2230 case QCD::K_star_P:
2231 return (I_2s(q2, 0) + I_2s(q2, 1))/2.;
2232 break;
2233 case QCD::PHI:
2234 return (I_2s(q2, 0) + I_2s(q2, 1) - ys * h_2s(q2, 0))/2.;
2235 break;
2236 default:
2237 std::stringstream out;
2238 out << vectorM;
2239 throw std::runtime_error("MVll::getSigma2s : vector " + out.str() + " not implemented");
2240 }
2241 };
2242
2248 double getSigma3(double q2)
2249 {
2250 switch(vectorM){
2251 case QCD::K_star:
2252 return (I_3(q2, 0) + I_3(q2, 1))/2.;
2253 break;
2254 case QCD::K_star_P:
2255 return (I_3(q2, 0) + I_3(q2, 1))/2.;
2256 break;
2257 case QCD::PHI:
2258 return (I_3(q2, 0) + I_3(q2, 1) - ys * h_3(q2, 0))/2.;
2259 break;
2260 default:
2261 std::stringstream out;
2262 out << vectorM;
2263 throw std::runtime_error("MVll::getSigma3 : vector " + out.str() + " not implemented");
2264 }
2265 };
2266
2272 double getSigma4(double q2)
2273 {
2274 switch(vectorM){
2275 case QCD::K_star:
2276 return (I_4(q2, 0) + I_4(q2, 1))/2.;
2277 break;
2278 case QCD::K_star_P:
2279 return (I_4(q2, 0) + I_4(q2, 1))/2.;
2280 break;
2281 case QCD::PHI:
2282 return (I_4(q2, 0) + I_4(q2, 1) - ys * h_4(q2, 0))/2.;
2283 break;
2284 default:
2285 std::stringstream out;
2286 out << vectorM;
2287 throw std::runtime_error("MVll::getSigma4 : vector " + out.str() + " not implemented");
2288 }
2289 };
2290
2296 double getSigma5(double q2)
2297 {
2298 return (I_5(q2, 0) + I_5(q2, 1))/2.;
2299 };
2300
2306 double getSigma6s(double q2)
2307 {
2308 return (I_6s(q2, 0) + I_6s(q2, 1))/2.;
2309 };
2310
2316 double getSigma6c(double q2)
2317 {
2318 return (I_6c(q2, 0) + I_6c(q2, 1))/2.;
2319 };
2320
2326 double getSigma7(double q2)
2327 {
2328 switch(vectorM){
2329 case QCD::K_star:
2330 return (I_7(q2, 0) + I_7(q2, 1))/2.;
2331 break;
2332 case QCD::K_star_P:
2333 return (I_7(q2, 0) + I_7(q2, 1))/2.;
2334 break;
2335 case QCD::PHI:
2336 return (I_7(q2, 0) + I_7(q2, 1) - ys * h_7(q2, 0))/2.;
2337 break;
2338 default:
2339 std::stringstream out;
2340 out << vectorM;
2341 throw std::runtime_error("MVll::getSigma7 : vector " + out.str() + " not implemented");
2342 }
2343 };
2344
2350 double getSigma8(double q2)
2351 {
2352 return (I_8(q2, 0) + I_8(q2, 1))/2.;
2353 };
2354
2360 double getSigma9(double q2)
2361 {
2362 return (I_9(q2, 0) + I_9(q2, 1))/2.;
2363 };
2364
2370 double getDelta1c(double q2)
2371 {
2372 return (I_1c(q2, 0) - I_1c(q2, 1)) / 2.;
2373 };
2374
2380 double getDelta1s(double q2)
2381 {
2382 return (I_1s(q2, 0) - I_1s(q2, 1)) / 2.;
2383 };
2384
2390 double getDelta2c(double q2)
2391 {
2392 return (I_2c(q2, 0) - I_2c(q2, 1)) / 2.;
2393 };
2394
2400 double getDelta2s(double q2)
2401 {
2402 return (I_2s(q2, 0) - I_2s(q2, 1))/2.;
2403 };
2404
2410 double getDelta3(double q2)
2411 {
2412 return (I_3(q2, 0) - I_3(q2, 1))/2.;
2413 };
2414
2420 double getDelta4(double q2)
2421 {
2422 return (I_4(q2, 0) - I_4(q2, 1))/2.;
2423 };
2424
2430 double getDelta5(double q2)
2431 {
2432 switch(vectorM){
2433 case QCD::K_star:
2434 return (I_5(q2, 0) - I_5(q2, 1))/2.;
2435 break;
2436 case QCD::K_star_P:
2437 return (I_5(q2, 0) - I_5(q2, 1))/2.;
2438 break;
2439 case QCD::PHI:
2440 return (1. - ys*ys)/(1. + xs*xs) * (I_5(q2, 0) - I_5(q2, 1) - xs * s_5(q2, 0))/2.;
2441 break;
2442 default:
2443 std::stringstream out;
2444 out << vectorM;
2445 throw std::runtime_error("MVll::getDelta5 : vector " + out.str() + " not implemented");
2446 }
2447 };
2448
2454 double getDelta6s(double q2)
2455 {
2456 switch(vectorM){
2457 case QCD::K_star:
2458 return (I_6s(q2, 0) - I_6s(q2, 1))/2.;
2459 break;
2460 case QCD::K_star_P:
2461 return (I_6s(q2, 0) - I_6s(q2, 1))/2.;
2462 break;
2463 case QCD::PHI:
2464 return (1. - ys*ys)/(1. + xs*xs) * (I_6s(q2, 0) - I_6s(q2, 1) - xs * s_6s(q2, 0))/2.;
2465 break;
2466 default:
2467 std::stringstream out;
2468 out << vectorM;
2469 throw std::runtime_error("MVll::getDelta6s : vector " + out.str() + " not implemented");
2470 }
2471 };
2472
2478 double getDelta6c(double q2)
2479 {
2480 switch(vectorM){
2481 case QCD::K_star:
2482 return (I_6c(q2, 0) - I_6c(q2, 1))/2.;
2483 break;
2484 case QCD::K_star_P:
2485 return (I_6c(q2, 0) - I_6c(q2, 1))/2.;
2486 break;
2487 case QCD::PHI:
2488 return (1. - ys*ys)/(1. + xs*xs) * (I_6c(q2, 0) - I_6c(q2, 1) - xs * s_6c(q2, 0))/2.;
2489 break;
2490 default:
2491 std::stringstream out;
2492 out << vectorM;
2493 throw std::runtime_error("MVll::getDelta6c : vector " + out.str() + " not implemented");
2494 }
2495 };
2496
2502 double getDelta7(double q2)
2503 {
2504 return (I_7(q2, 0) - I_7(q2, 1))/2.;
2505 };
2506
2512 double getDelta8(double q2)
2513 {
2514 switch(vectorM){
2515 case QCD::K_star:
2516 return (I_8(q2, 0) - I_8(q2, 1))/2.;
2517 break;
2518 case QCD::K_star_P:
2519 return (I_8(q2, 0) - I_8(q2, 1))/2.;
2520 break;
2521 case QCD::PHI:
2522 return (1. - ys*ys)/(1. + xs*xs) * (I_8(q2, 0) - I_8(q2, 1) - xs * s_8(q2, 0))/2.;
2523 break;
2524 default:
2525 std::stringstream out;
2526 out << vectorM;
2527 throw std::runtime_error("MVll::getDelta8 : vector " + out.str() + " not implemented");
2528 }
2529 };
2530
2536 double getDelta9(double q2)
2537 {
2538 switch(vectorM){
2539 case QCD::K_star:
2540 return (I_9(q2, 0) - I_9(q2, 1))/2.;
2541 break;
2542 case QCD::K_star_P:
2543 return (I_9(q2, 0) - I_9(q2, 1))/2.;
2544 break;
2545 case QCD::PHI:
2546 return (1. - ys*ys)/(1. + xs*xs) * (I_9(q2, 0) - I_9(q2, 1) - xs * s_9(q2, 0))/2.;
2547 break;
2548 default:
2549 std::stringstream out;
2550 out << vectorM;
2551 throw std::runtime_error("MVll::getDelta9 : vector " + out.str() + " not implemented");
2552 }
2553 };
2554
2560 double SigmaTree(double q2);
2561
2566 double getintegratedSigmaTree();
2567
2568 gslpp::complex A_Seidel(double q2, double mb2);
2569
2570 gslpp::complex B_Seidel(double q2, double mb2);
2571
2572 gslpp::complex C_Seidel(double q2);
2573
2579 gslpp::complex deltaC7_QCDF(double q2, bool conjugate, bool spline = false);
2580
2586 gslpp::complex deltaC9_QCDF(double q2, bool conjugate, bool spline = false);
2587
2593 gslpp::complex Cq34(bool conjugate);
2594
2599 gslpp::complex T_para_minus_WA(bool conjugate);
2600
2605 gslpp::complex T_perp_WA_1();
2606
2612 gslpp::complex T_perp_WA_2(bool conjugate);
2613
2619 gslpp::complex T_perp_plus_O8(double q2, double u);
2620
2626 gslpp::complex T_para_minus_O8(double q2, double u);
2627
2634 gslpp::complex t_perp(double q2, double u, double m2);
2635
2642 gslpp::complex t_para(double q2, double u, double m2);
2643
2644 gslpp::complex I1(double q2, double u, double m2);
2645
2646 gslpp::complex B0diff(double q2, double u, double m2);
2647
2648 gslpp::complex B0(double s, double m2);
2649
2650 gslpp::complex h_func(double s, double m2);
2651
2658 gslpp::complex T_perp_plus_QSS(double q2, double u, bool conjugate);
2659
2666 gslpp::complex T_para_plus_QSS(double q2, double u, bool conjugate);
2667
2674 gslpp::complex T_para_minus_QSS(double q2, double u, bool conjugate);
2675
2681 double phi_V(double u);
2682
2683 gslpp::complex lambda_B_minus(double q2);
2684
2692 double T_perp_real(double q2, double u, bool conjugate);
2693
2701 double T_perp_imag(double q2, double u, bool conjugate);
2702
2710 double T_para_real(double q2, double u, bool conjugate);
2711
2719 double T_para_imag(double q2, double u, bool conjugate);
2720
2727 double T_perp_real(double q2, bool conjugate);
2728
2735 double T_perp_imag(double q2, bool conjugate);
2736
2743 double T_para_real(double q2, bool conjugate);
2744
2751 double T_para_imag(double q2, bool conjugate);
2752
2753 double QCDF_fit_func(double* x, double* p);
2754
2755 void fit_QCDF_func();
2756
2757 void spline_QCDF_func();
2758
2759 gslpp::complex T_minus(double q2, bool conjugate);
2760
2761 gslpp::complex T_0(double q2, bool conjugate);
2762
2772 double FF_fit(double q2, double a_0, double a_1, double a_2, double MR2);
2773
2774};
2775
2776#endif /* MVLL_H */
2777
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:1883
gslpp::complex getQCDf_1(double q2)
Definition: MVll.h:643
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:1776
bool FixedWCbtos
Definition: MVll.h:862
std::vector< std::string > mvllParameters
Definition: MVll.h:857
const StandardModel & mySM
Definition: MVll.h:853
double xs
Definition: MVll.h:885
double mu_h
Definition: MVll.h:877
double getDelta_C9_zExp_p()
The non-pertubative ccbar contributions to the helicity amplitudes.
Definition: MVll.h:539
bool zExpansion
Definition: MVll.h:861
double getgtilde_2_re(double q2)
The real part of .
Definition: MVll.h:727
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:2040
double getQCDfC9_2(double q2, double cutoff)
Definition: MVll.h:670
void spline_QCDF_func()
Definition: MVll.cpp:2202
double getQCDfC9p_1(double cutoff)
Definition: MVll.h:682
gslpp::complex H_m_nunu(double q2, bool bar, QCD::lepton lep)
The helicity amplitude for the invisible decay .
Definition: MVll.cpp:2646
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:1922
gslpp::complex B_Seidel(double q2, double mb2)
Definition: MVll.cpp:1738
double geth_0_re(double q2)
The real part of .
Definition: MVll.h:773
bool MVll_DM_flag
Definition: MVll.h:864
double geth_p_im(double q2)
The imaginary part of .
Definition: MVll.h:812
gslpp::complex H_A_p(double q2, bool bar)
The helicity amplitude .
Definition: MVll.cpp:2590
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:1974
double ale
Definition: MVll.h:871
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:2077
gslpp::complex T_perp_WA_1()
QCDF Correction from various BFS paper (hep-ph/0412400). Weak Annihilation.
Definition: MVll.cpp:1888
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:1824
double Mb
Definition: MVll.h:875
std::unique_ptr< F_2 > myF_2
Definition: MVll.h:859
gslpp::complex Cq34(bool conjugate)
QCDF Correction from various BFS paper (hep-ph/0412400). Part of Weak Annihilation.
Definition: MVll.cpp:1873
double QCDF_fit_func(double *x, double *p)
Definition: MVll.cpp:2133
double mPsi2S2
Definition: MVll.h:866
double MM
Definition: MVll.h:873
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:2051
double getgtilde_2_im(double q2)
The immaginary part of .
Definition: MVll.h:738
double geth_m_im(double q2)
The imaginary part of .
Definition: MVll.h:841
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:1996
gslpp::complex getQCDf_2(double q2)
Definition: MVll.h:650
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:1906
gslpp::complex C_Seidel(double q2)
Definition: MVll.cpp:1770
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:2602
double mD2
Definition: MVll.h:867
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:858
double integrateDelta(int i, double q_min, double q_max)
The integral of from to .
Definition: MVll.cpp:3040
double width
Definition: MVll.h:883
double getgtilde_3_im(double q2)
The imaginary part of .
Definition: MVll.h:761
double alpha_s_mub
Definition: MVll.h:888
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:2566
QCD::meson meson
Definition: MVll.h:855
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:2089
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:2138
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:2064
bool dispersion
Definition: MVll.h:860
gslpp::complex h_func(double s, double m2)
Definition: MVll.cpp:1960
double GF
Definition: MVll.h:870
gslpp::complex T_minus(double q2, bool conjugate)
Definition: MVll.cpp:2271
double getSigma(int i, double q_2)
The value of from to .
Definition: MVll.cpp:2992
gslpp::complex geth_p_0()
.
Definition: MVll.h:792
int etaV
Definition: MVll.h:887
gslpp::complex H_V_p(double q2, bool bar)
The helicity amplitude .
Definition: MVll.cpp:2572
double getgtilde_3_re(double q2)
The real part of .
Definition: MVll.h:749
gslpp::complex getQCDf_3(double q2)
Definition: MVll.h:657
gslpp::complex lambda_B_minus(double q2)
Definition: MVll.cpp:2045
double getwidth()
The width of the meson M.
Definition: MVll.h:363
gslpp::complex T_0(double q2, bool conjugate)
Definition: MVll.cpp:2293
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:881
double get_unitarity_bound_T2_T0()
The unitarity constraints on form factors and .
Definition: MVll.h:504
double mPsi2S
Definition: MVll.h:866
gslpp::complex h_lambda(int hel, double q2)
The non-pertubative ccbar contributions to the helicity amplitudes.
Definition: MVll.cpp:2509
gslpp::complex exp_Phase[3]
Definition: MVll.h:868
double getDelta_C9_zExp_0()
The non-pertubative ccbar contributions to the helicity amplitudes.
Definition: MVll.h:529
double mJpsi
Definition: MVll.h:865
double MV
Definition: MVll.h:874
gslpp::complex geth_m_0()
.
Definition: MVll.h:821
double geth_p_re(double q2)
The real part of .
Definition: MVll.h:802
double getintegratedSigmaTree()
The integral of from 0 to .
Definition: MVll.cpp:3174
double integrateSigmaTree(double q_min, double q_max)
The integral of from to (arxiv/2301.06990)
Definition: MVll.cpp:3139
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:2018
double geth_0_im(double q2)
The imaginary part of .
Definition: MVll.h:783
double mc_pole
Definition: MVll.h:880
double getTp(double q2)
The form factor .
Definition: MVll.h:432
double angmomV
Definition: MVll.h:886
gslpp::complex T_perp_WA_2(bool conjugate)
QCDF Correction from various BFS paper (hep-ph/0412400). Weak Annihilation.
Definition: MVll.cpp:1893
gslpp::complex H_p_nunu(double q2, bool bar, QCD::lepton lep)
The helicity amplitude for the invisible decay .
Definition: MVll.cpp:2632
double FF_fit(double q2, double a_0, double a_1, double a_2, double MR2)
The fit function from , .
Definition: MVll.cpp:1501
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:1913
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:856
gslpp::complex H_V_m(double q2, bool bar)
The helicity amplitude .
Definition: MVll.cpp:2578
gslpp::complex A_Seidel(double q2, double mb2)
Definition: MVll.cpp:1724
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:2660
double spectator_charge
Definition: MVll.h:882
double Mlep
Definition: MVll.h:872
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:1946
double Delta_C9_zExp(int hel)
The non-pertubative ccbar contributions to the helicity amplitudes.
Definition: MVll.cpp:2539
double geth_m_re(double q2)
The real part of .
Definition: MVll.h:831
gslpp::complex H_A_0(double q2, bool bar)
The helicity amplitude .
Definition: MVll.cpp:2584
double getQCDfC9p_2(double cutoff)
Definition: MVll.h:688
double SigmaTree(double q2)
Definition: MVll.cpp:3169
double getQCDfC9p_3(double cutoff)
Definition: MVll.h:694
gslpp::complex H_0_nunu(double q2, bool bar, QCD::lepton lep)
The helicity amplitude for the invisible decay .
Definition: MVll.cpp:2618
QCD::lepton lep
Definition: MVll.h:854
gslpp::complex I1(double q2, double u, double m2)
Definition: MVll.cpp:1929
double getgtilde_1_im(double q2)
The immaginary part of .
Definition: MVll.h:716
double getQCDfC9_1(double q2, double cutoff)
Definition: MVll.h:664
gslpp::complex B0(double s, double m2)
Definition: MVll.cpp:1954
gslpp::complex H_A_m(double q2, bool bar)
The helicity amplitude .
Definition: MVll.cpp:2596
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:1898
gslpp::complex H_P(double q2, bool bar)
The helicity amplitude .
Definition: MVll.cpp:2610
double Mc
Definition: MVll.h:878
double integrateSigma(int i, double q_min, double q_max)
The integral of from to .
Definition: MVll.cpp:2865
bool NeutrinoTree_flag
Definition: MVll.h:863
double mb_pole
Definition: MVll.h:879
double getQCDfC9_3(double q2, double cutoff)
Definition: MVll.h:676
double beta(double q2)
The factor used in the angular coefficients .
Definition: MVll.cpp:2692
double mu_b
Definition: MVll.h:876
double ys
Definition: MVll.h:884
double mJ2
Definition: MVll.h:865
double getgtilde_1_re(double q2)
The real part of .
Definition: MVll.h:705
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 .