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
466 gslpp::complex h_lambda(int hel, double q2);
467
473 double Delta_C9_zExp(int hel);
474
480 {
481 updateParameters();
482 return Delta_C9_zExp(0);
483 };
484
490 {
491 updateParameters();
492 return Delta_C9_zExp(1);
493 };
494
500 {
501 updateParameters();
502 return Delta_C9_zExp(2);
503 };
504
510 gslpp::complex H_V_0(double q2, bool bar);
511
517 gslpp::complex H_V_p(double q2, bool bar);
518
524 gslpp::complex H_V_m(double q2, bool bar);
525
531 gslpp::complex H_A_0(double q2, bool bar);
532
538 gslpp::complex H_A_p(double q2, bool bar);
539
545 gslpp::complex H_A_m(double q2, bool bar);
546
552 gslpp::complex H_S(double q2, bool bar);
553
559 gslpp::complex H_P(double q2, bool bar);
560
567 gslpp::complex AmpMVpsi_zExpansion(double mpsi, int tran);
568
569 gslpp::complex getQCDf_1(double q2)
570 {
571 updateParameters();
572// 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)));
573 return 0.;
574 }
575
576 gslpp::complex getQCDf_2(double q2)
577 {
578 updateParameters();
579// 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)));
580 return 0.;
581 }
582
583 gslpp::complex getQCDf_3(double q2)
584 {
585 updateParameters();
586// 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))));
587 return 0.;
588 }
589
590 double getQCDfC9_1(double q2, double cutoff)
591 {
592 updateParameters();
593 return (getQCDf_1(q2) - getQCDf_1(cutoff) * cutoff/q2).abs();
594 }
595
596 double getQCDfC9_2(double q2, double cutoff)
597 {
598 updateParameters();
599 return (getQCDf_2(q2) - getQCDf_2(cutoff) * cutoff/q2).abs();
600 }
601
602 double getQCDfC9_3(double q2, double cutoff)
603 {
604 updateParameters();
605 return (getQCDf_3(q2) - getQCDf_3(cutoff) * cutoff/q2).abs();
606 }
607
608 double getQCDfC9p_1(double cutoff)
609 {
610 updateParameters();
611 return (getQCDf_1(cutoff) * cutoff).abs();
612 }
613
614 double getQCDfC9p_2(double cutoff)
615 {
616 updateParameters();
617 return (getQCDf_2(cutoff) * cutoff).abs();
618 }
619
620 double getQCDfC9p_3(double cutoff)
621 {
622 updateParameters();
623 return (getQCDf_3(cutoff) * cutoff).abs();
624 }
625
631 double getgtilde_1_re(double q2)
632 {
633 updateParameters();
634 return C2_inv * (gtilde_1_pre/(sqrt(lambda(q2)) * V(q2)) * (h_lambda(2,q2)-h_lambda(1,q2))).real()/q2;
635 }
636
642 double getgtilde_1_im(double q2)
643 {
644 updateParameters();
645 return C2_inv * (gtilde_1_pre/(sqrt(lambda(q2)) * V(q2)) * (h_lambda(2,q2)-h_lambda(1,q2))).imag()/q2;
646 }
647
653 double getgtilde_2_re(double q2)
654 {
655 updateParameters();
656 return C2_inv * (gtilde_2_pre/A_1(q2) * (h_lambda(1,q2)+h_lambda(2,q2))).real()/q2;
657 }
658
664 double getgtilde_2_im(double q2)
665 {
666 updateParameters();
667 return C2_inv * (gtilde_2_pre/A_1(q2) * (h_lambda(1,q2)+h_lambda(2,q2))).imag()/q2;
668 }
669
675 double getgtilde_3_re(double q2)
676 {
677 updateParameters();
678 return C2_inv * (gtilde_3_pre/(lambda(q2) * A_2(q2)) * (sqrt(q2)*h_lambda(0,q2)/q2-
679 (MM2mMV2 - q2)/(4.*MV) * (h_lambda(1,q2)+h_lambda(2,q2))/q2)).real();
680 }
681
687 double getgtilde_3_im(double q2)
688 {
689 updateParameters();
690 return C2_inv * (gtilde_3_pre/(lambda(q2) * A_2(q2)) * (sqrt(q2)*h_lambda(0,q2)/q2-
691 (MM2mMV2 - q2)/(4.*MV) * (h_lambda(1,q2)+h_lambda(2,q2))/q2)).imag();
692 }
693
699 double geth_0_re(double q2)
700 {
701 return (sixteenM_PI2MM2 * h_lambda(0,q2)/q2).real();
702 }
703
709 double geth_0_im(double q2)
710 {
711 return (sixteenM_PI2MM2 * h_lambda(0,q2)/q2).imag();
712 }
713
718 gslpp::complex geth_p_0()
719 {
720 return h_lambda(1,0.);
721 }
722
728 double geth_p_re(double q2)
729 {
730 return (sixteenM_PI2MM2 * h_lambda(1,q2)/q2).real();
731 }
732
738 double geth_p_im(double q2)
739 {
740 return (sixteenM_PI2MM2 * h_lambda(1,q2)/q2).imag();
741 }
742
747 gslpp::complex geth_m_0()
748 {
749 return h_lambda(2,0.);
750 }
751
757 double geth_m_re(double q2)
758 {
759 return (sixteenM_PI2MM2 * h_lambda(2,q2)/q2).real();
760 }
761
767 double geth_m_im(double q2)
768 {
769 return (sixteenM_PI2MM2 * h_lambda(2,q2)/q2).imag();
770 }
771
776 std::vector<std::string> initializeMVllParameters();
777
778private:
783 std::vector<std::string> mvllParameters;
784 std::unique_ptr<F_1> myF_1;
785 std::unique_ptr<F_2> myF_2;
791 double mJpsi, mJ2;
793 double mD2;
794 gslpp::complex exp_Phase[3];
795
796 double GF;
797 double ale;
798 double Mlep;
799 double MM;
800 double MV;
801 double Mb;
802 double mu_b;
803 double mu_h;
804 double Mc;
805 double mb_pole;
806 double mc_pole;
807 double Ms;
809 double width;
810 double ys;
811 double xs;
812 double angmomV;
813 int etaV;
814 double alpha_s_mub;
815 double fB;
816 double fpara;
817 double fperp;
818
819 double MW;
820 gslpp::complex lambda_t;
821 gslpp::complex lambda_u;
822 double b;
823 gslpp::complex h_0[3];
824 gslpp::complex h_1[3];
825 gslpp::complex h_2[3];
826 gslpp::complex SU3_breaking;
827
828 gslpp::complex beta_0[7];
829 gslpp::complex beta_1[7];
830 gslpp::complex beta_2[7];
832 double t_p;
833 double t_m;
834 double t_0;
835 double z_0;
836 double s_p;
837 double s_0;
838 double Q2;
839 double chiOPE;
840 double twoalphaBtoKst;
841 double rho_0;
842 double rho_1;
843 double rho_2;
844 double rho_3;
845 double rho_4;
846 double rho_5;
847 double onemrho_0_2;
848 double onemrho_1_2;
849 double onemrho_2_2;
850 double onemrho_3_2;
851 double onemrho_4_2;
852 double onemrho_5_2;
853 double DeltaC9;
854 double DeltaC10;
855 double MMpMV;
856 double MMpMV2;
857 double MMmMV;
858 double MMmMV2;
859 double rV;
860 double Chi1minus;
861 double Chi1plus;
862 double Chi0plus;
863 double Chi0minus;
864 double ChiTT;
865 double ChiBB;
866 double MM2;
867 double MM4;
868 double MV2;
869 double MV4;
870 double MMMV;
871 double MM2mMV2;
872 double MM2pMV2;
873 double fourMV;
874 double onepMMoMV;
875 double MM_MMpMV;
876 double twoMM2;
877 double twoMV2;
878 double twoMM_mbpms;
879 double fourMM2;
880 double Mlep2;
881 double twoMlepMb;
882 double MboMW;
883 double MsoMb;
884 double M_PI2osix;
885 double N_QCDF;
886 double twoMM;
887 double ninetysixM_PI3MM3;
888 double M_PI2;
889 double sixteenM_PI2;
890 double sixteenM_PI2MM2;
891 double twoMboMM;
892 gslpp::complex H_0_pre;
893 double mu_b2;
894 double Mc2;
895 double Mb2;
896 double fourMc2;
897 double fourMb2;
898 double logMc;
899 double logMb;
900 gslpp::complex H_0_WC;
901 gslpp::complex H_c_WC;
902 gslpp::complex H_b_WC;
903 double fournineth;
904 double half;
905 double twothird;
906 double sqrt3;
907 gslpp::complex ihalfMPI;
908 double twoMM3;
909 double gtilde_1_pre;
910 double gtilde_2_pre;
911 double gtilde_3_pre;
912 double C2_inv;
913 double S_L_pre;
914 gslpp::complex NN;
915 gslpp::complex NN_conjugate;
916 double CF;
917 double deltaT_0;
918 double deltaT_1par;
919 double deltaT_1perp;
921 bool h_pole;
922
923 gslpp::complex ubar;
924 gslpp::complex arg1;
925 gslpp::complex B01;
926 gslpp::complex B00;
927 gslpp::complex xp;
928 gslpp::complex xm;
929 gslpp::complex yp;
930 gslpp::complex ym;
931 gslpp::complex L1xp;
932 gslpp::complex L1xm;
933 gslpp::complex L1yp;
934 gslpp::complex L1ym;
935 gslpp::complex F87_0;
936 gslpp::complex F87_1;
937 gslpp::complex F87_2;
938 gslpp::complex F87_3;
939 double F89_0;
940 double F89_1;
941 double F89_2;
942 double F89_3;
944 double a_0V;
945 double a_1V;
946 double a_2V;
947 double MRV_2;
948 double a_0A0;
949 double a_1A0;
950 double a_2A0;
951 double MRA0_2;
952 double a_0A1;
953 double a_1A1;
954 double a_2A1;
955 double MRA1_2;
956 double a_0A12;
957 double a_1A12;
958 double a_2A12;
959 double MRA12_2;
960 double a_0T1;
961 double a_1T1;
962 double a_2T1;
963 double MRT1_2;
964 double a_0T2;
965 double a_1T2;
966 double a_2T2;
967 double MRT2_2;
968 double a_0T23;
969 double a_1T23;
970 double a_2T23;
971 double MRT23_2;
973 double a_0f;
974 double a_1f;
975 double a_2f;
976 double MRf_2;
977 double a_0g;
978 double a_1g;
979 double a_2g;
980 double MRg_2;
981 double a_0F1;
982 double a_1F1;
983 double a_2F1;
984 double MRF1_2;
985 double a_0F2;
986 double a_1F2;
987 double a_2F2;
988 double MRF2_2;
989 double a_0T0;
990 double a_1T0;
991 double a_2T0;
992 double MRT0_2;
994 //additional variables for B to K nu nu
995 double GF4;
996 double MM3;
997 double fM2;
998 double fV2;
999
1000 double mtau;
1001 double mtau2;
1002 double Gammatau;
1003 double VusVub_abs2;
1004
1005 gslpp::vector<gslpp::complex> ** allcoeff;
1006 gslpp::vector<gslpp::complex> ** allcoeffh;
1007 gslpp::vector<gslpp::complex> ** allcoeffprime;
1009 gslpp::vector<gslpp::complex> ** allcoeff_noSM;
1011 gslpp::vector<gslpp::complex> ** allcoeff_nu;
1013 gslpp::vector<gslpp::complex> ** allcoeff_noSM_nu;
1015 gslpp::complex C_1;
1016 gslpp::complex C_1L_bar;
1017 gslpp::complex C_1Lh_bar;
1018 gslpp::complex C_2;
1019 gslpp::complex C_2L_bar;
1020 gslpp::complex C_2Lh_bar;
1021 gslpp::complex C_3;
1022 gslpp::complex C_4;
1023 gslpp::complex C_5;
1024 gslpp::complex C_6;
1025 gslpp::complex C_7;
1026 gslpp::complex C_8;
1027 gslpp::complex C_8L;
1028 gslpp::complex C_8Lh;
1029 gslpp::complex C_9;
1030 gslpp::complex C_10;
1031 gslpp::complex C_S;
1032 gslpp::complex C_P;
1034 gslpp::complex C_7p;
1035 gslpp::complex C_9p;
1036 gslpp::complex C_10p;
1037 gslpp::complex C_Sp;
1038 gslpp::complex C_Pp;
1040 gslpp::complex C_L_nunu;
1041 gslpp::complex C_R_nunu;
1043 gslpp::complex C_L_nunu_e;
1044 gslpp::complex C_R_nunu_e;
1045 gslpp::complex C_L_nunu_mu;
1046 gslpp::complex C_R_nunu_mu;
1047 gslpp::complex C_L_nunu_tau;
1048 gslpp::complex C_R_nunu_tau;
1050 std::vector<double> Re_T_perp;
1051 std::vector<double> Im_T_perp;
1052 std::vector<double> Re_T_para;
1053 std::vector<double> Im_T_para;
1055 gsl_interp_accel *acc_Re_T_perp;
1056 gsl_interp_accel *acc_Im_T_perp;
1057 gsl_interp_accel *acc_Re_T_para;
1058 gsl_interp_accel *acc_Im_T_para;
1059
1060 gsl_spline *spline_Re_T_perp;
1061 gsl_spline *spline_Im_T_perp;
1062 gsl_spline *spline_Re_T_para;
1063 gsl_spline *spline_Im_T_para;
1064
1065 gsl_interp_accel *acc_Re_deltaC7_QCDF;
1066 gsl_interp_accel *acc_Im_deltaC7_QCDF;
1067 gsl_interp_accel *acc_Re_deltaC9_QCDF;
1068 gsl_interp_accel *acc_Im_deltaC9_QCDF;
1069
1070 gsl_spline *spline_Re_deltaC7_QCDF;
1071 gsl_spline *spline_Im_deltaC7_QCDF;
1072 gsl_spline *spline_Re_deltaC9_QCDF;
1073 gsl_spline *spline_Im_deltaC9_QCDF;
1074
1075#if COMPUTECP
1076 std::vector<double> Re_T_perp_conj;
1077 std::vector<double> Im_T_perp_conj;
1078 std::vector<double> Re_T_para_conj;
1079 std::vector<double> Im_T_para_conj;
1081 gsl_interp_accel *acc_Re_T_perp_conj;
1082 gsl_interp_accel *acc_Im_T_perp_conj;
1083 gsl_interp_accel *acc_Re_T_para_conj;
1084 gsl_interp_accel *acc_Im_T_para_conj;
1085
1086 gsl_interp_accel *acc_Re_deltaC7_QCDF_conj;
1087 gsl_interp_accel *acc_Im_deltaC7_QCDF_conj;
1088 gsl_interp_accel *acc_Re_deltaC9_QCDF_conj;
1089 gsl_interp_accel *acc_Im_deltaC9_QCDF_conj;
1090
1091 gsl_spline *spline_Re_T_perp_conj;
1092 gsl_spline *spline_Im_T_perp_conj;
1093 gsl_spline *spline_Re_T_para_conj;
1094 gsl_spline *spline_Im_T_para_conj;
1095
1096 gsl_spline *spline_Re_deltaC7_QCDF_conj;
1097 gsl_spline *spline_Im_deltaC7_QCDF_conj;
1098 gsl_spline *spline_Re_deltaC9_QCDF_conj;
1099 gsl_spline *spline_Im_deltaC9_QCDF_conj;
1100#endif
1101
1102 std::vector<double> myq2;
1104 TFitResultPtr Re_T_perp_res;
1105 TFitResultPtr Im_T_perp_res;
1106 TFitResultPtr Re_T_para_res;
1107 TFitResultPtr Im_T_para_res;
1109 TFitResultPtr Re_T_perp_res_conj;
1110 TFitResultPtr Im_T_perp_res_conj;
1111 TFitResultPtr Re_T_para_res_conj;
1112 TFitResultPtr Im_T_para_res_conj;
1114 TGraph gr1;
1115 TGraph gr2;
1117 TF1 QCDFfit;
1119 TF1 reffit;
1120 TF1 imffit;
1122 double avaSigma;
1123 double avaDelta;
1124 double avaSigmaTree;
1126 double errSigma;
1127 double errDelta;
1128 double errSigmaTree;
1130 gsl_function FS;
1131 gsl_function FD;
1133 gsl_integration_cquad_workspace * w_sigma;
1134 gsl_integration_cquad_workspace * w_delta;
1135 gsl_integration_cquad_workspace * w_sigmaTree;
1137 gsl_error_handler_t * old_handler;
1139 std::map<std::pair<double, double>, gslpp::complex > cacheI1;
1141 std::map<std::pair<double, double>, double > cacheSigma0;
1142 std::map<std::pair<double, double>, double > cacheSigma1;
1143 std::map<std::pair<double, double>, double > cacheSigma2;
1144 std::map<std::pair<double, double>, double > cacheSigma3;
1145 std::map<std::pair<double, double>, double > cacheSigma4;
1146 std::map<std::pair<double, double>, double > cacheSigma5;
1147 std::map<std::pair<double, double>, double > cacheSigma6;
1148 std::map<std::pair<double, double>, double > cacheSigma7;
1149 std::map<std::pair<double, double>, double > cacheSigma9;
1150 std::map<std::pair<double, double>, double > cacheSigma10;
1151 std::map<std::pair<double, double>, double > cacheSigma11;
1153 std::map<std::pair<double, double>, double > cacheDelta0;
1154 std::map<std::pair<double, double>, double > cacheDelta1;
1155 std::map<std::pair<double, double>, double > cacheDelta2;
1156 std::map<std::pair<double, double>, double > cacheDelta3;
1157 std::map<std::pair<double, double>, double > cacheDelta6;
1158 std::map<std::pair<double, double>, double > cacheDelta7;
1159 std::map<std::pair<double, double>, double > cacheDelta8;
1160 std::map<std::pair<double, double>, double > cacheDelta10;
1161 std::map<std::pair<double, double>, double > cacheDelta11;
1163 std::map<std::pair<double, double>, double > cacheSigmaTree;
1165 unsigned int N_updated;
1166 gslpp::vector<double> N_cache;
1167 gslpp::complex Nc_cache;
1169 unsigned int V_updated;
1170 gslpp::vector<double> V_cache;
1172 unsigned int A0_updated;
1173 gslpp::vector<double> A0_cache;
1175 unsigned int A1_updated;
1176 gslpp::vector<double> A1_cache;
1178 unsigned int T1_updated;
1179 gslpp::vector<double> T1_cache;
1181 unsigned int T2_updated;
1182 gslpp::vector<double> T2_cache;
1184 unsigned int k2_updated;
1185 gslpp::vector<double> k2_cache;
1187 unsigned int z_updated;
1189 unsigned int lambda_updated;
1191 unsigned int beta_updated;
1192 double beta_cache;
1194 unsigned int F_updated;
1196 unsigned int VL1_updated;
1197 unsigned int VL2_updated;
1199 unsigned int TL1_updated;
1200 unsigned int TL2_updated;
1202 unsigned int VR1_updated;
1203 unsigned int VR2_updated;
1205 unsigned int TR1_updated;
1206 unsigned int TR2_updated;
1208 unsigned int VL0_updated;
1209 gslpp::vector<double> VL0_cache;
1211 unsigned int TL0_updated;
1212 gslpp::vector<double> TL0_cache;
1214 unsigned int VR0_updated;
1216 unsigned int TR0_updated;
1218 unsigned int Mb_Ms_updated;
1220 unsigned int SL_updated;
1221 gslpp::vector<double> SL_cache;
1223 unsigned int SR_updated;
1225 unsigned int C_1_updated;
1226 gslpp::complex C_1_cache;
1228 unsigned int C_2_updated;
1229 gslpp::complex C_2_cache;
1231 unsigned int C_3_updated;
1232 gslpp::complex C_3_cache;
1234 unsigned int C_4_updated;
1235 gslpp::complex C_4_cache;
1237 unsigned int C_5_updated;
1238 gslpp::complex C_5_cache;
1240 unsigned int C_6_updated;
1241 gslpp::complex C_6_cache;
1243 unsigned int C_7_updated;
1244 gslpp::complex C_7_cache;
1246 unsigned int C_9_updated;
1247 gslpp::complex C_9_cache;
1249 unsigned int C_10_updated;
1250 gslpp::complex C_10_cache;
1252 unsigned int C_7p_updated;
1253 gslpp::complex C_7p_cache;
1255 unsigned int C_9p_updated;
1256 gslpp::complex C_9p_cache;
1258 unsigned int C_10p_updated;
1259 gslpp::complex C_10p_cache;
1261 unsigned int C_S_updated;
1262 gslpp::complex C_S_cache;
1264 unsigned int C_P_updated;
1265 gslpp::complex C_P_cache;
1267 unsigned int C_Sp_updated;
1268 gslpp::complex C_Sp_cache;
1270 unsigned int C_Pp_updated;
1271 gslpp::complex C_Pp_cache;
1273 unsigned int C_2Lh_updated;
1274 gslpp::complex C_2Lh_cache;
1276 unsigned int C_8Lh_updated;
1277 gslpp::complex C_8Lh_cache;
1279 unsigned int C_L_nunu_updated;
1280 gslpp::complex C_L_nunu_cache;
1282 unsigned int C_R_nunu_updated;
1283 gslpp::complex C_R_nunu_cache;
1285 unsigned int Yupdated;
1286 gslpp::vector<double> Ycache;
1288 gslpp::complex h0Ccache[4];
1289 gslpp::complex h1Ccache[4];
1290 gslpp::complex h2Ccache[4];
1292 gslpp::complex beta0Ccache[8];
1293 gslpp::complex beta1Ccache[8];
1294 gslpp::complex beta2Ccache[8];
1296 unsigned int h0_updated;
1297 unsigned int h1_updated;
1298 unsigned int h2_updated;
1300 unsigned int H_V0updated;
1301 gslpp::vector<double> H_V0cache;
1303 unsigned int H_V1updated;
1304 gslpp::vector<double> H_V1cache;
1306 unsigned int H_V2updated;
1307 gslpp::vector<double> H_V2cache;
1309 unsigned int H_A0updated;
1310 unsigned int H_A1updated;
1311 unsigned int H_A2updated;
1313 unsigned int H_Supdated;
1314 gslpp::vector<double> H_Scache;
1316 unsigned int H_Pupdated;
1317 gslpp::vector<double> H_Pcache;
1319 unsigned int I0_updated;
1320 unsigned int I1_updated;
1321 unsigned int I2_updated;
1322 unsigned int I3_updated;
1323 unsigned int I4_updated;
1324 unsigned int I5_updated;
1325 unsigned int I6_updated;
1326 unsigned int I7_updated;
1327 unsigned int I8_updated;
1328 unsigned int I9_updated;
1329 unsigned int I10_updated;
1330 unsigned int I11_updated;
1332 unsigned int Itree_updated;
1333 gslpp::vector<double> Itree_cache;
1335 std::map<std::pair<double, double>, unsigned int > I1Cached;
1337 std::map<std::pair<double, double>, unsigned int > sigma0Cached;
1338 std::map<std::pair<double, double>, unsigned int > sigma1Cached;
1339 std::map<std::pair<double, double>, unsigned int > sigma2Cached;
1340 std::map<std::pair<double, double>, unsigned int > sigma3Cached;
1341 std::map<std::pair<double, double>, unsigned int > sigma4Cached;
1342 std::map<std::pair<double, double>, unsigned int > sigma5Cached;
1343 std::map<std::pair<double, double>, unsigned int > sigma6Cached;
1344 std::map<std::pair<double, double>, unsigned int > sigma7Cached;
1345 std::map<std::pair<double, double>, unsigned int > sigma9Cached;
1346 std::map<std::pair<double, double>, unsigned int > sigma10Cached;
1347 std::map<std::pair<double, double>, unsigned int > sigma11Cached;
1349 std::map<std::pair<double, double>, unsigned int > delta0Cached;
1350 std::map<std::pair<double, double>, unsigned int > delta1Cached;
1351 std::map<std::pair<double, double>, unsigned int > delta2Cached;
1352 std::map<std::pair<double, double>, unsigned int > delta3Cached;
1353 std::map<std::pair<double, double>, unsigned int > delta6Cached;
1354 std::map<std::pair<double, double>, unsigned int > delta7Cached;
1355 std::map<std::pair<double, double>, unsigned int > delta8Cached;
1356 std::map<std::pair<double, double>, unsigned int > delta10Cached;
1357 std::map<std::pair<double, double>, unsigned int > delta11Cached;
1359 std::map<std::pair<double, double>, unsigned int > sigmaTreeCached;
1361 std::map<double, unsigned int> deltaTparpCached;
1362 std::map<double, unsigned int> deltaTparmCached;
1363 std::map<double, unsigned int> deltaTperpCached;
1365 std::map<double, gslpp::complex> cacheDeltaTparp;
1366 std::map<double, gslpp::complex> cacheDeltaTparm;
1367 std::map<double, gslpp::complex> cacheDeltaTperp;
1369 unsigned int deltaTparpupdated;
1370 unsigned int deltaTparmupdated;
1371 unsigned int deltaTperpupdated;
1373 unsigned int T_updated;
1374 gslpp::vector<double> T_cache;
1379 void updateParameters();
1380
1384 void checkCache();
1385
1391 double z(double q2);
1392
1398 double z_DM(double q2);
1399
1406 double phi_f(double q2, double MRf_2);
1407
1414 double phi_g(double q2, double MRg_2);
1415
1422 double phi_F1(double q2, double MRF1_2);
1423
1430 double phi_F2(double q2, double MRF2_2);
1431
1438 double phi_T0(double q2, double MRT0_2);
1439
1446 double phi_T1(double q2, double MRT1_2);
1447
1454 double phi_T2(double q2, double MRT2_2);
1455
1465 double f_DM(double q2, double a_0f, double a_1f, double a_2f, double MRf_2);
1466
1476 double g_DM(double q2, double a_0g, double a_1g, double a_2g, double MRg_2);
1477
1487 double F1_DM(double q2, double a_0F1, double a_1F1, double a_2F1, double MRF1_2);
1488
1498 double F2_DM(double q2, double a_0F2, double a_1F2, double a_2F2, double MRF2_2);
1499
1509 double T0_DM(double q2, double a_0T0, double a_1T0, double a_2T0, double MRT0_2);
1510
1520 double T1_DM(double q2, double a_0T1, double a_1T1, double a_2T1, double MRT1_2);
1521
1531 double T2_DM(double q2, double a_0T2, double a_1T2, double a_2T2, double MRT2_2);
1532
1538 double V(double q2);
1539
1540
1546 double A_0(double q2);
1547
1548
1554 double A_1(double q2);
1555
1561 double A_2(double q2);
1562
1568 double T_1(double q2);
1569
1570
1576 double T_2(double q2);
1577
1583 double V_0t(double q2);
1584
1590 double V_p(double q2);
1591
1597 double V_m(double q2);
1598
1604 double T_0t(double q2);
1605
1611 double T_p(double q2);
1612
1618 double T_m(double q2);
1619
1625 double S_L(double q2);
1626
1632 gslpp::complex H_0(double q2);
1633
1641 gslpp::complex H(double q2, double m2, double mu2);
1642
1648 gslpp::complex Y(double q2);
1649
1650 gslpp::complex funct_g(double q2);
1651
1652 gslpp::complex DeltaC9_KD(double q2, int com);
1653
1659 gslpp::complex zh(double q2);
1660
1666 gslpp::complex P(double q2);
1667
1673 gslpp::complex Phi_1(double q2);
1674
1680 gslpp::complex Phi_1_st(double q2);
1681
1687 gslpp::complex Phi_2(double q2);
1688
1694 gslpp::complex Phi_2_st(double q2);
1695
1701 gslpp::complex Phi_3(double q2);
1702
1708 gslpp::complex Phi_3_st(double q2);
1709
1715 gslpp::complex Phi_4(double q2);
1716
1722 gslpp::complex Phi_4_st(double q2);
1723
1729 gslpp::complex Phi_5(double q2);
1730
1736 gslpp::complex Phi_5_st(double q2);
1737
1743 gslpp::complex Phi_6(double q2);
1744
1750 gslpp::complex Phi_6_st(double q2);
1751
1756 gslpp::complex p0();
1757
1763 gslpp::complex p1(double q2);
1764
1770 gslpp::complex p2(double q2);
1771
1777 gslpp::complex p3(double q2);
1778
1784 gslpp::complex p4(double q2);
1785
1791 gslpp::complex p5(double q2);
1792
1798 gslpp::complex p6(double q2);
1799
1805 gslpp::complex phi_1(double q2);
1806
1812 gslpp::complex phi_2(double q2);
1813
1819 gslpp::complex phi_3(double q2);
1820
1826 gslpp::complex phi_4(double q2);
1827
1834 gslpp::complex DeltaC9_zExpansion(double q2, int tran);
1835
1841 double k2 (double q2);
1842
1848 double beta2 (double q2);
1849
1855 double lambda(double q2);
1856
1863 double F(double q2, double b_i);
1864
1871 double I_1c(double q2, bool bar);
1872
1879 double I_1s(double q2, bool bar);
1880
1887 double I_2c(double q2, bool bar);
1888
1895 double I_2s(double q2, bool bar);
1896
1903 double I_3(double q2, bool bar);
1904
1911 double I_4(double q2, bool bar);
1912
1919 double I_5(double q2, bool bar);
1920
1927 double I_6c(double q2, bool bar);
1928
1935 double I_6s(double q2, bool bar);
1936
1943 double I_7(double q2, bool bar);
1944
1951 double I_8(double q2, bool bar);
1952
1959 double I_9(double q2, bool bar);
1960
1967 double h_1s(double q2, bool bar);
1968
1975 double h_1c(double q2, bool bar);
1976
1983 double h_2s(double q2, bool bar);
1984
1991 double h_2c(double q2, bool bar);
1992
1999 double h_3(double q2, bool bar);
2000
2007 double h_4(double q2, bool bar);
2008
2015 double h_7(double q2, bool bar);
2016
2023 double s_5(double q2, bool bar);
2024
2031 double s_6s(double q2, bool bar);
2032
2039 double s_6c(double q2, bool bar);
2040
2047 double s_8(double q2, bool bar);
2048
2055 double s_9(double q2, bool bar);
2056
2062 double getSigma1c(double q2)
2063 {
2064 switch(vectorM){
2065 case QCD::K_star:
2066 return (I_1c(q2, 0) + I_1c(q2, 1))/2.;
2067 break;
2068 case QCD::K_star_P:
2069 return (I_1c(q2, 0) + I_1c(q2, 1))/2.;
2070 break;
2071 case QCD::PHI:
2072 return (I_1c(q2, 0) + I_1c(q2, 1) - ys * h_1c(q2, 0) )/2.;
2073 break;
2074 default:
2075 std::stringstream out;
2076 out << vectorM;
2077 throw std::runtime_error("MVll::getSigma1c : vector " + out.str() + " not implemented");
2078 }
2079 };
2080
2086 double getSigma1s(double q2)
2087 {
2088 switch(vectorM){
2089 case QCD::K_star:
2090 return (I_1s(q2, 0) + I_1s(q2, 1))/2.;
2091 break;
2092 case QCD::K_star_P:
2093 return (I_1s(q2, 0) + I_1s(q2, 1))/2.;
2094 break;
2095 case QCD::PHI:
2096 return (I_1s(q2, 0) + I_1s(q2, 1) - ys * h_1s(q2, 0))/2.;
2097 break;
2098 default:
2099 std::stringstream out;
2100 out << vectorM;
2101 throw std::runtime_error("MVll::getSigma1s : vector " + out.str() + " not implemented");
2102 }
2103 };
2104
2110 double getSigma2c(double q2)
2111 {
2112 switch(vectorM){
2113 case QCD::K_star:
2114 return (I_2c(q2, 0) + I_2c(q2, 1))/2.;
2115 break;
2116 case QCD::K_star_P:
2117 return (I_2c(q2, 0) + I_2c(q2, 1))/2.;
2118 break;
2119 case QCD::PHI:
2120 return (I_2c(q2, 0) + I_2c(q2, 1) - ys * h_2c(q2, 0))/2.;
2121 break;
2122 default:
2123 std::stringstream out;
2124 out << vectorM;
2125 throw std::runtime_error("MVll::getSigma2c : vector " + out.str() + " not implemented");
2126 }
2127 };
2128
2134 double getSigma2s(double q2)
2135 {
2136 switch(vectorM){
2137 case QCD::K_star:
2138 return (I_2s(q2, 0) + I_2s(q2, 1))/2.;
2139 break;
2140 case QCD::K_star_P:
2141 return (I_2s(q2, 0) + I_2s(q2, 1))/2.;
2142 break;
2143 case QCD::PHI:
2144 return (I_2s(q2, 0) + I_2s(q2, 1) - ys * h_2s(q2, 0))/2.;
2145 break;
2146 default:
2147 std::stringstream out;
2148 out << vectorM;
2149 throw std::runtime_error("MVll::getSigma2s : vector " + out.str() + " not implemented");
2150 }
2151 };
2152
2158 double getSigma3(double q2)
2159 {
2160 switch(vectorM){
2161 case QCD::K_star:
2162 return (I_3(q2, 0) + I_3(q2, 1))/2.;
2163 break;
2164 case QCD::K_star_P:
2165 return (I_3(q2, 0) + I_3(q2, 1))/2.;
2166 break;
2167 case QCD::PHI:
2168 return (I_3(q2, 0) + I_3(q2, 1) - ys * h_3(q2, 0))/2.;
2169 break;
2170 default:
2171 std::stringstream out;
2172 out << vectorM;
2173 throw std::runtime_error("MVll::getSigma3 : vector " + out.str() + " not implemented");
2174 }
2175 };
2176
2182 double getSigma4(double q2)
2183 {
2184 switch(vectorM){
2185 case QCD::K_star:
2186 return (I_4(q2, 0) + I_4(q2, 1))/2.;
2187 break;
2188 case QCD::K_star_P:
2189 return (I_4(q2, 0) + I_4(q2, 1))/2.;
2190 break;
2191 case QCD::PHI:
2192 return (I_4(q2, 0) + I_4(q2, 1) - ys * h_4(q2, 0))/2.;
2193 break;
2194 default:
2195 std::stringstream out;
2196 out << vectorM;
2197 throw std::runtime_error("MVll::getSigma4 : vector " + out.str() + " not implemented");
2198 }
2199 };
2200
2206 double getSigma5(double q2)
2207 {
2208 return (I_5(q2, 0) + I_5(q2, 1))/2.;
2209 };
2210
2216 double getSigma6s(double q2)
2217 {
2218 return (I_6s(q2, 0) + I_6s(q2, 1))/2.;
2219 };
2220
2226 double getSigma6c(double q2)
2227 {
2228 return (I_6c(q2, 0) + I_6c(q2, 1))/2.;
2229 };
2230
2236 double getSigma7(double q2)
2237 {
2238 switch(vectorM){
2239 case QCD::K_star:
2240 return (I_7(q2, 0) + I_7(q2, 1))/2.;
2241 break;
2242 case QCD::K_star_P:
2243 return (I_7(q2, 0) + I_7(q2, 1))/2.;
2244 break;
2245 case QCD::PHI:
2246 return (I_7(q2, 0) + I_7(q2, 1) - ys * h_7(q2, 0))/2.;
2247 break;
2248 default:
2249 std::stringstream out;
2250 out << vectorM;
2251 throw std::runtime_error("MVll::getSigma7 : vector " + out.str() + " not implemented");
2252 }
2253 };
2254
2260 double getSigma8(double q2)
2261 {
2262 return (I_8(q2, 0) + I_8(q2, 1))/2.;
2263 };
2264
2270 double getSigma9(double q2)
2271 {
2272 return (I_9(q2, 0) + I_9(q2, 1))/2.;
2273 };
2274
2280 double getDelta1c(double q2)
2281 {
2282 return (I_1c(q2, 0) - I_1c(q2, 1)) / 2.;
2283 };
2284
2290 double getDelta1s(double q2)
2291 {
2292 return (I_1s(q2, 0) - I_1s(q2, 1)) / 2.;
2293 };
2294
2300 double getDelta2c(double q2)
2301 {
2302 return (I_2c(q2, 0) - I_2c(q2, 1)) / 2.;
2303 };
2304
2310 double getDelta2s(double q2)
2311 {
2312 return (I_2s(q2, 0) - I_2s(q2, 1))/2.;
2313 };
2314
2320 double getDelta3(double q2)
2321 {
2322 return (I_3(q2, 0) - I_3(q2, 1))/2.;
2323 };
2324
2330 double getDelta4(double q2)
2331 {
2332 return (I_4(q2, 0) - I_4(q2, 1))/2.;
2333 };
2334
2340 double getDelta5(double q2)
2341 {
2342 switch(vectorM){
2343 case QCD::K_star:
2344 return (I_5(q2, 0) - I_5(q2, 1))/2.;
2345 break;
2346 case QCD::K_star_P:
2347 return (I_5(q2, 0) - I_5(q2, 1))/2.;
2348 break;
2349 case QCD::PHI:
2350 return (1. - ys*ys)/(1. + xs*xs) * (I_5(q2, 0) - I_5(q2, 1) - xs * s_5(q2, 0))/2.;
2351 break;
2352 default:
2353 std::stringstream out;
2354 out << vectorM;
2355 throw std::runtime_error("MVll::getDelta5 : vector " + out.str() + " not implemented");
2356 }
2357 };
2358
2364 double getDelta6s(double q2)
2365 {
2366 switch(vectorM){
2367 case QCD::K_star:
2368 return (I_6s(q2, 0) - I_6s(q2, 1))/2.;
2369 break;
2370 case QCD::K_star_P:
2371 return (I_6s(q2, 0) - I_6s(q2, 1))/2.;
2372 break;
2373 case QCD::PHI:
2374 return (1. - ys*ys)/(1. + xs*xs) * (I_6s(q2, 0) - I_6s(q2, 1) - xs * s_6s(q2, 0))/2.;
2375 break;
2376 default:
2377 std::stringstream out;
2378 out << vectorM;
2379 throw std::runtime_error("MVll::getDelta6s : vector " + out.str() + " not implemented");
2380 }
2381 };
2382
2388 double getDelta6c(double q2)
2389 {
2390 switch(vectorM){
2391 case QCD::K_star:
2392 return (I_6c(q2, 0) - I_6c(q2, 1))/2.;
2393 break;
2394 case QCD::K_star_P:
2395 return (I_6c(q2, 0) - I_6c(q2, 1))/2.;
2396 break;
2397 case QCD::PHI:
2398 return (1. - ys*ys)/(1. + xs*xs) * (I_6c(q2, 0) - I_6c(q2, 1) - xs * s_6c(q2, 0))/2.;
2399 break;
2400 default:
2401 std::stringstream out;
2402 out << vectorM;
2403 throw std::runtime_error("MVll::getDelta6c : vector " + out.str() + " not implemented");
2404 }
2405 };
2406
2412 double getDelta7(double q2)
2413 {
2414 return (I_7(q2, 0) - I_7(q2, 1))/2.;
2415 };
2416
2422 double getDelta8(double q2)
2423 {
2424 switch(vectorM){
2425 case QCD::K_star:
2426 return (I_8(q2, 0) - I_8(q2, 1))/2.;
2427 break;
2428 case QCD::K_star_P:
2429 return (I_8(q2, 0) - I_8(q2, 1))/2.;
2430 break;
2431 case QCD::PHI:
2432 return (1. - ys*ys)/(1. + xs*xs) * (I_8(q2, 0) - I_8(q2, 1) - xs * s_8(q2, 0))/2.;
2433 break;
2434 default:
2435 std::stringstream out;
2436 out << vectorM;
2437 throw std::runtime_error("MVll::getDelta8 : vector " + out.str() + " not implemented");
2438 }
2439 };
2440
2446 double getDelta9(double q2)
2447 {
2448 switch(vectorM){
2449 case QCD::K_star:
2450 return (I_9(q2, 0) - I_9(q2, 1))/2.;
2451 break;
2452 case QCD::K_star_P:
2453 return (I_9(q2, 0) - I_9(q2, 1))/2.;
2454 break;
2455 case QCD::PHI:
2456 return (1. - ys*ys)/(1. + xs*xs) * (I_9(q2, 0) - I_9(q2, 1) - xs * s_9(q2, 0))/2.;
2457 break;
2458 default:
2459 std::stringstream out;
2460 out << vectorM;
2461 throw std::runtime_error("MVll::getDelta9 : vector " + out.str() + " not implemented");
2462 }
2463 };
2464
2470 double SigmaTree(double q2);
2471
2476 double getintegratedSigmaTree();
2477
2478 gslpp::complex A_Seidel(double q2, double mb2);
2479
2480 gslpp::complex B_Seidel(double q2, double mb2);
2481
2482 gslpp::complex C_Seidel(double q2);
2483
2489 gslpp::complex deltaC7_QCDF(double q2, bool conjugate, bool spline = false);
2490
2496 gslpp::complex deltaC9_QCDF(double q2, bool conjugate, bool spline = false);
2497
2503 gslpp::complex Cq34(bool conjugate);
2504
2509 gslpp::complex T_para_minus_WA(bool conjugate);
2510
2515 gslpp::complex T_perp_WA_1();
2516
2522 gslpp::complex T_perp_WA_2(bool conjugate);
2523
2529 gslpp::complex T_perp_plus_O8(double q2, double u);
2530
2536 gslpp::complex T_para_minus_O8(double q2, double u);
2537
2544 gslpp::complex t_perp(double q2, double u, double m2);
2545
2552 gslpp::complex t_para(double q2, double u, double m2);
2553
2554 gslpp::complex I1(double q2, double u, double m2);
2555
2556 gslpp::complex B0diff(double q2, double u, double m2);
2557
2558 gslpp::complex B0(double s, double m2);
2559
2560 gslpp::complex h_func(double s, double m2);
2561
2568 gslpp::complex T_perp_plus_QSS(double q2, double u, bool conjugate);
2569
2576 gslpp::complex T_para_plus_QSS(double q2, double u, bool conjugate);
2577
2584 gslpp::complex T_para_minus_QSS(double q2, double u, bool conjugate);
2585
2591 double phi_V(double u);
2592
2593 gslpp::complex lambda_B_minus(double q2);
2594
2602 double T_perp_real(double q2, double u, bool conjugate);
2603
2611 double T_perp_imag(double q2, double u, bool conjugate);
2612
2620 double T_para_real(double q2, double u, bool conjugate);
2621
2629 double T_para_imag(double q2, double u, bool conjugate);
2630
2637 double T_perp_real(double q2, bool conjugate);
2638
2645 double T_perp_imag(double q2, bool conjugate);
2646
2653 double T_para_real(double q2, bool conjugate);
2654
2661 double T_para_imag(double q2, bool conjugate);
2662
2663 double QCDF_fit_func(double* x, double* p);
2664
2665 void fit_QCDF_func();
2666
2667 void spline_QCDF_func();
2668
2669 gslpp::complex T_minus(double q2, bool conjugate);
2670
2671 gslpp::complex T_0(double q2, bool conjugate);
2672
2682 double FF_fit(double q2, double a_0, double a_1, double a_2, double MR2);
2683
2684};
2685
2686#endif /* MVLL_H */
2687
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:1778
gslpp::complex getQCDf_1(double q2)
Definition: MVll.h:569
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:1671
bool FixedWCbtos
Definition: MVll.h:788
std::vector< std::string > mvllParameters
Definition: MVll.h:783
const StandardModel & mySM
Definition: MVll.h:779
double xs
Definition: MVll.h:811
double mu_h
Definition: MVll.h:803
double getDelta_C9_zExp_p()
The non-pertubative ccbar contributions to the helicity amplitudes.
Definition: MVll.h:489
bool zExpansion
Definition: MVll.h:787
double getgtilde_2_re(double q2)
The real part of .
Definition: MVll.h:653
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:1935
double getQCDfC9_2(double q2, double cutoff)
Definition: MVll.h:596
void spline_QCDF_func()
Definition: MVll.cpp:2097
double getQCDfC9p_1(double cutoff)
Definition: MVll.h:608
double getDelta_C9_zExp_m()
The non-pertubative ccbar contributions to the helicity amplitudes.
Definition: MVll.h:499
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:1817
gslpp::complex B_Seidel(double q2, double mb2)
Definition: MVll.cpp:1633
double geth_0_re(double q2)
The real part of .
Definition: MVll.h:699
bool MVll_DM_flag
Definition: MVll.h:790
double geth_p_im(double q2)
The imaginary part of .
Definition: MVll.h:738
gslpp::complex H_A_p(double q2, bool bar)
The helicity amplitude .
Definition: MVll.cpp:2502
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:1869
double ale
Definition: MVll.h:797
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:1972
gslpp::complex T_perp_WA_1()
QCDF Correction from various BFS paper (hep-ph/0412400). Weak Annihilation.
Definition: MVll.cpp:1783
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:1719
double Mb
Definition: MVll.h:801
std::unique_ptr< F_2 > myF_2
Definition: MVll.h:785
gslpp::complex Cq34(bool conjugate)
QCDF Correction from various BFS paper (hep-ph/0412400). Part of Weak Annihilation.
Definition: MVll.cpp:1768
double QCDF_fit_func(double *x, double *p)
Definition: MVll.cpp:2028
double mPsi2S2
Definition: MVll.h:792
double MM
Definition: MVll.h:799
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:1946
double getgtilde_2_im(double q2)
The immaginary part of .
Definition: MVll.h:664
double geth_m_im(double q2)
The imaginary part of .
Definition: MVll.h:767
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:1891
gslpp::complex getQCDf_2(double q2)
Definition: MVll.h:576
double getVm(double q2)
The form factor .
Definition: MVll.h:411
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:1801
gslpp::complex C_Seidel(double q2)
Definition: MVll.cpp:1665
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:2522
double mD2
Definition: MVll.h:793
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:784
double integrateDelta(int i, double q_min, double q_max)
The integral of from to .
Definition: MVll.cpp:2890
double width
Definition: MVll.h:809
double getgtilde_3_im(double q2)
The imaginary part of .
Definition: MVll.h:687
double alpha_s_mub
Definition: MVll.h:814
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:2461
QCD::meson meson
Definition: MVll.h:781
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:1984
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:2033
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:1959
bool dispersion
Definition: MVll.h:786
gslpp::complex h_func(double s, double m2)
Definition: MVll.cpp:1855
double GF
Definition: MVll.h:796
gslpp::complex T_minus(double q2, bool conjugate)
Definition: MVll.cpp:2166
double getSigma(int i, double q_2)
The value of from to .
Definition: MVll.cpp:2845
gslpp::complex geth_p_0()
.
Definition: MVll.h:718
int etaV
Definition: MVll.h:813
gslpp::complex H_V_p(double q2, bool bar)
The helicity amplitude .
Definition: MVll.cpp:2472
double getgtilde_3_re(double q2)
The real part of .
Definition: MVll.h:675
gslpp::complex getQCDf_3(double q2)
Definition: MVll.h:583
gslpp::complex lambda_B_minus(double q2)
Definition: MVll.cpp:1940
double getwidth()
The width of the meson M.
Definition: MVll.h:363
gslpp::complex T_0(double q2, bool conjugate)
Definition: MVll.cpp:2188
double getV0(double q2)
The form factor .
Definition: MVll.h:389
double Ms
Definition: MVll.h:807
double mPsi2S
Definition: MVll.h:792
gslpp::complex h_lambda(int hel, double q2)
The non-pertubative ccbar contributions to the helicity amplitudes.
Definition: MVll.cpp:2404
gslpp::complex exp_Phase[3]
Definition: MVll.h:794
double getDelta_C9_zExp_0()
The non-pertubative ccbar contributions to the helicity amplitudes.
Definition: MVll.h:479
double mJpsi
Definition: MVll.h:791
double MV
Definition: MVll.h:800
gslpp::complex geth_m_0()
.
Definition: MVll.h:747
double geth_p_re(double q2)
The real part of .
Definition: MVll.h:728
double getintegratedSigmaTree()
The integral of from 0 to .
Definition: MVll.cpp:3024
double integrateSigmaTree(double q_min, double q_max)
The integral of from to (arxiv/2301.06990)
Definition: MVll.cpp:2989
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:1913
double geth_0_im(double q2)
The imaginary part of .
Definition: MVll.h:709
double mc_pole
Definition: MVll.h:806
double getTp(double q2)
The form factor .
Definition: MVll.h:432
double angmomV
Definition: MVll.h:812
gslpp::complex T_perp_WA_2(bool conjugate)
QCDF Correction from various BFS paper (hep-ph/0412400). Weak Annihilation.
Definition: MVll.cpp:1788
double FF_fit(double q2, double a_0, double a_1, double a_2, double MR2)
The fit function from , .
Definition: MVll.cpp:1403
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:1808
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:782
gslpp::complex H_V_m(double q2, bool bar)
The helicity amplitude .
Definition: MVll.cpp:2482
gslpp::complex A_Seidel(double q2, double mb2)
Definition: MVll.cpp:1619
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:2538
double spectator_charge
Definition: MVll.h:808
double Mlep
Definition: MVll.h:798
gslpp::complex B0diff(double q2, double u, double m2)
Definition: MVll.cpp:1841
double Delta_C9_zExp(int hel)
The non-pertubative ccbar contributions to the helicity amplitudes.
Definition: MVll.cpp:2434
double geth_m_re(double q2)
The real part of .
Definition: MVll.h:757
gslpp::complex H_A_0(double q2, bool bar)
The helicity amplitude .
Definition: MVll.cpp:2492
double getQCDfC9p_2(double cutoff)
Definition: MVll.h:614
double SigmaTree(double q2)
Definition: MVll.cpp:3019
double getQCDfC9p_3(double cutoff)
Definition: MVll.h:620
QCD::lepton lep
Definition: MVll.h:780
gslpp::complex I1(double q2, double u, double m2)
Definition: MVll.cpp:1824
double getgtilde_1_im(double q2)
The immaginary part of .
Definition: MVll.h:642
double getQCDfC9_1(double q2, double cutoff)
Definition: MVll.h:590
gslpp::complex B0(double s, double m2)
Definition: MVll.cpp:1849
gslpp::complex H_A_m(double q2, bool bar)
The helicity amplitude .
Definition: MVll.cpp:2512
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:1793
gslpp::complex H_P(double q2, bool bar)
The helicity amplitude .
Definition: MVll.cpp:2530
double Mc
Definition: MVll.h:804
double integrateSigma(int i, double q_min, double q_max)
The integral of from to .
Definition: MVll.cpp:2727
bool NeutrinoTree_flag
Definition: MVll.h:789
double mb_pole
Definition: MVll.h:805
double getQCDfC9_3(double q2, double cutoff)
Definition: MVll.h:602
double beta(double q2)
The factor used in the angular coefficients .
Definition: MVll.cpp:2570
double mu_b
Definition: MVll.h:802
double ys
Definition: MVll.h:810
double mJ2
Definition: MVll.h:791
double getgtilde_1_re(double q2)
The real part of .
Definition: MVll.h:631
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.