a Code for the Combination of Indirect and Direct Constraints on High Energy Physics Models Logo
StandardModel.h
Go to the documentation of this file.
1 /*
2  * Copyright (C) 2012 HEPfit Collaboration
3  *
4  *
5  * For the licensing terms see doc/COPYING.
6  */
7 
8 #ifndef STANDARDMODEL_H
9 #define STANDARDMODEL_H
10 
11 #include "Flavour.h"
12 #include "QCD.h"
13 #include "CKM.h"
14 #include "PMNS.h"
15 #include "StandardModelMatching.h"
16 #include "Matching.h"
17 #include <gsl/gsl_integration.h>
18 
19 class EWSMcache;
20 class EWSMOneLoopEW;
21 class EWSMTwoLoopQCD;
22 class EWSMTwoLoopEW;
23 class EWSMThreeLoopQCD;
25 class EWSMThreeLoopEW;
27 class LeptonFlavour;
28 /* BEGIN: REMOVE FROM THE PACKAGE */
30 /* END: REMOVE FROM THE PACKAGE */
31 
32 
477 class StandardModel : public QCD {
478 public:
479 
480 
481  // Radiative Corrections for the LEP-II observables
482  enum LEP2RCs {
483  Weak = 0,
485  ISR,
489  };
490 
495  enum orders_EW {
496  EW1 = 0,
499  EW2,
501  EW3,
503  };
504 
505 
506 
507 // static const int NSMvars = 38; ///< The number of the model parameters in %StandardModel. !!! PMNS INCLUDED
508  static const int NSMvars = 26;
509 
512  static std::string SMvars[NSMvars];
513 
514  static const double GeVminus2_to_nb;
515 
520  static const double Mw_error;
521 
525  StandardModel();
526 
530  virtual ~StandardModel();
531 
532 
534  // Initialization
535 
542  virtual bool InitializeModel();
543 
544 
546  // Model parameters
547 
554  virtual bool Init(const std::map<std::string, double>& DPars);
555 
563  virtual bool PreUpdate();
564 
572  virtual bool Update(const std::map<std::string, double>& DPars);
573 
583  virtual bool PostUpdate();
584 
585  int getIterationNo() const
586  {
587  return iterationNo;
588  }
589 
597  virtual bool CheckParameters(const std::map<std::string, double>& DPars);
598 
599 
601  // Flags
602 
609  virtual bool setFlag(const std::string name, const bool value);
610 
617  virtual bool setFlagStr(const std::string name, const std::string value);
618 
623  virtual bool CheckFlags() const;
624 
635  {
637  }
638 
648  {
650  }
651 
653  {
654  this->FlagNoApproximateGammaZ = FlagNoApproximateGammaZ;
655  }
656 
662  std::string getFlagMw() const
663  {
664  return FlagMw;
665  }
666 
672  std::string getFlagRhoZ() const
673  {
674  return FlagRhoZ;
675  }
676 
682  std::string getFlagKappaZ() const
683  {
684  return FlagKappaZ;
685  }
686 
699  {
700  this->FlagCacheInStandardModel = FlagCacheInStandardModel;
701  }
702 
703 
705  // get and set methods for class members
706 
713  {
714  return leptons[p];
715  }
716 
721  double getMz() const
722  {
723  return Mz;
724  }
725 
730  double getAlsMz() const
731  {
732  return AlsMz;
733  }
734 
739  double getGF() const
740  {
741  return GF;
742  }
743 
748  double getAle() const
749  {
750  return ale;
751  }
752 
759  double getDAle5Mz() const
760  {
761  return dAle5Mz;
762  }
763 
768  virtual double getMHl() const
769  {
770  return mHl;
771  }
772 
778  double getDelMw() const
779  {
780  return delMw;
781  }
782 
789  double getDelSin2th_l() const
790  {
791  return delSin2th_l;
792  }
793 
800  double getDelSin2th_q() const
801  {
802  return delSin2th_q;
803  }
804 
811  double getDelSin2th_b() const
812  {
813  return delSin2th_b;
814  }
815 
821  double getDelGammaZ() const
822  {
823  return delGammaZ;
824  }
825 
831  double getDelSigma0H() const
832  {
833  return delsigma0H;
834  }
835 
841  double getDelR0l() const
842  {
843  return delR0l;
844  }
845 
851  double getDelR0c() const
852  {
853  return delR0c;
854  }
855 
861  double getDelR0b() const
862  {
863  return delR0b;
864  }
865 
870  gslpp::matrix<gslpp::complex> getVCKM() const // why don't we return a const reference?
871  {
872  return myCKM.getCKM();
873  }
874 
879  CKM getCKM() const
880  {
881  return myCKM;
882  }
883 
889  {
890  return myPMNS.getPMNS();
891  }
892 
899  {
900  return Yu;
901  }
902 
909  {
910  return Yd;
911  }
912 
919  {
920  return Yn;
921  }
922 
929  {
930  return Ye;
931  }
932 
938  double getMuw() const
939  {
940  return muw;
941  }
942 
943  virtual StandardModel getTrueSM() const
944  {
945  throw std::runtime_error("StandardModel::getTrueSM() must be overridden by the NP extension.");
946  }
947 
953  {
954  return SMM.getObj();
955  }
956 
962  {
963  return myEWSMcache;
964  }
965 
971  {
972  return myOneLoopEW;
973  }
974 
980  {
981  return myApproximateFormulae;
982  }
983 
984  /* BEGIN: REMOVE FROM THE PACKAGE */
990  {
991  return myTwoFermionsLEP2;
992  }
993  /* END: REMOVE FROM THE PACKAGE */
994 
996  {
997  return myThreeLoopEW;
998  }
999 
1001  {
1002  return myThreeLoopEW2QCD;
1003  }
1004 
1006  {
1007  return myThreeLoopQCD;
1008  }
1009 
1011  {
1012  return myTwoLoopEW;
1013  }
1014 
1016  {
1017  return myTwoLoopQCD;
1018  }
1019 
1020  const Flavour& getFlavour() const
1021  {
1022  return SMFlavour;
1023  }
1024 
1026  {
1027  return myLeptonFlavour;
1028  }
1030  // QED coupling
1031 
1043  double ale_OS(const double mu, orders order = FULLNLO) const;
1044 
1051  double Beta_s(int nm, unsigned int nf) const;
1052 
1059  double Beta_e(int nm, unsigned int nf) const;
1060 
1071 #ifdef __clang__
1072 #pragma clang diagnostic push
1073 #pragma clang diagnostic ignored "-Woverloaded-virtual"
1074 #endif
1075  double Als(double mu, orders order = FULLNLO, bool qed_flag = false, bool Nf_thr = true) const;
1076  double AlsByOrder(double mu, orders order = FULLNLO, bool qed_flag = false, bool Nf_thr = true) const;
1077 #ifdef __clang__
1078 #pragma clang diagnostic pop
1079 #endif
1080 
1090  double Ale(double mu, orders order, bool Nf_thr = true) const;
1091 
1098  double DeltaAlphaLepton(const double s) const;
1099 
1112  double DeltaAlphaL5q() const;
1113 
1120  double DeltaAlphaTop(const double s) const;
1121 
1134  double DeltaAlpha() const;
1135 
1145  double alphaMz() const;
1146 
1153  double Alstilde5(const double mu) const;
1154 
1158  virtual double getCCC1() const { return 0.; };
1159 
1163  virtual double getCCC2() const { return 0.; };
1164 
1168  virtual double getCCC3() const { return 0.; };
1169 
1173  virtual double getCCC4() const { return 0.; };
1174 
1178  virtual double getCCC5() const { return 0.; };
1179 
1180 
1181 
1183  // Higgs VEV
1184 
1193  virtual double v() const;
1194 
1195 
1197  // The W-boson mass
1198 
1203  virtual double Mw_tree() const;
1204 
1220  double s02() const;
1221 
1234  double c02() const;
1235 
1262  virtual double Mw() const;
1263 
1273  virtual double cW2(const double Mw_i) const;
1274  virtual double cW2() const;
1275 
1285  virtual double sW2(const double Mw_i) const;
1286  double sW2() const;
1287 
1309  virtual double DeltaR() const;
1310 
1319  void ComputeDeltaRho(const double Mw_i, double DeltaRho[orders_EW_size]) const;
1320 
1329  void ComputeDeltaR_rem(const double Mw_i, double DeltaR_rem[orders_EW_size]) const;
1330 
1331 
1333  // The W and Z masses in the complex-pole/fixed-width scheme
1334 
1365  double Mzbar() const;
1366 
1387  double MwbarFromMw(const double Mw) const;
1388 
1412  double MwFromMwbar(const double Mwbar) const;
1413 
1431  virtual double DeltaRbar() const;
1432 
1433 
1435  // The W-boson decay width
1436 
1445  virtual double rho_GammaW(const Particle fi, const Particle fj) const;
1446 
1473  virtual double GammaW(const Particle fi, const Particle fj) const;
1474 
1479  virtual double GammaW() const;
1480 
1481 
1483  // EWPO at Z-pole
1484 
1503  virtual double A_f(const Particle f) const;
1504 
1510  virtual double AFB(const Particle f) const;
1511 
1531  virtual double sin2thetaEff(const Particle f) const;
1532 
1558  virtual double GammaZ(const Particle f) const;
1559 
1570  virtual double Gamma_inv() const;
1571 
1586  virtual double Gamma_had() const;
1587 
1601  virtual double Gamma_Z() const;
1602 
1616  virtual double sigma0_had() const;
1617 
1632  virtual double R0_f(const Particle f) const;
1633 
1642  virtual double R_inv() const;
1643 
1653  virtual double N_nu() const;
1654 
1655 
1657  // Zff effective couplings
1658 
1668  virtual gslpp::complex gV_f(const Particle f) const;
1669 
1679  virtual gslpp::complex gA_f(const Particle f) const;
1680 
1696  virtual gslpp::complex rhoZ_f(const Particle f) const;
1697 
1725  virtual gslpp::complex kappaZ_f(const Particle f) const;
1726 
1747  virtual gslpp::complex deltaRhoZ_f(const Particle f) const;
1748 
1773  virtual gslpp::complex deltaKappaZ_f(const Particle f) const;
1774 
1775 
1777  // Epsilon parameters for EWPO
1778 
1790  virtual double epsilon1() const;
1791 
1812  virtual double epsilon2() const;
1813 
1831  virtual double epsilon3() const;
1832 
1850  virtual double epsilonb() const;
1851 
1852 
1854  // For EWPO caches
1855 
1863  static const int NumSMParamsForEWPO = 33;
1864 
1879  bool checkSMparamsForEWPO();
1880 
1882  // Several Higgs-related quantities used in Higgs coupling analysis
1883 
1897  double computeSigmaggH(const double sqrt_s) const
1898  {
1899  if (sqrt_s == 7.0) {
1900  return 16.83; // in pb for Mh=125.1 GeV
1901  } else if (sqrt_s == 8.0) {
1902  return 21.40; // in pb for Mh=125.1 GeV
1903  } else if (sqrt_s == 13.0) {
1904  return 48.61; // in pb for Mh=125.09 GeV
1905  } else if (sqrt_s == 14.0) {
1906  return 54.72; // in pb for Mh=125.09 GeV
1907  } else if (sqrt_s == 27.0) {
1908  return 146.65; // in pb for Mh=125.09 GeV
1909  } else if (sqrt_s == 100.0) {
1910  return 740.3; // in pb for Mh=125. GeV
1911  } else if (sqrt_s == 1.96) {
1912  return 0.9493; // in pb for Mh=125 GeV
1913  } else
1914  throw std::runtime_error("Bad argument in StandardModel::computeSigmaggH()");
1915  }
1916 
1924  double computeSigmaggH_tt(const double sqrt_s) const
1925  {
1926  if (sqrt_s == 7.0) {
1927  return 16.69; // in pb for Mh=125.09 GeV
1928  } else if (sqrt_s == 8.0) {
1929  return 21.20; // in pb for Mh=125.09 GeV
1930  } else if (sqrt_s == 13.0) {
1931  return 47.94; // in pb for Mh=125.09 GeV
1932  } else if (sqrt_s == 14.0) {
1933  return 53.93; // in pb for Mh=125.09 GeV
1934  } else if (sqrt_s == 27.0) {
1935  return computeSigmaggH(sqrt_s) / computeSigmaggH(14.) * computeSigmaggH_tt(14.); // in the absence of this value we rescale the LHC result at 14 TeV
1936  } else if (sqrt_s == 100.0) {
1937  return computeSigmaggH(sqrt_s) / computeSigmaggH(14.) * computeSigmaggH_tt(14.); // in the absence of this value we rescale the LHC result at 14 TeV
1938  } else
1939  throw std::runtime_error("Bad argument in StandardModel::computeSigmaggH_tt()");
1940  }
1941 
1949  double computeSigmaggH_bb(const double sqrt_s) const
1950  {
1951  if (sqrt_s == 7.0) {
1952  return 0.04; // in pb for Mh=125.09 GeV
1953  } else if (sqrt_s == 8.0) {
1954  return 0.05; // in pb for Mh=125.09 GeV
1955  } else if (sqrt_s == 13.0) {
1956  return 0.10; // in pb for Mh=125.09 GeV
1957  } else if (sqrt_s == 14.0) {
1958  return 0.11; // in pb for Mh=125.09 GeV
1959  } else if (sqrt_s == 27.0) {
1960  return computeSigmaggH(sqrt_s) / computeSigmaggH(14.) * computeSigmaggH_bb(14.); // in the absence of this value we rescale the LHC result at 14 TeV
1961  } else if (sqrt_s == 100.0) {
1962  return computeSigmaggH(sqrt_s) / computeSigmaggH(14.) * computeSigmaggH_bb(14.); // in the absence of this value we rescale the LHC result at 14 TeV
1963  } else
1964  throw std::runtime_error("Bad argument in StandardModel::computeSigmaggH_bb()");
1965  }
1966 
1974  double computeSigmaggH_tb(const double sqrt_s) const
1975  {
1976  if (sqrt_s == 7.0) {
1977  return -0.66; // in pb for Mh=125.09 GeV
1978  } else if (sqrt_s == 8.0) {
1979  return -0.82; // in pb for Mh=125.09 GeV
1980  } else if (sqrt_s == 13.0) {
1981  return -1.73; // in pb for Mh=125.09 GeV
1982  } else if (sqrt_s == 14.0) {
1983  return -1.92; // in pb for Mh=125.09 GeV
1984  } else if (sqrt_s == 27.0) {
1985  return computeSigmaggH(sqrt_s) / computeSigmaggH(14.) * computeSigmaggH_tb(14.); // in the absence of this value we rescale the LHC result at 14 TeV
1986  } else if (sqrt_s == 100.0) {
1987  return computeSigmaggH(sqrt_s) / computeSigmaggH(14.) * computeSigmaggH_tb(14.); // in the absence of this value we rescale the LHC result at 14 TeV
1988  } else
1989  throw std::runtime_error("Bad argument in StandardModel::computeSigmaggH_tb()");
1990  }
1991 
2003  double computeSigmaVBF(const double sqrt_s) const
2004  {
2005  if (sqrt_s == 7.0) {
2006  return 1.241; // in pb for Mh=125.09 GeV
2007  } else if (sqrt_s == 8.0) {
2008  return 1.601; // in pb for Mh=125.09 GeV
2009  } else if (sqrt_s == 13.0) {
2010  return 3.766; // in pb for Mh=125.09 GeV
2011  } else if (sqrt_s == 14.0) {
2012  return 4.260; // in pb for Mh=125.09 GeV
2013  } else if (sqrt_s == 27.0) {
2014  return 11.838; // in pb for Mh=125.09 GeV
2015  } else if (sqrt_s == 100.0) {
2016  return 82.0; // in pb for Mh=125. GeV
2017  } else if (sqrt_s == 1.96) {
2018  return 0.0653; // in pb for Mh=125 GeV
2019  } else
2020  throw std::runtime_error("Bad argument in StandardModel::computeSigmaVBF()");
2021  }
2022 
2031  double computeSigmaWF(const double sqrt_s) const
2032  {
2033  if (sqrt_s == 7.0) {
2034  return 0.946; // in pb for Mh=125 GeV
2035  } else if (sqrt_s == 8.0) {
2036  return 1.220; // in pb for Mh=125 GeV
2037  } else if (sqrt_s == 13.0) {
2038  return 2.882; // in pb for Mh=125 GeV
2039  } else if (sqrt_s == 14.0) {
2040  return 3.260; // in pb for Mh=125 GeV
2041  } else if (sqrt_s == 27.0) {
2042  return computeSigmaVBF(sqrt_s) / computeSigmaVBF(14.) * computeSigmaWF(14.); // in the absence of this value we rescale the LHC result at 14 TeV
2043  } else if (sqrt_s == 100.0) {
2044  return computeSigmaVBF(sqrt_s) / computeSigmaVBF(14.) * computeSigmaWF(14.); // in the absence of this value we rescale the LHC result at 14 TeV
2045  } else if (sqrt_s == 1.96) {
2046  return computeSigmaVBF(sqrt_s) / computeSigmaVBF(7.) * computeSigmaWF(7.); // in the absence of individual cross sections for TeVatron we rescale the LHC ones
2047  } else
2048  throw std::runtime_error("Bad argument in StandardModel::computeSigmaWF()");
2049  }
2050 
2059  double computeSigmaZF(const double sqrt_s) const
2060  {
2061  if (sqrt_s == 7.0) {
2062  return 0.333; // in pb for Mh=125 GeV
2063  } else if (sqrt_s == 8.0) {
2064  return 0.432; // in pb for Mh=125 GeV
2065  } else if (sqrt_s == 13.0) {
2066  return 1.049; // in pb for Mh=125 GeV
2067  } else if (sqrt_s == 14.0) {
2068  return 1.191; // in pb for Mh=125 GeV
2069  } else if (sqrt_s == 27.0) {
2070  return computeSigmaVBF(sqrt_s) / computeSigmaVBF(14.) * computeSigmaZF(14.); // in the absence of this value we rescale the LHC result at 14 TeV
2071  } else if (sqrt_s == 100.0) {
2072  return computeSigmaVBF(sqrt_s) / computeSigmaVBF(14.) * computeSigmaZF(14.); // in the absence of this value we rescale the LHC result at 14 TeV
2073  } else if (sqrt_s == 1.96) {
2074  return computeSigmaVBF(sqrt_s) / computeSigmaVBF(7.) * computeSigmaZF(7.); // in the absence of individual cross sections for TeVatron we rescale the LHC ones
2075  } else
2076  throw std::runtime_error("Bad argument in StandardModel::computeSigmaZF()");
2077  }
2078 
2086  double computeSigmaZWF(const double sqrt_s) const
2087  {
2088  return 0.;
2089  }
2090 
2102  double computeSigmaWH(const double sqrt_s) const
2103  {
2104  if (sqrt_s == 7.0) {
2105  return 0.577; // in pb for Mh=125.1 GeV
2106  //return 0.5688; // in pb for Mh=125.6 GeV
2107  } else if (sqrt_s == 8.0) {
2108  return 0.7027; // in pb for Mh=125.1 GeV
2109  //return 0.6931; // in pb for Mh=125.6 GeV
2110  } else if (sqrt_s == 13.0) {
2111  return 1.358; // in pb for Mh=125.09 GeV
2112  } else if (sqrt_s == 14.0) {
2113  return 1.498; // in pb for Mh=125.09 GeV
2114  } else if (sqrt_s == 27.0) {
2115  return 3.397; // in pb for Mh=125.09 GeV
2116  } else if (sqrt_s == 100.0) {
2117  return 15.9; // in pb for Mh=125. GeV
2118  } else if (sqrt_s == 1.96) {
2119  return 0.1295; // in pb for Mh=125 GeV
2120  } else
2121  throw std::runtime_error("Bad argument in StandardModel::computeSigmaWH()");
2122  }
2123 
2135  double computeSigmaZH(const double sqrt_s) const
2136  {
2137  if (sqrt_s == 7.0) {
2138  return 0.3341; // in pb for Mh=125.1 GeV
2139  //return 0.3299; // in pb for Mh=125.6 GeV
2140  } else if (sqrt_s == 8.0) {
2141  return 0.4142; // in pb for Mh=125.1 GeV
2142  //return 0.4091; // in pb for Mh=125.6 GeV
2143  } else if (sqrt_s == 13.0) {
2144  return 0.880; // in pb for Mh=125.09 GeV
2145  } else if (sqrt_s == 14.0) {
2146  return 0.981; // in pb for Mh=125.09 GeV
2147  } else if (sqrt_s == 27.0) {
2148  return 2.463; // in pb for Mh=125.09 GeV
2149  } else if (sqrt_s == 100.0) {
2150  return 11.26; // in pb for Mh=125. GeV
2151  } else if (sqrt_s == 1.96) {
2152  return 0.0785; // in pb for Mh=125 GeV
2153  } else
2154  throw std::runtime_error("Bad argument in StandardModel::computeSigmaZH()");
2155  }
2156 
2171  double computeSigmattH(const double sqrt_s) const
2172  {
2173  if (sqrt_s == 7.0) {
2174  return 0.0861; // in pb for Mh=125.1 GeV
2175  //return 0.0851; // in pb for Mh=125.6 GeV
2176  } else if (sqrt_s == 8.0) {
2177  return 0.129; // in pb for Mh=125.1 GeV
2178  //return 0.1274; // in pb for Mh=125.6 GeV
2179  } else if (sqrt_s == 13.0) {
2180  return 0.5060; // in pb for Mh=125.1 GeV
2181  } else if (sqrt_s == 14.0) {
2182  return 0.6128; // in pb for Mh=125.09 GeV
2183  } else if (sqrt_s == 27.0) {
2184  return 2.86; // in pb for Mh=125.09 GeV
2185  } else if (sqrt_s == 100.0) {
2186  return 37.9; // in pb for Mh=125. GeV
2187  } else if (sqrt_s == 1.96) {
2188  return 0.0043; // in pb for Mh=125 GeV
2189  } else
2190  throw std::runtime_error("Bad argument in StandardModel::computeSigmattH()");
2191  }
2192 
2199  double computeBrHtogg() const
2200  {
2201  return 8.179e-2; // Mh=125.1 GeV
2202  }
2203 
2210  double computeBrHtoWW() const
2211  {
2212  //return 2.23e-1; // Mh=125.5 GeV
2213  return 2.154e-1; // Mh=125.1 GeV
2214  }
2215 
2222  double computeBrHtoZZ() const
2223  {
2224  return 2.643e-2; // Mh=125.1 GeV
2225  //return 2.79e-2; // Mh=125.6 GeV
2226  }
2227 
2233  double computeBrHtoZZinv() const
2234  {
2235  return 1.06e-3;
2236  }
2237 
2244  double computeBrHtoZga() const
2245  {
2246  return 1.541e-3; // Mh=125.1 GeV
2247  //return 1.59e-3; // Mh=125.6 GeV
2248  }
2249 
2256  double computeBrHtogaga() const
2257  {
2258  return 2.27e-3; // Mh=125.1 GeV
2259  }
2260 
2267  double computeBrHtomumu() const
2268  {
2269  return 2.17e-4; // Mh=125.1 GeV
2270  }
2271 
2278  double computeBrHtotautau() const
2279  {
2280  return 6.256e-2; // Mh=125.1 GeV
2281  //return 6.22e-2; // Mh=125.6 GeV
2282  }
2283 
2290  double computeBrHtocc() const
2291  {
2292  return 2.883e-2; // Mh=125.1 GeV
2293  //return 2.86e-2; // Mh=125.6 GeV
2294  }
2295 
2302  double computeBrHtoss() const
2303  {
2304  return 4.0e-4;
2305  }
2306 
2313  double computeBrHtobb() const
2314  {
2315  return 5.807e-1; // Mh=125.1 GeV
2316  //return 5.67e-1; // Mh=125.6 GeV
2317  }
2318 
2325  double computeGammaHTotal() const
2326  {
2327  return 4.101e-3; // Mh=125.1 GeV
2328  //return 4.15e-3; // Mh=125.6 GeV
2329  }
2330 
2336  double computeGammaHgg_tt() const
2337  {
2338  return 380.8; // in keV for Mh=125 GeV
2339  //return 389.6; // in keV for Mh=126 GeV
2340  }
2341 
2347  double computeGammaHgg_bb() const
2348  {
2349  return 3.96; // in keV for Mh=125 GeV
2350  //return 3.95; // in keV for Mh=126 GeV
2351  }
2352 
2358  double computeGammaHgg_tb() const
2359  {
2360  return -42.1; // in keV for Mh=125 GeV
2361  //return -42.7; // in keV for Mh=126 GeV
2362  }
2363 
2369  double computeGammaHZga_tt() const
2370  {
2371  return 21.74; // in eV for Mh=125 GeV
2372  //return 23.51; // in eV for Mh=126 GeV
2373  }
2374 
2380  double computeGammaHZga_WW() const
2381  {
2382  return 7005.6; // in eV for Mh=125 GeV
2383  //return 7648.4; // in eV for Mh=126 GeV
2384  }
2385 
2391  double computeGammaHZga_tW() const
2392  {
2393  return -780.4; // in eV for Mh=125 GeV
2394  //return -848.1; // in eV for Mh=126 GeV
2395  }
2396 
2402  double computeGammaHgaga_tt() const
2403  {
2404  return 662.84; // in eV for Mh=125 GeV
2405  //return 680.39; // in eV for Mh=126 GeV
2406  }
2407 
2413  double computeGammaHgaga_WW() const
2414  {
2415  return 14731.86; // in eV for Mh=125 GeV
2416  //return 15221.98; // in eV for Mh=126 GeV
2417  }
2418 
2424  double computeGammaHgaga_tW() const
2425  {
2426  return -6249.93; // in eV for Mh=125 GeV
2427  //return -6436.35; // in eV for Mh=126 GeV
2428  }
2429 
2434  virtual double getCBd() const
2435  {
2436  return 1.;
2437  }
2438 
2443  virtual double getCBs() const
2444  {
2445  return 1.;
2446  }
2447 
2452  virtual double getCDMK() const
2453  {
2454  return 1.;
2455  }
2456 
2461  virtual double getCepsK() const
2462  {
2463  return 1.;
2464  }
2465 
2470  virtual double getPhiBs() const
2471  {
2472  return 0.;
2473  }
2474 
2479  virtual double getPhiBd() const
2480  {
2481  return 0.;
2482  }
2483 
2484 /* BEGIN: REMOVE FROM THE PACKAGE */
2486  //LEP2 Observables
2487 
2488  virtual double LEP2sigmaMu(const double s) const;
2489  virtual double LEP2sigmaTau(const double s) const;
2490  virtual double LEP2sigmaHadron(const double s) const;
2491  virtual double LEP2sigmaCharm(const double s) const;
2492  virtual double LEP2sigmaBottom(const double s) const;
2493 
2494  virtual double LEP2AFBmu(const double s) const;
2495  virtual double LEP2AFBtau(const double s) const;
2496  virtual double LEP2AFBcharm(const double s) const;
2497  virtual double LEP2AFBbottom(const double s) const;
2498  virtual double LEP2Rcharm(const double s) const;
2499  virtual double LEP2Rbottom(const double s) const;
2500 /* END: REMOVE FROM THE PACKAGE */
2501 
2502 bool setFlagSigmaForAFB(const bool flagSigmaForAFB_i)
2503 {
2504  bSigmaForAFB = flagSigmaForAFB_i;
2505  return true;
2506 }
2507 
2508 bool setFlagSigmaForR(const bool flagSigmaForR_i)
2509 {
2510  bSigmaForR = flagSigmaForR_i;
2511  return true;
2512 }
2513 
2514 virtual double getmq(const QCD::quark q, const double mu) const
2515 {
2516  return m_q(q, mu, FULLNLO);
2517 }
2519 protected:
2520 
2526  virtual void setParameter(const std::string name, const double& value);
2527 
2531  virtual void computeCKM();
2532 
2538  virtual void computeYukawas();
2539 
2543 // gslpp::matrix<gslpp::complex> VCKM; ///< The %CKM matrix.
2544 // gslpp::matrix<gslpp::complex> UPMNS; ///< The %PMNS matrix.
2549 
2551 
2552  // model parameters
2553  double AlsMz;
2554  double Mz;
2555  double GF;
2556  double ale;
2557  double dAle5Mz;
2558  double mHl;
2559  double delMw;
2560  double delSin2th_l;
2561  double delSin2th_q;
2562  double delSin2th_b;
2563  double delGammaZ;
2564  double delsigma0H;
2565  double delR0l;
2566  double delR0c;
2567  double delR0b;
2568  double lambda;
2569  double A;
2570  double rhob;
2571  double etab;
2572  double Vus;
2573  double Vcb;
2574  double Vub;
2575  double gamma;
2576  double muw;
2578 
2579 
2581  // For EWPO
2582 
2590 
2601  double SchemeToDouble(const std::string scheme) const
2602  {
2603  if (scheme.compare("NORESUM") == 0)
2604  return 0.0;
2605  else if (scheme.compare("OMSI") == 0)
2606  return 1.0;
2607  else if (scheme.compare("INTERMEDIATE") == 0)
2608  return 2.0;
2609  else if (scheme.compare("OMSII") == 0)
2610  return 3.0;
2611  else if (scheme.compare("APPROXIMATEFORMULA") == 0)
2612  return 4.0;
2613  else
2614  throw std::runtime_error("EWSM::SchemeToDouble: bad scheme");
2615  }
2616 
2622  bool checkEWPOscheme(const std::string scheme) const
2623  {
2624  if (scheme.compare("NORESUM") == 0
2625  || scheme.compare("OMSI") == 0
2626  || scheme.compare("INTERMEDIATE") == 0
2627  || scheme.compare("OMSII") == 0
2628  || scheme.compare("APPROXIMATEFORMULA") == 0)
2629  return true;
2630  else
2631  return false;
2632  }
2633 
2634 
2635  double m_q(const QCD::quark q, const double mu, const orders order=FULLNLO) const
2636  {
2637  switch(q) {
2638  case QCD::UP:
2639  case QCD::DOWN:
2640  case QCD::STRANGE:
2641  return Mrun(mu, getQuarks(q).getMass_scale(),
2642  getQuarks(q).getMass(), order);
2643  case QCD::CHARM:
2644  case QCD::BOTTOM:
2645  return Mrun(mu, getQuarks(q).getMass(), order);
2646  case QCD::TOP:
2647  return getMtpole(); // the pole mass
2648  default:
2649  throw std::runtime_error("Error in StandardModel::m_q()");
2650  }
2651  }
2652 
2685  double resumMw(const double Mw_i, const double DeltaRho[orders_EW_size],
2686  const double DeltaR_rem[orders_EW_size]) const;
2687 
2714  double resumRhoZ(const double DeltaRho[orders_EW_size],
2715  const double deltaRho_rem[orders_EW_size],
2716  const double DeltaRbar_rem, const bool bool_Zbb) const;
2717 
2744  double resumKappaZ(const double DeltaRho[orders_EW_size],
2745  const double deltaKappa_rem[orders_EW_size],
2746  const double DeltaRbar_rem, const bool bool_Zbb) const;
2747 
2770  double taub() const;
2771 
2780  double Delta_EWQCD(const QCD::quark q) const;
2781 
2791  double RVq(const QCD::quark q) const;
2792 
2802  double RAq(const QCD::quark q) const;
2803 
2817  double RVh() const;
2818 
2819  bool requireCKM;
2820  bool requireYe;
2821  bool requireYn;
2822 
2824 
2825  /* BEGIN: REMOVE FROM THE PACKAGE */
2827  //Migrated from LEP2ThObservables.h
2828 
2829 
2830 
2832  mutable bool bSigmaForAFB;
2833  mutable bool bSigmaForR;
2834 
2835 
2836 
2837  double sigma_NoISR_l(const QCD::lepton l_flavor, const double s) const;
2838  double sigma_NoISR_q(const QCD::quark q_flavor, const double s) const;
2839  double AFB_NoISR_l(const QCD::lepton l_flavor, const double s) const;
2840  double AFB_NoISR_q(const QCD::quark q_flavor, const double s) const;
2841 
2842  double Integrand_sigmaWithISR_l(double x, const QCD::lepton l_flavor, const double s) const;
2843 
2844  double getIntegrand_sigmaWithISR_mu130(double x) const;
2845  double getIntegrand_sigmaWithISR_mu136(double x) const;
2846  double getIntegrand_sigmaWithISR_mu161(double x) const;
2847  double getIntegrand_sigmaWithISR_mu172(double x) const;
2848  double getIntegrand_sigmaWithISR_mu183(double x) const;
2849  double getIntegrand_sigmaWithISR_mu189(double x) const;
2850  double getIntegrand_sigmaWithISR_mu192(double x) const;
2851  double getIntegrand_sigmaWithISR_mu196(double x) const;
2852  double getIntegrand_sigmaWithISR_mu200(double x) const;
2853  double getIntegrand_sigmaWithISR_mu202(double x) const;
2854  double getIntegrand_sigmaWithISR_mu205(double x) const;
2855  double getIntegrand_sigmaWithISR_mu207(double x) const;
2856 
2857 
2858  double getIntegrand_sigmaWithISR_tau130(double x) const;
2859  double getIntegrand_sigmaWithISR_tau136(double x) const;
2860  double getIntegrand_sigmaWithISR_tau161(double x) const;
2861  double getIntegrand_sigmaWithISR_tau172(double x) const;
2862  double getIntegrand_sigmaWithISR_tau183(double x) const;
2863  double getIntegrand_sigmaWithISR_tau189(double x) const;
2864  double getIntegrand_sigmaWithISR_tau192(double x) const;
2865  double getIntegrand_sigmaWithISR_tau196(double x) const;
2866  double getIntegrand_sigmaWithISR_tau200(double x) const;
2867  double getIntegrand_sigmaWithISR_tau202(double x) const;
2868  double getIntegrand_sigmaWithISR_tau205(double x) const;
2869  double getIntegrand_sigmaWithISR_tau207(double x) const;
2870 
2871 
2872  double Integrand_sigmaWithISR_q(double x, const QCD::quark q_flavor, const double s) const;
2873 
2874  double getIntegrand_sigmaWithISR_up130(double x) const;
2875  double getIntegrand_sigmaWithISR_up133(double x) const;
2876  double getIntegrand_sigmaWithISR_up136(double x) const;
2877  double getIntegrand_sigmaWithISR_up161(double x) const;
2878  double getIntegrand_sigmaWithISR_up167(double x) const;
2879  double getIntegrand_sigmaWithISR_up172(double x) const;
2880  double getIntegrand_sigmaWithISR_up183(double x) const;
2881  double getIntegrand_sigmaWithISR_up189(double x) const;
2882  double getIntegrand_sigmaWithISR_up192(double x) const;
2883  double getIntegrand_sigmaWithISR_up196(double x) const;
2884  double getIntegrand_sigmaWithISR_up200(double x) const;
2885  double getIntegrand_sigmaWithISR_up202(double x) const;
2886  double getIntegrand_sigmaWithISR_up205(double x) const;
2887  double getIntegrand_sigmaWithISR_up207(double x) const;
2888 
2889  double getIntegrand_sigmaWithISR_down130(double x) const;
2890  double getIntegrand_sigmaWithISR_down133(double x) const;
2891  double getIntegrand_sigmaWithISR_down136(double x) const;
2892  double getIntegrand_sigmaWithISR_down161(double x) const;
2893  double getIntegrand_sigmaWithISR_down167(double x) const;
2894  double getIntegrand_sigmaWithISR_down172(double x) const;
2895  double getIntegrand_sigmaWithISR_down183(double x) const;
2896  double getIntegrand_sigmaWithISR_down189(double x) const;
2897  double getIntegrand_sigmaWithISR_down192(double x) const;
2898  double getIntegrand_sigmaWithISR_down196(double x) const;
2899  double getIntegrand_sigmaWithISR_down200(double x) const;
2900  double getIntegrand_sigmaWithISR_down202(double x) const;
2901  double getIntegrand_sigmaWithISR_down205(double x) const;
2902  double getIntegrand_sigmaWithISR_down207(double x) const;
2903 
2904  double getIntegrand_sigmaWithISR_charm130(double x) const;
2905  double getIntegrand_sigmaWithISR_charm133(double x) const;
2906  double getIntegrand_sigmaWithISR_charm136(double x) const;
2907  double getIntegrand_sigmaWithISR_charm161(double x) const;
2908  double getIntegrand_sigmaWithISR_charm167(double x) const;
2909  double getIntegrand_sigmaWithISR_charm172(double x) const;
2910  double getIntegrand_sigmaWithISR_charm183(double x) const;
2911  double getIntegrand_sigmaWithISR_charm189(double x) const;
2912  double getIntegrand_sigmaWithISR_charm192(double x) const;
2913  double getIntegrand_sigmaWithISR_charm196(double x) const;
2914  double getIntegrand_sigmaWithISR_charm200(double x) const;
2915  double getIntegrand_sigmaWithISR_charm202(double x) const;
2916  double getIntegrand_sigmaWithISR_charm205(double x) const;
2917  double getIntegrand_sigmaWithISR_charm207(double x) const;
2918 
2919  double getIntegrand_sigmaWithISR_strange130(double x) const;
2920  double getIntegrand_sigmaWithISR_strange133(double x) const;
2921  double getIntegrand_sigmaWithISR_strange136(double x) const;
2922  double getIntegrand_sigmaWithISR_strange161(double x) const;
2923  double getIntegrand_sigmaWithISR_strange167(double x) const;
2924  double getIntegrand_sigmaWithISR_strange172(double x) const;
2925  double getIntegrand_sigmaWithISR_strange183(double x) const;
2926  double getIntegrand_sigmaWithISR_strange189(double x) const;
2927  double getIntegrand_sigmaWithISR_strange192(double x) const;
2928  double getIntegrand_sigmaWithISR_strange196(double x) const;
2929  double getIntegrand_sigmaWithISR_strange200(double x) const;
2930  double getIntegrand_sigmaWithISR_strange202(double x) const;
2931  double getIntegrand_sigmaWithISR_strange205(double x) const;
2932  double getIntegrand_sigmaWithISR_strange207(double x) const;
2933 
2934  double getIntegrand_sigmaWithISR_bottom130(double x) const;
2935  double getIntegrand_sigmaWithISR_bottom133(double x) const;
2936  double getIntegrand_sigmaWithISR_bottom136(double x) const;
2937  double getIntegrand_sigmaWithISR_bottom161(double x) const;
2938  double getIntegrand_sigmaWithISR_bottom167(double x) const;
2939  double getIntegrand_sigmaWithISR_bottom172(double x) const;
2940  double getIntegrand_sigmaWithISR_bottom183(double x) const;
2941  double getIntegrand_sigmaWithISR_bottom189(double x) const;
2942  double getIntegrand_sigmaWithISR_bottom192(double x) const;
2943  double getIntegrand_sigmaWithISR_bottom196(double x) const;
2944  double getIntegrand_sigmaWithISR_bottom200(double x) const;
2945  double getIntegrand_sigmaWithISR_bottom202(double x) const;
2946  double getIntegrand_sigmaWithISR_bottom205(double x) const;
2947  double getIntegrand_sigmaWithISR_bottom207(double x) const;
2948 
2949  double Integrand_dsigmaBox_l(double cosTheta, const QCD::lepton l_flavor, const double s) const;
2950 
2951  double getIntegrand_dsigmaBox_mu130(double x) const;
2952  double getIntegrand_dsigmaBox_mu133(double x) const;
2953  double getIntegrand_dsigmaBox_mu136(double x) const;
2954  double getIntegrand_dsigmaBox_mu161(double x) const;
2955  double getIntegrand_dsigmaBox_mu167(double x) const;
2956  double getIntegrand_dsigmaBox_mu172(double x) const;
2957  double getIntegrand_dsigmaBox_mu183(double x) const;
2958  double getIntegrand_dsigmaBox_mu189(double x) const;
2959  double getIntegrand_dsigmaBox_mu192(double x) const;
2960  double getIntegrand_dsigmaBox_mu196(double x) const;
2961  double getIntegrand_dsigmaBox_mu200(double x) const;
2962  double getIntegrand_dsigmaBox_mu202(double x) const;
2963  double getIntegrand_dsigmaBox_mu205(double x) const;
2964  double getIntegrand_dsigmaBox_mu207(double x) const;
2965 
2966  double getIntegrand_dsigmaBox_tau130(double x) const;
2967  double getIntegrand_dsigmaBox_tau133(double x) const;
2968  double getIntegrand_dsigmaBox_tau136(double x) const;
2969  double getIntegrand_dsigmaBox_tau161(double x) const;
2970  double getIntegrand_dsigmaBox_tau167(double x) const;
2971  double getIntegrand_dsigmaBox_tau172(double x) const;
2972  double getIntegrand_dsigmaBox_tau183(double x) const;
2973  double getIntegrand_dsigmaBox_tau189(double x) const;
2974  double getIntegrand_dsigmaBox_tau192(double x) const;
2975  double getIntegrand_dsigmaBox_tau196(double x) const;
2976  double getIntegrand_dsigmaBox_tau200(double x) const;
2977  double getIntegrand_dsigmaBox_tau202(double x) const;
2978  double getIntegrand_dsigmaBox_tau205(double x) const;
2979  double getIntegrand_dsigmaBox_tau207(double x) const;
2980 
2981 
2982  double Integrand_dsigmaBox_q(double cosTheta, const QCD::quark q_flavor, const double s) const;
2983 
2984  double getIntegrand_dsigmaBox_up130(double x) const;
2985  double getIntegrand_dsigmaBox_up133(double x) const;
2986  double getIntegrand_dsigmaBox_up136(double x) const;
2987  double getIntegrand_dsigmaBox_up161(double x) const;
2988  double getIntegrand_dsigmaBox_up167(double x) const;
2989  double getIntegrand_dsigmaBox_up172(double x) const;
2990  double getIntegrand_dsigmaBox_up183(double x) const;
2991  double getIntegrand_dsigmaBox_up189(double x) const;
2992  double getIntegrand_dsigmaBox_up192(double x) const;
2993  double getIntegrand_dsigmaBox_up196(double x) const;
2994  double getIntegrand_dsigmaBox_up200(double x) const;
2995  double getIntegrand_dsigmaBox_up202(double x) const;
2996  double getIntegrand_dsigmaBox_up205(double x) const;
2997  double getIntegrand_dsigmaBox_up207(double x) const;
2998 
2999  double getIntegrand_dsigmaBox_down130(double x) const;
3000  double getIntegrand_dsigmaBox_down133(double x) const;
3001  double getIntegrand_dsigmaBox_down136(double x) const;
3002  double getIntegrand_dsigmaBox_down161(double x) const;
3003  double getIntegrand_dsigmaBox_down167(double x) const;
3004  double getIntegrand_dsigmaBox_down172(double x) const;
3005  double getIntegrand_dsigmaBox_down183(double x) const;
3006  double getIntegrand_dsigmaBox_down189(double x) const;
3007  double getIntegrand_dsigmaBox_down192(double x) const;
3008  double getIntegrand_dsigmaBox_down196(double x) const;
3009  double getIntegrand_dsigmaBox_down200(double x) const;
3010  double getIntegrand_dsigmaBox_down202(double x) const;
3011  double getIntegrand_dsigmaBox_down205(double x) const;
3012  double getIntegrand_dsigmaBox_down207(double x) const;
3013 
3014  double getIntegrand_dsigmaBox_charm130(double x) const;
3015  double getIntegrand_dsigmaBox_charm133(double x) const;
3016  double getIntegrand_dsigmaBox_charm136(double x) const;
3017  double getIntegrand_dsigmaBox_charm161(double x) const;
3018  double getIntegrand_dsigmaBox_charm167(double x) const;
3019  double getIntegrand_dsigmaBox_charm172(double x) const;
3020  double getIntegrand_dsigmaBox_charm183(double x) const;
3021  double getIntegrand_dsigmaBox_charm189(double x) const;
3022  double getIntegrand_dsigmaBox_charm192(double x) const;
3023  double getIntegrand_dsigmaBox_charm196(double x) const;
3024  double getIntegrand_dsigmaBox_charm200(double x) const;
3025  double getIntegrand_dsigmaBox_charm202(double x) const;
3026  double getIntegrand_dsigmaBox_charm205(double x) const;
3027  double getIntegrand_dsigmaBox_charm207(double x) const;
3028 
3029  double getIntegrand_dsigmaBox_strange130(double x) const;
3030  double getIntegrand_dsigmaBox_strange133(double x) const;
3031  double getIntegrand_dsigmaBox_strange136(double x) const;
3032  double getIntegrand_dsigmaBox_strange161(double x) const;
3033  double getIntegrand_dsigmaBox_strange167(double x) const;
3034  double getIntegrand_dsigmaBox_strange172(double x) const;
3035  double getIntegrand_dsigmaBox_strange183(double x) const;
3036  double getIntegrand_dsigmaBox_strange189(double x) const;
3037  double getIntegrand_dsigmaBox_strange192(double x) const;
3038  double getIntegrand_dsigmaBox_strange196(double x) const;
3039  double getIntegrand_dsigmaBox_strange200(double x) const;
3040  double getIntegrand_dsigmaBox_strange202(double x) const;
3041  double getIntegrand_dsigmaBox_strange205(double x) const;
3042  double getIntegrand_dsigmaBox_strange207(double x) const;
3043 
3044  double getIntegrand_dsigmaBox_bottom130(double x) const;
3045  double getIntegrand_dsigmaBox_bottom133(double x) const;
3046  double getIntegrand_dsigmaBox_bottom136(double x) const;
3047  double getIntegrand_dsigmaBox_bottom161(double x) const;
3048  double getIntegrand_dsigmaBox_bottom167(double x) const;
3049  double getIntegrand_dsigmaBox_bottom172(double x) const;
3050  double getIntegrand_dsigmaBox_bottom183(double x) const;
3051  double getIntegrand_dsigmaBox_bottom189(double x) const;
3052  double getIntegrand_dsigmaBox_bottom192(double x) const;
3053  double getIntegrand_dsigmaBox_bottom196(double x) const;
3054  double getIntegrand_dsigmaBox_bottom200(double x) const;
3055  double getIntegrand_dsigmaBox_bottom202(double x) const;
3056  double getIntegrand_dsigmaBox_bottom205(double x) const;
3057  double getIntegrand_dsigmaBox_bottom207(double x) const;
3058 
3059 
3060  double Integrand_AFBnumeratorWithISR_l(double x, const QCD::lepton l_flavor, const double s) const;
3061 
3062  double getIntegrand_AFBnumeratorWithISR_mu130(double x) const;
3063  double getIntegrand_AFBnumeratorWithISR_mu136(double x) const;
3064  double getIntegrand_AFBnumeratorWithISR_mu161(double x) const;
3065  double getIntegrand_AFBnumeratorWithISR_mu172(double x) const;
3066  double getIntegrand_AFBnumeratorWithISR_mu183(double x) const;
3067  double getIntegrand_AFBnumeratorWithISR_mu189(double x) const;
3068  double getIntegrand_AFBnumeratorWithISR_mu192(double x) const;
3069  double getIntegrand_AFBnumeratorWithISR_mu196(double x) const;
3070  double getIntegrand_AFBnumeratorWithISR_mu200(double x) const;
3071  double getIntegrand_AFBnumeratorWithISR_mu202(double x) const;
3072  double getIntegrand_AFBnumeratorWithISR_mu205(double x) const;
3073  double getIntegrand_AFBnumeratorWithISR_mu207(double x) const;
3074 
3075  double getIntegrand_AFBnumeratorWithISR_tau130(double x) const;
3076  double getIntegrand_AFBnumeratorWithISR_tau136(double x) const;
3077  double getIntegrand_AFBnumeratorWithISR_tau161(double x) const;
3078  double getIntegrand_AFBnumeratorWithISR_tau172(double x) const;
3079  double getIntegrand_AFBnumeratorWithISR_tau183(double x) const;
3080  double getIntegrand_AFBnumeratorWithISR_tau189(double x) const;
3081  double getIntegrand_AFBnumeratorWithISR_tau192(double x) const;
3082  double getIntegrand_AFBnumeratorWithISR_tau196(double x) const;
3083  double getIntegrand_AFBnumeratorWithISR_tau200(double x) const;
3084  double getIntegrand_AFBnumeratorWithISR_tau202(double x) const;
3085  double getIntegrand_AFBnumeratorWithISR_tau205(double x) const;
3086  double getIntegrand_AFBnumeratorWithISR_tau207(double x) const;
3087 
3088  double Integrand_AFBnumeratorWithISR_q(double x, const QCD::quark q_flavor, const double s) const;
3089 
3090  double getIntegrand_AFBnumeratorWithISR_charm133(double x) const;
3091  double getIntegrand_AFBnumeratorWithISR_charm167(double x) const;
3092  double getIntegrand_AFBnumeratorWithISR_charm172(double x) const;
3093  double getIntegrand_AFBnumeratorWithISR_charm183(double x) const;
3094  double getIntegrand_AFBnumeratorWithISR_charm189(double x) const;
3095  double getIntegrand_AFBnumeratorWithISR_charm192(double x) const;
3096  double getIntegrand_AFBnumeratorWithISR_charm196(double x) const;
3097  double getIntegrand_AFBnumeratorWithISR_charm200(double x) const;
3098  double getIntegrand_AFBnumeratorWithISR_charm202(double x) const;
3099  double getIntegrand_AFBnumeratorWithISR_charm205(double x) const;
3100  double getIntegrand_AFBnumeratorWithISR_charm207(double x) const;
3101 
3102  double getIntegrand_AFBnumeratorWithISR_bottom133(double x) const;
3103  double getIntegrand_AFBnumeratorWithISR_bottom167(double x) const;
3104  double getIntegrand_AFBnumeratorWithISR_bottom172(double x) const;
3105  double getIntegrand_AFBnumeratorWithISR_bottom183(double x) const;
3106  double getIntegrand_AFBnumeratorWithISR_bottom189(double x) const;
3107  double getIntegrand_AFBnumeratorWithISR_bottom192(double x) const;
3108  double getIntegrand_AFBnumeratorWithISR_bottom196(double x) const;
3109  double getIntegrand_AFBnumeratorWithISR_bottom200(double x) const;
3110  double getIntegrand_AFBnumeratorWithISR_bottom202(double x) const;
3111  double getIntegrand_AFBnumeratorWithISR_bottom205(double x) const;
3112  double getIntegrand_AFBnumeratorWithISR_bottom207(double x) const;
3113 
3114  /* END: REMOVE FROM THE PACKAGE */
3116 private:
3126  /* BEGIN: REMOVE FROM THE PACKAGE */
3128  /* END: REMOVE FROM THE PACKAGE */
3129 
3132  std::string FlagMw;
3133  std::string FlagRhoZ;
3134  std::string FlagKappaZ;
3136 
3137  bool FlagSMAux;
3138 
3140  // Caches for EWPO
3141 
3144  mutable double DeltaAlphaLepton_cache;
3145  mutable double DeltaAlpha_cache;
3146  mutable double Mw_cache;
3147  mutable double GammaW_cache;
3151  mutable bool useDeltaAlpha_cache;
3152  mutable bool useMw_cache;
3153  mutable bool useGammaW_cache;
3154  mutable bool useRhoZ_f_cache[12];
3155  mutable bool useKappaZ_f_cache[12];
3156 
3157 /* BEGIN: REMOVE FROM THE PACKAGE */
3158  // caches for the SM prediction of LEP2 Obs
3159 // mutable double SMparams_cache[NumSMParamsForEWPO+3];
3160  mutable double SMresult_cache;
3161 // mutable bool flag_cache[NUMofLEP2RCs];
3162 // mutable double ml_cache, mq_cache, mqForHad_cache[6];
3163 
3164  mutable double average;
3165  mutable double error;
3166  mutable gsl_function f_GSL;
3167  gsl_integration_workspace * w_GSL1;
3168 /* END: REMOVE FROM THE PACKAGE */
3169 
3171 
3172  double AlsWithInit(double mu, double alsi, double mu_i, orders order, bool qed_flag) const;
3173  double AleWithInit(double mu, double alsi, double mu_i, orders order) const;
3174  static const int CacheSize = 5;
3175  mutable double als_cache[11][CacheSize];
3176  mutable double ale_cache[10][CacheSize];
3178 
3179 
3180 };
3181 
3182 #endif /* STANDARDMODEL_H */
StandardModel::getIntegrand_sigmaWithISR_down136
double getIntegrand_sigmaWithISR_down136(double x) const
Definition: StandardModel.cpp:6628
StandardModel::delR0c
double delR0c
The theoretical uncertainty in , denoted as .
Definition: StandardModel.h:2566
StandardModel::getIntegrand_dsigmaBox_charm192
double getIntegrand_dsigmaBox_charm192(double x) const
Definition: StandardModel.cpp:7362
StandardModel::getIntegrand_dsigmaBox_bottom189
double getIntegrand_dsigmaBox_bottom189(double x) const
Definition: StandardModel.cpp:7534
StandardModel::m_q
double m_q(const QCD::quark q, const double mu, const orders order=FULLNLO) const
Definition: StandardModel.h:2635
StandardModel::Integrand_sigmaWithISR_l
double Integrand_sigmaWithISR_l(double x, const QCD::lepton l_flavor, const double s) const
Definition: StandardModel.cpp:6344
StandardModel::getIntegrand_sigmaWithISR_strange200
double getIntegrand_sigmaWithISR_strange200(double x) const
Definition: StandardModel.cpp:6852
StandardModel::getMyOneLoopEW
EWSMOneLoopEW * getMyOneLoopEW() const
A get method to retrieve the member pointer of type EWSMOneLoopEW,.
Definition: StandardModel.h:970
StandardModel::getIntegrand_dsigmaBox_mu172
double getIntegrand_dsigmaBox_mu172(double x) const
Definition: StandardModel.cpp:6992
StandardModel::getIntegrand_dsigmaBox_down207
double getIntegrand_dsigmaBox_down207(double x) const
Definition: StandardModel.cpp:7304
StandardModel::getIntegrand_sigmaWithISR_bottom202
double getIntegrand_sigmaWithISR_bottom202(double x) const
Definition: StandardModel.cpp:6946
StandardModel::getIntegrand_sigmaWithISR_strange189
double getIntegrand_sigmaWithISR_strange189(double x) const
Definition: StandardModel.cpp:6834
StandardModel::Mw_cache
double Mw_cache
A cache of the value of .
Definition: StandardModel.h:3146
StandardModel::getIntegrand_dsigmaBox_mu192
double getIntegrand_dsigmaBox_mu192(double x) const
Definition: StandardModel.cpp:7010
StandardModel::FlagWithoutNonUniversalVC
bool FlagWithoutNonUniversalVC
A boolean for the model flag WithoutNonUniversalVC.
Definition: StandardModel.h:3130
StandardModel::~StandardModel
virtual ~StandardModel()
The default destructor.
Definition: StandardModel.cpp:141
StandardModel::getIntegrand_sigmaWithISR_charm207
double getIntegrand_sigmaWithISR_charm207(double x) const
Definition: StandardModel.cpp:6782
StandardModel::Gamma_inv
virtual double Gamma_inv() const
The invisible partial decay width of the boson, .
Definition: StandardModel.cpp:1303
StandardModel::EW2
Two-loop of .
Definition: StandardModel.h:499
StandardModel::rhoZ_f
virtual gslpp::complex rhoZ_f(const Particle f) const
The effective leptonic neutral-current coupling in the SM.
Definition: StandardModel.cpp:1579
StandardModel::taub
double taub() const
Top-mass corrections to the vertex, denoted by .
Definition: StandardModel.cpp:2104
StandardModel::setParameter
virtual void setParameter(const std::string name, const double &value)
A method to set the value of a parameter of StandardModel.
Definition: StandardModel.cpp:257
StandardModel::getIntegrand_AFBnumeratorWithISR_tau189
double getIntegrand_AFBnumeratorWithISR_tau189(double x) const
Definition: StandardModel.cpp:7703
StandardModel::getIntegrand_dsigmaBox_up133
double getIntegrand_dsigmaBox_up133(double x) const
Definition: StandardModel.cpp:7145
StandardModel::Integrand_sigmaWithISR_q
double Integrand_sigmaWithISR_q(double x, const QCD::quark q_flavor, const double s) const
Definition: StandardModel.cpp:6504
StandardModel::v
virtual double v() const
The Higgs vacuum expectation value.
Definition: StandardModel.cpp:943
StandardModel::getIterationNo
int getIterationNo() const
Definition: StandardModel.h:585
StandardModel::NSMvars
static const int NSMvars
The number of the model parameters in StandardModel.
Definition: StandardModel.h:508
StandardModel::getIntegrand_sigmaWithISR_tau207
double getIntegrand_sigmaWithISR_tau207(double x) const
Definition: StandardModel.cpp:6498
StandardModel::delSin2th_q
double delSin2th_q
The theoretical uncertainty in , denoted as .
Definition: StandardModel.h:2561
StandardModel::getIntegrand_sigmaWithISR_charm130
double getIntegrand_sigmaWithISR_charm130(double x) const
Definition: StandardModel.cpp:6704
StandardModel::getIntegrand_dsigmaBox_mu183
double getIntegrand_dsigmaBox_mu183(double x) const
Definition: StandardModel.cpp:6998
StandardModel::getIntegrand_AFBnumeratorWithISR_bottom189
double getIntegrand_AFBnumeratorWithISR_bottom189(double x) const
Definition: StandardModel.cpp:7855
StandardModel::getIntegrand_dsigmaBox_mu205
double getIntegrand_dsigmaBox_mu205(double x) const
Definition: StandardModel.cpp:7034
StandardModel::getIntegrand_sigmaWithISR_bottom200
double getIntegrand_sigmaWithISR_bottom200(double x) const
Definition: StandardModel.cpp:6940
StandardModel::useRhoZ_f_cache
bool useRhoZ_f_cache[12]
Definition: StandardModel.h:3154
QCD::BOTTOM
Definition: QCD.h:329
EWSMTwoLoopQCD
A class for two-loop corrections to the EW precision observables.
Definition: EWSMTwoLoopQCD.h:55
StandardModel::average
double average
Definition: StandardModel.h:3164
StandardModel::A
double A
The CKM parameter in the Wolfenstein parameterization.
Definition: StandardModel.h:2569
StandardModel::getIntegrand_sigmaWithISR_mu205
double getIntegrand_sigmaWithISR_mu205(double x) const
Definition: StandardModel.cpp:6419
StandardModel::DeltaRbar
virtual double DeltaRbar() const
The SM prediction for derived from that for the -boson mass.
Definition: StandardModel.cpp:1146
StandardModel::getIntegrand_dsigmaBox_tau183
double getIntegrand_dsigmaBox_tau183(double x) const
Definition: StandardModel.cpp:7074
QCD
A class for parameters related to QCD, hadrons and quarks.
Definition: QCD.h:304
Particle
A class for particles.
Definition: Particle.h:26
StandardModel::computeSigmaWH
double computeSigmaWH(const double sqrt_s) const
The WH production cross section in the Standard Model.
Definition: StandardModel.h:2102
StandardModel::ISR
Definition: StandardModel.h:485
StandardModel::rhob
double rhob
The CKM parameter in the Wolfenstein parameterization.
Definition: StandardModel.h:2570
StandardModel::setFlagNoApproximateGammaZ
void setFlagNoApproximateGammaZ(bool FlagNoApproximateGammaZ)
Definition: StandardModel.h:652
StandardModel::computeGammaHZga_tt
double computeGammaHZga_tt() const
The top loop contribution to in the Standard Model.
Definition: StandardModel.h:2369
StandardModel::gamma
double gamma
used as an input for FlagWolfenstein = FALSE
Definition: StandardModel.h:2575
StandardModel::getIntegrand_dsigmaBox_strange136
double getIntegrand_dsigmaBox_strange136(double x) const
Definition: StandardModel.cpp:7414
StandardModel::LEP2sigmaTau
virtual double LEP2sigmaTau(const double s) const
Definition: StandardModel.cpp:2674
StandardModel::getMyLeptonFlavour
LeptonFlavour * getMyLeptonFlavour() const
Definition: StandardModel.h:1025
StandardModel::getIntegrand_sigmaWithISR_charm183
double getIntegrand_sigmaWithISR_charm183(double x) const
Definition: StandardModel.cpp:6740
StandardModel::getCCC5
virtual double getCCC5() const
A virtual implementation for the RealWeakEFTCC class.
Definition: StandardModel.h:1178
StandardModel::realorder
orders realorder
Definition: StandardModel.h:3177
StandardModel::getIntegrand_sigmaWithISR_mu196
double getIntegrand_sigmaWithISR_mu196(double x) const
Definition: StandardModel.cpp:6401
StandardModel::getIntegrand_dsigmaBox_mu167
double getIntegrand_dsigmaBox_mu167(double x) const
StandardModel::getIntegrand_AFBnumeratorWithISR_mu202
double getIntegrand_AFBnumeratorWithISR_mu202(double x) const
Definition: StandardModel.cpp:7654
StandardModel::getIntegrand_dsigmaBox_down172
double getIntegrand_dsigmaBox_down172(double x) const
Definition: StandardModel.cpp:7256
StandardModel::getIntegrand_dsigmaBox_strange196
double getIntegrand_dsigmaBox_strange196(double x) const
Definition: StandardModel.cpp:7458
StandardModel::Ye
gslpp::matrix< gslpp::complex > Ye
The Yukawa matrix of the charged leptons.
Definition: StandardModel.h:2548
StandardModel::getIntegrand_sigmaWithISR_strange202
double getIntegrand_sigmaWithISR_strange202(double x) const
Definition: StandardModel.cpp:6858
StandardModel::computeBrHtotautau
double computeBrHtotautau() const
The Br in the Standard Model.
Definition: StandardModel.h:2278
StandardModel::rho_GammaW
virtual double rho_GammaW(const Particle fi, const Particle fj) const
EW radiative corrections to the width of , denoted as .
Definition: StandardModel.cpp:1158
StandardModel::getIntegrand_AFBnumeratorWithISR_mu136
double getIntegrand_AFBnumeratorWithISR_mu136(double x) const
Definition: StandardModel.cpp:7606
StandardModel::getIntegrand_sigmaWithISR_up172
double getIntegrand_sigmaWithISR_up172(double x) const
Definition: StandardModel.cpp:6558
StandardModel::myTwoFermionsLEP2
EWSMTwoFermionsLEP2 * myTwoFermionsLEP2
A pointer to an object of type EWSMTwoFermionsLEP2.
Definition: StandardModel.h:3127
StandardModel::getIntegrand_sigmaWithISR_tau183
double getIntegrand_sigmaWithISR_tau183(double x) const
Definition: StandardModel.cpp:6456
StandardModel::getDelSigma0H
double getDelSigma0H() const
A get method to retrieve the theoretical uncertainty in , denoted as .
Definition: StandardModel.h:831
StandardModel::getIntegrand_sigmaWithISR_charm172
double getIntegrand_sigmaWithISR_charm172(double x) const
Definition: StandardModel.cpp:6734
StandardModel::Yu
gslpp::matrix< gslpp::complex > Yu
The Yukawa matrix of the up-type quarks.
Definition: StandardModel.h:2545
EWSMThreeLoopEW2QCD
A class for three-loop corrections to the EW precision observables.
Definition: EWSMThreeLoopEW2QCD.h:34
StandardModel::getIntegrand_sigmaWithISR_mu200
double getIntegrand_sigmaWithISR_mu200(double x) const
Definition: StandardModel.cpp:6407
StandardModel::NumSMParamsForEWPO
static const int NumSMParamsForEWPO
The number of the SM parameters that are relevant to the EW precision observables.
Definition: StandardModel.h:1863
StandardModel::getIntegrand_sigmaWithISR_tau200
double getIntegrand_sigmaWithISR_tau200(double x) const
Definition: StandardModel.cpp:6480
StandardModel::ale_cache
double ale_cache[10][CacheSize]
Cache for .
Definition: StandardModel.h:3176
StandardModel::QCDFSR
Definition: StandardModel.h:487
StandardModel::getFlagRhoZ
std::string getFlagRhoZ() const
A method to retrieve the model flag RhoZ.
Definition: StandardModel.h:672
StandardModel::EW1
One-loop of .
Definition: StandardModel.h:496
StandardModel::getCCC1
virtual double getCCC1() const
A virtual implementation for the RealWeakEFTCC class.
Definition: StandardModel.h:1158
StandardModel::getIntegrand_dsigmaBox_charm207
double getIntegrand_dsigmaBox_charm207(double x) const
Definition: StandardModel.cpp:7392
StandardModel::getMyEWSMcache
EWSMcache * getMyEWSMcache() const
A get method to retrieve the member pointer of type EWSMcache.
Definition: StandardModel.h:961
StandardModel::MwbarFromMw
double MwbarFromMw(const double Mw) const
A method to convert the -boson mass in the experimental/running-width scheme to that in the complex-p...
Definition: StandardModel.cpp:1128
StandardModel::getIntegrand_sigmaWithISR_down167
double getIntegrand_sigmaWithISR_down167(double x) const
Definition: StandardModel.cpp:6640
StandardModel::getIntegrand_AFBnumeratorWithISR_bottom196
double getIntegrand_AFBnumeratorWithISR_bottom196(double x) const
Definition: StandardModel.cpp:7867
StandardModel::getIntegrand_dsigmaBox_tau172
double getIntegrand_dsigmaBox_tau172(double x) const
Definition: StandardModel.cpp:7068
StandardModel::getIntegrand_sigmaWithISR_mu192
double getIntegrand_sigmaWithISR_mu192(double x) const
Definition: StandardModel.cpp:6395
StandardModel::getIntegrand_dsigmaBox_bottom205
double getIntegrand_dsigmaBox_bottom205(double x) const
Definition: StandardModel.cpp:7564
StandardModel::getIntegrand_sigmaWithISR_strange183
double getIntegrand_sigmaWithISR_strange183(double x) const
Definition: StandardModel.cpp:6828
StandardModel::delR0l
double delR0l
The theoretical uncertainty in , denoted as .
Definition: StandardModel.h:2565
StandardModel::getMyTwoLoopQCD
EWSMTwoLoopQCD * getMyTwoLoopQCD() const
Definition: StandardModel.h:1015
StandardModel::Integrand_AFBnumeratorWithISR_l
double Integrand_AFBnumeratorWithISR_l(double x, const QCD::lepton l_flavor, const double s) const
Definition: StandardModel.cpp:7588
EWSMThreeLoopQCD
A class for three-loop corrections to the EW precision observables.
Definition: EWSMThreeLoopQCD.h:33
StandardModel::Delta_EWQCD
double Delta_EWQCD(const QCD::quark q) const
The non-factorizable EW-QCD corrections to the partial widths for , denoted as .
Definition: StandardModel.cpp:2124
StandardModel::getAlsMz
double getAlsMz() const
A get method to access the value of .
Definition: StandardModel.h:730
StandardModel::getIntegrand_sigmaWithISR_charm192
double getIntegrand_sigmaWithISR_charm192(double x) const
Definition: StandardModel.cpp:6752
Matching::getObj
T & getObj()
Definition: Matching.h:14
StandardModel::getDelR0l
double getDelR0l() const
A get method to retrieve the theoretical uncertainty in , denoted as .
Definition: StandardModel.h:841
StandardModel::getIntegrand_dsigmaBox_bottom200
double getIntegrand_dsigmaBox_bottom200(double x) const
Definition: StandardModel.cpp:7552
StandardModel::getIntegrand_dsigmaBox_bottom192
double getIntegrand_dsigmaBox_bottom192(double x) const
Definition: StandardModel.cpp:7540
StandardModel::getIntegrand_sigmaWithISR_mu202
double getIntegrand_sigmaWithISR_mu202(double x) const
Definition: StandardModel.cpp:6413
StandardModel::computeSigmaWF
double computeSigmaWF(const double sqrt_s) const
The W fusion contribution to higgs-production cross section in the Standard Model.
Definition: StandardModel.h:2031
StandardModel::sigma_NoISR_l
double sigma_NoISR_l(const QCD::lepton l_flavor, const double s) const
Definition: StandardModel.cpp:6298
StandardModel::getIntegrand_dsigmaBox_up196
double getIntegrand_dsigmaBox_up196(double x) const
Definition: StandardModel.cpp:7193
StandardModel::FlagSMAux
bool FlagSMAux
A boolean for the model flag SMAux.
Definition: StandardModel.h:3137
StandardModel::getYu
gslpp::matrix< gslpp::complex > getYu() const
A get method to retrieve the Yukawa matrix of the up-type quarks, .
Definition: StandardModel.h:898
StandardModel::getIntegrand_AFBnumeratorWithISR_mu172
double getIntegrand_AFBnumeratorWithISR_mu172(double x) const
Definition: StandardModel.cpp:7618
StandardModel::getIntegrand_dsigmaBox_up136
double getIntegrand_dsigmaBox_up136(double x) const
Definition: StandardModel.cpp:7151
QCD::UP
Definition: QCD.h:324
StandardModel::getIntegrand_dsigmaBox_bottom133
double getIntegrand_dsigmaBox_bottom133(double x) const
Definition: StandardModel.cpp:7498
StandardModel::getIntegrand_dsigmaBox_up161
double getIntegrand_dsigmaBox_up161(double x) const
Definition: StandardModel.cpp:7157
StandardModel::GF
double GF
The Fermi constant in .
Definition: StandardModel.h:2555
StandardModel::AFB_NoISR_q
double AFB_NoISR_q(const QCD::quark q_flavor, const double s) const
Definition: StandardModel.cpp:6333
StandardModel::A_f
virtual double A_f(const Particle f) const
The left-right asymmetry in at the -pole, .
Definition: StandardModel.cpp:1209
StandardModel::getIntegrand_AFBnumeratorWithISR_mu192
double getIntegrand_AFBnumeratorWithISR_mu192(double x) const
Definition: StandardModel.cpp:7636
StandardModel::sigma0_had
virtual double sigma0_had() const
The hadronic cross section for at the -pole, .
Definition: StandardModel.cpp:1370
StandardModel::getMyThreeLoopQCD
EWSMThreeLoopQCD * getMyThreeLoopQCD() const
Definition: StandardModel.h:1005
StandardModel::getIntegrand_AFBnumeratorWithISR_charm183
double getIntegrand_AFBnumeratorWithISR_charm183(double x) const
Definition: StandardModel.cpp:7780
StandardModel::computeGammaHgaga_tW
double computeGammaHgaga_tW() const
The mixed loop contribution to in the Standard Model.
Definition: StandardModel.h:2424
StandardModel::getIntegrand_sigmaWithISR_tau202
double getIntegrand_sigmaWithISR_tau202(double x) const
Definition: StandardModel.cpp:6486
StandardModel::getIntegrand_AFBnumeratorWithISR_charm172
double getIntegrand_AFBnumeratorWithISR_charm172(double x) const
Definition: StandardModel.cpp:7774
StandardModel::getIntegrand_sigmaWithISR_up167
double getIntegrand_sigmaWithISR_up167(double x) const
Definition: StandardModel.cpp:6552
StandardModel::getIntegrand_dsigmaBox_charm133
double getIntegrand_dsigmaBox_charm133(double x) const
Definition: StandardModel.cpp:7320
StandardModel::getIntegrand_sigmaWithISR_up189
double getIntegrand_sigmaWithISR_up189(double x) const
Definition: StandardModel.cpp:6570
StandardModel::getmq
virtual double getmq(const QCD::quark q, const double mu) const
Definition: StandardModel.h:2514
StandardModel::delMw
double delMw
The theoretical uncertainty in , denoted as , in GeV.
Definition: StandardModel.h:2559
StandardModel::getIntegrand_dsigmaBox_bottom130
double getIntegrand_dsigmaBox_bottom130(double x) const
Definition: StandardModel.cpp:7492
StandardModel::SchemeToDouble
double SchemeToDouble(const std::string scheme) const
A method to convert a given scheme name in string form into a floating-point number with double preci...
Definition: StandardModel.h:2601
StandardModel::computeGammaHgg_tb
double computeGammaHgg_tb() const
The top-bottom interference contribution to in the Standard Model.
Definition: StandardModel.h:2358
StandardModel::getIntegrand_AFBnumeratorWithISR_mu161
double getIntegrand_AFBnumeratorWithISR_mu161(double x) const
Definition: StandardModel.cpp:7612
StandardModel::Beta_e
double Beta_e(int nm, unsigned int nf) const
QED beta function coefficients - eq. (36) hep-ph/0512066.
Definition: StandardModel.cpp:582
StandardModel::CheckParameters
virtual bool CheckParameters(const std::map< std::string, double > &DPars)
A method to check if all the mandatory parameters for StandardModel have been provided in model initi...
Definition: StandardModel.cpp:339
StandardModel::getIntegrand_dsigmaBox_charm136
double getIntegrand_dsigmaBox_charm136(double x) const
Definition: StandardModel.cpp:7326
StandardModel::alphaMz
double alphaMz() const
The electromagnetic coupling at the -mass scale, .
Definition: StandardModel.cpp:893
StandardModel::getIntegrand_dsigmaBox_down130
double getIntegrand_dsigmaBox_down130(double x) const
Definition: StandardModel.cpp:7226
StandardModel::LEP2AFBcharm
virtual double LEP2AFBcharm(const double s) const
Definition: StandardModel.cpp:4882
QCD::CHARM
Definition: QCD.h:326
StandardModel::DeltaAlphaL5q
double DeltaAlphaL5q() const
The sum of the leptonic and the five-flavour hadronic corrections to the electromagnetic coupling at...
Definition: StandardModel.cpp:856
StandardModel::bSigmaForAFB
bool bSigmaForAFB
Definition: StandardModel.h:2832
StandardModel::getIntegrand_sigmaWithISR_down172
double getIntegrand_sigmaWithISR_down172(double x) const
Definition: StandardModel.cpp:6646
StandardModel::computeBrHtobb
double computeBrHtobb() const
The Br in the Standard Model.
Definition: StandardModel.h:2313
StandardModel::getYe
gslpp::matrix< gslpp::complex > getYe() const
A get method to retrieve the Yukawa matrix of the charged leptons, .
Definition: StandardModel.h:928
StandardModel::getCCC2
virtual double getCCC2() const
A virtual implementation for the RealWeakEFTCC class.
Definition: StandardModel.h:1163
gslpp::complex
A class for defining operations on and functions of complex numbers.
Definition: gslpp_complex.h:35
StandardModel::LEP2AFBtau
virtual double LEP2AFBtau(const double s) const
Definition: StandardModel.cpp:5784
PMNS.h
StandardModel::getIntegrand_dsigmaBox_charm172
double getIntegrand_dsigmaBox_charm172(double x) const
Definition: StandardModel.cpp:7344
StandardModel::getIntegrand_sigmaWithISR_mu183
double getIntegrand_sigmaWithISR_mu183(double x) const
Definition: StandardModel.cpp:6383
StandardModel::ComputeDeltaR_rem
void ComputeDeltaR_rem(const double Mw_i, double DeltaR_rem[orders_EW_size]) const
A method to collect computed via subclasses.
Definition: StandardModel.cpp:1079
StandardModel::mHl
double mHl
The Higgs mass in GeV.
Definition: StandardModel.h:2558
StandardModel::epsilonb
virtual double epsilonb() const
The SM contribution to the epsilon parameter .
Definition: StandardModel.cpp:1815
StandardModel::getIntegrand_dsigmaBox_tau189
double getIntegrand_dsigmaBox_tau189(double x) const
Definition: StandardModel.cpp:7080
StandardModel::getPhiBs
virtual double getPhiBs() const
Half the relative phase of the $B_s$ mixing amplitude w.r.t. the Standard Model one.
Definition: StandardModel.h:2470
StandardModel::s23
double s23
Definition: StandardModel.h:2577
StandardModel::f_GSL
gsl_function f_GSL
Definition: StandardModel.h:3166
StandardModel::Integrand_AFBnumeratorWithISR_q
double Integrand_AFBnumeratorWithISR_q(double x, const QCD::quark q_flavor, const double s) const
Definition: StandardModel.cpp:7747
StandardModel::getIntegrand_dsigmaBox_charm183
double getIntegrand_dsigmaBox_charm183(double x) const
Definition: StandardModel.cpp:7350
StandardModel::GammaW
virtual double GammaW() const
The total width of the boson, .
Definition: StandardModel.cpp:1190
StandardModel::LEP2Rbottom
virtual double LEP2Rbottom(const double s) const
Definition: StandardModel.cpp:6274
gslpp::matrix< gslpp::complex >
StandardModel::getIntegrand_dsigmaBox_strange202
double getIntegrand_dsigmaBox_strange202(double x) const
Definition: StandardModel.cpp:7470
CKM.h
StandardModel::flag_order
bool flag_order[orders_EW_size]
An array of internal flags controlling the inclusions of higher-order corrections.
Definition: StandardModel.h:2589
StandardModel::getIntegrand_sigmaWithISR_charm202
double getIntegrand_sigmaWithISR_charm202(double x) const
Definition: StandardModel.cpp:6770
StandardModel::getIntegrand_sigmaWithISR_up192
double getIntegrand_sigmaWithISR_up192(double x) const
Definition: StandardModel.cpp:6576
StandardModel::getIntegrand_sigmaWithISR_up196
double getIntegrand_sigmaWithISR_up196(double x) const
Definition: StandardModel.cpp:6582
StandardModel::getMyThreeLoopEW
EWSMThreeLoopEW * getMyThreeLoopEW() const
Definition: StandardModel.h:995
StandardModel::getIntegrand_dsigmaBox_down200
double getIntegrand_dsigmaBox_down200(double x) const
Definition: StandardModel.cpp:7286
StandardModel::GammaW_cache
double GammaW_cache
A cache of the value of .
Definition: StandardModel.h:3147
StandardModel::getIntegrand_dsigmaBox_down133
double getIntegrand_dsigmaBox_down133(double x) const
Definition: StandardModel.cpp:7232
StandardModel::getCepsK
virtual double getCepsK() const
The ratio of the imaginary part of the $K$ mixing amplitude over the Standard Model value.
Definition: StandardModel.h:2461
StandardModel::getIntegrand_dsigmaBox_charm189
double getIntegrand_dsigmaBox_charm189(double x) const
Definition: StandardModel.cpp:7356
StandardModel
A model class for the Standard Model.
Definition: StandardModel.h:477
StandardModel::DeltaAlphaLepton_cache
double DeltaAlphaLepton_cache
A cache of the value of .
Definition: StandardModel.h:3144
EWSMOneLoopEW
A class for one-loop corrections to the EW precision observables.
Definition: EWSMOneLoopEW.h:105
StandardModel::DeltaAlpha_cache
double DeltaAlpha_cache
A cache of the value of .
Definition: StandardModel.h:3145
StandardModel::FlagRhoZ
std::string FlagRhoZ
A string for the model flag RhoZ.
Definition: StandardModel.h:3133
StandardModel::DeltaAlpha
double DeltaAlpha() const
The total corrections to the electromagnetic coupling at the -mass scale, denoted as .
Definition: StandardModel.cpp:881
StandardModel::SMM
Matching< StandardModelMatching, StandardModel > SMM
An object of type Matching.
Definition: StandardModel.h:2550
StandardModel::getIntegrand_sigmaWithISR_tau161
double getIntegrand_sigmaWithISR_tau161(double x) const
Definition: StandardModel.cpp:6444
StandardModel::getDelMw
double getDelMw() const
A get method to retrieve the theoretical uncertainty in , denoted as .
Definition: StandardModel.h:778
StandardModel::getIntegrand_sigmaWithISR_tau205
double getIntegrand_sigmaWithISR_tau205(double x) const
Definition: StandardModel.cpp:6492
StandardModel::DeltaR
virtual double DeltaR() const
The SM prediction for derived from that for the boson mass.
Definition: StandardModel.cpp:1036
StandardModel::getIntegrand_dsigmaBox_up172
double getIntegrand_dsigmaBox_up172(double x) const
Definition: StandardModel.cpp:7169
StandardModel::getIntegrand_dsigmaBox_down136
double getIntegrand_dsigmaBox_down136(double x) const
Definition: StandardModel.cpp:7238
StandardModel::getIntegrand_AFBnumeratorWithISR_mu207
double getIntegrand_AFBnumeratorWithISR_mu207(double x) const
Definition: StandardModel.cpp:7666
StandardModel::getIntegrand_sigmaWithISR_up183
double getIntegrand_sigmaWithISR_up183(double x) const
Definition: StandardModel.cpp:6564
StandardModel::cW2
virtual double cW2() const
Definition: StandardModel.cpp:1020
StandardModel::getIntegrand_dsigmaBox_strange189
double getIntegrand_dsigmaBox_strange189(double x) const
Definition: StandardModel.cpp:7446
StandardModel::myThreeLoopQCD
EWSMThreeLoopQCD * myThreeLoopQCD
A pointer to an object of type EWSMThreeLoopQCD.
Definition: StandardModel.h:3120
StandardModel::getIntegrand_dsigmaBox_up167
double getIntegrand_dsigmaBox_up167(double x) const
Definition: StandardModel.cpp:7163
StandardModel::getIntegrand_dsigmaBox_charm196
double getIntegrand_dsigmaBox_charm196(double x) const
Definition: StandardModel.cpp:7368
StandardModel::getIntegrand_sigmaWithISR_strange133
double getIntegrand_sigmaWithISR_strange133(double x) const
Definition: StandardModel.cpp:6798
StandardModel::LEP2AFBbottom
virtual double LEP2AFBbottom(const double s) const
Definition: StandardModel.cpp:4470
StandardModel::ale
double ale
The fine-structure constant .
Definition: StandardModel.h:2556
StandardModel::getIntegrand_sigmaWithISR_charm161
double getIntegrand_sigmaWithISR_charm161(double x) const
Definition: StandardModel.cpp:6722
StandardModel::Alstilde5
double Alstilde5(const double mu) const
The value of at any scale with the number of flavours and full EW corrections.
Definition: StandardModel.cpp:899
StandardModel::setFlag
virtual bool setFlag(const std::string name, const bool value)
A method to set a flag of StandardModel.
Definition: StandardModel.cpp:404
StandardModel::myTwoLoopEW
EWSMTwoLoopEW * myTwoLoopEW
A pointer to an object of type EWSMTwoLoopEW.
Definition: StandardModel.h:3121
StandardModel::DeltaAlphaTop
double DeltaAlphaTop(const double s) const
Top-quark contribution to the electromagnetic coupling , denoted as .
Definition: StandardModel.cpp:862
StandardModel::getIntegrand_dsigmaBox_bottom207
double getIntegrand_dsigmaBox_bottom207(double x) const
Definition: StandardModel.cpp:7570
StandardModel::getIntegrand_dsigmaBox_charm205
double getIntegrand_dsigmaBox_charm205(double x) const
Definition: StandardModel.cpp:7386
StandardModel::AlsByOrder
double AlsByOrder(double mu, orders order=FULLNLO, bool qed_flag=false, bool Nf_thr=true) const
Definition: StandardModel.cpp:623
StandardModel::getIntegrand_sigmaWithISR_bottom172
double getIntegrand_sigmaWithISR_bottom172(double x) const
Definition: StandardModel.cpp:6910
StandardModel::useDeltaAlphaLepton_cache
bool useDeltaAlphaLepton_cache
Definition: StandardModel.h:3150
StandardModel::DeltaAlphaLepton
double DeltaAlphaLepton(const double s) const
Leptonic contribution to the electromagnetic coupling , denoted as .
Definition: StandardModel.cpp:828
StandardModel::IsFlagWithoutNonUniversalVC
bool IsFlagWithoutNonUniversalVC() const
A method to retrieve the model flag WithoutNonUniversalVC.
Definition: StandardModel.h:634
StandardModel::getIntegrand_AFBnumeratorWithISR_mu196
double getIntegrand_AFBnumeratorWithISR_mu196(double x) const
Definition: StandardModel.cpp:7642
EWSMTwoFermionsLEP2
A class for the form factors , and in the processes at LEP-II.
Definition: EWSMTwoFermionsLEP2.h:25
StandardModel::RVh
double RVh() const
The singlet vector corrections to the hadronic -boson width, denoted as .
Definition: StandardModel.cpp:2399
StandardModel::getIntegrand_dsigmaBox_charm167
double getIntegrand_dsigmaBox_charm167(double x) const
Definition: StandardModel.cpp:7338
StandardModel::getIntegrand_dsigmaBox_tau205
double getIntegrand_dsigmaBox_tau205(double x) const
Definition: StandardModel.cpp:7110
StandardModel::getIntegrand_dsigmaBox_bottom196
double getIntegrand_dsigmaBox_bottom196(double x) const
Definition: StandardModel.cpp:7546
StandardModel::epsilon2
virtual double epsilon2() const
The SM contribution to the epsilon parameter .
Definition: StandardModel.cpp:1793
LeptonFlavour
The parent class in LeptonFlavour for calculating all the Wilson coefficients for various Lepton Flav...
Definition: LeptonFlavour.h:26
StandardModel::StandardModel
StandardModel()
The default constructor.
Definition: StandardModel.cpp:40
StandardModel::getIntegrand_sigmaWithISR_charm196
double getIntegrand_sigmaWithISR_charm196(double x) const
Definition: StandardModel.cpp:6758
StandardModel::getIntegrand_dsigmaBox_strange207
double getIntegrand_dsigmaBox_strange207(double x) const
Definition: StandardModel.cpp:7482
StandardModel::getIntegrand_dsigmaBox_tau202
double getIntegrand_dsigmaBox_tau202(double x) const
Definition: StandardModel.cpp:7104
StandardModel::computeSigmaZF
double computeSigmaZF(const double sqrt_s) const
The Z fusion contribution to higgs-production cross section in the Standard Model.
Definition: StandardModel.h:2059
Matching.h
StandardModel::getIntegrand_AFBnumeratorWithISR_mu130
double getIntegrand_AFBnumeratorWithISR_mu130(double x) const
Definition: StandardModel.cpp:7600
StandardModel::delGammaZ
double delGammaZ
The theoretical uncertainty in , denoted as , in GeV.
Definition: StandardModel.h:2563
StandardModel::getIntegrand_AFBnumeratorWithISR_mu205
double getIntegrand_AFBnumeratorWithISR_mu205(double x) const
Definition: StandardModel.cpp:7660
StandardModel::computeBrHtoZZ
double computeBrHtoZZ() const
The Br in the Standard Model.
Definition: StandardModel.h:2222
StandardModel::computeGammaHgg_bb
double computeGammaHgg_bb() const
The bottom loop contribution to in the Standard Model.
Definition: StandardModel.h:2347
StandardModel::rhoZ_f_cache
gslpp::complex rhoZ_f_cache[12]
A cache of the value of .
Definition: StandardModel.h:3148
StandardModel::getIntegrand_AFBnumeratorWithISR_charm205
double getIntegrand_AFBnumeratorWithISR_charm205(double x) const
Definition: StandardModel.cpp:7816
StandardModel::c02
double c02() const
The square of the cosine of the weak mixing angle defined without weak radiative corrections.
Definition: StandardModel.cpp:965
StandardModel::getIntegrand_dsigmaBox_tau167
double getIntegrand_dsigmaBox_tau167(double x) const
StandardModel::getIntegrand_dsigmaBox_tau133
double getIntegrand_dsigmaBox_tau133(double x) const
StandardModel::Integrand_dsigmaBox_l
double Integrand_dsigmaBox_l(double cosTheta, const QCD::lepton l_flavor, const double s) const
Definition: StandardModel.cpp:6968
StandardModel::getIntegrand_sigmaWithISR_tau172
double getIntegrand_sigmaWithISR_tau172(double x) const
Definition: StandardModel.cpp:6450
StandardModel::SMresult_cache
double SMresult_cache
Definition: StandardModel.h:3160
StandardModel::Init
virtual bool Init(const std::map< std::string, double > &DPars)
A method to initialize the model parameters.
Definition: StandardModel.cpp:185
StandardModel::w_GSL1
gsl_integration_workspace * w_GSL1
Definition: StandardModel.h:3167
StandardModel::getIntegrand_dsigmaBox_up130
double getIntegrand_dsigmaBox_up130(double x) const
Definition: StandardModel.cpp:7139
StandardModel::computeGammaHgaga_tt
double computeGammaHgaga_tt() const
The top loop contribution to in the Standard Model.
Definition: StandardModel.h:2402
StandardModelMatching
A class for the matching in the Standard Model.
Definition: StandardModelMatching.h:26
StandardModel::ale_OS
double ale_OS(const double mu, orders order=FULLNLO) const
The running electromagnetic coupling in the on-shell scheme.
Definition: StandardModel.cpp:533
StandardModel::dAle5Mz
double dAle5Mz
The five-flavour hadronic contribution to the electromagnetic coupling, .
Definition: StandardModel.h:2557
StandardModel::Vub
double Vub
used as an input for FlagWolfenstein = FALSE
Definition: StandardModel.h:2574
StandardModel::getIntegrand_dsigmaBox_up192
double getIntegrand_dsigmaBox_up192(double x) const
Definition: StandardModel.cpp:7187
StandardModel::getIntegrand_sigmaWithISR_bottom183
double getIntegrand_sigmaWithISR_bottom183(double x) const
Definition: StandardModel.cpp:6916
StandardModel::useMw_cache
bool useMw_cache
Definition: StandardModel.h:3152
StandardModel::AlsMz
double AlsMz
The strong coupling constant at the Z-boson mass, .
Definition: StandardModel.h:2553
StandardModel::getIntegrand_dsigmaBox_strange130
double getIntegrand_dsigmaBox_strange130(double x) const
Definition: StandardModel.cpp:7402
StandardModel::getIntegrand_sigmaWithISR_up200
double getIntegrand_sigmaWithISR_up200(double x) const
Definition: StandardModel.cpp:6588
StandardModel::getIntegrand_sigmaWithISR_strange136
double getIntegrand_sigmaWithISR_strange136(double x) const
Definition: StandardModel.cpp:6804
StandardModel::sW2
double sW2() const
Definition: StandardModel.cpp:1031
StandardModel::getIntegrand_AFBnumeratorWithISR_mu189
double getIntegrand_AFBnumeratorWithISR_mu189(double x) const
Definition: StandardModel.cpp:7630
CKM::getCKM
gslpp::matrix< gslpp::complex > getCKM() const
A member for returning the CKM matrix.
Definition: CKM.h:49
StandardModel::getCCC4
virtual double getCCC4() const
A virtual implementation for the RealWeakEFTCC class.
Definition: StandardModel.h:1173
StandardModel::getIntegrand_dsigmaBox_down189
double getIntegrand_dsigmaBox_down189(double x) const
Definition: StandardModel.cpp:7268
StandardModel::EW1QCD1
Two-loop of .
Definition: StandardModel.h:497
Flavour.h
StandardModel::PreUpdate
virtual bool PreUpdate()
The pre-update method for StandardModel.
Definition: StandardModel.cpp:198
StandardModel::Gamma_had
virtual double Gamma_had() const
The hadronic decay width of the boson, .
Definition: StandardModel.cpp:1309
StandardModel::resumRhoZ
double resumRhoZ(const double DeltaRho[orders_EW_size], const double deltaRho_rem[orders_EW_size], const double DeltaRbar_rem, const bool bool_Zbb) const
A method to compute the real part of the effective coupling from , and .
Definition: StandardModel.cpp:1948
StandardModel::computeSigmaggH_bb
double computeSigmaggH_bb(const double sqrt_s) const
The square of the bottom-quark contribution to the ggH cross section in the Standard Model.
Definition: StandardModel.h:1949
StandardModel::getIntegrand_dsigmaBox_bottom136
double getIntegrand_dsigmaBox_bottom136(double x) const
Definition: StandardModel.cpp:7504
StandardModel::getIntegrand_dsigmaBox_tau136
double getIntegrand_dsigmaBox_tau136(double x) const
Definition: StandardModel.cpp:7056
QCD::TOP
Definition: QCD.h:328
StandardModel::getIntegrand_dsigmaBox_up200
double getIntegrand_dsigmaBox_up200(double x) const
Definition: StandardModel.cpp:7199
StandardModel::getIntegrand_sigmaWithISR_up130
double getIntegrand_sigmaWithISR_up130(double x) const
Definition: StandardModel.cpp:6528
StandardModel::setFlagSigmaForR
bool setFlagSigmaForR(const bool flagSigmaForR_i)
Definition: StandardModel.h:2508
StandardModel::getIntegrand_dsigmaBox_bottom167
double getIntegrand_dsigmaBox_bottom167(double x) const
Definition: StandardModel.cpp:7516
StandardModel::computeGammaHgaga_WW
double computeGammaHgaga_WW() const
The loop contribution to in the Standard Model.
Definition: StandardModel.h:2413
StandardModel::computeGammaHTotal
double computeGammaHTotal() const
The Higgs total width in the Standard Model.
Definition: StandardModel.h:2325
StandardModel::getIntegrand_sigmaWithISR_bottom192
double getIntegrand_sigmaWithISR_bottom192(double x) const
Definition: StandardModel.cpp:6928
StandardModel::getIntegrand_AFBnumeratorWithISR_charm202
double getIntegrand_AFBnumeratorWithISR_charm202(double x) const
Definition: StandardModel.cpp:7810
StandardModel::RAq
double RAq(const QCD::quark q) const
The radiator factor associated with the final-state QED and QCD corrections to the the axial-vector-c...
Definition: StandardModel.cpp:2262
StandardModel::getIntegrand_sigmaWithISR_down189
double getIntegrand_sigmaWithISR_down189(double x) const
Definition: StandardModel.cpp:6658
StandardModel::getIntegrand_AFBnumeratorWithISR_charm192
double getIntegrand_AFBnumeratorWithISR_charm192(double x) const
Definition: StandardModel.cpp:7792
StandardModel::getIntegrand_sigmaWithISR_charm133
double getIntegrand_sigmaWithISR_charm133(double x) const
Definition: StandardModel.cpp:6710
StandardModel::getGF
double getGF() const
A get method to retrieve the Fermi constant .
Definition: StandardModel.h:739
StandardModel::SMvars
static std::string SMvars[NSMvars]
A string array containing the labels of the model parameters in StandardModel.
Definition: StandardModel.h:512
StandardModel::getIntegrand_sigmaWithISR_mu172
double getIntegrand_sigmaWithISR_mu172(double x) const
Definition: StandardModel.cpp:6377
StandardModel::myLeptonFlavour
LeptonFlavour * myLeptonFlavour
A pointer to an object of the type LeptonFlavour.
Definition: StandardModel.h:3125
StandardModel::getIntegrand_sigmaWithISR_mu207
double getIntegrand_sigmaWithISR_mu207(double x) const
Definition: StandardModel.cpp:6425
StandardModel::getIntegrand_dsigmaBox_tau161
double getIntegrand_dsigmaBox_tau161(double x) const
Definition: StandardModel.cpp:7062
StandardModel::getIntegrand_sigmaWithISR_down192
double getIntegrand_sigmaWithISR_down192(double x) const
Definition: StandardModel.cpp:6664
StandardModel::getIntegrand_dsigmaBox_up189
double getIntegrand_dsigmaBox_up189(double x) const
Definition: StandardModel.cpp:7181
StandardModel::getDAle5Mz
double getDAle5Mz() const
A get method to retrieve the five-flavour hadronic contribution to the electromagnetic coupling,...
Definition: StandardModel.h:759
StandardModel::EW3
Three-loop of .
Definition: StandardModel.h:501
StandardModel::setFlagSigmaForAFB
bool setFlagSigmaForAFB(const bool flagSigmaForAFB_i)
Definition: StandardModel.h:2502
EWSMcache
A class for cache variables used in computing radiative corrections to the EW precision observables.
Definition: EWSMcache.h:40
StandardModel::getTrueSM
virtual StandardModel getTrueSM() const
Definition: StandardModel.h:943
StandardModel::FlagMw
std::string FlagMw
A string for the model flag Mw.
Definition: StandardModel.h:3132
StandardModel::computeBrHtoZga
double computeBrHtoZga() const
The Br in the Standard Model.
Definition: StandardModel.h:2244
QCD::getMtpole
double getMtpole() const
A get method to access the pole mass of the top quark.
Definition: QCD.h:588
StandardModel::getIntegrand_sigmaWithISR_charm189
double getIntegrand_sigmaWithISR_charm189(double x) const
Definition: StandardModel.cpp:6746
StandardModel::getIntegrand_AFBnumeratorWithISR_tau130
double getIntegrand_AFBnumeratorWithISR_tau130(double x) const
Definition: StandardModel.cpp:7673
StandardModel::lambda
double lambda
The CKM parameter in the Wolfenstein parameterization.
Definition: StandardModel.h:2568
StandardModel::MwFromMwbar
double MwFromMwbar(const double Mwbar) const
A method to convert the -boson mass in the complex-pole/fixed-width scheme to that in the experimenta...
Definition: StandardModel.cpp:1137
StandardModel::getIntegrand_AFBnumeratorWithISR_mu200
double getIntegrand_AFBnumeratorWithISR_mu200(double x) const
Definition: StandardModel.cpp:7648
StandardModel::getIntegrand_sigmaWithISR_strange192
double getIntegrand_sigmaWithISR_strange192(double x) const
Definition: StandardModel.cpp:6840
StandardModel::getIntegrand_AFBnumeratorWithISR_bottom207
double getIntegrand_AFBnumeratorWithISR_bottom207(double x) const
Definition: StandardModel.cpp:7891
StandardModel::getIntegrand_dsigmaBox_strange133
double getIntegrand_dsigmaBox_strange133(double x) const
Definition: StandardModel.cpp:7408
EWSMTwoLoopEW
A class for two-loop corrections to the EW precision observables.
Definition: EWSMTwoLoopEW.h:57
EWSMApproximateFormulae
A class for approximate formulae of the EW precision observables.
Definition: EWSMApproximateFormulae.h:33
StandardModel::useDeltaAlpha_cache
bool useDeltaAlpha_cache
Definition: StandardModel.h:3151
StandardModel::getUPMNS
gslpp::matrix< gslpp::complex > getUPMNS() const
A get method to retrieve the object of the PMNS matrix.
Definition: StandardModel.h:888
StandardModel::getIntegrand_dsigmaBox_bottom183
double getIntegrand_dsigmaBox_bottom183(double x) const
Definition: StandardModel.cpp:7528
StandardModel::computeBrHtogaga
double computeBrHtogaga() const
The Br in the Standard Model.
Definition: StandardModel.h:2256
StandardModel::getIntegrand_sigmaWithISR_down161
double getIntegrand_sigmaWithISR_down161(double x) const
Definition: StandardModel.cpp:6634
StandardModel::computeSigmaggH
double computeSigmaggH(const double sqrt_s) const
The ggH cross section in the Standard Model.
Definition: StandardModel.h:1897
StandardModel::getIntegrand_dsigmaBox_strange161
double getIntegrand_dsigmaBox_strange161(double x) const
Definition: StandardModel.cpp:7420
StandardModel::getIntegrand_AFBnumeratorWithISR_tau207
double getIntegrand_AFBnumeratorWithISR_tau207(double x) const
Definition: StandardModel.cpp:7739
StandardModel::PostUpdate
virtual bool PostUpdate()
The post-update method for StandardModel.
Definition: StandardModel.cpp:225
StandardModel::getCBs
virtual double getCBs() const
The ratio of the absolute value of the $B_s$ mixing amplitude over the Standard Model value.
Definition: StandardModel.h:2443
StandardModel::getYn
gslpp::matrix< gslpp::complex > getYn() const
A get method to retrieve the Yukawa matrix of the neutrinos, .
Definition: StandardModel.h:918
StandardModel::IsFlagNoApproximateGammaZ
bool IsFlagNoApproximateGammaZ() const
A method to retrieve the model flag NoApproximateGammaZ.
Definition: StandardModel.h:647
StandardModel::getIntegrand_sigmaWithISR_strange161
double getIntegrand_sigmaWithISR_strange161(double x) const
Definition: StandardModel.cpp:6810
StandardModel::myOneLoopEW
EWSMOneLoopEW * myOneLoopEW
A pointer to an object of type EWSMOneLoopEW.
Definition: StandardModel.h:3118
StandardModel::setFlagStr
virtual bool setFlagStr(const std::string name, const std::string value)
A method to set a flag of StandardModel.
Definition: StandardModel.cpp:444
StandardModel::LEP2sigmaCharm
virtual double LEP2sigmaCharm(const double s) const
Definition: StandardModel.cpp:2924
StandardModel::Mw_error
static const double Mw_error
The target accuracy of the iterative calculation of the -boson mass in units of GeV.
Definition: StandardModel.h:520
StandardModel::getMyTwoLoopEW
EWSMTwoLoopEW * getMyTwoLoopEW() const
Definition: StandardModel.h:1010
StandardModel::Ale
double Ale(double mu, orders order, bool Nf_thr=true) const
The running electromagnetic coupling in the scheme.
Definition: StandardModel.cpp:732
StandardModel::epsilon1
virtual double epsilon1() const
The SM contribution to the epsilon parameter .
Definition: StandardModel.cpp:1785
StandardModel::getIntegrand_dsigmaBox_charm200
double getIntegrand_dsigmaBox_charm200(double x) const
Definition: StandardModel.cpp:7374
StandardModel::computeGammaHZga_tW
double computeGammaHZga_tW() const
The mixed loop contribution to in the Standard Model.
Definition: StandardModel.h:2391
StandardModel::getIntegrand_sigmaWithISR_charm200
double getIntegrand_sigmaWithISR_charm200(double x) const
Definition: StandardModel.cpp:6764
StandardModel::getIntegrand_sigmaWithISR_down133
double getIntegrand_sigmaWithISR_down133(double x) const
Definition: StandardModel.cpp:6621
StandardModel::getIntegrand_sigmaWithISR_bottom205
double getIntegrand_sigmaWithISR_bottom205(double x) const
Definition: StandardModel.cpp:6952
StandardModel::s13
double s13
Definition: StandardModel.h:2577
StandardModel::Als
double Als(double mu, orders order=FULLNLO, bool qed_flag=false, bool Nf_thr=true) const
The running QCD coupling in the scheme including QED corrections.
Definition: StandardModel.cpp:602
StandardModel::computeSigmattH
double computeSigmattH(const double sqrt_s) const
The ttH production cross section in the Standard Model.
Definition: StandardModel.h:2171
StandardModel::checkEWPOscheme
bool checkEWPOscheme(const std::string scheme) const
A method to check if a given scheme name in string form is valid.
Definition: StandardModel.h:2622
QCD::getQuarks
Particle getQuarks(const QCD::quark q) const
A get method to access a quark as an object of the type Particle.
Definition: QCD.h:534
StandardModel::getIntegrand_dsigmaBox_strange167
double getIntegrand_dsigmaBox_strange167(double x) const
Definition: StandardModel.cpp:7426
StandardModel::getIntegrand_dsigmaBox_tau130
double getIntegrand_dsigmaBox_tau130(double x) const
Definition: StandardModel.cpp:7050
StandardModel::getIntegrand_dsigmaBox_up202
double getIntegrand_dsigmaBox_up202(double x) const
Definition: StandardModel.cpp:7205
StandardModel::getDelSin2th_l
double getDelSin2th_l() const
A get method to retrieve the theoretical uncertainty in , denoted as .
Definition: StandardModel.h:789
StandardModel::deltaKappaZ_f
virtual gslpp::complex deltaKappaZ_f(const Particle f) const
Flavour non-universal vertex corrections to , denoted by .
Definition: StandardModel.cpp:1755
StandardModel::getIntegrand_sigmaWithISR_strange130
double getIntegrand_sigmaWithISR_strange130(double x) const
Definition: StandardModel.cpp:6792
StandardModel::getIntegrand_sigmaWithISR_tau189
double getIntegrand_sigmaWithISR_tau189(double x) const
Definition: StandardModel.cpp:6462
StandardModel::getIntegrand_sigmaWithISR_bottom167
double getIntegrand_sigmaWithISR_bottom167(double x) const
Definition: StandardModel.cpp:6904
StandardModel::getIntegrand_dsigmaBox_mu200
double getIntegrand_dsigmaBox_mu200(double x) const
Definition: StandardModel.cpp:7022
StandardModel::getVCKM
gslpp::matrix< gslpp::complex > getVCKM() const
A get method to retrieve the CKM matrix.
Definition: StandardModel.h:870
StandardModel::FlagKappaZ
std::string FlagKappaZ
A string for the model flag KappaZ.
Definition: StandardModel.h:3134
StandardModel::getMyApproximateFormulae
EWSMApproximateFormulae * getMyApproximateFormulae() const
A get method to retrieve the member pointer of type EWSMApproximateFormulae.
Definition: StandardModel.h:979
StandardModel::getIntegrand_AFBnumeratorWithISR_tau200
double getIntegrand_AFBnumeratorWithISR_tau200(double x) const
Definition: StandardModel.cpp:7721
StandardModel::getIntegrand_AFBnumeratorWithISR_charm189
double getIntegrand_AFBnumeratorWithISR_charm189(double x) const
Definition: StandardModel.cpp:7786
StandardModel::getIntegrand_dsigmaBox_mu136
double getIntegrand_dsigmaBox_mu136(double x) const
Definition: StandardModel.cpp:6980
StandardModel::getIntegrand_AFBnumeratorWithISR_tau183
double getIntegrand_AFBnumeratorWithISR_tau183(double x) const
Definition: StandardModel.cpp:7697
QCD::quark
quark
An enum type for quarks.
Definition: QCD.h:323
StandardModel::getIntegrand_AFBnumeratorWithISR_bottom192
double getIntegrand_AFBnumeratorWithISR_bottom192(double x) const
Definition: StandardModel.cpp:7861
StandardModel::CacheSize
static const int CacheSize
Defines the depth of the cache.
Definition: StandardModel.h:3174
StandardModel::getIntegrand_sigmaWithISR_down183
double getIntegrand_sigmaWithISR_down183(double x) const
Definition: StandardModel.cpp:6652
StandardModel::getIntegrand_dsigmaBox_down183
double getIntegrand_dsigmaBox_down183(double x) const
Definition: StandardModel.cpp:7262
StandardModel::getDelSin2th_q
double getDelSin2th_q() const
A get method to retrieve the theoretical uncertainty in , denoted as .
Definition: StandardModel.h:800
StandardModel::getFlavour
const Flavour & getFlavour() const
Definition: StandardModel.h:1020
StandardModel::R0_f
virtual double R0_f(const Particle f) const
The ratio .
Definition: StandardModel.cpp:1395
StandardModel::getIntegrand_AFBnumeratorWithISR_bottom172
double getIntegrand_AFBnumeratorWithISR_bottom172(double x) const
Definition: StandardModel.cpp:7843
StandardModel::getIntegrand_sigmaWithISR_tau196
double getIntegrand_sigmaWithISR_tau196(double x) const
Definition: StandardModel.cpp:6474
StandardModel::AleWithInit
double AleWithInit(double mu, double alsi, double mu_i, orders order) const
Definition: StandardModel.cpp:805
QCD.h
StandardModel::myCKM
CKM myCKM
An object of type CKM.
Definition: StandardModel.h:2541
StandardModel::getPhiBd
virtual double getPhiBd() const
Half the relative phase of the $B_d$ mixing amplitude w.r.t. the Standard Model one.
Definition: StandardModel.h:2479
StandardModel::resumKappaZ
double resumKappaZ(const double DeltaRho[orders_EW_size], const double deltaKappa_rem[orders_EW_size], const double DeltaRbar_rem, const bool bool_Zbb) const
A method to compute the real part of the effetvive coupling from , and .
Definition: StandardModel.cpp:2027
StandardModel::getIntegrand_dsigmaBox_up205
double getIntegrand_dsigmaBox_up205(double x) const
Definition: StandardModel.cpp:7211
StandardModel::getIntegrand_dsigmaBox_mu202
double getIntegrand_dsigmaBox_mu202(double x) const
Definition: StandardModel.cpp:7028
StandardModel::getIntegrand_AFBnumeratorWithISR_tau172
double getIntegrand_AFBnumeratorWithISR_tau172(double x) const
Definition: StandardModel.cpp:7691
StandardModel::epsilon3
virtual double epsilon3() const
The SM contribution to the epsilon parameter .
Definition: StandardModel.cpp:1805
PMNS::getPMNS
gslpp::matrix< gslpp::complex > getPMNS() const
A member for returning the PMNS matrix.
Definition: PMNS.h:42
StandardModel::Beta_s
double Beta_s(int nm, unsigned int nf) const
QCD beta function coefficients including QED corrections - eq. (36) hep-ph/0512066.
Definition: StandardModel.cpp:554
StandardModel::getIntegrand_AFBnumeratorWithISR_tau192
double getIntegrand_AFBnumeratorWithISR_tau192(double x) const
Definition: StandardModel.cpp:7709
StandardModel::LEP2RCs
LEP2RCs
Definition: StandardModel.h:482
StandardModel::getIntegrand_dsigmaBox_tau196
double getIntegrand_dsigmaBox_tau196(double x) const
Definition: StandardModel.cpp:7092
StandardModel::myThreeLoopEW
EWSMThreeLoopEW * myThreeLoopEW
A pointer to an object of type EWSMThreeLoopEW.
Definition: StandardModel.h:3123
StandardModel::getIntegrand_sigmaWithISR_strange205
double getIntegrand_sigmaWithISR_strange205(double x) const
Definition: StandardModel.cpp:6864
StandardModel::InitializeModel
virtual bool InitializeModel()
A method to initialize the model.
Definition: StandardModel.cpp:163
orders
orders
An enum type for orders in QCD.
Definition: OrderScheme.h:31
Matching< StandardModelMatching, StandardModel >
StandardModel::getIntegrand_dsigmaBox_up207
double getIntegrand_dsigmaBox_up207(double x) const
Definition: StandardModel.cpp:7217
StandardModelMatching.h
StandardModel::myEWSMcache
EWSMcache * myEWSMcache
A pointer to an object of type EWSMcache.
Definition: StandardModel.h:3117
StandardModel::getIntegrand_dsigmaBox_down161
double getIntegrand_dsigmaBox_down161(double x) const
Definition: StandardModel.cpp:7244
StandardModel::s12
double s12
Definition: StandardModel.h:2577
StandardModel::getIntegrand_AFBnumeratorWithISR_tau136
double getIntegrand_AFBnumeratorWithISR_tau136(double x) const
Definition: StandardModel.cpp:7679
StandardModel::getIntegrand_sigmaWithISR_up136
double getIntegrand_sigmaWithISR_up136(double x) const
Definition: StandardModel.cpp:6540
StandardModel::getIntegrand_dsigmaBox_strange200
double getIntegrand_dsigmaBox_strange200(double x) const
Definition: StandardModel.cpp:7464
StandardModel::resumMw
double resumMw(const double Mw_i, const double DeltaRho[orders_EW_size], const double DeltaR_rem[orders_EW_size]) const
A method to compute the -boson mass from and .
Definition: StandardModel.cpp:1861
StandardModel::computeGammaHZga_WW
double computeGammaHZga_WW() const
The loop contribution to in the Standard Model. Currently it returns the value of tab 41 in ref....
Definition: StandardModel.h:2380
Flavour
The parent class in Flavour for calculating all the Wilson coefficients for various Flavor Violating ...
Definition: Flavour.h:33
StandardModel::getIntegrand_sigmaWithISR_up202
double getIntegrand_sigmaWithISR_up202(double x) const
Definition: StandardModel.cpp:6594
StandardModel::computeYukawas
virtual void computeYukawas()
The method to compute the Yukawa matrices.
Definition: StandardModel.cpp:371
StandardModel::computeSigmaZH
double computeSigmaZH(const double sqrt_s) const
The ZH production cross section in the Standard Model.
Definition: StandardModel.h:2135
StandardModel::getMz
double getMz() const
A get method to access the mass of the boson .
Definition: StandardModel.h:721
StandardModel::getIntegrand_sigmaWithISR_mu130
double getIntegrand_sigmaWithISR_mu130(double x) const
Definition: StandardModel.cpp:6359
StandardModel::AFB
virtual double AFB(const Particle f) const
Definition: StandardModel.cpp:1216
StandardModel::computeBrHtocc
double computeBrHtocc() const
The Br in the Standard Model.
Definition: StandardModel.h:2290
StandardModel::getIntegrand_dsigmaBox_down167
double getIntegrand_dsigmaBox_down167(double x) const
Definition: StandardModel.cpp:7250
StandardModel::computeCKM
virtual void computeCKM()
The method to compute the CKM matrix.
Definition: StandardModel.cpp:351
StandardModel::getIntegrand_sigmaWithISR_strange172
double getIntegrand_sigmaWithISR_strange172(double x) const
Definition: StandardModel.cpp:6822
StandardModel::getIntegrand_dsigmaBox_mu189
double getIntegrand_dsigmaBox_mu189(double x) const
Definition: StandardModel.cpp:7004
StandardModel::iterationNo
int iterationNo
Definition: StandardModel.h:3170
StandardModel::getIntegrand_dsigmaBox_mu161
double getIntegrand_dsigmaBox_mu161(double x) const
Definition: StandardModel.cpp:6986
StandardModel::getIntegrand_sigmaWithISR_down200
double getIntegrand_sigmaWithISR_down200(double x) const
Definition: StandardModel.cpp:6676
StandardModel::Update
virtual bool Update(const std::map< std::string, double > &DPars)
The update method for StandardModel.
Definition: StandardModel.cpp:209
StandardModel::Weak
Definition: StandardModel.h:483
StandardModel::RVq
double RVq(const QCD::quark q) const
The radiator factor associated with the final-state QED and QCD corrections to the the vector-current...
Definition: StandardModel.cpp:2142
EWSMThreeLoopEW
A class for three-loop corrections to the EW precision observables.
Definition: EWSMThreeLoopEW.h:35
StandardModel::getIntegrand_sigmaWithISR_bottom196
double getIntegrand_sigmaWithISR_bottom196(double x) const
Definition: StandardModel.cpp:6934
StandardModel::alpha31
double alpha31
Definition: StandardModel.h:2577
StandardModel::gV_f
virtual gslpp::complex gV_f(const Particle f) const
The effective leptonic neutral-current vector coupling in the SM.
Definition: StandardModel.cpp:1568
StandardModel::getIntegrand_sigmaWithISR_up205
double getIntegrand_sigmaWithISR_up205(double x) const
Definition: StandardModel.cpp:6600
StandardModel::Mw_tree
virtual double Mw_tree() const
The tree-level mass of the boson, .
Definition: StandardModel.cpp:951
StandardModel::s02
double s02() const
The square of the sine of the weak mixing angle defined without weak radiative corrections.
Definition: StandardModel.cpp:956
StandardModel::Yn
gslpp::matrix< gslpp::complex > Yn
The Yukawa matrix of the neutrinos.
Definition: StandardModel.h:2547
StandardModel::computeBrHtomumu
double computeBrHtomumu() const
The Br in the Standard Model.
Definition: StandardModel.h:2267
StandardModel::requireCKM
bool requireCKM
An internal flag to control whether the CKM matrix has to be recomputed.
Definition: StandardModel.h:2819
StandardModel::getIntegrand_AFBnumeratorWithISR_bottom167
double getIntegrand_AFBnumeratorWithISR_bottom167(double x) const
Definition: StandardModel.cpp:7837
StandardModel::computeSigmaVBF
double computeSigmaVBF(const double sqrt_s) const
The VBF cross section in the Standard Model.
Definition: StandardModel.h:2003
StandardModel::getIntegrand_dsigmaBox_strange172
double getIntegrand_dsigmaBox_strange172(double x) const
Definition: StandardModel.cpp:7434
StandardModel::getIntegrand_AFBnumeratorWithISR_tau161
double getIntegrand_AFBnumeratorWithISR_tau161(double x) const
Definition: StandardModel.cpp:7685
StandardModel::getIntegrand_dsigmaBox_mu196
double getIntegrand_dsigmaBox_mu196(double x) const
Definition: StandardModel.cpp:7016
StandardModel::getDelR0b
double getDelR0b() const
A get method to retrieve the theoretical uncertainty in , denoted as .
Definition: StandardModel.h:861
StandardModel::useKappaZ_f_cache
bool useKappaZ_f_cache[12]
Definition: StandardModel.h:3155
StandardModel::getFlagKappaZ
std::string getFlagKappaZ() const
A method to retrieve the model flag KappaZ.
Definition: StandardModel.h:682
StandardModel::getIntegrand_dsigmaBox_mu130
double getIntegrand_dsigmaBox_mu130(double x) const
Definition: StandardModel.cpp:6974
StandardModel::getIntegrand_dsigmaBox_bottom202
double getIntegrand_dsigmaBox_bottom202(double x) const
Definition: StandardModel.cpp:7558
Mw
An observable class for the -boson mass.
Definition: Mw.h:22
StandardModel::R_inv
virtual double R_inv() const
The ratio of the invisible and leptonic (electron) decay widths of the boson, .
Definition: StandardModel.cpp:1541
StandardModel::computeBrHtoWW
double computeBrHtoWW() const
The Br in the Standard Model.
Definition: StandardModel.h:2210
StandardModel::getMatching
virtual StandardModelMatching & getMatching() const
A get method to access the member reference of type StandardModelMatching.
Definition: StandardModel.h:952
StandardModel::computeBrHtoss
double computeBrHtoss() const
The Br in the Standard Model.
Definition: StandardModel.h:2302
StandardModel::getDelR0c
double getDelR0c() const
A get method to retrieve the theoretical uncertainty in , denoted as .
Definition: StandardModel.h:851
StandardModel::computeSigmaggH_tt
double computeSigmaggH_tt(const double sqrt_s) const
The square of the top-quark contribution to the ggH cross section in the Standard Model.
Definition: StandardModel.h:1924
StandardModel::LEP2sigmaHadron
virtual double LEP2sigmaHadron(const double s) const
Definition: StandardModel.cpp:3346
QCD::Mrun
double Mrun(const double mu, const double m, const orders order=FULLNNLO) const
Computes a running quark mass from .
Definition: QCD.cpp:1064
StandardModel::SMFlavour
Flavour SMFlavour
An object of type Flavour.
Definition: StandardModel.h:2823
StandardModel::getIntegrand_AFBnumeratorWithISR_charm200
double getIntegrand_AFBnumeratorWithISR_charm200(double x) const
Definition: StandardModel.cpp:7804
QCD::STRANGE
Definition: QCD.h:327
StandardModel::sin2thetaEff
virtual double sin2thetaEff(const Particle f) const
The effective weak mixing angle for at the the -mass scale.
Definition: StandardModel.cpp:1221
StandardModel::AlsWithInit
double AlsWithInit(double mu, double alsi, double mu_i, orders order, bool qed_flag) const
Definition: StandardModel.cpp:689
StandardModel::getIntegrand_sigmaWithISR_bottom133
double getIntegrand_sigmaWithISR_bottom133(double x) const
Definition: StandardModel.cpp:6886
StandardModel::getIntegrand_AFBnumeratorWithISR_bottom133
double getIntegrand_AFBnumeratorWithISR_bottom133(double x) const
Definition: StandardModel.cpp:7831
StandardModel::bSigmaForR
bool bSigmaForR
Definition: StandardModel.h:2833
StandardModel::getIntegrand_AFBnumeratorWithISR_charm133
double getIntegrand_AFBnumeratorWithISR_charm133(double x) const
Definition: StandardModel.cpp:7762
StandardModel::getIntegrand_dsigmaBox_charm202
double getIntegrand_dsigmaBox_charm202(double x) const
Definition: StandardModel.cpp:7380
StandardModel::getIntegrand_sigmaWithISR_bottom207
double getIntegrand_sigmaWithISR_bottom207(double x) const
Definition: StandardModel.cpp:6958
StandardModel::als_cache
double als_cache[11][CacheSize]
Cache for .
Definition: StandardModel.h:3175
StandardModel::getIntegrand_sigmaWithISR_bottom136
double getIntegrand_sigmaWithISR_bottom136(double x) const
Definition: StandardModel.cpp:6892
StandardModel::getIntegrand_sigmaWithISR_tau130
double getIntegrand_sigmaWithISR_tau130(double x) const
Definition: StandardModel.cpp:6432
StandardModel::getIntegrand_AFBnumeratorWithISR_tau196
double getIntegrand_AFBnumeratorWithISR_tau196(double x) const
Definition: StandardModel.cpp:7715
StandardModel::getFlagMw
std::string getFlagMw() const
A method to retrieve the model flag Mw.
Definition: StandardModel.h:662
StandardModel::getIntegrand_sigmaWithISR_down196
double getIntegrand_sigmaWithISR_down196(double x) const
Definition: StandardModel.cpp:6670
StandardModel::leptons
Particle leptons[6]
An array of Particle objects for the leptons.
Definition: StandardModel.h:2540
StandardModel::getIntegrand_sigmaWithISR_charm205
double getIntegrand_sigmaWithISR_charm205(double x) const
Definition: StandardModel.cpp:6776
StandardModel::getIntegrand_dsigmaBox_strange183
double getIntegrand_dsigmaBox_strange183(double x) const
Definition: StandardModel.cpp:7440
StandardModel::flagLEP2
bool flagLEP2[NUMofLEP2RCs]
Definition: StandardModel.h:2831
StandardModel::getIntegrand_sigmaWithISR_down207
double getIntegrand_sigmaWithISR_down207(double x) const
Definition: StandardModel.cpp:6694
StandardModel::getIntegrand_dsigmaBox_strange192
double getIntegrand_dsigmaBox_strange192(double x) const
Definition: StandardModel.cpp:7452
StandardModel::getMyTwoFermionsLEP2
EWSMTwoFermionsLEP2 * getMyTwoFermionsLEP2() const
A get method to retrieve the member pointer of type EWSMTwoFermionsLEP2.
Definition: StandardModel.h:989
StandardModel::getIntegrand_AFBnumeratorWithISR_charm167
double getIntegrand_AFBnumeratorWithISR_charm167(double x) const
Definition: StandardModel.cpp:7768
StandardModel::GammaZ
virtual double GammaZ(const Particle f) const
The partial decay width, .
Definition: StandardModel.cpp:1227
StandardModel::getIntegrand_sigmaWithISR_bottom161
double getIntegrand_sigmaWithISR_bottom161(double x) const
Definition: StandardModel.cpp:6898
StandardModel::getIntegrand_sigmaWithISR_bottom189
double getIntegrand_sigmaWithISR_bottom189(double x) const
Definition: StandardModel.cpp:6922
StandardModel::etab
double etab
The CKM parameter in the Wolfenstein parameterization.
Definition: StandardModel.h:2571
StandardModel::getCCC3
virtual double getCCC3() const
A virtual implementation for the RealWeakEFTCC class.
Definition: StandardModel.h:1168
CKM
A class for the CKM matrix elements.
Definition: CKM.h:23
StandardModel::getIntegrand_sigmaWithISR_bottom130
double getIntegrand_sigmaWithISR_bottom130(double x) const
Definition: StandardModel.cpp:6880
StandardModel::LEP2AFBmu
virtual double LEP2AFBmu(const double s) const
Definition: StandardModel.cpp:5294
StandardModel::requireYe
bool requireYe
An internal flag to control whether the charged-lepton Yukawa matrix has to be recomputed.
Definition: StandardModel.h:2820
StandardModel::computeBrHtogg
double computeBrHtogg() const
The Br in the Standard Model.
Definition: StandardModel.h:2199
Model::name
std::string name
The name of the model.
Definition: Model.h:275
StandardModel::Mz
double Mz
The mass of the boson in GeV.
Definition: StandardModel.h:2554
StandardModel::getIntegrand_AFBnumeratorWithISR_bottom205
double getIntegrand_AFBnumeratorWithISR_bottom205(double x) const
Definition: StandardModel.cpp:7885
StandardModel::myPMNS
PMNS myPMNS
Definition: StandardModel.h:2542
StandardModel::kappaZ_f_cache
gslpp::complex kappaZ_f_cache[12]
A cache of the value of .
Definition: StandardModel.h:3149
StandardModel::NUMofLEP2RCs
Definition: StandardModel.h:488
StandardModel::getIntegrand_dsigmaBox_up183
double getIntegrand_dsigmaBox_up183(double x) const
Definition: StandardModel.cpp:7175
StandardModel::N_nu
virtual double N_nu() const
The number of neutrinos obtained indirectly from the measurements at the Z pole, .
Definition: StandardModel.cpp:1547
StandardModel::getIntegrand_dsigmaBox_down205
double getIntegrand_dsigmaBox_down205(double x) const
Definition: StandardModel.cpp:7298
StandardModel::AFB_NoISR_l
double AFB_NoISR_l(const QCD::lepton l_flavor, const double s) const
Definition: StandardModel.cpp:6325
StandardModel::computeBrHtoZZinv
double computeBrHtoZZinv() const
The Br in the Standard Model.
Definition: StandardModel.h:2233
StandardModel::getIntegrand_sigmaWithISR_charm136
double getIntegrand_sigmaWithISR_charm136(double x) const
Definition: StandardModel.cpp:6716
StandardModel::EW2QCD1
Three-loop of .
Definition: StandardModel.h:500
StandardModel::getIntegrand_sigmaWithISR_strange196
double getIntegrand_sigmaWithISR_strange196(double x) const
Definition: StandardModel.cpp:6846
StandardModel::getIntegrand_AFBnumeratorWithISR_bottom200
double getIntegrand_AFBnumeratorWithISR_bottom200(double x) const
Definition: StandardModel.cpp:7873
StandardModel::computeSigmaZWF
double computeSigmaZWF(const double sqrt_s) const
The Z W interference fusion contribution to higgs-production cross section in the Standard Model.
Definition: StandardModel.h:2086
StandardModel::FlagCacheInStandardModel
bool FlagCacheInStandardModel
A flag for caching (true by default).
Definition: StandardModel.h:3142
StandardModel::ComputeDeltaRho
void ComputeDeltaRho(const double Mw_i, double DeltaRho[orders_EW_size]) const
A method to collect computed via subclasses.
Definition: StandardModel.cpp:1050
StandardModel::getMHl
virtual double getMHl() const
A get method to retrieve the Higgs mass .
Definition: StandardModel.h:768
StandardModel::getIntegrand_AFBnumeratorWithISR_charm196
double getIntegrand_AFBnumeratorWithISR_charm196(double x) const
Definition: StandardModel.cpp:7798
StandardModel::getIntegrand_sigmaWithISR_up207
double getIntegrand_sigmaWithISR_up207(double x) const
Definition: StandardModel.cpp:6606
StandardModel::Mw
virtual double Mw() const
The SM prediction for the -boson mass in the on-shell scheme, .
Definition: StandardModel.cpp:970
StandardModel::myApproximateFormulae
EWSMApproximateFormulae * myApproximateFormulae
A pointer to an object of type EWSMApproximateFormulae.
Definition: StandardModel.h:3124
StandardModel::Mzbar
double Mzbar() const
The -boson mass in the complex-pole/fixed-width scheme.
Definition: StandardModel.cpp:1111
StandardModel::getIntegrand_sigmaWithISR_down130
double getIntegrand_sigmaWithISR_down130(double x) const
Definition: StandardModel.cpp:6615
StandardModel::getIntegrand_dsigmaBox_strange205
double getIntegrand_dsigmaBox_strange205(double x) const
Definition: StandardModel.cpp:7476
StandardModel::getIntegrand_sigmaWithISR_mu189
double getIntegrand_sigmaWithISR_mu189(double x) const
Definition: StandardModel.cpp:6389
StandardModel::Integrand_dsigmaBox_q
double Integrand_dsigmaBox_q(double cosTheta, const QCD::quark q_flavor, const double s) const
Definition: StandardModel.cpp:7127
StandardModel::LEP2sigmaBottom
virtual double LEP2sigmaBottom(const double s) const
Definition: StandardModel.cpp:3135
StandardModel::gA_f
virtual gslpp::complex gA_f(const Particle f) const
The effective leptonic neutral-current axial-vector coupling in the SM.
Definition: StandardModel.cpp:1574
StandardModel::requireYn
bool requireYn
An internal flag to control whether the neutrino Yukawa matrix has to be recomputed.
Definition: StandardModel.h:2821
StandardModel::LEP2sigmaMu
virtual double LEP2sigmaMu(const double s) const
Definition: StandardModel.cpp:2425
StandardModel::delR0b
double delR0b
The theoretical uncertainty in , denoted as .
Definition: StandardModel.h:2567
StandardModel::getIntegrand_dsigmaBox_mu207
double getIntegrand_dsigmaBox_mu207(double x) const
Definition: StandardModel.cpp:7040
StandardModel::sigma_NoISR_q
double sigma_NoISR_q(const QCD::quark q_flavor, const double s) const
Definition: StandardModel.cpp:6310
StandardModel::LEP2Rcharm
virtual double LEP2Rcharm(const double s) const
Definition: StandardModel.cpp:6286
StandardModel::getIntegrand_AFBnumeratorWithISR_charm207
double getIntegrand_AFBnumeratorWithISR_charm207(double x) const
Definition: StandardModel.cpp:7822
StandardModel::Vus
double Vus
used as an input for FlagWolfenstein = FALSE
Definition: StandardModel.h:2572
StandardModel::GeVminus2_to_nb
static const double GeVminus2_to_nb
Definition: StandardModel.h:514
StandardModel::myThreeLoopEW2QCD
EWSMThreeLoopEW2QCD * myThreeLoopEW2QCD
A pointer to an object of type EWSMThreeLoopEW2QCD.
Definition: StandardModel.h:3122
StandardModel::useGammaW_cache
bool useGammaW_cache
Definition: StandardModel.h:3153
StandardModel::delSin2th_l
double delSin2th_l
The theoretical uncertainty in , denoted as .
Definition: StandardModel.h:2560
StandardModel::getIntegrand_dsigmaBox_charm161
double getIntegrand_dsigmaBox_charm161(double x) const
Definition: StandardModel.cpp:7332
StandardModel::getIntegrand_dsigmaBox_mu133
double getIntegrand_dsigmaBox_mu133(double x) const
StandardModel::getDelSin2th_b
double getDelSin2th_b() const
A get method to retrieve the theoretical uncertainty in , denoted as .
Definition: StandardModel.h:811
StandardModel::getIntegrand_AFBnumeratorWithISR_mu183
double getIntegrand_AFBnumeratorWithISR_mu183(double x) const
Definition: StandardModel.cpp:7624
StandardModel::getDelGammaZ
double getDelGammaZ() const
A get method to retrieve the theoretical uncertainty in , denoted as .
Definition: StandardModel.h:821
StandardModel::QEDFSR
Definition: StandardModel.h:486
StandardModel::getIntegrand_dsigmaBox_down196
double getIntegrand_dsigmaBox_down196(double x) const
Definition: StandardModel.cpp:7280
StandardModel::SMparamsForEWPO_cache
double SMparamsForEWPO_cache[NumSMParamsForEWPO]
Definition: StandardModel.h:3143
StandardModel::getIntegrand_dsigmaBox_charm130
double getIntegrand_dsigmaBox_charm130(double x) const
Definition: StandardModel.cpp:7314
StandardModel::getIntegrand_sigmaWithISR_tau192
double getIntegrand_sigmaWithISR_tau192(double x) const
Definition: StandardModel.cpp:6468
StandardModel::getIntegrand_AFBnumeratorWithISR_bottom202
double getIntegrand_AFBnumeratorWithISR_bottom202(double x) const
Definition: StandardModel.cpp:7879
StandardModel::Gamma_Z
virtual double Gamma_Z() const
The total decay width of the boson, .
Definition: StandardModel.cpp:1344
StandardModel::getCBd
virtual double getCBd() const
The ratio of the absolute value of the $B_d$ mixing amplitude over the Standard Model value.
Definition: StandardModel.h:2434
StandardModel::deltaRhoZ_f
virtual gslpp::complex deltaRhoZ_f(const Particle f) const
Flavour non-universal vertex corrections to , denoted by .
Definition: StandardModel.cpp:1730
PMNS
A class for the PMNS matrix elements.
Definition: PMNS.h:23
StandardModel::getIntegrand_sigmaWithISR_charm167
double getIntegrand_sigmaWithISR_charm167(double x) const
Definition: StandardModel.cpp:6728
StandardModel::delsigma0H
double delsigma0H
The theoretical uncertainty in , denoted as in nb.
Definition: StandardModel.h:2564
StandardModel::setFlagCacheInStandardModel
void setFlagCacheInStandardModel(bool FlagCacheInStandardModel)
A set method to change the model flag CacheInStandardModel of StandardModel.
Definition: StandardModel.h:698
QCD::DOWN
Definition: QCD.h:325
StandardModel::getIntegrand_dsigmaBox_tau192
double getIntegrand_dsigmaBox_tau192(double x) const
Definition: StandardModel.cpp:7086
StandardModel::CheckFlags
virtual bool CheckFlags() const
A method to check the sanity of the set of model flags.
Definition: StandardModel.cpp:475
StandardModel::kappaZ_f
virtual gslpp::complex kappaZ_f(const Particle f) const
The effective leptonic neutral-current coupling in the SM.
Definition: StandardModel.cpp:1644
StandardModel::getIntegrand_sigmaWithISR_strange207
double getIntegrand_sigmaWithISR_strange207(double x) const
Definition: StandardModel.cpp:6870
StandardModel::getIntegrand_dsigmaBox_bottom161
double getIntegrand_dsigmaBox_bottom161(double x) const
Definition: StandardModel.cpp:7510
StandardModel::getIntegrand_dsigmaBox_down192
double getIntegrand_dsigmaBox_down192(double x) const
Definition: StandardModel.cpp:7274
StandardModel::muw
double muw
A matching scale around the weak scale in GeV.
Definition: StandardModel.h:2576
StandardModel::getIntegrand_AFBnumeratorWithISR_tau202
double getIntegrand_AFBnumeratorWithISR_tau202(double x) const
Definition: StandardModel.cpp:7727
StandardModel::getIntegrand_sigmaWithISR_mu136
double getIntegrand_sigmaWithISR_mu136(double x) const
Definition: StandardModel.cpp:6365
StandardModel::getIntegrand_dsigmaBox_down202
double getIntegrand_dsigmaBox_down202(double x) const
Definition: StandardModel.cpp:7292
StandardModel::computeGammaHgg_tt
double computeGammaHgg_tt() const
The top loop contribution to in the Standard Model.
Definition: StandardModel.h:2336
StandardModel::Yd
gslpp::matrix< gslpp::complex > Yd
The Yukawa matrix of the down-type quarks.
Definition: StandardModel.h:2546
StandardModel::EW1QCD2
Three-loop of .
Definition: StandardModel.h:498
StandardModel::getIntegrand_sigmaWithISR_down202
double getIntegrand_sigmaWithISR_down202(double x) const
Definition: StandardModel.cpp:6682
FULLNLO
Definition: OrderScheme.h:37
StandardModel::error
double error
Definition: StandardModel.h:3165
StandardModel::getAle
double getAle() const
A get method to retrieve the fine-structure constant .
Definition: StandardModel.h:748
StandardModel::FlagWolfenstein
bool FlagWolfenstein
A boolean for the model flag Wolfenstein.
Definition: StandardModel.h:3135
StandardModel::getIntegrand_sigmaWithISR_up133
double getIntegrand_sigmaWithISR_up133(double x) const
Definition: StandardModel.cpp:6534
StandardModel::getCDMK
virtual double getCDMK() const
The ratio of the real part of the $K$ mixing amplitude over the Standard Model value.
Definition: StandardModel.h:2452
StandardModel::getYd
gslpp::matrix< gslpp::complex > getYd() const
A get method to retrieve the Yukawa matrix of the down-type quarks, .
Definition: StandardModel.h:908
StandardModel::checkSMparamsForEWPO
bool checkSMparamsForEWPO()
A method to check whether the parameters relevant to the EWPO are updated.
Definition: StandardModel.cpp:484
StandardModel::WeakBox
Definition: StandardModel.h:484
StandardModel::getIntegrand_dsigmaBox_tau200
double getIntegrand_dsigmaBox_tau200(double x) const
Definition: StandardModel.cpp:7098
StandardModel::getIntegrand_sigmaWithISR_tau136
double getIntegrand_sigmaWithISR_tau136(double x) const
Definition: StandardModel.cpp:6438
StandardModel::orders_EW_size
The size of this enum.
Definition: StandardModel.h:502
StandardModel::myTwoLoopQCD
EWSMTwoLoopQCD * myTwoLoopQCD
A pointer to an object of type EWSMTwoLoopQCD.
Definition: StandardModel.h:3119
StandardModel::getCKM
CKM getCKM() const
A get method to retrieve the member object of type CKM.
Definition: StandardModel.h:879
StandardModel::Vcb
double Vcb
used as an input for FlagWolfenstein = FALSE
Definition: StandardModel.h:2573
StandardModel::getMuw
double getMuw() const
A get method to retrieve the matching scale around the weak scale.
Definition: StandardModel.h:938
StandardModel::getIntegrand_sigmaWithISR_down205
double getIntegrand_sigmaWithISR_down205(double x) const
Definition: StandardModel.cpp:6688
StandardModel::getIntegrand_AFBnumeratorWithISR_bottom183
double getIntegrand_AFBnumeratorWithISR_bottom183(double x) const
Definition: StandardModel.cpp:7849
StandardModel::delSin2th_b
double delSin2th_b
The theoretical uncertainty in , denoted as .
Definition: StandardModel.h:2562
StandardModel::getMyThreeLoopEW2QCD
EWSMThreeLoopEW2QCD * getMyThreeLoopEW2QCD() const
Definition: StandardModel.h:1000
StandardModel::getIntegrand_sigmaWithISR_mu161
double getIntegrand_sigmaWithISR_mu161(double x) const
Definition: StandardModel.cpp:6371
StandardModel::getIntegrand_dsigmaBox_tau207
double getIntegrand_dsigmaBox_tau207(double x) const
Definition: StandardModel.cpp:7116
StandardModel::getIntegrand_sigmaWithISR_strange167
double getIntegrand_sigmaWithISR_strange167(double x) const
Definition: StandardModel.cpp:6816
QCD::lepton
lepton
An enum type for leptons.
Definition: QCD.h:310
StandardModel::computeSigmaggH_tb
double computeSigmaggH_tb(const double sqrt_s) const
The top-bottom interference contribution to the ggH cross section in the Standard Model.
Definition: StandardModel.h:1974
StandardModel::delta
double delta
Definition: StandardModel.h:2577
StandardModel::getIntegrand_AFBnumeratorWithISR_tau205
double getIntegrand_AFBnumeratorWithISR_tau205(double x) const
Definition: StandardModel.cpp:7733
StandardModel::getIntegrand_sigmaWithISR_up161
double getIntegrand_sigmaWithISR_up161(double x) const
Definition: StandardModel.cpp:6546
StandardModel::alpha21
double alpha21
Definition: StandardModel.h:2577
StandardModel::FlagNoApproximateGammaZ
bool FlagNoApproximateGammaZ
A boolean for the model flag NoApproximateGammaZ.
Definition: StandardModel.h:3131
StandardModel::getIntegrand_dsigmaBox_bottom172
double getIntegrand_dsigmaBox_bottom172(double x) const
Definition: StandardModel.cpp:7522
StandardModel::getLeptons
Particle getLeptons(const QCD::lepton p) const
A get method to retrieve the member object of a lepton.
Definition: StandardModel.h:712
StandardModel::orders_EW
orders_EW
An enumerated type representing perturbative orders of radiative corrections to EW precision observab...
Definition: StandardModel.h:495