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 
11 class StandardModel;
12 class F_1;
13 class 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 true
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 
308 class MVll {
309 public:
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 
355  double getwidth(){
356  updateParameters();
357  return width;
358  }
359 
364  double getMlep(){
365  updateParameters();
366  return Mlep;
367  }
368 
374  double beta (double q2);
375 
381  double getV0(double q2)
382  {
383  updateParameters();
384  return (2. * MM * sqrt(q2))/sqrt(lambda(q2)) * V_0t(q2);
385  };
386 
392  double getVp(double q2)
393  {
394  updateParameters();
395  return V_p(q2);
396  };
397 
403  double getVm(double q2)
404  {
405  return V_m(q2);
406  };
407 
413  double getT0(double q2)
414  {
415  updateParameters();
416  return twoMM3/sqrt(q2 * lambda(q2)) * T_0t(q2);
417  };
418 
424  double getTp(double q2)
425  {
426  updateParameters();
427  return T_p(q2);
428  };
429 
435  double getTm(double q2)
436  {
437  updateParameters();
438  return T_m(q2);
439  };
440 
446  double getS(double q2)
447  {
448  updateParameters();
449  return S_L_pre/sqrt(lambda(q2)) * S_L(q2);
450  };
451 
458  gslpp::complex h_lambda(int hel, double q2);
459 
465  gslpp::complex H_V_0(double q2, bool bar);
466 
472  gslpp::complex H_V_p(double q2, bool bar);
473 
479  gslpp::complex H_V_m(double q2, bool bar);
480 
486  gslpp::complex H_A_0(double q2, bool bar);
487 
493  gslpp::complex H_A_p(double q2, bool bar);
494 
500  gslpp::complex H_A_m(double q2, bool bar);
501 
507  gslpp::complex H_S(double q2, bool bar);
508 
514  gslpp::complex H_P(double q2, bool bar);
515 
517  {
518  updateParameters();
519 // 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)));
520  return 0.;
521  }
522 
524  {
525  updateParameters();
526 // 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)));
527  return 0.;
528  }
529 
531  {
532  updateParameters();
533 // 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))));
534  return 0.;
535  }
536 
537  double getQCDfC9_1(double q2, double cutoff)
538  {
539  updateParameters();
540  return (getQCDf_1(q2) - getQCDf_1(cutoff) * cutoff/q2).abs();
541  }
542 
543  double getQCDfC9_2(double q2, double cutoff)
544  {
545  updateParameters();
546  return (getQCDf_2(q2) - getQCDf_2(cutoff) * cutoff/q2).abs();
547  }
548 
549  double getQCDfC9_3(double q2, double cutoff)
550  {
551  updateParameters();
552  return (getQCDf_3(q2) - getQCDf_3(cutoff) * cutoff/q2).abs();
553  }
554 
555  double getQCDfC9p_1(double cutoff)
556  {
557  updateParameters();
558  return (getQCDf_1(cutoff) * cutoff).abs();
559  }
560 
561  double getQCDfC9p_2(double cutoff)
562  {
563  updateParameters();
564  return (getQCDf_2(cutoff) * cutoff).abs();
565  }
566 
567  double getQCDfC9p_3(double cutoff)
568  {
569  updateParameters();
570  return (getQCDf_3(cutoff) * cutoff).abs();
571  }
572 
578  double getgtilde_1_re(double q2)
579  {
580  updateParameters();
581  return C2_inv * (gtilde_1_pre/(sqrt(lambda(q2)) * V(q2)) * (h_lambda(2,q2)-h_lambda(1,q2))).real()/q2;
582  }
583 
589  double getgtilde_1_im(double q2)
590  {
591  updateParameters();
592  return C2_inv * (gtilde_1_pre/(sqrt(lambda(q2)) * V(q2)) * (h_lambda(2,q2)-h_lambda(1,q2))).imag()/q2;
593  }
594 
600  double getgtilde_2_re(double q2)
601  {
602  updateParameters();
603  return C2_inv * (gtilde_2_pre/A_1(q2) * (h_lambda(1,q2)+h_lambda(2,q2))).real()/q2;
604  }
605 
611  double getgtilde_2_im(double q2)
612  {
613  updateParameters();
614  return C2_inv * (gtilde_2_pre/A_1(q2) * (h_lambda(1,q2)+h_lambda(2,q2))).imag()/q2;
615  }
616 
622  double getgtilde_3_re(double q2)
623  {
624  updateParameters();
625  return C2_inv * (gtilde_3_pre/(lambda(q2) * A_2(q2)) * (sqrt(q2)*h_lambda(0,q2)/q2-
626  (MM2mMV2 - q2)/(4.*MV) * (h_lambda(1,q2)+h_lambda(2,q2))/q2)).real();
627  }
628 
634  double getgtilde_3_im(double q2)
635  {
636  updateParameters();
637  return C2_inv * (gtilde_3_pre/(lambda(q2) * A_2(q2)) * (sqrt(q2)*h_lambda(0,q2)/q2-
638  (MM2mMV2 - q2)/(4.*MV) * (h_lambda(1,q2)+h_lambda(2,q2))/q2)).imag();
639  }
640 
646  double geth_0_re(double q2)
647  {
648  return (sixteenM_PI2MM2 * h_lambda(0,q2)/q2).real();
649  }
650 
656  double geth_0_im(double q2)
657  {
658  return (sixteenM_PI2MM2 * h_lambda(0,q2)/q2).imag();
659  }
660 
666  {
667  return h_lambda(1,0.);
668  }
669 
675  double geth_p_re(double q2)
676  {
677  return (sixteenM_PI2MM2 * h_lambda(1,q2)/q2).real();
678  }
679 
685  double geth_p_im(double q2)
686  {
687  return (sixteenM_PI2MM2 * h_lambda(1,q2)/q2).imag();
688  }
689 
695  {
696  return h_lambda(2,0.);
697  }
698 
704  double geth_m_re(double q2)
705  {
706  return (sixteenM_PI2MM2 * h_lambda(2,q2)/q2).real();
707  }
708 
714  double geth_m_im(double q2)
715  {
716  return (sixteenM_PI2MM2 * h_lambda(2,q2)/q2).imag();
717  }
718 
723  std::vector<std::string> initializeMVllParameters();
724 
725 private:
730  std::vector<std::string> mvllParameters;
731  std::unique_ptr<F_1> myF_1;
732  std::unique_ptr<F_2> myF_2;
735  double mJ2;
737 
738  double GF;
739  double ale;
740  double Mlep;
741  double MM;
742  double MV;
743  double Mb;
744  double mu_b;
745  double mu_h;
746  double Mc;
747  double mb_pole;
748  double mc_pole;
749  double Ms;
751  double width;
752  double ys;
753  double xs;
754  double angmomV;
755  int etaV;
756  double alpha_s_mub;
757  double fB;
758  double fpara;
759  double fperp;
760 
761  double MW;
762  gslpp::complex lambda_t;
763  gslpp::complex lambda_u;
764  double b;
765  gslpp::complex h_0[3];
766  gslpp::complex h_1[3];
767  gslpp::complex h_2[3];
768  gslpp::complex SU3_breaking;
769 
770  double t_p;
771  double t_m;
772  double t_0;
773  double z_0;
774  double MMpMV;
775  double MMpMV2;
776  double MMmMV;
777  double MMmMV2;
778  double MM2;
779  double MM4;
780  double MV2;
781  double MV4;
782  double MMMV;
783  double MM2mMV2;
784  double MM2pMV2;
785  double fourMV;
786  double onepMMoMV;
787  double MM_MMpMV;
788  double twoMM2;
789  double twoMV2;
790  double twoMM_mbpms;
791  double fourMM2;
792  double Mlep2;
793  double twoMlepMb;
794  double MboMW;
795  double MsoMb;
796  double M_PI2osix;
797  double N_QCDF;
798  double twoMM;
799  double ninetysixM_PI3MM3;
800  double M_PI2;
801  double sixteenM_PI2;
802  double sixteenM_PI2MM2;
803  double twoMboMM;
804  gslpp::complex H_0_pre;
805  double mu_b2;
806  double Mc2;
807  double Mb2;
808  double fourMc2;
809  double fourMb2;
810  double logMc;
811  double logMb;
812  gslpp::complex H_0_WC;
813  gslpp::complex H_c_WC;
814  gslpp::complex H_b_WC;
815  double fournineth;
816  double half;
817  double twothird;
818  gslpp::complex ihalfMPI;
819  double twoMM3;
820  double gtilde_1_pre;
821  double gtilde_2_pre;
822  double gtilde_3_pre;
823  double C2_inv;
824  double S_L_pre;
825  gslpp::complex NN;
826  gslpp::complex NN_conjugate;
827  double CF;
828  double deltaT_0;
829  double deltaT_1par;
830  double deltaT_1perp;
832  bool h_pole;
833 
834  gslpp::complex ubar;
835  gslpp::complex arg1;
836  gslpp::complex B01;
837  gslpp::complex B00;
838  gslpp::complex xp;
839  gslpp::complex xm;
840  gslpp::complex yp;
841  gslpp::complex ym;
842  gslpp::complex L1xp;
843  gslpp::complex L1xm;
844  gslpp::complex L1yp;
845  gslpp::complex L1ym;
846  gslpp::complex F87_0;
847  gslpp::complex F87_1;
848  gslpp::complex F87_2;
849  gslpp::complex F87_3;
850  double F89_0;
851  double F89_1;
852  double F89_2;
853  double F89_3;
855  double a_0V;
856  double a_1V;
857  double a_2V;
858  double MRV_2;
859  double a_0A0;
860  double a_1A0;
861  double a_2A0;
862  double MRA0_2;
863  double a_0A1;
864  double a_1A1;
865  double a_2A1;
866  double MRA1_2;
867  double a_0A12;
868  double a_1A12;
869  double a_2A12;
870  double MRA12_2;
871  double a_0T1;
872  double a_1T1;
873  double a_2T1;
874  double MRT1_2;
875  double a_0T2;
876  double a_1T2;
877  double a_2T2;
878  double MRT2_2;
879  double a_0T23;
880  double a_1T23;
881  double a_2T23;
882  double MRT23_2;
884  gslpp::vector<gslpp::complex> ** allcoeff;
885  gslpp::vector<gslpp::complex> ** allcoeffh;
886  gslpp::vector<gslpp::complex> ** allcoeffprime;
888  gslpp::vector<gslpp::complex> ** allcoeff_noSM;
890  gslpp::complex C_1;
891  gslpp::complex C_1L_bar;
892  gslpp::complex C_1Lh_bar;
893  gslpp::complex C_2;
894  gslpp::complex C_2L_bar;
895  gslpp::complex C_2Lh_bar;
898  gslpp::complex C_5;
899  gslpp::complex C_6;
900  gslpp::complex C_7;
901  gslpp::complex C_8;
902  gslpp::complex C_8L;
903  gslpp::complex C_8Lh;
904  gslpp::complex C_9;
905  gslpp::complex C_10;
906  gslpp::complex C_S;
907  gslpp::complex C_P;
909  gslpp::complex C_7p;
910  gslpp::complex C_9p;
911  gslpp::complex C_10p;
912  gslpp::complex C_Sp;
913  gslpp::complex C_Pp;
916  std::vector<double> Re_T_perp;
917  std::vector<double> Im_T_perp;
918  std::vector<double> Re_T_para;
919  std::vector<double> Im_T_para;
921  gsl_interp_accel *acc_Re_T_perp;
922  gsl_interp_accel *acc_Im_T_perp;
923  gsl_interp_accel *acc_Re_T_para;
924  gsl_interp_accel *acc_Im_T_para;
925 
926  gsl_spline *spline_Re_T_perp;
927  gsl_spline *spline_Im_T_perp;
928  gsl_spline *spline_Re_T_para;
929  gsl_spline *spline_Im_T_para;
930 
931  gsl_interp_accel *acc_Re_deltaC7_QCDF;
932  gsl_interp_accel *acc_Im_deltaC7_QCDF;
933  gsl_interp_accel *acc_Re_deltaC9_QCDF;
934  gsl_interp_accel *acc_Im_deltaC9_QCDF;
935 
936  gsl_spline *spline_Re_deltaC7_QCDF;
937  gsl_spline *spline_Im_deltaC7_QCDF;
938  gsl_spline *spline_Re_deltaC9_QCDF;
939  gsl_spline *spline_Im_deltaC9_QCDF;
940 
941 #if COMPUTECP
942  std::vector<double> Re_T_perp_conj;
943  std::vector<double> Im_T_perp_conj;
944  std::vector<double> Re_T_para_conj;
945  std::vector<double> Im_T_para_conj;
947  gsl_interp_accel *acc_Re_T_perp_conj;
948  gsl_interp_accel *acc_Im_T_perp_conj;
949  gsl_interp_accel *acc_Re_T_para_conj;
950  gsl_interp_accel *acc_Im_T_para_conj;
951 
952  gsl_interp_accel *acc_Re_deltaC7_QCDF_conj;
953  gsl_interp_accel *acc_Im_deltaC7_QCDF_conj;
954  gsl_interp_accel *acc_Re_deltaC9_QCDF_conj;
955  gsl_interp_accel *acc_Im_deltaC9_QCDF_conj;
956 
957  gsl_spline *spline_Re_T_perp_conj;
958  gsl_spline *spline_Im_T_perp_conj;
959  gsl_spline *spline_Re_T_para_conj;
960  gsl_spline *spline_Im_T_para_conj;
961 
962  gsl_spline *spline_Re_deltaC7_QCDF_conj;
963  gsl_spline *spline_Im_deltaC7_QCDF_conj;
964  gsl_spline *spline_Re_deltaC9_QCDF_conj;
965  gsl_spline *spline_Im_deltaC9_QCDF_conj;
966 #endif
967 
968  std::vector<double> myq2;
970  TFitResultPtr Re_T_perp_res;
971  TFitResultPtr Im_T_perp_res;
972  TFitResultPtr Re_T_para_res;
973  TFitResultPtr Im_T_para_res;
975  TFitResultPtr Re_T_perp_res_conj;
976  TFitResultPtr Im_T_perp_res_conj;
977  TFitResultPtr Re_T_para_res_conj;
978  TFitResultPtr Im_T_para_res_conj;
980  TGraph gr1;
981  TGraph gr2;
983  TF1 QCDFfit;
985  TF1 reffit;
986  TF1 imffit;
988  double avaSigma;
989  double avaDelta;
991  double errSigma;
992  double errDelta;
994  gsl_function FS;
995  gsl_function FD;
997  gsl_integration_cquad_workspace * w_sigma;
998  gsl_integration_cquad_workspace * w_delta;
1000  gsl_error_handler_t * old_handler;
1002  std::map<std::pair<double, double>, gslpp::complex > cacheI1;
1004  std::map<std::pair<double, double>, double > cacheSigma0;
1005  std::map<std::pair<double, double>, double > cacheSigma1;
1006  std::map<std::pair<double, double>, double > cacheSigma2;
1007  std::map<std::pair<double, double>, double > cacheSigma3;
1008  std::map<std::pair<double, double>, double > cacheSigma4;
1009  std::map<std::pair<double, double>, double > cacheSigma5;
1010  std::map<std::pair<double, double>, double > cacheSigma6;
1011  std::map<std::pair<double, double>, double > cacheSigma7;
1012  std::map<std::pair<double, double>, double > cacheSigma9;
1013  std::map<std::pair<double, double>, double > cacheSigma10;
1014  std::map<std::pair<double, double>, double > cacheSigma11;
1016  std::map<std::pair<double, double>, double > cacheDelta0;
1017  std::map<std::pair<double, double>, double > cacheDelta1;
1018  std::map<std::pair<double, double>, double > cacheDelta2;
1019  std::map<std::pair<double, double>, double > cacheDelta3;
1020  std::map<std::pair<double, double>, double > cacheDelta6;
1021  std::map<std::pair<double, double>, double > cacheDelta7;
1022  std::map<std::pair<double, double>, double > cacheDelta8;
1023  std::map<std::pair<double, double>, double > cacheDelta10;
1024  std::map<std::pair<double, double>, double > cacheDelta11;
1026  unsigned int N_updated;
1027  gslpp::vector<double> N_cache;
1028  gslpp::complex Nc_cache;
1030  unsigned int V_updated;
1031  gslpp::vector<double> V_cache;
1033  unsigned int A0_updated;
1034  gslpp::vector<double> A0_cache;
1036  unsigned int A1_updated;
1037  gslpp::vector<double> A1_cache;
1039  unsigned int T1_updated;
1040  gslpp::vector<double> T1_cache;
1042  unsigned int T2_updated;
1043  gslpp::vector<double> T2_cache;
1045  unsigned int k2_updated;
1046  gslpp::vector<double> k2_cache;
1048  unsigned int z_updated;
1050  unsigned int lambda_updated;
1052  unsigned int beta_updated;
1053  double beta_cache;
1055  unsigned int F_updated;
1057  unsigned int VL1_updated;
1058  unsigned int VL2_updated;
1060  unsigned int TL1_updated;
1061  unsigned int TL2_updated;
1063  unsigned int VR1_updated;
1064  unsigned int VR2_updated;
1066  unsigned int TR1_updated;
1067  unsigned int TR2_updated;
1069  unsigned int VL0_updated;
1070  gslpp::vector<double> VL0_cache;
1072  unsigned int TL0_updated;
1073  gslpp::vector<double> TL0_cache;
1075  unsigned int VR0_updated;
1077  unsigned int TR0_updated;
1079  unsigned int Mb_Ms_updated;
1081  unsigned int SL_updated;
1082  gslpp::vector<double> SL_cache;
1084  unsigned int SR_updated;
1086  unsigned int C_1_updated;
1087  gslpp::complex C_1_cache;
1089  unsigned int C_2_updated;
1090  gslpp::complex C_2_cache;
1092  unsigned int C_3_updated;
1093  gslpp::complex C_3_cache;
1095  unsigned int C_4_updated;
1096  gslpp::complex C_4_cache;
1098  unsigned int C_5_updated;
1099  gslpp::complex C_5_cache;
1101  unsigned int C_6_updated;
1102  gslpp::complex C_6_cache;
1104  unsigned int C_7_updated;
1105  gslpp::complex C_7_cache;
1107  unsigned int C_9_updated;
1108  gslpp::complex C_9_cache;
1110  unsigned int C_10_updated;
1111  gslpp::complex C_10_cache;
1113  unsigned int C_7p_updated;
1114  gslpp::complex C_7p_cache;
1116  unsigned int C_9p_updated;
1117  gslpp::complex C_9p_cache;
1119  unsigned int C_10p_updated;
1120  gslpp::complex C_10p_cache;
1122  unsigned int C_S_updated;
1123  gslpp::complex C_S_cache;
1125  unsigned int C_P_updated;
1126  gslpp::complex C_P_cache;
1128  unsigned int C_Sp_updated;
1129  gslpp::complex C_Sp_cache;
1131  unsigned int C_Pp_updated;
1132  gslpp::complex C_Pp_cache;
1134  unsigned int C_2Lh_updated;
1135  gslpp::complex C_2Lh_cache;
1137  unsigned int C_8Lh_updated;
1138  gslpp::complex C_8Lh_cache;
1140  unsigned int Yupdated;
1141  gslpp::vector<double> Ycache;
1143  gslpp::complex h0Ccache[4];
1144  gslpp::complex h1Ccache[4];
1145  gslpp::complex h2Ccache[4];
1147  unsigned int h0_updated;
1148  unsigned int h1_updated;
1149  unsigned int h2_updated;
1151  unsigned int H_V0updated;
1152  gslpp::vector<double> H_V0cache;
1154  unsigned int H_V1updated;
1155  gslpp::vector<double> H_V1cache;
1157  unsigned int H_V2updated;
1158  gslpp::vector<double> H_V2cache;
1160  unsigned int H_A0updated;
1161  unsigned int H_A1updated;
1162  unsigned int H_A2updated;
1164  unsigned int H_Supdated;
1165  gslpp::vector<double> H_Scache;
1167  unsigned int H_Pupdated;
1168  gslpp::vector<double> H_Pcache;
1170  unsigned int I0_updated;
1171  unsigned int I1_updated;
1172  unsigned int I2_updated;
1173  unsigned int I3_updated;
1174  unsigned int I4_updated;
1175  unsigned int I5_updated;
1176  unsigned int I6_updated;
1177  unsigned int I7_updated;
1178  unsigned int I8_updated;
1179  unsigned int I9_updated;
1180  unsigned int I10_updated;
1181  unsigned int I11_updated;
1183  std::map<std::pair<double, double>, unsigned int > I1Cached;
1185  std::map<std::pair<double, double>, unsigned int > sigma0Cached;
1186  std::map<std::pair<double, double>, unsigned int > sigma1Cached;
1187  std::map<std::pair<double, double>, unsigned int > sigma2Cached;
1188  std::map<std::pair<double, double>, unsigned int > sigma3Cached;
1189  std::map<std::pair<double, double>, unsigned int > sigma4Cached;
1190  std::map<std::pair<double, double>, unsigned int > sigma5Cached;
1191  std::map<std::pair<double, double>, unsigned int > sigma6Cached;
1192  std::map<std::pair<double, double>, unsigned int > sigma7Cached;
1193  std::map<std::pair<double, double>, unsigned int > sigma9Cached;
1194  std::map<std::pair<double, double>, unsigned int > sigma10Cached;
1195  std::map<std::pair<double, double>, unsigned int > sigma11Cached;
1197  std::map<std::pair<double, double>, unsigned int > delta0Cached;
1198  std::map<std::pair<double, double>, unsigned int > delta1Cached;
1199  std::map<std::pair<double, double>, unsigned int > delta2Cached;
1200  std::map<std::pair<double, double>, unsigned int > delta3Cached;
1201  std::map<std::pair<double, double>, unsigned int > delta6Cached;
1202  std::map<std::pair<double, double>, unsigned int > delta7Cached;
1203  std::map<std::pair<double, double>, unsigned int > delta8Cached;
1204  std::map<std::pair<double, double>, unsigned int > delta10Cached;
1205  std::map<std::pair<double, double>, unsigned int > delta11Cached;
1207  std::map<double, unsigned int> deltaTparpCached;
1208  std::map<double, unsigned int> deltaTparmCached;
1209  std::map<double, unsigned int> deltaTperpCached;
1211  std::map<double, gslpp::complex> cacheDeltaTparp;
1212  std::map<double, gslpp::complex> cacheDeltaTparm;
1213  std::map<double, gslpp::complex> cacheDeltaTperp;
1215  unsigned int deltaTparpupdated;
1216  unsigned int deltaTparmupdated;
1217  unsigned int deltaTperpupdated;
1219  unsigned int T_updated;
1220  gslpp::vector<double> T_cache;
1225  void updateParameters();
1226 
1230  void checkCache();
1231 
1237  double z(double q2);
1238 
1244  double V(double q2);
1245 
1246 
1252  double A_0(double q2);
1253 
1254 
1260  double A_1(double q2);
1261 
1267  double A_2(double q2);
1268 
1274  double T_1(double q2);
1275 
1276 
1282  double T_2(double q2);
1283 
1289  double V_0t(double q2);
1290 
1296  double V_p(double q2);
1297 
1303  double V_m(double q2);
1304 
1310  double T_0t(double q2);
1311 
1317  double T_p(double q2);
1318 
1324  double T_m(double q2);
1325 
1331  double S_L(double q2);
1332 
1338  gslpp::complex H_0(double q2);
1339 
1347  gslpp::complex H(double q2, double m2, double mu2);
1348 
1354  gslpp::complex Y(double q2);
1355 
1356  gslpp::complex funct_g(double q2);
1357 
1358  gslpp::complex DeltaC9_KD(double q2, int com);
1359 
1365  double k2 (double q2);
1366 
1372  double beta2 (double q2);
1373 
1379  double lambda(double q2);
1380 
1387  double F(double q2, double b_i);
1388 
1395  double I_1c(double q2, bool bar);
1396 
1403  double I_1s(double q2, bool bar);
1404 
1411  double I_2c(double q2, bool bar);
1412 
1419  double I_2s(double q2, bool bar);
1420 
1427  double I_3(double q2, bool bar);
1428 
1435  double I_4(double q2, bool bar);
1436 
1443  double I_5(double q2, bool bar);
1444 
1451  double I_6c(double q2, bool bar);
1452 
1459  double I_6s(double q2, bool bar);
1460 
1467  double I_7(double q2, bool bar);
1468 
1475  double I_8(double q2, bool bar);
1476 
1483  double I_9(double q2, bool bar);
1484 
1491  double h_1s(double q2, bool bar);
1492 
1499  double h_1c(double q2, bool bar);
1500 
1507  double h_2s(double q2, bool bar);
1508 
1515  double h_2c(double q2, bool bar);
1516 
1523  double h_3(double q2, bool bar);
1524 
1531  double h_4(double q2, bool bar);
1532 
1539  double h_7(double q2, bool bar);
1540 
1547  double s_5(double q2, bool bar);
1548 
1555  double s_6s(double q2, bool bar);
1556 
1563  double s_6c(double q2, bool bar);
1564 
1571  double s_8(double q2, bool bar);
1572 
1579  double s_9(double q2, bool bar);
1580 
1586  double getSigma1c(double q2)
1587  {
1588  switch(vectorM){
1589  case QCD::K_star:
1590  return (I_1c(q2, 0) + I_1c(q2, 1))/2.;
1591  break;
1592  case QCD::K_star_P:
1593  return (I_1c(q2, 0) + I_1c(q2, 1))/2.;
1594  break;
1595  case QCD::PHI:
1596  return (I_1c(q2, 0) + I_1c(q2, 1) - ys * h_1c(q2, 0) )/2.;
1597  break;
1598  default:
1599  std::stringstream out;
1600  out << vectorM;
1601  throw std::runtime_error("MVll::getSigma1c : vector " + out.str() + " not implemented");
1602  }
1603  };
1604 
1610  double getSigma1s(double q2)
1611  {
1612  switch(vectorM){
1613  case QCD::K_star:
1614  return (I_1s(q2, 0) + I_1s(q2, 1))/2.;
1615  break;
1616  case QCD::K_star_P:
1617  return (I_1s(q2, 0) + I_1s(q2, 1))/2.;
1618  break;
1619  case QCD::PHI:
1620  return (I_1s(q2, 0) + I_1s(q2, 1) - ys * h_1s(q2, 0))/2.;
1621  break;
1622  default:
1623  std::stringstream out;
1624  out << vectorM;
1625  throw std::runtime_error("MVll::getSigma1s : vector " + out.str() + " not implemented");
1626  }
1627  };
1628 
1634  double getSigma2c(double q2)
1635  {
1636  switch(vectorM){
1637  case QCD::K_star:
1638  return (I_2c(q2, 0) + I_2c(q2, 1))/2.;
1639  break;
1640  case QCD::K_star_P:
1641  return (I_2c(q2, 0) + I_2c(q2, 1))/2.;
1642  break;
1643  case QCD::PHI:
1644  return (I_2c(q2, 0) + I_2c(q2, 1) - ys * h_2c(q2, 0))/2.;
1645  break;
1646  default:
1647  std::stringstream out;
1648  out << vectorM;
1649  throw std::runtime_error("MVll::getSigma2c : vector " + out.str() + " not implemented");
1650  }
1651  };
1652 
1658  double getSigma2s(double q2)
1659  {
1660  switch(vectorM){
1661  case QCD::K_star:
1662  return (I_2s(q2, 0) + I_2s(q2, 1))/2.;
1663  break;
1664  case QCD::K_star_P:
1665  return (I_2s(q2, 0) + I_2s(q2, 1))/2.;
1666  break;
1667  case QCD::PHI:
1668  return (I_2s(q2, 0) + I_2s(q2, 1) - ys * h_2s(q2, 0))/2.;
1669  break;
1670  default:
1671  std::stringstream out;
1672  out << vectorM;
1673  throw std::runtime_error("MVll::getSigma2s : vector " + out.str() + " not implemented");
1674  }
1675  };
1676 
1682  double getSigma3(double q2)
1683  {
1684  switch(vectorM){
1685  case QCD::K_star:
1686  return (I_3(q2, 0) + I_3(q2, 1))/2.;
1687  break;
1688  case QCD::K_star_P:
1689  return (I_3(q2, 0) + I_3(q2, 1))/2.;
1690  break;
1691  case QCD::PHI:
1692  return (I_3(q2, 0) + I_3(q2, 1) - ys * h_3(q2, 0))/2.;
1693  break;
1694  default:
1695  std::stringstream out;
1696  out << vectorM;
1697  throw std::runtime_error("MVll::getSigma3 : vector " + out.str() + " not implemented");
1698  }
1699  };
1700 
1706  double getSigma4(double q2)
1707  {
1708  switch(vectorM){
1709  case QCD::K_star:
1710  return (I_4(q2, 0) + I_4(q2, 1))/2.;
1711  break;
1712  case QCD::K_star_P:
1713  return (I_4(q2, 0) + I_4(q2, 1))/2.;
1714  break;
1715  case QCD::PHI:
1716  return (I_4(q2, 0) + I_4(q2, 1) - ys * h_4(q2, 0))/2.;
1717  break;
1718  default:
1719  std::stringstream out;
1720  out << vectorM;
1721  throw std::runtime_error("MVll::getSigma4 : vector " + out.str() + " not implemented");
1722  }
1723  };
1724 
1730  double getSigma5(double q2)
1731  {
1732  return (I_5(q2, 0) + I_5(q2, 1))/2.;
1733  };
1734 
1740  double getSigma6s(double q2)
1741  {
1742  return (I_6s(q2, 0) + I_6s(q2, 1))/2.;
1743  };
1744 
1750  double getSigma6c(double q2)
1751  {
1752  return (I_6c(q2, 0) + I_6c(q2, 1))/2.;
1753  };
1754 
1760  double getSigma7(double q2)
1761  {
1762  switch(vectorM){
1763  case QCD::K_star:
1764  return (I_7(q2, 0) + I_7(q2, 1))/2.;
1765  break;
1766  case QCD::K_star_P:
1767  return (I_7(q2, 0) + I_7(q2, 1))/2.;
1768  break;
1769  case QCD::PHI:
1770  return (I_7(q2, 0) + I_7(q2, 1) - ys * h_7(q2, 0))/2.;
1771  break;
1772  default:
1773  std::stringstream out;
1774  out << vectorM;
1775  throw std::runtime_error("MVll::getSigma7 : vector " + out.str() + " not implemented");
1776  }
1777  };
1778 
1784  double getSigma8(double q2)
1785  {
1786  return (I_8(q2, 0) + I_8(q2, 1))/2.;
1787  };
1788 
1794  double getSigma9(double q2)
1795  {
1796  return (I_9(q2, 0) + I_9(q2, 1))/2.;
1797  };
1798 
1804  double getDelta1c(double q2)
1805  {
1806  return (I_1c(q2, 0) - I_1c(q2, 1)) / 2.;
1807  };
1808 
1814  double getDelta1s(double q2)
1815  {
1816  return (I_1s(q2, 0) - I_1s(q2, 1)) / 2.;
1817  };
1818 
1824  double getDelta2c(double q2)
1825  {
1826  return (I_2c(q2, 0) - I_2c(q2, 1)) / 2.;
1827  };
1828 
1834  double getDelta2s(double q2)
1835  {
1836  return (I_2s(q2, 0) - I_2s(q2, 1))/2.;
1837  };
1838 
1844  double getDelta3(double q2)
1845  {
1846  return (I_3(q2, 0) - I_3(q2, 1))/2.;
1847  };
1848 
1854  double getDelta4(double q2)
1855  {
1856  return (I_4(q2, 0) - I_4(q2, 1))/2.;
1857  };
1858 
1864  double getDelta5(double q2)
1865  {
1866  switch(vectorM){
1867  case QCD::K_star:
1868  return (I_5(q2, 0) - I_5(q2, 1))/2.;
1869  break;
1870  case QCD::K_star_P:
1871  return (I_5(q2, 0) - I_5(q2, 1))/2.;
1872  break;
1873  case QCD::PHI:
1874  return (1. - ys*ys)/(1. + xs*xs) * (I_5(q2, 0) - I_5(q2, 1) - xs * s_5(q2, 0))/2.;
1875  break;
1876  default:
1877  std::stringstream out;
1878  out << vectorM;
1879  throw std::runtime_error("MVll::getDelta5 : vector " + out.str() + " not implemented");
1880  }
1881  };
1882 
1888  double getDelta6s(double q2)
1889  {
1890  switch(vectorM){
1891  case QCD::K_star:
1892  return (I_6s(q2, 0) - I_6s(q2, 1))/2.;
1893  break;
1894  case QCD::K_star_P:
1895  return (I_6s(q2, 0) - I_6s(q2, 1))/2.;
1896  break;
1897  case QCD::PHI:
1898  return (1. - ys*ys)/(1. + xs*xs) * (I_6s(q2, 0) - I_6s(q2, 1) - xs * s_6s(q2, 0))/2.;
1899  break;
1900  default:
1901  std::stringstream out;
1902  out << vectorM;
1903  throw std::runtime_error("MVll::getDelta6s : vector " + out.str() + " not implemented");
1904  }
1905  };
1906 
1912  double getDelta6c(double q2)
1913  {
1914  switch(vectorM){
1915  case QCD::K_star:
1916  return (I_6c(q2, 0) - I_6c(q2, 1))/2.;
1917  break;
1918  case QCD::K_star_P:
1919  return (I_6c(q2, 0) - I_6c(q2, 1))/2.;
1920  break;
1921  case QCD::PHI:
1922  return (1. - ys*ys)/(1. + xs*xs) * (I_6c(q2, 0) - I_6c(q2, 1) - xs * s_6c(q2, 0))/2.;
1923  break;
1924  default:
1925  std::stringstream out;
1926  out << vectorM;
1927  throw std::runtime_error("MVll::getDelta6c : vector " + out.str() + " not implemented");
1928  }
1929  };
1930 
1936  double getDelta7(double q2)
1937  {
1938  return (I_7(q2, 0) - I_7(q2, 1))/2.;
1939  };
1940 
1946  double getDelta8(double q2)
1947  {
1948  switch(vectorM){
1949  case QCD::K_star:
1950  return (I_8(q2, 0) - I_8(q2, 1))/2.;
1951  break;
1952  case QCD::K_star_P:
1953  return (I_8(q2, 0) - I_8(q2, 1))/2.;
1954  break;
1955  case QCD::PHI:
1956  return (1. - ys*ys)/(1. + xs*xs) * (I_8(q2, 0) - I_8(q2, 1) - xs * s_8(q2, 0))/2.;
1957  break;
1958  default:
1959  std::stringstream out;
1960  out << vectorM;
1961  throw std::runtime_error("MVll::getDelta8 : vector " + out.str() + " not implemented");
1962  }
1963  };
1964 
1970  double getDelta9(double q2)
1971  {
1972  switch(vectorM){
1973  case QCD::K_star:
1974  return (I_9(q2, 0) - I_9(q2, 1))/2.;
1975  break;
1976  case QCD::K_star_P:
1977  return (I_9(q2, 0) - I_9(q2, 1))/2.;
1978  break;
1979  case QCD::PHI:
1980  return (1. - ys*ys)/(1. + xs*xs) * (I_9(q2, 0) - I_9(q2, 1) - xs * s_9(q2, 0))/2.;
1981  break;
1982  default:
1983  std::stringstream out;
1984  out << vectorM;
1985  throw std::runtime_error("MVll::getDelta9 : vector " + out.str() + " not implemented");
1986  }
1987  };
1988 
1989  gslpp::complex A_Seidel(double q2, double mb2);
1990 
1991  gslpp::complex B_Seidel(double q2, double mb2);
1992 
1993  gslpp::complex C_Seidel(double q2);
1994 
2000  gslpp::complex deltaC7_QCDF(double q2, bool conjugate, bool spline = false);
2001 
2007  gslpp::complex deltaC9_QCDF(double q2, bool conjugate, bool spline = false);
2008 
2014  gslpp::complex Cq34(bool conjugate);
2015 
2020  gslpp::complex T_para_minus_WA(bool conjugate);
2021 
2026  gslpp::complex T_perp_WA_1();
2027 
2033  gslpp::complex T_perp_WA_2(bool conjugate);
2034 
2040  gslpp::complex T_perp_plus_O8(double q2, double u);
2041 
2047  gslpp::complex T_para_minus_O8(double q2, double u);
2048 
2055  gslpp::complex t_perp(double q2, double u, double m2);
2056 
2063  gslpp::complex t_para(double q2, double u, double m2);
2064 
2065  gslpp::complex I1(double q2, double u, double m2);
2066 
2067  gslpp::complex B0diff(double q2, double u, double m2);
2068 
2069  gslpp::complex B0(double s, double m2);
2070 
2071  gslpp::complex h_func(double s, double m2);
2072 
2079  gslpp::complex T_perp_plus_QSS(double q2, double u, bool conjugate);
2080 
2087  gslpp::complex T_para_plus_QSS(double q2, double u, bool conjugate);
2088 
2095  gslpp::complex T_para_minus_QSS(double q2, double u, bool conjugate);
2096 
2102  double phi_V(double u);
2103 
2104  gslpp::complex lambda_B_minus(double q2);
2105 
2113  double T_perp_real(double q2, double u, bool conjugate);
2114 
2122  double T_perp_imag(double q2, double u, bool conjugate);
2123 
2131  double T_para_real(double q2, double u, bool conjugate);
2132 
2140  double T_para_imag(double q2, double u, bool conjugate);
2141 
2148  double T_perp_real(double q2, bool conjugate);
2149 
2156  double T_perp_imag(double q2, bool conjugate);
2157 
2164  double T_para_real(double q2, bool conjugate);
2165 
2172  double T_para_imag(double q2, bool conjugate);
2173 
2174  double QCDF_fit_func(double* x, double* p);
2175 
2176  void fit_QCDF_func();
2177 
2178  void spline_QCDF_func();
2179 
2180  gslpp::complex T_minus(double q2, bool conjugate);
2181 
2182  gslpp::complex T_0(double q2, bool conjugate);
2183 
2193  double FF_fit(double q2, double a_0, double a_1, double a_2, double MR2);
2194 
2195 };
2196 
2197 #endif /* MVLL_H */
2198 
MVll::H_A_0
gslpp::complex H_A_0(double q2, bool bar)
The helicity amplitude .
Definition: MVll.cpp:1644
MVll::H_A_m
gslpp::complex H_A_m(double q2, bool bar)
The helicity amplitude .
Definition: MVll.cpp:1656
MVll::geth_m_im
double geth_m_im(double q2)
The imaginary part of .
Definition: MVll.h:714
MVll::GF
double GF
Definition: MVll.h:738
MVll::getQCDfC9p_3
double getQCDfC9p_3(double cutoff)
Definition: MVll.h:567
MVll::vectorM
QCD::meson vectorM
Definition: MVll.h:729
MVll::initializeMVllParameters
std::vector< std::string > initializeMVllParameters()
A method for initializing the parameters necessary for MVll.
Definition: MVll.cpp:149
MVll::geth_p_re
double geth_p_re(double q2)
The real part of .
Definition: MVll.h:675
MVll::MVll
MVll(const StandardModel &SM_i, QCD::meson meson_i, QCD::meson vector_i, QCD::lepton lep_i)
Constructor.
Definition: MVll.cpp:21
MVll::MV
double MV
Definition: MVll.h:742
MVll::H_P
gslpp::complex H_P(double q2, bool bar)
The helicity amplitude .
Definition: MVll.cpp:1668
MVll::getQCDfC9_2
double getQCDfC9_2(double q2, double cutoff)
Definition: MVll.h:543
MVll::h_lambda
gslpp::complex h_lambda(int hel, double q2)
The non-pertubative ccbar contributions to the helicity amplitudes.
Definition: MVll.cpp:1602
MVll::getgtilde_1_re
double getgtilde_1_re(double q2)
The real part of .
Definition: MVll.h:578
MVll::getQCDf_2
gslpp::complex getQCDf_2(double q2)
Definition: MVll.h:523
MVll::getQCDfC9p_1
double getQCDfC9p_1(double cutoff)
Definition: MVll.h:555
MVll::getQCDfC9_1
double getQCDfC9_1(double q2, double cutoff)
Definition: MVll.h:537
MVll::ale
double ale
Definition: MVll.h:739
MVll::H_V_p
gslpp::complex H_V_p(double q2, bool bar)
The helicity amplitude .
Definition: MVll.cpp:1632
MVll::getQCDf_1
gslpp::complex getQCDf_1(double q2)
Definition: MVll.h:516
MVll::exp_Phase
gslpp::complex exp_Phase[3]
Definition: MVll.h:736
MVll::xs
double xs
Definition: MVll.h:753
MVll::getSigma
double getSigma(int i, double q_2)
The value of from to .
Definition: MVll.cpp:1956
gslpp::complex
A class for defining operations on and functions of complex numbers.
Definition: gslpp_complex.h:35
MVll::getgtilde_2_im
double getgtilde_2_im(double q2)
The immaginary part of .
Definition: MVll.h:611
h_0
A class for the correction in .
Definition: MVllObservables.h:1587
MVll::getQCDfC9p_2
double getQCDfC9p_2(double cutoff)
Definition: MVll.h:561
StandardModel
A model class for the Standard Model.
Definition: StandardModel.h:474
MVll::mc_pole
double mc_pole
Definition: MVll.h:748
MVll::mb_pole
double mb_pole
Definition: MVll.h:747
MVll::H_V_0
gslpp::complex H_V_0(double q2, bool bar)
The helicity amplitude .
Definition: MVll.cpp:1625
MVll::lep
QCD::lepton lep
Definition: MVll.h:727
MVll::getMlep
double getMlep()
The mass of the lepton l.
Definition: MVll.h:364
MVll::mu_h
double mu_h
Definition: MVll.h:745
MVll::mySM
const StandardModel & mySM
Definition: MVll.h:726
MVll::myF_1
std::unique_ptr< F_1 > myF_1
Definition: MVll.h:731
MVll::H_A_p
gslpp::complex H_A_p(double q2, bool bar)
The helicity amplitude .
Definition: MVll.cpp:1650
MVll::etaV
int etaV
Definition: MVll.h:755
MVll::H_V_m
gslpp::complex H_V_m(double q2, bool bar)
The helicity amplitude .
Definition: MVll.cpp:1638
MVll::getTp
double getTp(double q2)
The form factor .
Definition: MVll.h:424
MVll::ys
double ys
Definition: MVll.h:752
QCD::K_star
Definition: QCD.h:348
QCD::meson
meson
An enum type for mesons.
Definition: QCD.h:336
MVll::Mlep
double Mlep
Definition: MVll.h:740
F_2
Definition: F_2.h:15
MVll::mJ2
double mJ2
Definition: MVll.h:735
gslpp::sqrt
complex sqrt(const complex &z)
Definition: gslpp_complex.cpp:385
MVll::getwidth
double getwidth()
The width of the meson M.
Definition: MVll.h:355
MVll::meson
QCD::meson meson
Definition: MVll.h:728
MVll::beta
double beta(double q2)
The factor used in the angular coefficients .
Definition: MVll.cpp:1681
MVll::~MVll
virtual ~MVll()
Destructor.
Definition: MVll.cpp:145
MVll::getgtilde_1_im
double getgtilde_1_im(double q2)
The immaginary part of .
Definition: MVll.h:589
MVll::getgtilde_3_im
double getgtilde_3_im(double q2)
The imaginary part of .
Definition: MVll.h:634
MVll::getT0
double getT0(double q2)
The form factor .
Definition: MVll.h:413
MVll::H_S
gslpp::complex H_S(double q2, bool bar)
The helicity amplitude .
Definition: MVll.cpp:1662
MVll::geth_p_0
gslpp::complex geth_p_0()
.
Definition: MVll.h:665
gslpp::vector< double >
A class for constructing and defining operations on real vectors.
Definition: gslpp_vector_double.h:33
MVll::geth_p_im
double geth_p_im(double q2)
The imaginary part of .
Definition: MVll.h:685
MVll::getgtilde_3_re
double getgtilde_3_re(double q2)
The real part of .
Definition: MVll.h:622
MVll::integrateSigma
double integrateSigma(int i, double q_min, double q_max)
The integral of from to .
Definition: MVll.cpp:1838
MVll::getTm
double getTm(double q2)
The form factor .
Definition: MVll.h:435
MVll::geth_0_im
double geth_0_im(double q2)
The imaginary part of .
Definition: MVll.h:656
C_4
Definition: Doxygen/examples-src/myModel/src/myObservables.h:76
F_1
Definition: F_1.h:15
MVll::FixedWCbtos
bool FixedWCbtos
Definition: MVll.h:734
QCD::K_star_P
Definition: QCD.h:349
MVll::geth_m_0
gslpp::complex geth_m_0()
.
Definition: MVll.h:694
MVll::width
double width
Definition: MVll.h:751
QCD::PHI
Definition: QCD.h:347
MVll::Mb
double Mb
Definition: MVll.h:743
MVll::myF_2
std::unique_ptr< F_2 > myF_2
Definition: MVll.h:732
MVll::angmomV
double angmomV
Definition: MVll.h:754
MVll::getVp
double getVp(double q2)
The form factor .
Definition: MVll.h:392
MVll
A class for the decay.
Definition: MVll.h:308
MVll::geth_m_re
double geth_m_re(double q2)
The real part of .
Definition: MVll.h:704
MVll::getV0
double getV0(double q2)
The form factor .
Definition: MVll.h:381
MVll::getQCDfC9_3
double getQCDfC9_3(double q2, double cutoff)
Definition: MVll.h:549
MVll::getS
double getS(double q2)
The form factor .
Definition: MVll.h:446
MVll::mu_b
double mu_b
Definition: MVll.h:744
MVll::spectator_charge
double spectator_charge
Definition: MVll.h:750
MVll::integrateDelta
double integrateDelta(int i, double q_min, double q_max)
The integral of from to .
Definition: MVll.cpp:2001
C_3
Definition: Doxygen/examples-src/myModel/src/myObservables.h:59
MVll::geth_0_re
double geth_0_re(double q2)
The real part of .
Definition: MVll.h:646
MVll::getgtilde_2_re
double getgtilde_2_re(double q2)
The real part of .
Definition: MVll.h:600
MVll::MM
double MM
Definition: MVll.h:741
MVll::dispersion
bool dispersion
Definition: MVll.h:733
MVll::Ms
double Ms
Definition: MVll.h:749
gslpp::vector< gslpp::complex >
MVll::getQCDf_3
gslpp::complex getQCDf_3(double q2)
Definition: MVll.h:530
MVll::mvllParameters
std::vector< std::string > mvllParameters
Definition: MVll.h:730
MVll::Mc
double Mc
Definition: MVll.h:746
QCD::lepton
lepton
An enum type for leptons.
Definition: QCD.h:310
MVll::getVm
double getVm(double q2)
The form factor .
Definition: MVll.h:403