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 t_p;
907 double t_m;
908 double t_0;
909 double z_0;
910 double s_p;
911 double s_0;
912 double Q2;
913 double chiOPE;
914 double twoalphaBtoKst;
915 double rho_0;
916 double rho_1;
917 double rho_2;
918 double rho_3;
919 double rho_4;
920 double rho_5;
921 double onemrho_0_2;
922 double onemrho_1_2;
923 double onemrho_2_2;
924 double onemrho_3_2;
925 double onemrho_4_2;
926 double onemrho_5_2;
927 double DeltaC9;
928 double DeltaC10;
929 double MMpMV;
930 double MMpMV2;
931 double MMmMV;
932 double MMmMV2;
933 double rV;
934 double Chi1minus;
935 double Chi1plus;
936 double Chi0plus;
937 double Chi0minus;
938 double ChiTT;
939 double ChiBB;
940 double MM2;
941 double MM4;
942 double MV2;
943 double MV4;
944 double MMMV;
945 double MM2mMV2;
946 double MM2pMV2;
947 double fourMV;
948 double onepMMoMV;
949 double MM_MMpMV;
950 double twoMM2;
951 double twoMV2;
952 double twoMM_mbpms;
953 double fourMM2;
954 double Mlep2;
955 double twoMlepMb;
956 double MboMW;
957 double MsoMb;
958 double M_PI2osix;
959 double N_QCDF;
960 double twoMM;
961 double ninetysixM_PI3MM3;
962 double M_PI2;
963 double sixteenM_PI2;
964 double sixteenM_PI2MM2;
965 double twoMboMM;
966 gslpp::complex H_0_pre;
967 double mu_b2;
968 double Mc2;
969 double Mb2;
970 double fourMc2;
971 double fourMb2;
972 double logMc;
973 double logMb;
974 gslpp::complex H_0_WC;
975 gslpp::complex H_c_WC;
976 gslpp::complex H_b_WC;
977 double fournineth;
978 double half;
979 double twothird;
980 double sqrt3;
981 gslpp::complex ihalfMPI;
982 double twoMM3;
983 double gtilde_1_pre;
984 double gtilde_2_pre;
985 double gtilde_3_pre;
986 double C2_inv;
987 double S_L_pre;
988 gslpp::complex NN;
989 gslpp::complex NN_conjugate;
990 double CF;
991 double deltaT_0;
992 double deltaT_1par;
993 double deltaT_1perp;
995 bool h_pole;
996
997 gslpp::complex ubar;
998 gslpp::complex arg1;
999 gslpp::complex B01;
1000 gslpp::complex B00;
1001 gslpp::complex xp;
1002 gslpp::complex xm;
1003 gslpp::complex yp;
1004 gslpp::complex ym;
1005 gslpp::complex L1xp;
1006 gslpp::complex L1xm;
1007 gslpp::complex L1yp;
1008 gslpp::complex L1ym;
1009 gslpp::complex F87_0;
1010 gslpp::complex F87_1;
1011 gslpp::complex F87_2;
1012 gslpp::complex F87_3;
1013 double F89_0;
1014 double F89_1;
1015 double F89_2;
1016 double F89_3;
1018 double a_0V;
1019 double a_1V;
1020 double a_2V;
1021 double MRV_2;
1022 double a_0A0;
1023 double a_1A0;
1024 double a_2A0;
1025 double MRA0_2;
1026 double a_0A1;
1027 double a_1A1;
1028 double a_2A1;
1029 double MRA1_2;
1030 double a_0A12;
1031 double a_1A12;
1032 double a_2A12;
1033 double MRA12_2;
1034 double a_0T1;
1035 double a_1T1;
1036 double a_2T1;
1037 double MRT1_2;
1038 double a_0T2;
1039 double a_1T2;
1040 double a_2T2;
1041 double MRT2_2;
1042 double a_0T23;
1043 double a_1T23;
1044 double a_2T23;
1045 double MRT23_2;
1047 double a_0f;
1048 double a_1f;
1049 double a_2f;
1050 double MRf_2;
1051 double a_0g;
1052 double a_1g;
1053 double a_2g;
1054 double MRg_2;
1055 double a_0F1;
1056 double a_1F1;
1057 double a_2F1;
1058 double MRF1_2;
1059 double a_0F2;
1060 double a_1F2;
1061 double a_2F2;
1062 double MRF2_2;
1063 double a_0T0;
1064 double a_1T0;
1065 double a_2T0;
1066 double MRT0_2;
1068 double unitarity_bound_f_F1;
1069 double unitarity_bound_g;
1070 double unitarity_bound_F2;
1071 double unitarity_bound_T1;
1072 double unitarity_bound_T2_T0;
1074 //additional variables for B to K nu nu
1075 double GF4;
1076 double MM3;
1077 double fM2;
1078 double fV2;
1079
1080 double mtau;
1081 double mtau2;
1082 double Gammatau;
1083 double VusVub_abs2;
1084
1085 gslpp::vector<gslpp::complex> ** allcoeff;
1086 gslpp::vector<gslpp::complex> ** allcoeffh;
1087 gslpp::vector<gslpp::complex> ** allcoeffprime;
1089 gslpp::vector<gslpp::complex> ** allcoeff_noSM;
1091 gslpp::vector<gslpp::complex> ** allcoeff_nu;
1093 gslpp::vector<gslpp::complex> ** allcoeff_noSM_nu;
1095 gslpp::complex C_1;
1096 gslpp::complex C_1L_bar;
1097 gslpp::complex C_1Lh_bar;
1098 gslpp::complex C_2;
1099 gslpp::complex C_2L_bar;
1100 gslpp::complex C_2Lh_bar;
1101 gslpp::complex C_3;
1102 gslpp::complex C_4;
1103 gslpp::complex C_5;
1104 gslpp::complex C_6;
1105 gslpp::complex C_7;
1106 gslpp::complex C_8;
1107 gslpp::complex C_8L;
1108 gslpp::complex C_8Lh;
1109 gslpp::complex C_9;
1110 gslpp::complex C_10;
1111 gslpp::complex C_S;
1112 gslpp::complex C_P;
1114 gslpp::complex C_7p;
1115 gslpp::complex C_9p;
1116 gslpp::complex C_10p;
1117 gslpp::complex C_Sp;
1118 gslpp::complex C_Pp;
1120 gslpp::complex C_L_nunu_e;
1121 gslpp::complex C_R_nunu_e;
1122 gslpp::complex C_L_nunu_mu;
1123 gslpp::complex C_R_nunu_mu;
1124 gslpp::complex C_L_nunu_tau;
1125 gslpp::complex C_R_nunu_tau;
1127 std::vector<double> Re_T_perp;
1128 std::vector<double> Im_T_perp;
1129 std::vector<double> Re_T_para;
1130 std::vector<double> Im_T_para;
1132 gsl_interp_accel *acc_Re_T_perp;
1133 gsl_interp_accel *acc_Im_T_perp;
1134 gsl_interp_accel *acc_Re_T_para;
1135 gsl_interp_accel *acc_Im_T_para;
1136
1137 gsl_spline *spline_Re_T_perp;
1138 gsl_spline *spline_Im_T_perp;
1139 gsl_spline *spline_Re_T_para;
1140 gsl_spline *spline_Im_T_para;
1141
1142 gsl_interp_accel *acc_Re_deltaC7_QCDF;
1143 gsl_interp_accel *acc_Im_deltaC7_QCDF;
1144 gsl_interp_accel *acc_Re_deltaC9_QCDF;
1145 gsl_interp_accel *acc_Im_deltaC9_QCDF;
1146
1147 gsl_spline *spline_Re_deltaC7_QCDF;
1148 gsl_spline *spline_Im_deltaC7_QCDF;
1149 gsl_spline *spline_Re_deltaC9_QCDF;
1150 gsl_spline *spline_Im_deltaC9_QCDF;
1151
1152#if COMPUTECP
1153 std::vector<double> Re_T_perp_conj;
1154 std::vector<double> Im_T_perp_conj;
1155 std::vector<double> Re_T_para_conj;
1156 std::vector<double> Im_T_para_conj;
1158 gsl_interp_accel *acc_Re_T_perp_conj;
1159 gsl_interp_accel *acc_Im_T_perp_conj;
1160 gsl_interp_accel *acc_Re_T_para_conj;
1161 gsl_interp_accel *acc_Im_T_para_conj;
1162
1163 gsl_interp_accel *acc_Re_deltaC7_QCDF_conj;
1164 gsl_interp_accel *acc_Im_deltaC7_QCDF_conj;
1165 gsl_interp_accel *acc_Re_deltaC9_QCDF_conj;
1166 gsl_interp_accel *acc_Im_deltaC9_QCDF_conj;
1167
1168 gsl_spline *spline_Re_T_perp_conj;
1169 gsl_spline *spline_Im_T_perp_conj;
1170 gsl_spline *spline_Re_T_para_conj;
1171 gsl_spline *spline_Im_T_para_conj;
1172
1173 gsl_spline *spline_Re_deltaC7_QCDF_conj;
1174 gsl_spline *spline_Im_deltaC7_QCDF_conj;
1175 gsl_spline *spline_Re_deltaC9_QCDF_conj;
1176 gsl_spline *spline_Im_deltaC9_QCDF_conj;
1177#endif
1178
1179 std::vector<double> myq2;
1181 TFitResultPtr Re_T_perp_res;
1182 TFitResultPtr Im_T_perp_res;
1183 TFitResultPtr Re_T_para_res;
1184 TFitResultPtr Im_T_para_res;
1186 TFitResultPtr Re_T_perp_res_conj;
1187 TFitResultPtr Im_T_perp_res_conj;
1188 TFitResultPtr Re_T_para_res_conj;
1189 TFitResultPtr Im_T_para_res_conj;
1191 TGraph gr1;
1192 TGraph gr2;
1194 TF1 QCDFfit;
1196 TF1 reffit;
1197 TF1 imffit;
1199 double avaSigma;
1200 double avaDelta;
1201 double avaSigmaTree;
1203 double errSigma;
1204 double errDelta;
1205 double errSigmaTree;
1207 gsl_function FS;
1208 gsl_function FD;
1210 gsl_integration_cquad_workspace * w_sigma;
1211 gsl_integration_cquad_workspace * w_delta;
1212 gsl_integration_cquad_workspace * w_sigmaTree;
1214 gsl_error_handler_t * old_handler;
1216 std::map<std::pair<double, double>, gslpp::complex > cacheI1;
1218 std::map<std::pair<double, double>, double > cacheSigma0;
1219 std::map<std::pair<double, double>, double > cacheSigma1;
1220 std::map<std::pair<double, double>, double > cacheSigma2;
1221 std::map<std::pair<double, double>, double > cacheSigma3;
1222 std::map<std::pair<double, double>, double > cacheSigma4;
1223 std::map<std::pair<double, double>, double > cacheSigma5;
1224 std::map<std::pair<double, double>, double > cacheSigma6;
1225 std::map<std::pair<double, double>, double > cacheSigma7;
1226 std::map<std::pair<double, double>, double > cacheSigma8;
1227 std::map<std::pair<double, double>, double > cacheSigma9;
1228 std::map<std::pair<double, double>, double > cacheSigma10;
1229 std::map<std::pair<double, double>, double > cacheSigma11;
1231 std::map<std::pair<double, double>, double > cacheDelta0;
1232 std::map<std::pair<double, double>, double > cacheDelta1;
1233 std::map<std::pair<double, double>, double > cacheDelta2;
1234 std::map<std::pair<double, double>, double > cacheDelta3;
1235 std::map<std::pair<double, double>, double > cacheDelta6;
1236 std::map<std::pair<double, double>, double > cacheDelta7;
1237 std::map<std::pair<double, double>, double > cacheDelta8;
1238 std::map<std::pair<double, double>, double > cacheDelta10;
1239 std::map<std::pair<double, double>, double > cacheDelta11;
1241 std::map<std::pair<double, double>, double > cacheSigmaTree;
1243 unsigned int N_updated;
1244 gslpp::vector<double> N_cache;
1245 gslpp::complex Nc_cache;
1247 unsigned int V_updated;
1248 gslpp::vector<double> V_cache;
1250 unsigned int A0_updated;
1251 gslpp::vector<double> A0_cache;
1253 unsigned int A1_updated;
1254 gslpp::vector<double> A1_cache;
1256 unsigned int T1_updated;
1257 gslpp::vector<double> T1_cache;
1259 unsigned int T2_updated;
1260 gslpp::vector<double> T2_cache;
1262 unsigned int k2_updated;
1263 gslpp::vector<double> k2_cache;
1265 unsigned int z_updated;
1267 unsigned int lambda_updated;
1269 unsigned int beta_updated;
1270 double beta_cache;
1272 unsigned int F_updated;
1274 unsigned int VL1_updated;
1275 unsigned int VL2_updated;
1277 unsigned int TL1_updated;
1278 unsigned int TL2_updated;
1280 unsigned int VR1_updated;
1281 unsigned int VR2_updated;
1283 unsigned int TR1_updated;
1284 unsigned int TR2_updated;
1286 unsigned int VL0_updated;
1287 gslpp::vector<double> VL0_cache;
1289 unsigned int TL0_updated;
1290 gslpp::vector<double> TL0_cache;
1292 unsigned int VR0_updated;
1294 unsigned int TR0_updated;
1296 unsigned int Mb_Ms_updated;
1298 unsigned int SL_updated;
1299 gslpp::vector<double> SL_cache;
1301 unsigned int SR_updated;
1303 unsigned int C_1_updated;
1304 gslpp::complex C_1_cache;
1306 unsigned int C_2_updated;
1307 gslpp::complex C_2_cache;
1309 unsigned int C_3_updated;
1310 gslpp::complex C_3_cache;
1312 unsigned int C_4_updated;
1313 gslpp::complex C_4_cache;
1315 unsigned int C_5_updated;
1316 gslpp::complex C_5_cache;
1318 unsigned int C_6_updated;
1319 gslpp::complex C_6_cache;
1321 unsigned int C_7_updated;
1322 gslpp::complex C_7_cache;
1324 unsigned int C_9_updated;
1325 gslpp::complex C_9_cache;
1327 unsigned int C_10_updated;
1328 gslpp::complex C_10_cache;
1330 unsigned int C_7p_updated;
1331 gslpp::complex C_7p_cache;
1333 unsigned int C_9p_updated;
1334 gslpp::complex C_9p_cache;
1336 unsigned int C_10p_updated;
1337 gslpp::complex C_10p_cache;
1339 unsigned int C_S_updated;
1340 gslpp::complex C_S_cache;
1342 unsigned int C_P_updated;
1343 gslpp::complex C_P_cache;
1345 unsigned int C_Sp_updated;
1346 gslpp::complex C_Sp_cache;
1348 unsigned int C_Pp_updated;
1349 gslpp::complex C_Pp_cache;
1351 unsigned int C_2Lh_updated;
1352 gslpp::complex C_2Lh_cache;
1354 unsigned int C_8Lh_updated;
1355 gslpp::complex C_8Lh_cache;
1357 unsigned int C_L_nunu_e_updated;
1358 unsigned int C_L_nunu_mu_updated;
1359 unsigned int C_L_nunu_tau_updated;
1360 gslpp::complex C_L_nunu_e_cache;
1361 gslpp::complex C_L_nunu_mu_cache;
1362 gslpp::complex C_L_nunu_tau_cache;
1364 unsigned int C_R_nunu_e_updated;
1365 unsigned int C_R_nunu_mu_updated;
1366 unsigned int C_R_nunu_tau_updated;
1367 gslpp::complex C_R_nunu_e_cache;
1368 gslpp::complex C_R_nunu_mu_cache;
1369 gslpp::complex C_R_nunu_tau_cache;
1371 unsigned int Yupdated;
1372 gslpp::vector<double> Ycache;
1374 gslpp::complex h0Ccache[4];
1375 gslpp::complex h1Ccache[4];
1376 gslpp::complex h2Ccache[4];
1378 gslpp::complex beta0Ccache[8];
1379 gslpp::complex beta1Ccache[8];
1380 gslpp::complex beta2Ccache[8];
1382 unsigned int h0_updated;
1383 unsigned int h1_updated;
1384 unsigned int h2_updated;
1386 unsigned int H_V0updated;
1387 gslpp::vector<double> H_V0cache;
1389 unsigned int H_V1updated;
1390 gslpp::vector<double> H_V1cache;
1392 unsigned int H_V2updated;
1393 gslpp::vector<double> H_V2cache;
1395 unsigned int H_A0updated;
1396 unsigned int H_A1updated;
1397 unsigned int H_A2updated;
1399 unsigned int H_Supdated;
1400 gslpp::vector<double> H_Scache;
1402 unsigned int H_Pupdated;
1403 gslpp::vector<double> H_Pcache;
1405 unsigned int I0_updated;
1406 unsigned int I1_updated;
1407 unsigned int I2_updated;
1408 unsigned int I3_updated;
1409 unsigned int I4_updated;
1410 unsigned int I5_updated;
1411 unsigned int I6_updated;
1412 unsigned int I7_updated;
1413 unsigned int I8_updated;
1414 unsigned int I9_updated;
1415 unsigned int I10_updated;
1416 unsigned int I11_updated;
1418 unsigned int Itree_updated;
1419 gslpp::vector<double> Itree_cache;
1421 std::map<std::pair<double, double>, unsigned int > I1Cached;
1423 std::map<std::pair<double, double>, unsigned int > sigma0Cached;
1424 std::map<std::pair<double, double>, unsigned int > sigma1Cached;
1425 std::map<std::pair<double, double>, unsigned int > sigma2Cached;
1426 std::map<std::pair<double, double>, unsigned int > sigma3Cached;
1427 std::map<std::pair<double, double>, unsigned int > sigma4Cached;
1428 std::map<std::pair<double, double>, unsigned int > sigma5Cached;
1429 std::map<std::pair<double, double>, unsigned int > sigma6Cached;
1430 std::map<std::pair<double, double>, unsigned int > sigma7Cached;
1431 std::map<std::pair<double, double>, unsigned int > sigma8Cached;
1432 std::map<std::pair<double, double>, unsigned int > sigma9Cached;
1433 std::map<std::pair<double, double>, unsigned int > sigma10Cached;
1434 std::map<std::pair<double, double>, unsigned int > sigma11Cached;
1436 std::map<std::pair<double, double>, unsigned int > delta0Cached;
1437 std::map<std::pair<double, double>, unsigned int > delta1Cached;
1438 std::map<std::pair<double, double>, unsigned int > delta2Cached;
1439 std::map<std::pair<double, double>, unsigned int > delta3Cached;
1440 std::map<std::pair<double, double>, unsigned int > delta6Cached;
1441 std::map<std::pair<double, double>, unsigned int > delta7Cached;
1442 std::map<std::pair<double, double>, unsigned int > delta8Cached;
1443 std::map<std::pair<double, double>, unsigned int > delta10Cached;
1444 std::map<std::pair<double, double>, unsigned int > delta11Cached;
1446 std::map<std::pair<double, double>, unsigned int > sigmaTreeCached;
1448 std::map<double, unsigned int> deltaTparpCached;
1449 std::map<double, unsigned int> deltaTparmCached;
1450 std::map<double, unsigned int> deltaTperpCached;
1452 std::map<double, gslpp::complex> cacheDeltaTparp;
1453 std::map<double, gslpp::complex> cacheDeltaTparm;
1454 std::map<double, gslpp::complex> cacheDeltaTperp;
1456 unsigned int deltaTparpupdated;
1457 unsigned int deltaTparmupdated;
1458 unsigned int deltaTperpupdated;
1460 unsigned int T_updated;
1461 gslpp::vector<double> T_cache;
1466 void updateParameters();
1467
1471 void checkCache();
1472
1478 double z(double q2);
1479
1485 double z_DM(double q2);
1486
1493 double phi_f(double q2, double MRf_2);
1494
1501 double phi_g(double q2, double MRg_2);
1502
1509 double phi_F1(double q2, double MRF1_2);
1510
1517 double phi_F2(double q2, double MRF2_2);
1518
1525 double phi_T0(double q2, double MRT0_2);
1526
1533 double phi_T1(double q2, double MRT1_2);
1534
1541 double phi_T2(double q2, double MRT2_2);
1542
1552 double f_DM(double q2, double a_0f, double a_1f, double a_2f, double MRf_2);
1553
1563 double g_DM(double q2, double a_0g, double a_1g, double a_2g, double MRg_2);
1564
1574 double F1_DM(double q2, double a_0F1, double a_1F1, double a_2F1, double MRF1_2);
1575
1585 double F2_DM(double q2, double a_0F2, double a_1F2, double a_2F2, double MRF2_2);
1586
1596 double T0_DM(double q2, double a_0T0, double a_1T0, double a_2T0, double MRT0_2);
1597
1607 double T1_DM(double q2, double a_0T1, double a_1T1, double a_2T1, double MRT1_2);
1608
1618 double T2_DM(double q2, double a_0T2, double a_1T2, double a_2T2, double MRT2_2);
1619
1625 double V(double q2);
1626
1627
1633 double A_0(double q2);
1634
1635
1641 double A_1(double q2);
1642
1648 double A_2(double q2);
1649
1655 double T_1(double q2);
1656
1657
1663 double T_2(double q2);
1664
1670 double V_0t(double q2);
1671
1677 double V_p(double q2);
1678
1684 double V_m(double q2);
1685
1691 double T_0t(double q2);
1692
1698 double T_p(double q2);
1699
1705 double T_m(double q2);
1706
1712 double S_L(double q2);
1713
1719 gslpp::complex H_0(double q2);
1720
1728 gslpp::complex H(double q2, double m2, double mu2);
1729
1735 gslpp::complex Y(double q2);
1736
1737 gslpp::complex funct_g(double q2);
1738
1739 gslpp::complex DeltaC9_KD(double q2, int com);
1740
1746 gslpp::complex zh(double q2);
1747
1753 gslpp::complex P(double q2);
1754
1760 gslpp::complex Phi_1(double q2);
1761
1767 gslpp::complex Phi_1_st(double q2);
1768
1774 gslpp::complex Phi_2(double q2);
1775
1781 gslpp::complex Phi_2_st(double q2);
1782
1788 gslpp::complex Phi_3(double q2);
1789
1795 gslpp::complex Phi_3_st(double q2);
1796
1802 gslpp::complex Phi_4(double q2);
1803
1809 gslpp::complex Phi_4_st(double q2);
1810
1816 gslpp::complex Phi_5(double q2);
1817
1823 gslpp::complex Phi_5_st(double q2);
1824
1830 gslpp::complex Phi_6(double q2);
1831
1837 gslpp::complex Phi_6_st(double q2);
1838
1843 gslpp::complex p0();
1844
1850 gslpp::complex p1(double q2);
1851
1857 gslpp::complex p2(double q2);
1858
1864 gslpp::complex p3(double q2);
1865
1871 gslpp::complex p4(double q2);
1872
1878 gslpp::complex p5(double q2);
1879
1885 gslpp::complex p6(double q2);
1886
1892 gslpp::complex phi_1(double q2);
1893
1899 gslpp::complex phi_2(double q2);
1900
1906 gslpp::complex phi_3(double q2);
1907
1913 gslpp::complex phi_4(double q2);
1914
1921 gslpp::complex DeltaC9_zExpansion(double q2, int tran);
1922
1928 double k2 (double q2);
1929
1935 double beta2 (double q2);
1936
1942 double lambda(double q2);
1943
1950 double F(double q2, double b_i);
1951
1958 double I_1c(double q2, bool bar);
1959
1966 double I_1s(double q2, bool bar);
1967
1974 double I_2c(double q2, bool bar);
1975
1982 double I_2s(double q2, bool bar);
1983
1990 double I_3(double q2, bool bar);
1991
1998 double I_4(double q2, bool bar);
1999
2006 double I_5(double q2, bool bar);
2007
2014 double I_6c(double q2, bool bar);
2015
2022 double I_6s(double q2, bool bar);
2023
2030 double I_7(double q2, bool bar);
2031
2038 double I_8(double q2, bool bar);
2039
2046 double I_9(double q2, bool bar);
2047
2054 double h_1s(double q2, bool bar);
2055
2062 double h_1c(double q2, bool bar);
2063
2070 double h_2s(double q2, bool bar);
2071
2078 double h_2c(double q2, bool bar);
2079
2086 double h_3(double q2, bool bar);
2087
2094 double h_4(double q2, bool bar);
2095
2102 double h_7(double q2, bool bar);
2103
2110 double s_5(double q2, bool bar);
2111
2118 double s_6s(double q2, bool bar);
2119
2126 double s_6c(double q2, bool bar);
2127
2134 double s_8(double q2, bool bar);
2135
2142 double s_9(double q2, bool bar);
2143
2149 double getSigma1c(double q2)
2150 {
2151 switch(vectorM){
2152 case QCD::K_star:
2153 return (I_1c(q2, 0) + I_1c(q2, 1))/2.;
2154 break;
2155 case QCD::K_star_P:
2156 return (I_1c(q2, 0) + I_1c(q2, 1))/2.;
2157 break;
2158 case QCD::PHI:
2159 return (I_1c(q2, 0) + I_1c(q2, 1) - ys * h_1c(q2, 0) )/2.;
2160 break;
2161 default:
2162 std::stringstream out;
2163 out << vectorM;
2164 throw std::runtime_error("MVll::getSigma1c : vector " + out.str() + " not implemented");
2165 }
2166 };
2167
2173 double getSigma1s(double q2)
2174 {
2175 switch(vectorM){
2176 case QCD::K_star:
2177 return (I_1s(q2, 0) + I_1s(q2, 1))/2.;
2178 break;
2179 case QCD::K_star_P:
2180 return (I_1s(q2, 0) + I_1s(q2, 1))/2.;
2181 break;
2182 case QCD::PHI:
2183 return (I_1s(q2, 0) + I_1s(q2, 1) - ys * h_1s(q2, 0))/2.;
2184 break;
2185 default:
2186 std::stringstream out;
2187 out << vectorM;
2188 throw std::runtime_error("MVll::getSigma1s : vector " + out.str() + " not implemented");
2189 }
2190 };
2191
2197 double getSigma2c(double q2)
2198 {
2199 switch(vectorM){
2200 case QCD::K_star:
2201 return (I_2c(q2, 0) + I_2c(q2, 1))/2.;
2202 break;
2203 case QCD::K_star_P:
2204 return (I_2c(q2, 0) + I_2c(q2, 1))/2.;
2205 break;
2206 case QCD::PHI:
2207 return (I_2c(q2, 0) + I_2c(q2, 1) - ys * h_2c(q2, 0))/2.;
2208 break;
2209 default:
2210 std::stringstream out;
2211 out << vectorM;
2212 throw std::runtime_error("MVll::getSigma2c : vector " + out.str() + " not implemented");
2213 }
2214 };
2215
2221 double getSigma2s(double q2)
2222 {
2223 switch(vectorM){
2224 case QCD::K_star:
2225 return (I_2s(q2, 0) + I_2s(q2, 1))/2.;
2226 break;
2227 case QCD::K_star_P:
2228 return (I_2s(q2, 0) + I_2s(q2, 1))/2.;
2229 break;
2230 case QCD::PHI:
2231 return (I_2s(q2, 0) + I_2s(q2, 1) - ys * h_2s(q2, 0))/2.;
2232 break;
2233 default:
2234 std::stringstream out;
2235 out << vectorM;
2236 throw std::runtime_error("MVll::getSigma2s : vector " + out.str() + " not implemented");
2237 }
2238 };
2239
2245 double getSigma3(double q2)
2246 {
2247 switch(vectorM){
2248 case QCD::K_star:
2249 return (I_3(q2, 0) + I_3(q2, 1))/2.;
2250 break;
2251 case QCD::K_star_P:
2252 return (I_3(q2, 0) + I_3(q2, 1))/2.;
2253 break;
2254 case QCD::PHI:
2255 return (I_3(q2, 0) + I_3(q2, 1) - ys * h_3(q2, 0))/2.;
2256 break;
2257 default:
2258 std::stringstream out;
2259 out << vectorM;
2260 throw std::runtime_error("MVll::getSigma3 : vector " + out.str() + " not implemented");
2261 }
2262 };
2263
2269 double getSigma4(double q2)
2270 {
2271 switch(vectorM){
2272 case QCD::K_star:
2273 return (I_4(q2, 0) + I_4(q2, 1))/2.;
2274 break;
2275 case QCD::K_star_P:
2276 return (I_4(q2, 0) + I_4(q2, 1))/2.;
2277 break;
2278 case QCD::PHI:
2279 return (I_4(q2, 0) + I_4(q2, 1) - ys * h_4(q2, 0))/2.;
2280 break;
2281 default:
2282 std::stringstream out;
2283 out << vectorM;
2284 throw std::runtime_error("MVll::getSigma4 : vector " + out.str() + " not implemented");
2285 }
2286 };
2287
2293 double getSigma5(double q2)
2294 {
2295 return (I_5(q2, 0) + I_5(q2, 1))/2.;
2296 };
2297
2303 double getSigma6s(double q2)
2304 {
2305 return (I_6s(q2, 0) + I_6s(q2, 1))/2.;
2306 };
2307
2313 double getSigma6c(double q2)
2314 {
2315 return (I_6c(q2, 0) + I_6c(q2, 1))/2.;
2316 };
2317
2323 double getSigma7(double q2)
2324 {
2325 switch(vectorM){
2326 case QCD::K_star:
2327 return (I_7(q2, 0) + I_7(q2, 1))/2.;
2328 break;
2329 case QCD::K_star_P:
2330 return (I_7(q2, 0) + I_7(q2, 1))/2.;
2331 break;
2332 case QCD::PHI:
2333 return (I_7(q2, 0) + I_7(q2, 1) - ys * h_7(q2, 0))/2.;
2334 break;
2335 default:
2336 std::stringstream out;
2337 out << vectorM;
2338 throw std::runtime_error("MVll::getSigma7 : vector " + out.str() + " not implemented");
2339 }
2340 };
2341
2347 double getSigma8(double q2)
2348 {
2349 return (I_8(q2, 0) + I_8(q2, 1))/2.;
2350 };
2351
2357 double getSigma9(double q2)
2358 {
2359 return (I_9(q2, 0) + I_9(q2, 1))/2.;
2360 };
2361
2367 double getDelta1c(double q2)
2368 {
2369 return (I_1c(q2, 0) - I_1c(q2, 1)) / 2.;
2370 };
2371
2377 double getDelta1s(double q2)
2378 {
2379 return (I_1s(q2, 0) - I_1s(q2, 1)) / 2.;
2380 };
2381
2387 double getDelta2c(double q2)
2388 {
2389 return (I_2c(q2, 0) - I_2c(q2, 1)) / 2.;
2390 };
2391
2397 double getDelta2s(double q2)
2398 {
2399 return (I_2s(q2, 0) - I_2s(q2, 1))/2.;
2400 };
2401
2407 double getDelta3(double q2)
2408 {
2409 return (I_3(q2, 0) - I_3(q2, 1))/2.;
2410 };
2411
2417 double getDelta4(double q2)
2418 {
2419 return (I_4(q2, 0) - I_4(q2, 1))/2.;
2420 };
2421
2427 double getDelta5(double q2)
2428 {
2429 switch(vectorM){
2430 case QCD::K_star:
2431 return (I_5(q2, 0) - I_5(q2, 1))/2.;
2432 break;
2433 case QCD::K_star_P:
2434 return (I_5(q2, 0) - I_5(q2, 1))/2.;
2435 break;
2436 case QCD::PHI:
2437 return (1. - ys*ys)/(1. + xs*xs) * (I_5(q2, 0) - I_5(q2, 1) - xs * s_5(q2, 0))/2.;
2438 break;
2439 default:
2440 std::stringstream out;
2441 out << vectorM;
2442 throw std::runtime_error("MVll::getDelta5 : vector " + out.str() + " not implemented");
2443 }
2444 };
2445
2451 double getDelta6s(double q2)
2452 {
2453 switch(vectorM){
2454 case QCD::K_star:
2455 return (I_6s(q2, 0) - I_6s(q2, 1))/2.;
2456 break;
2457 case QCD::K_star_P:
2458 return (I_6s(q2, 0) - I_6s(q2, 1))/2.;
2459 break;
2460 case QCD::PHI:
2461 return (1. - ys*ys)/(1. + xs*xs) * (I_6s(q2, 0) - I_6s(q2, 1) - xs * s_6s(q2, 0))/2.;
2462 break;
2463 default:
2464 std::stringstream out;
2465 out << vectorM;
2466 throw std::runtime_error("MVll::getDelta6s : vector " + out.str() + " not implemented");
2467 }
2468 };
2469
2475 double getDelta6c(double q2)
2476 {
2477 switch(vectorM){
2478 case QCD::K_star:
2479 return (I_6c(q2, 0) - I_6c(q2, 1))/2.;
2480 break;
2481 case QCD::K_star_P:
2482 return (I_6c(q2, 0) - I_6c(q2, 1))/2.;
2483 break;
2484 case QCD::PHI:
2485 return (1. - ys*ys)/(1. + xs*xs) * (I_6c(q2, 0) - I_6c(q2, 1) - xs * s_6c(q2, 0))/2.;
2486 break;
2487 default:
2488 std::stringstream out;
2489 out << vectorM;
2490 throw std::runtime_error("MVll::getDelta6c : vector " + out.str() + " not implemented");
2491 }
2492 };
2493
2499 double getDelta7(double q2)
2500 {
2501 return (I_7(q2, 0) - I_7(q2, 1))/2.;
2502 };
2503
2509 double getDelta8(double q2)
2510 {
2511 switch(vectorM){
2512 case QCD::K_star:
2513 return (I_8(q2, 0) - I_8(q2, 1))/2.;
2514 break;
2515 case QCD::K_star_P:
2516 return (I_8(q2, 0) - I_8(q2, 1))/2.;
2517 break;
2518 case QCD::PHI:
2519 return (1. - ys*ys)/(1. + xs*xs) * (I_8(q2, 0) - I_8(q2, 1) - xs * s_8(q2, 0))/2.;
2520 break;
2521 default:
2522 std::stringstream out;
2523 out << vectorM;
2524 throw std::runtime_error("MVll::getDelta8 : vector " + out.str() + " not implemented");
2525 }
2526 };
2527
2533 double getDelta9(double q2)
2534 {
2535 switch(vectorM){
2536 case QCD::K_star:
2537 return (I_9(q2, 0) - I_9(q2, 1))/2.;
2538 break;
2539 case QCD::K_star_P:
2540 return (I_9(q2, 0) - I_9(q2, 1))/2.;
2541 break;
2542 case QCD::PHI:
2543 return (1. - ys*ys)/(1. + xs*xs) * (I_9(q2, 0) - I_9(q2, 1) - xs * s_9(q2, 0))/2.;
2544 break;
2545 default:
2546 std::stringstream out;
2547 out << vectorM;
2548 throw std::runtime_error("MVll::getDelta9 : vector " + out.str() + " not implemented");
2549 }
2550 };
2551
2557 double SigmaTree(double q2);
2558
2563 double getintegratedSigmaTree();
2564
2565 gslpp::complex A_Seidel(double q2, double mb2);
2566
2567 gslpp::complex B_Seidel(double q2, double mb2);
2568
2569 gslpp::complex C_Seidel(double q2);
2570
2576 gslpp::complex deltaC7_QCDF(double q2, bool conjugate, bool spline = false);
2577
2583 gslpp::complex deltaC9_QCDF(double q2, bool conjugate, bool spline = false);
2584
2590 gslpp::complex Cq34(bool conjugate);
2591
2596 gslpp::complex T_para_minus_WA(bool conjugate);
2597
2602 gslpp::complex T_perp_WA_1();
2603
2609 gslpp::complex T_perp_WA_2(bool conjugate);
2610
2616 gslpp::complex T_perp_plus_O8(double q2, double u);
2617
2623 gslpp::complex T_para_minus_O8(double q2, double u);
2624
2631 gslpp::complex t_perp(double q2, double u, double m2);
2632
2639 gslpp::complex t_para(double q2, double u, double m2);
2640
2641 gslpp::complex I1(double q2, double u, double m2);
2642
2643 gslpp::complex B0diff(double q2, double u, double m2);
2644
2645 gslpp::complex B0(double s, double m2);
2646
2647 gslpp::complex h_func(double s, double m2);
2648
2655 gslpp::complex T_perp_plus_QSS(double q2, double u, bool conjugate);
2656
2663 gslpp::complex T_para_plus_QSS(double q2, double u, bool conjugate);
2664
2671 gslpp::complex T_para_minus_QSS(double q2, double u, bool conjugate);
2672
2678 double phi_V(double u);
2679
2680 gslpp::complex lambda_B_minus(double q2);
2681
2689 double T_perp_real(double q2, double u, bool conjugate);
2690
2698 double T_perp_imag(double q2, double u, bool conjugate);
2699
2707 double T_para_real(double q2, double u, bool conjugate);
2708
2716 double T_para_imag(double q2, double u, bool conjugate);
2717
2724 double T_perp_real(double q2, bool conjugate);
2725
2732 double T_perp_imag(double q2, bool conjugate);
2733
2740 double T_para_real(double q2, bool conjugate);
2741
2748 double T_para_imag(double q2, bool conjugate);
2749
2750 double QCDF_fit_func(double* x, double* p);
2751
2752 void fit_QCDF_func();
2753
2754 void spline_QCDF_func();
2755
2756 gslpp::complex T_minus(double q2, bool conjugate);
2757
2758 gslpp::complex T_0(double q2, bool conjugate);
2759
2769 double FF_fit(double q2, double a_0, double a_1, double a_2, double MR2);
2770
2771};
2772
2773#endif /* MVLL_H */
2774
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:1862
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:1755
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:2019
double getQCDfC9_2(double q2, double cutoff)
Definition: MVll.h:670
void spline_QCDF_func()
Definition: MVll.cpp:2181
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:2625
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:1901
gslpp::complex B_Seidel(double q2, double mb2)
Definition: MVll.cpp:1717
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:2569
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:1953
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:2056
gslpp::complex T_perp_WA_1()
QCDF Correction from various BFS paper (hep-ph/0412400). Weak Annihilation.
Definition: MVll.cpp:1867
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:1803
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:1852
double QCDF_fit_func(double *x, double *p)
Definition: MVll.cpp:2112
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:2030
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:1975
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:1885
gslpp::complex C_Seidel(double q2)
Definition: MVll.cpp:1749
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: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:3015
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:2545
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:2068
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:2117
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:2043
bool dispersion
Definition: MVll.h:860
gslpp::complex h_func(double s, double m2)
Definition: MVll.cpp:1939
double GF
Definition: MVll.h:870
gslpp::complex T_minus(double q2, bool conjugate)
Definition: MVll.cpp:2250
double getSigma(int i, double q_2)
The value of from to .
Definition: MVll.cpp:2967
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:2551
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:2024
double getwidth()
The width of the meson M.
Definition: MVll.h:363
gslpp::complex T_0(double q2, bool conjugate)
Definition: MVll.cpp:2272
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:2488
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:3149
double integrateSigmaTree(double q_min, double q_max)
The integral of from to (arxiv/2301.06990)
Definition: MVll.cpp:3114
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:1997
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:1872
gslpp::complex H_p_nunu(double q2, bool bar, QCD::lepton lep)
The helicity amplitude for the invisible decay .
Definition: MVll.cpp:2611
double FF_fit(double q2, double a_0, double a_1, double a_2, double MR2)
The fit function from , .
Definition: MVll.cpp:1480
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:1892
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:2557
gslpp::complex A_Seidel(double q2, double mb2)
Definition: MVll.cpp:1703
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:2639
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:1925
double Delta_C9_zExp(int hel)
The non-pertubative ccbar contributions to the helicity amplitudes.
Definition: MVll.cpp:2518
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:2563
double getQCDfC9p_2(double cutoff)
Definition: MVll.h:688
double SigmaTree(double q2)
Definition: MVll.cpp:3144
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:2597
QCD::lepton lep
Definition: MVll.h:854
gslpp::complex I1(double q2, double u, double m2)
Definition: MVll.cpp:1908
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:1933
gslpp::complex H_A_m(double q2, bool bar)
The helicity amplitude .
Definition: MVll.cpp:2575
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:1877
gslpp::complex H_P(double q2, bool bar)
The helicity amplitude .
Definition: MVll.cpp:2589
double Mc
Definition: MVll.h:878
double integrateSigma(int i, double q_min, double q_max)
The integral of from to .
Definition: MVll.cpp:2840
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:2671
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 .