a Code for the Combination of Indirect and Direct Constraints on High Energy Physics Models Logo
MPll.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 MPLL_H
9 #define MPLL_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 MPllSWITCH 8.2
22 #define LATTICE true
23 #define GSL_INTERP_DIM_DC 10
24 #define SPLINE true
25 
26 #define NFPOLARBASIS_MPLL true
27 
172 class MPll{
173 public:
181  MPll(const StandardModel& SM_i, QCD::meson meson_i, QCD::meson pseudoscalar_i, QCD::lepton lep_i);
182 
186  virtual ~MPll();
187 
188  gslpp::complex funct_g(double q2);
189 
190  gslpp::complex DeltaC9_KD(double q2);
191 
198  gslpp::complex h_lambda(double q2);
199 
205  gslpp::complex H_V(double q2);
206 
207 
213  gslpp::complex H_A(double q2);
214 
215 
221  gslpp::complex H_S(double q2);
222 
223 
229  gslpp::complex H_P(double q2);
230 
231 
239  double integrateSigma(int i, double q_min, double q_max);
240 
247  double getSigma(int i, double q_2);
248 
256  double integrateDelta(int i, double q_min, double q_max);
257 
262  double getwidth(){
263  updateParameters();
264  return width;
265  }
266 
271  std::vector<std::string> initializeMPllParameters();
272 
273 
274 private:
279  std::vector<std::string> mpllParameters;
280  std::unique_ptr<F_1> myF_1;
281  std::unique_ptr<F_2> myF_2;
284  double mJ2;
285 
286  double GF;
287  double ale;
288  double Mlep;
289  double MM;
290  double MP;
291  double Mb;
292  double mu_b;
293  double mu_h;
294  double Mc;
295  double mb_pole;
296  double mc_pole;
297  double Ms;
299  double width;
300  double alpha_s_mub;
301  double angmomP;
302  int etaP;
304  double MW;
305  gslpp::complex lambda_t;
307  gslpp::complex h_1;
308  gslpp::complex h_2;
309  gslpp::complex r_1;
310  gslpp::complex r_2;
311  gslpp::complex Delta_C9;
312  gslpp::complex exp_Phase;
313 
314  /*LATTICE fit parameters*/
315 #if LATTICE
316  double b_0_fplus;
317  double b_1_fplus;
318  double b_2_fplus;
319  double m_fit2_fplus_lat;
320  double b_0_fT;
321  double b_1_fT;
322  double b_2_fT;
323  double m_fit2_fT_lat;
324  double b_0_f0;
325  double b_1_f0;
326  double b_2_f0;
327  double m_fit2_f0_lat;
328 #else
329  /*LCSR fit parameters*/
330  double r_1_fplus;
331  double r_2_fplus;
332  double m_fit2_fplus;
333  double r_1_fT;
334  double r_2_fT;
335  double m_fit2_fT;
336  double r_2_f0;
337  double m_fit2_f0;
338 #endif
339 
340 
341  gslpp::vector<gslpp::complex> ** allcoeff;
342  gslpp::vector<gslpp::complex> ** allcoeffh;
343  gslpp::vector<gslpp::complex> ** allcoeffprime;
345  gslpp::vector<gslpp::complex> ** allcoeff_noSM;
347  gslpp::complex C_1;
348  gslpp::complex C_1L_bar;
349  gslpp::complex C_1Lh_bar;
350  gslpp::complex C_2;
351  gslpp::complex C_2L_bar;
352  gslpp::complex C_2Lh_bar;
355  gslpp::complex C_5;
356  gslpp::complex C_6;
357  gslpp::complex C_7;
358  gslpp::complex C_8;
359  gslpp::complex C_8L;
360  gslpp::complex C_8Lh;
361  gslpp::complex C_9;
362  gslpp::complex C_10;
363  gslpp::complex C_S;
364  gslpp::complex C_P;
366  gslpp::complex C_7p;
367  gslpp::complex C_9p;
368  gslpp::complex C_10p;
369  gslpp::complex C_Sp;
370  gslpp::complex C_Pp;
372  gsl_interp_accel *acc_Re_deltaC7_QCDF;
373  gsl_interp_accel *acc_Im_deltaC7_QCDF;
374  gsl_interp_accel *acc_Re_deltaC9_QCDF;
375  gsl_interp_accel *acc_Im_deltaC9_QCDF;
376 
377  gsl_spline *spline_Re_deltaC7_QCDF;
378  gsl_spline *spline_Im_deltaC7_QCDF;
379  gsl_spline *spline_Re_deltaC9_QCDF;
380  gsl_spline *spline_Im_deltaC9_QCDF;
381 
382 
383  std::vector<double> ReDeltaC9;
384  std::vector<double> ImDeltaC9;
385  std::vector<double> myq2;
386  TFitResultPtr refres;
387  TFitResultPtr imfres;
389  TGraph gr1;
390  TGraph gr2;
392  TF1 reffit;
393  TF1 imffit;
395  double tmpq2;
397  gslpp::complex H_0_pre;
398  gslpp::complex H_0_WC;
399  gslpp::complex H_c_WC;
400  gslpp::complex H_b_WC;
402  gslpp::complex ihalfMPI;
403  double fournineth;
404  double half;
405  double twothird;
406  double Mc2;
407  double Mb2;
408  double logMc;
409  double logMb;
410  double mu_b2;
411  double fourMc2;
412  double fourMb2;
413  double Mlep2;
414  double NN;
415  double MM2;
416  double MM4;
417  double MP2;
418  double MP4;
419  double MM2mMP2;
420  double twoMP2;
421  double twoMM;
422  double twoMM2;
423  double twoMM2_MMpMP;
424  double twoMM_MbpMs;
425  double S_L_pre;
426  double fourMM2;
427  double twoMboMM;
428  double sixteenM_PI2;
429  double ninetysixM_PI3MM3;
430  double MboMW;
431  double MboMM;
432  double MsoMb;
433  double twoMlepMb;
434  double threeGegen0;
435  double threeGegen1otwo;
436  double M_PI2osix;
437  double twoMc2;
438  double sixMMoMb;
439  double CF;
440  double deltaT_0;
441  double deltaT_1par;
443  gslpp::complex ubar;
444  gslpp::complex arg1;
445  gslpp::complex B01;
446  gslpp::complex B00;
447  gslpp::complex xp;
448  gslpp::complex xm;
449  gslpp::complex yp;
450  gslpp::complex ym;
451  gslpp::complex L1xp;
452  gslpp::complex L1xm;
453  gslpp::complex L1yp;
454  gslpp::complex L1ym;
455  gslpp::complex F87_0;
456  gslpp::complex F87_1;
457  gslpp::complex F87_2;
458  gslpp::complex F87_3;
459  gslpp::complex F29_0;
460  gslpp::complex F29_L1;
461  gslpp::complex F29_1;
462  gslpp::complex F29_2;
463  gslpp::complex F29_3;
464  gslpp::complex F29_L1_1;
465  gslpp::complex F29_L1_2;
466  gslpp::complex F29_L1_3;
467  gslpp::complex F27_0;
468  gslpp::complex F27_1;
469  gslpp::complex F27_2;
470  gslpp::complex F27_3;
471  gslpp::complex F27_L1_1;
472  gslpp::complex F27_L1_2;
473  gslpp::complex F27_L1_3;
474  double F89_0;
475  double F89_1;
476  double F89_2;
477  double F89_3;
478  double Ee;
480 #if LATTICE
481  unsigned int fplus_lat_updated;
482  gslpp::vector<double> fplus_lat_cache;
484  unsigned int fT_lat_updated;
485  gslpp::vector<double> fT_lat_cache;
487  unsigned int f0_lat_updated;
488  gslpp::vector<double> f0_lat_cache;
489 #else
490  unsigned int fplus_updated;
491  gslpp::vector<double> fplus_cache;
493  unsigned int fT_updated;
494  gslpp::vector<double> fT_cache;
496  unsigned int f0_updated;
497  double f0_cache;
498 #endif
499 
500  unsigned int k2_updated;
501  gslpp::vector<double> k2_cache;
503  unsigned int beta_updated;
504  double beta_cache;
506  unsigned int lambda_updated;
508  unsigned int F_updated;
510  unsigned int VL_updated;
512  unsigned int TL_updated;
514  unsigned int SL_updated;
515  gslpp::vector<double> SL_cache;
517  unsigned int N_updated;
518  gslpp::vector<double> N_cache;
519  gslpp::complex Nc_cache;
521  unsigned int C_1_updated;
522  gslpp::complex C_1_cache;
524  unsigned int C_2_updated;
525  gslpp::complex C_2_cache;
527  unsigned int C_3_updated;
528  gslpp::complex C_3_cache;
530  unsigned int C_4_updated;
531  gslpp::complex C_4_cache;
533  unsigned int C_5_updated;
534  gslpp::complex C_5_cache;
536  unsigned int C_6_updated;
537  gslpp::complex C_6_cache;
539  unsigned int C_7_updated;
540  gslpp::complex C_7_cache;
542  unsigned int C_9_updated;
543  gslpp::complex C_9_cache;
545  unsigned int C_10_updated;
546  gslpp::complex C_10_cache;
548  unsigned int C_7p_updated;
549  gslpp::complex C_7p_cache;
551  unsigned int C_9p_updated;
552  gslpp::complex C_9p_cache;
554  unsigned int C_10p_updated;
555  gslpp::complex C_10p_cache;
557  unsigned int C_S_updated;
558  gslpp::complex C_S_cache;
560  unsigned int C_P_updated;
561  gslpp::complex C_P_cache;
563  unsigned int C_Sp_updated;
564  gslpp::complex C_Sp_cache;
566  unsigned int C_Pp_updated;
567  gslpp::complex C_Pp_cache;
569  unsigned int C_2Lh_updated;
570  gslpp::complex C_2Lh_cache;
572  unsigned int C_8Lh_updated;
573  gslpp::complex C_8Lh_cache;
575  unsigned int Yupdated;
576  gslpp::vector<double> Ycache;
578  unsigned int H_V0updated;
579  gslpp::vector<double> H_V0cache;
580  gslpp::complex H_V0Ccache[3];
581  gslpp::complex H_V0Ccache_dispersion[4];
583  unsigned int H_A0updated;
585  unsigned int H_Supdated;
586  gslpp::vector<double> H_Scache;
588  unsigned int H_P_updated;
589  gslpp::vector<double> H_P_cache;
591  unsigned int I0_updated;
592  unsigned int I2_updated;
593  unsigned int I8_updated;
595  std::map<std::pair<double, double>, unsigned int > I1Cached;
597  std::map<std::pair<double, double>, unsigned int > sigma0Cached;
598  std::map<std::pair<double, double>, unsigned int > sigma2Cached;
600  std::map<std::pair<double, double>, unsigned int > delta0Cached;
601  std::map<std::pair<double, double>, unsigned int > delta2Cached;
603  double avaSigma;
604  double avaDelta;
605  double avaDTPPR;
607  double errSigma;
608  double errDelta;
609  double errDTPPR;
611  gsl_function FS;
612  gsl_function FD;
613  gsl_function DTPPR;
615  gsl_integration_cquad_workspace * w_sigma;
616  gsl_integration_cquad_workspace * w_delta;
617  gsl_integration_cquad_workspace * w_DTPPR;
619  gsl_error_handler_t * old_handler;
621  std::map<std::pair<double, double>, gslpp::complex > cacheI1;
623  std::map<std::pair<double, double>, double > cacheSigma0;
624  std::map<std::pair<double, double>, double > cacheSigma2;
626  std::map<std::pair<double, double>, double > cacheDelta0;
627  std::map<std::pair<double, double>, double > cacheDelta2;
629  std::map<double, unsigned int> deltaTparpCached;
630  std::map<double, unsigned int> deltaTparmCached;
632  std::map<double, gslpp::complex> cacheDeltaTparp;
633  std::map<double, gslpp::complex> cacheDeltaTparm;
635  unsigned int deltaTparpupdated;
636  unsigned int deltaTparmupdated;
638  unsigned int T_updated;
639  gslpp::vector<double> T_cache;
646  void updateParameters();
647 
651  void checkCache();
652 
661  double LCSR_fit1(double q2, double r_1, double r_2, double m_fit2);
662 
663 
671  double LCSR_fit2(double q2, double r_2, double m_fit2);
672 
673 
679  double zeta(double q2);
680 
681 
691  double LATTICE_fit1(double q2, double b_0, double b_1, double b_2, double m_fit2);
692 
693 
703  double LATTICE_fit2(double q2, double b_0, double b_1, double b_2, double m_fit2);
704 
705 
711  double f_plus(double q2);
712 
713 
719  double f_T(double q2);
720 
721 
727  double f_0(double q2);
728 
734  gslpp::complex V_L(double q2);
735 
736 
742  gslpp::complex T_L(double q2);
743 
744 
750  double S_L(double q2);
751 
752 
758  gslpp::complex H_0(double q2);
759 
766  gslpp::complex H_c(double q2, double mu2);
767 
774  gslpp::complex H_b(double q2, double mu2);
775 
776 
782  gslpp::complex Y(double q2);
783 
784 
790  double k2 (double q2);
791 
792 
798  double beta (double q2);
799 
805  double beta2 (double q2);
806 
807 
813  double lambda(double q2);
814 
815 
821  double F(double q2);
822 
828  double I_1c(double q2);
829 
835  double I_2c(double q2);
836 
842  double I_6c(double q2);
843 
844 
851  double Delta(int i, double q2);
852 
858  double getSigma1c(double q2)
859  {
860  return I_1c(q2);
861  };
862 
868  double getSigma2c(double q2)
869  {
870  return I_2c(q2);
871  };
872 
878  double getSigma6c(double q2)
879  {
880  return I_6c(q2);
881  };
882 
888  double getDelta1c(double q2)
889  {
890  return Delta(0, q2);
891  };
892 
898  double getDelta2c(double q2)
899  {
900  return Delta(2, q2);
901  };
902 
909  gslpp::complex I1(double u, double q2);
910 
917  gslpp::complex Tparplus(double u, double q2);
918 
925  gslpp::complex Tparminus(double u, double q2);
926 
932  double Integrand_ReTparplus(double up);
933 
939  double Integrand_ImTparplus(double up);
940 
946  double Integrand_ReTparminus(double up);
947 
953  double Integrand_ImTparminus(double up);
954 
960  double Integrand_ReTpar_pm(double up);
961 
967  double Integrand_ImTpar_pm(double up);
968 
974  gslpp::complex F19(double q2);
975 
981  gslpp::complex F27(double q2);
982 
988  gslpp::complex F29(double q2);
989 
995  gslpp::complex F87(double q2);
996 
1002  double F89(double q2);
1003 
1009  gslpp::complex Cpar(double q2);
1010 
1016  gslpp::complex deltaTpar(double q2);
1017 
1024  double reDC9fit(double* x, double* p);
1025 
1032  double imDC9fit(double* x, double* p);
1033 
1037  void fit_DeltaC9_mumu();
1038 
1044  gslpp::complex fDeltaC9(double q2);
1045 
1051  gslpp::complex DeltaC9(double q2);
1052 
1058  gslpp::complex deltaC7_QCDF(double q2, bool spline = false);
1059 
1065  gslpp::complex deltaC9_QCDF(double q2, bool spline = false);
1066 
1067  void spline_QCDF_func();
1068 
1069 };
1070 
1071 #endif /* MPLL_H */
1072 
1073 
MPll::dispersion
bool dispersion
Definition: MPll.h:282
gslpp_special_functions::zeta
double zeta(int i)
Definition: gslpp_special_functions.cpp:20
MPll::H_P
gslpp::complex H_P(double q2)
The helicity amplitude .
Definition: MPll.cpp:1120
MPll::H_V
gslpp::complex H_V(double q2)
The helicity amplitude .
Definition: MPll.cpp:1105
MPll::mc_pole
double mc_pole
Definition: MPll.h:296
MPll::MP
double MP
Definition: MPll.h:290
MPll::pseudoscalar
QCD::meson pseudoscalar
Definition: MPll.h:278
MPll::meson
QCD::meson meson
Definition: MPll.h:277
MPll::~MPll
virtual ~MPll()
Destructor.
Definition: MPll.cpp:72
MPll::H_A
gslpp::complex H_A(double q2)
The helicity amplitude .
Definition: MPll.cpp:1110
MPll::Ms
double Ms
Definition: MPll.h:297
MPll::ale
double ale
Definition: MPll.h:287
gslpp::complex
A class for defining operations on and functions of complex numbers.
Definition: gslpp_complex.h:35
h_0
A class for the correction in .
Definition: MVllObservables.h:1587
MPll::myF_1
std::unique_ptr< F_1 > myF_1
Definition: MPll.h:280
MPll::spectator_charge
double spectator_charge
Definition: MPll.h:298
StandardModel
A model class for the Standard Model.
Definition: StandardModel.h:474
MPll::mpllParameters
std::vector< std::string > mpllParameters
Definition: MPll.h:279
MPll::integrateSigma
double integrateSigma(int i, double q_min, double q_max)
The integral of from to .
Definition: MPll.cpp:1174
QCD::meson
meson
An enum type for mesons.
Definition: QCD.h:336
MPll::integrateDelta
double integrateDelta(int i, double q_min, double q_max)
The integral of from to .
Definition: MPll.cpp:1232
F_2
Definition: F_2.h:15
MPll::lep
QCD::lepton lep
Definition: MPll.h:276
MPll::MPll
MPll(const StandardModel &SM_i, QCD::meson meson_i, QCD::meson pseudoscalar_i, QCD::lepton lep_i)
Constructor.
Definition: MPll.cpp:19
MPll::Mc
double Mc
Definition: MPll.h:294
MPll::H_S
gslpp::complex H_S(double q2)
The helicity amplitude .
Definition: MPll.cpp:1115
MPll::FixedWCbtos
bool FixedWCbtos
Definition: MPll.h:283
MPll::Mb
double Mb
Definition: MPll.h:291
MPll::MM
double MM
Definition: MPll.h:289
gslpp::vector< double >
A class for constructing and defining operations on real vectors.
Definition: gslpp_vector_double.h:33
MPll::GF
double GF
Definition: MPll.h:286
C_4
Definition: Doxygen/examples-src/myModel/src/myObservables.h:76
MPll::mu_b
double mu_b
Definition: MPll.h:292
MPll::Mlep
double Mlep
Definition: MPll.h:288
F_1
Definition: F_1.h:15
MPll::mJ2
double mJ2
Definition: MPll.h:284
MPll::mu_h
double mu_h
Definition: MPll.h:293
MPll::getSigma
double getSigma(int i, double q_2)
The value of from to .
Definition: MPll.cpp:1211
MPll::mySM
const StandardModel & mySM
Definition: MPll.h:275
MPll
A class for the decay.
Definition: MPll.h:172
MPll::h_lambda
gslpp::complex h_lambda(double q2)
The non-pertubative ccbar contributions to the helicity amplitudes.
Definition: MPll.cpp:1096
MPll::width
double width
Definition: MPll.h:299
MPll::funct_g
gslpp::complex funct_g(double q2)
Definition: MPll.cpp:1082
MPll::mb_pole
double mb_pole
Definition: MPll.h:295
MPll::initializeMPllParameters
std::vector< std::string > initializeMPllParameters()
A method for initializing the parameters necessary for MPll.
Definition: MPll.cpp:76
C_3
Definition: Doxygen/examples-src/myModel/src/myObservables.h:59
MPll::DeltaC9_KD
gslpp::complex DeltaC9_KD(double q2)
Definition: MPll.cpp:1090
MPll::myF_2
std::unique_ptr< F_2 > myF_2
Definition: MPll.h:281
MPll::getwidth
double getwidth()
The width of the meson M.
Definition: MPll.h:262
gslpp::vector< gslpp::complex >
QCD::lepton
lepton
An enum type for leptons.
Definition: QCD.h:310