a Code for the Combination of Indirect and Direct Constraints on High Energy Physics Models Logo
FlavourWilsonCoefficient_DF2 Class Reference

Model for NP contributions to \(\Delta F=2\) using modification to the Wilson coefficients. More...

#include <FlavourWilsonCoefficient_DF2.h>

+ Inheritance diagram for FlavourWilsonCoefficient_DF2:

Detailed Description

Model for NP contributions to \(\Delta F=2\) using modification to the Wilson coefficients.

Author
HEPfit Collaboration

Model parameters

The model parameters of FlavourWilsonCoefficient_DF2 model are summarized below:

Label LaTeX symbol Description
reC1s, reC2s, reC3s, reC4s, reC5s \( \mathcal{R}(C_{1,s}) \), \( \mathcal{R}(C_{2,s}) \), \( \mathcal{R}(C_{3,s}) \), \( \mathcal{R}(C_{4,s}) \), \( \mathcal{R}(C_{5,s}) \) The real parts of the rescaling of the Wilson coefficient \( C_{i,s}, (i,1,\ldots,5) \) that modified the SM prediction to \(K^0-\bar{K}^0\) mixing
WCscale_s \( \mu_s \) The scale of the Wilson coefficient \( C_{i,s}, (i,1,\ldots,5) \) that modified the SM prediction to \(K^0-\bar{K}^0\) mixing
reC1c, reC2c, reC3c, reC4c, reC5s \( \mathcal{R}(C_{1,c}) \), \( \mathcal{R}(C_{2,c}) \), \( \mathcal{R}(C_{3,c}) \), \( \mathcal{R}(C_{4,c}) \), \( \mathcal{R}(C_{5,c}) \) The real parts of the rescaling of the Wilson coefficient \( C_{i,c}, (i,1,\ldots,5) \) that modified the SM prediction to \(D^0-\bar{D}^0\) mixing
WCscale_c \( \mu_c \) The scale of the Wilson coefficient \( C_{i,c}, (i,1,\ldots,5) \) that modified the SM prediction to \(D^0-\bar{D}^0\) mixing
reC1bd, reC2bd, reC3bd, reC4bd, reC5bd \( \mathcal{R}(C_{1,bd}) \), \( \mathcal{R}(C_{2,bd}) \), \( \mathcal{R}(C_{3,bd}) \), \( \mathcal{R}(C_{4,bd}) \), \( \mathcal{R}(C_{5,bd}) \) The real parts of the rescaling of the Wilson coefficient \( C_{i,bd}, (i,1,\ldots,5) \) that modified the SM prediction to \(B^0-\bar{B}^0\) mixing
WCscale_bd \( \mu_{bd} \) The scale of the Wilson coefficient \( C_{i,bd}, (i,1,\ldots,5) \) that modified the SM prediction to \(B^0-\bar{B}^0\) mixing
reC1bs, reC2bs, reC3bs, reC4bs, reC5bs \( \mathcal{R}(C_{1,bs}) \), \( \mathcal{R}(C_{2,bs}) \), \( \mathcal{R}(C_{3,bs}) \), \( \mathcal{R}(C_{4,bs}) \), \( \mathcal{R}(C_{5,bs}) \) The real parts of the rescaling of the Wilson coefficient \( C_{i,bs}, (i,1,\ldots,5) \) that modified the SM prediction to \(B_s-\bar{B}_s\) mixing
WCscale_bs \( \mu_{bs} \) The scale of the Wilson coefficient \( C_{i,bs}, (i,1,\ldots,5) \) that modified the SM prediction to \(B_s-\bar{B}_s\) mixing

Definition at line 100 of file FlavourWilsonCoefficient_DF2.h.

Public Member Functions

virtual bool CheckParameters (const std::map< std::string, double > &DPars)
 A method to check if all the mandatory parameters for FlavourWilsonCoefficient_DF2 have been provided in model initialization. More...
 
 FlavourWilsonCoefficient_DF2 ()
 FlavourWilsonCoefficient constructor. More...
 
gslpp::vector< gslpp::complex > getC_bd () const
 
gslpp::vector< gslpp::complex > getC_bs () const
 
gslpp::vector< gslpp::complex > getC_c () const
 
gslpp::vector< gslpp::complex > getC_s () const
 
virtual FlavourWilsonCoefficient_DF2MatchinggetMatching () const
 A get method to access the member reference of type FlavourWilsonCoefficient_DF2Matching. More...
 
double getWCscale_bd () const
 
double getWCscale_bs () const
 
double getWCscale_c () const
 
double getWCscale_s () const
 
virtual bool InitializeModel ()
 A method to initialize the model. More...
 
virtual bool PostUpdate ()
 The post-update method for FlavourWilsonCoefficient_DF2. More...
 
- Public Member Functions inherited from StandardModel
virtual const double A_f (const Particle f) const
 The left-right asymmetry in \(e^+e^-\to Z\to \ell \bar{\ell}\) at the \(Z\)-pole, \(\mathcal{A}_\ell\). More...
 
virtual const double AFB (const Particle f) const
 
gslpp::complex AH_f (const double tau) const
 Fermionic loop function entering in the calculation of the effective \(Hgg\) and \(H\gamma\gamma\) couplings. More...
 
gslpp::complex AH_W (const double tau) const
 W loop function entering in the calculation of the effective \(H\gamma\gamma\) coupling. More...
 
gslpp::complex AHZga_f (const double tau, const double lambda) const
 Fermionic loop function entering in the calculation of the effective \(HZ\gamma\) coupling. More...
 
gslpp::complex AHZga_W (const double tau, const double lambda) const
 W loop function entering in the calculation of the effective \(HZ\gamma\) coupling. More...
 
const double Ale (double mu, orders order, bool Nf_thr=true) const
 The running electromagnetic coupling \(\alpha_e(\mu)\) in the \(\overline{MS}\) scheme. More...
 
const double ale_OS (const double mu, orders order=FULLNLO) const
 The running electromagnetic coupling \(\alpha(\mu)\) in the on-shell scheme. More...
 
virtual const double alphaMz () const
 The electromagnetic coupling at the \(Z\)-mass scale, \(\alpha(M_Z^2)=\alpha/(1-\Delta\alpha(M_Z^2))\). More...
 
virtual const double alrmoller (const double q2, const double y) const
 The computation of the parity violating asymmetry in Moller scattering. More...
 
const double Als (const double mu, const int Nf_in, const orders order=FULLNLO) const
 Computes the running strong coupling \(\alpha_s(\mu)\) with \(N_f\) active flavours in the \(\overline{\mathrm{MS}}\) scheme. In the cases of LO, NLO and FULLNLO, the coupling is computed with AlsWithInit(). On the other hand, in the cases of NNLO and FULLNNLO, the coupling is computed with AlsWithLambda(). More...
 
const double Als (const double mu, const orders order, const bool Nf_thr, const bool qed_flag) const
 The running QCD coupling \(\alpha(\mu)\) in the \(\overline{MS}\) scheme including QED corrections. More...
 
const double Als (const double mu, const orders order=FULLNLO, const bool Nf_thr=true) const
 
const double Alstilde5 (const double mu) const
 The value of \(\frac{\alpha_s^{\mathrm{FULLNLO}}}{4\pi}\) at any scale \(\mu\) with the number of flavours \(n_f = 4\) and full EW corrections. More...
 
virtual const double amuon () const
 The computation of the anomalous magnetic moment of the muon \(a_\mu=(g_\mu-2)/2\). More...
 
const double Beta_e (int nm, unsigned int nf) const
 QED beta function coefficients - eq. (36) hep-ph/0512066. More...
 
const double Beta_s (int nm, unsigned int nf) const
 QCD beta function coefficients including QED corrections - eq. (36) hep-ph/0512066. More...
 
virtual const double BrHtobb () const
 The Br \((H\to b \bar{b})\) in the Standard Model. More...
 
virtual const double BrHtocc () const
 The Br \((H\to c \bar{c})\) in the Standard Model. More...
 
virtual const double BrHtogaga () const
 The Br \((H\to \gamma \gamma)\) in the Standard Model. More...
 
virtual const double BrHtogg () const
 The Br \(\(H\to gg)\) in the Standard Model. More...
 
virtual const double BrHtomumu () const
 The Br \((H\to \mu^+ \mu^-)\) in the Standard Model. More...
 
virtual const double BrHtoss () const
 The Br \((H\to s \bar{s})\) in the Standard Model. More...
 
virtual const double BrHtotautau () const
 The Br \((H\to \tau^+ \tau^-)\) in the Standard Model. More...
 
virtual const double BrHtoWWstar () const
 The Br \((H\to W W^*)\) in the Standard Model. More...
 
virtual const double BrHtoZga () const
 The Br \((H\to Z \gamma)\) in the Standard Model. More...
 
virtual const double BrHtoZZstar () const
 The Br \((H\to Z Z^*)\) in the Standard Model. More...
 
virtual const double BrW (const Particle fi, const Particle fj) const
 The branching ratio of the \(W\) boson decaying into a SM fermion pair, \(Br(W\to f_i f_j)\). More...
 
const double c02 () const
 The square of the cosine of the weak mixing angle \(c_0^2\) defined without weak radiative corrections. More...
 
virtual bool CheckFlags () const
 A method to check the sanity of the set of model flags. More...
 
bool checkSMparamsForEWPO ()
 A method to check whether the parameters relevant to the EWPO are updated. More...
 
const double computeBrHto4f () const
 The Br \((H\to 4f)\) in the Standard Model. More...
 
const double computeBrHto4l2 () const
 The Br \((H\to 4l)\) \(l=e,\mu\) in the Standard Model. More...
 
const double computeBrHto4l3 () const
 The Br \((H\to 4l)\) \(l=e,\mu,\tau\) in the Standard Model. More...
 
const double computeBrHto4q () const
 The Br \((H\to 4q)\) in the Standard Model. More...
 
const double computeBrHto4v () const
 The Br \((H\to 4\nu)\) in the Standard Model. More...
 
const double computeBrHtobb () const
 The Br \((H\to bb)\) in the Standard Model. More...
 
const double computeBrHtocc () const
 The Br \((H\to cc)\) in the Standard Model. More...
 
const double computeBrHtoevmuv () const
 The Br \((H\to e \nu \mu \nu)\) in the Standard Model. More...
 
const double computeBrHtogaga () const
 The Br \((H\to\gamma\gamma)\) in the Standard Model. More...
 
const double computeBrHtogg () const
 The Br \((H\to gg)\) in the Standard Model. More...
 
const double computeBrHtollvv2 () const
 The Br \((H\to l^+ l^- \nu \nu)\) \(l=e,\mu\) in the Standard Model. More...
 
const double computeBrHtollvv3 () const
 The Br \((H\to l^+ l^- \nu \nu)\) \(l=e,\mu,\tau\) in the Standard Model. More...
 
const double computeBrHtomumu () const
 The Br \((H\to \mu\mu)\) in the Standard Model. More...
 
const double computeBrHtoss () const
 The Br \((H\to ss)\) in the Standard Model. More...
 
const double computeBrHtotautau () const
 The Br \((H\to \tau\tau)\) in the Standard Model. More...
 
const double computeBrHtoWW () const
 The Br \((H\to WW)\) in the Standard Model. More...
 
const double computeBrHtoZga () const
 The Br \((H\to Z\gamma)\) in the Standard Model. More...
 
const double computeBrHtoZZ () const
 The Br \((H\to ZZ)\) in the Standard Model. More...
 
void ComputeDeltaR_rem (const double Mw_i, double DeltaR_rem[orders_EW_size]) const
 A method to collect \(\Delta r_{\mathrm{rem}}\) computed via subclasses. More...
 
void ComputeDeltaRho (const double Mw_i, double DeltaRho[orders_EW_size]) const
 A method to collect \(\Delta\rho\) computed via subclasses. More...
 
const double computeGammaHgaga_tt () const
 The top loop contribution to \(H\to\gamma\gamma\) in the Standard Model. More...
 
const double computeGammaHgaga_tW () const
 The mixed \(t-W\) loop contribution to \(H\to\gamma\gamma\) in the Standard Model. More...
 
const double computeGammaHgaga_WW () const
 The \(W\) loop contribution to \(H\to\gamma\gamma\) in the Standard Model. More...
 
const double computeGammaHgg_bb () const
 The bottom loop contribution to \(H\to gg\) in the Standard Model. More...
 
const double computeGammaHgg_tb () const
 The top-bottom interference contribution to \(H\to gg\) in the Standard Model. More...
 
const double computeGammaHgg_tt () const
 The top loop contribution to \(H\to gg\) in the Standard Model. More...
 
const double computeGammaHTotal () const
 The Higgs total width in the Standard Model. More...
 
const double computeGammaHZga_tt () const
 The top loop contribution to \(H\to Z\gamma\) in the Standard Model. More...
 
const double computeGammaHZga_tW () const
 The mixed \(t-W\) loop contribution to \(H\to Z\gamma\) in the Standard Model. More...
 
const double computeGammaHZga_WW () const
 The \(W\) loop contribution to \(H\to Z\gamma\) in the Standard Model. Currently it returns the value of tab 41 in ref. [Heinemeyer:2013tqa]. More...
 
const double computeSigmabbH (const double sqrt_s) const
 The bbH production cross section in the Standard Model. More...
 
const double computeSigmaggH (const double sqrt_s) const
 The ggH cross section in the Standard Model. More...
 
const double computeSigmaggH_bb (const double sqrt_s) const
 The square of the bottom-quark contribution to the ggH cross section in the Standard Model. More...
 
const double computeSigmaggH_tb (const double sqrt_s) const
 The top-bottom interference contribution to the ggH cross section in the Standard Model. More...
 
const double computeSigmaggH_tt (const double sqrt_s) const
 The square of the top-quark contribution to the ggH cross section in the Standard Model. More...
 
const double computeSigmatHq (const double sqrt_s) const
 The tHq production cross section in the Standard Model. More...
 
const double computeSigmattH (const double sqrt_s) const
 The ttH production cross section in the Standard Model. More...
 
const double computeSigmaVBF (const double sqrt_s) const
 The VBF cross section in the Standard Model. More...
 
const double computeSigmaWF (const double sqrt_s) const
 The W fusion contribution \(\sigma_{WF}\) to higgs-production cross section in the Standard Model. More...
 
const double computeSigmaWH (const double sqrt_s) const
 The WH production cross section in the Standard Model. More...
 
const double computeSigmaZF (const double sqrt_s) const
 The Z fusion contribution \(\sigma_{ZF}\) to higgs-production cross section in the Standard Model. More...
 
const double computeSigmaZH (const double sqrt_s) const
 The ZH production cross section in the Standard Model. More...
 
const double computeSigmaZWF (const double sqrt_s) const
 The Z W interference fusion contribution \(\sigma_{ZWF}\) to higgs-production cross section in the Standard Model. More...
 
virtual const double cW2 () const
 
virtual const double cW2 (const double Mw_i) const
 The square of the cosine of the weak mixing angle in the on-shell scheme, denoted as \(c_W^2\). More...
 
virtual const double Dalpha5hMz () const
 The 5-quark contribution to the running of the em constant to the \(Z\) pole. \(\Delta\alpha_{had}^{(5)}(M_Z)\). More...
 
const double DeltaAlpha () const
 The total corrections to the electromagnetic coupling \(\alpha\) at the \(Z\)-mass scale, denoted as \(\Delta\alpha(M_Z^2)\). More...
 
const double DeltaAlphaL5q () const
 The sum of the leptonic and the five-flavour hadronic corrections to the electromagnetic coupling \(\alpha\) at the \(Z\)-mass scale, denoted as \(\Delta\alpha^{\ell+5q}(M_Z^2)\). More...
 
const double DeltaAlphaLepton (const double s) const
 Leptonic contribution to the electromagnetic coupling \(\alpha\), denoted as \(\Delta\alpha_{\mathrm{lept}}(s)\). More...
 
const double DeltaAlphaTop (const double s) const
 Top-quark contribution to the electromagnetic coupling \(\alpha\), denoted as \(\Delta\alpha_{\mathrm{top}}(s)\). More...
 
virtual const gslpp::complex deltaKappaZ_f (const Particle f) const
 Flavour non-universal vertex corrections to \(\kappa_Z^l\), denoted by \(\Delta\kappa_Z^l\). More...
 
virtual const double DeltaR () const
 The SM prediction for \(\Delta r\) derived from that for the \(W\) boson mass. More...
 
virtual const double DeltaRbar () const
 The SM prediction for \(\Delta \overline{r}\) derived from that for the \(W\)-boson mass. More...
 
virtual const gslpp::complex deltaRhoZ_f (const Particle f) const
 Flavour non-universal vertex corrections to \(\rho_Z^l\), denoted by \(\Delta\rho_Z^l\). More...
 
virtual const double eeffAFBbottom (const double pol_e, const double pol_p, const double s) const
 
virtual const double eeffAFBcharm (const double pol_e, const double pol_p, const double s) const
 
virtual const double eeffAFBe (const double pol_e, const double pol_p, const double s) const
 
virtual const double eeffAFBetsub (const double pol_e, const double pol_p, const double s) const
 
virtual const double eeffAFBmu (const double pol_e, const double pol_p, const double s) const
 
virtual const double eeffAFBstrange (const double pol_e, const double pol_p, const double s) const
 
virtual const double eeffAFBtau (const double pol_e, const double pol_p, const double s) const
 
virtual const double eeffRbottom (const double pol_e, const double pol_p, const double s) const
 
virtual const double eeffRcharm (const double pol_e, const double pol_p, const double s) const
 
virtual const double eeffRelectron (const double pol_e, const double pol_p, const double s) const
 
virtual const double eeffRelectrontsub (const double pol_e, const double pol_p, const double s) const
 
virtual const double eeffRmuon (const double pol_e, const double pol_p, const double s) const
 
virtual const double eeffRstrange (const double pol_e, const double pol_p, const double s) const
 
virtual const double eeffRtau (const double pol_e, const double pol_p, const double s) const
 
const double eeffsigma (const Particle f, const double pol_e, const double pol_p, const double s, const double cosmin, const double cosmax) const
 
virtual const double eeffsigmaBottom (const double pol_e, const double pol_p, const double s) const
 
virtual const double eeffsigmaCharm (const double pol_e, const double pol_p, const double s) const
 
virtual const double eeffsigmaE (const double pol_e, const double pol_p, const double s) const
 
const double eeffsigmaEbin (const double pol_e, const double pol_p, const double s, const double cosmin, const double cosmax) const
 
virtual const double eeffsigmaEtsub (const double pol_e, const double pol_p, const double s) const
 
virtual const double eeffsigmaHadron (const double pol_e, const double pol_p, const double s) const
 
virtual const double eeffsigmaMu (const double pol_e, const double pol_p, const double s) const
 
virtual const double eeffsigmaStrange (const double pol_e, const double pol_p, const double s) const
 
virtual const double eeffsigmaTau (const double pol_e, const double pol_p, const double s) const
 
virtual const double epsilon1 () const
 The SM contribution to the epsilon parameter \(\varepsilon_1\). More...
 
virtual const double epsilon2 () const
 The SM contribution to the epsilon parameter \(\varepsilon_2\). More...
 
virtual const double epsilon3 () const
 The SM contribution to the epsilon parameter \(\varepsilon_3\). More...
 
virtual const double epsilonb () const
 The SM contribution to the epsilon parameter \(\varepsilon_b\). More...
 
gslpp::complex f_triangle (const double tau) const
 Loop function entering in the calculation of the effective \(Hgg\) and \(H\gamma\gamma\) couplings. More...
 
gslpp::complex g_triangle (const double tau) const
 Loop function entering in the calculation of the effective \(HZ\gamma\) coupling. More...
 
virtual const gslpp::complex gA_f (const Particle f) const
 The effective leptonic neutral-current axial-vector coupling \(g_A^l\) in the SM. More...
 
virtual const double Gamma_had () const
 The hadronic decay width of the \(Z\) boson, \(\Gamma_{h}\). More...
 
virtual const double Gamma_inv () const
 The invisible partial decay width of the \(Z\) boson, \(\Gamma_{\mathrm{inv}}\). More...
 
virtual const double Gamma_muon () const
 The computation of the muon decay. More...
 
virtual const double Gamma_tau_l_nunu (const Particle l) const
 The computation of the leptonic tau decays. More...
 
virtual const double Gamma_Z () const
 The total decay width of the \(Z\) boson, \(\Gamma_Z\). More...
 
virtual const double GammaHtobb () const
 The \(\Gamma(H\to b \bar{b})\) in the Standard Model. More...
 
virtual const double GammaHtocc () const
 The \(\Gamma(H\to c \bar{c})\) in the Standard Model. More...
 
virtual const double GammaHtogaga () const
 The \(\Gamma(H\to \gamma \gamma)\) in the Standard Model. More...
 
virtual const double GammaHtogg () const
 The \(\Gamma(H\to gg)\) in the Standard Model. More...
 
virtual const double GammaHtomumu () const
 The \(\Gamma(H\to \mu^+ \mu^-)\) in the Standard Model. More...
 
virtual const double GammaHtoss () const
 The \(\Gamma(H\to s \bar{s})\) in the Standard Model. More...
 
virtual const double GammaHTot () const
 The total Higgs width \(\Gamma(H)\) in the Standard Model. More...
 
virtual const double GammaHtotautau () const
 The \(\Gamma(H\to \tau^+ \tau^-)\) in the Standard Model. More...
 
virtual const double GammaHtoWWstar () const
 The \(\Gamma(H\to W W^*)\) in the Standard Model. More...
 
virtual const double GammaHtoZga () const
 The \(\Gamma(H\to Z \gamma)\) in the Standard Model. More...
 
virtual const double GammaHtoZZstar () const
 The \(\Gamma(H\to Z Z^*)\) in the Standard Model. More...
 
virtual const double GammaW () const
 The total width of the \(W\) boson, \(\Gamma_W\). More...
 
virtual const double GammaW (const Particle fi, const Particle fj) const
 A partial decay width of the \(W\) boson decay into a SM fermion pair. More...
 
virtual const double GammaZ (const Particle f) const
 The \(Z\to \ell\bar{\ell}\) partial decay width, \(\Gamma_\ell\). More...
 
virtual const double gAnue () const
 The effective (muon) neutrino-electron axial-vector coupling: gAnue. More...
 
const double getAle () const
 A get method to retrieve the fine-structure constant \(\alpha\). More...
 
const double getAlsMz () const
 A get method to access the value of \(\alpha_s(M_Z)\). More...
 
virtual const double getCBd () const
 The ratio of the absolute value of the $B_d$ mixing amplitude over the Standard Model value. More...
 
virtual const double getCBs () const
 The ratio of the absolute value of the $B_s$ mixing amplitude over the Standard Model value. More...
 
virtual const double getCCC1 () const
 A virtual implementation for the RealWeakEFTCC class. More...
 
virtual const double getCCC2 () const
 A virtual implementation for the RealWeakEFTCC class. More...
 
virtual const double getCCC3 () const
 A virtual implementation for the RealWeakEFTCC class. More...
 
virtual const double getCCC4 () const
 A virtual implementation for the RealWeakEFTCC class. More...
 
virtual const double getCCC5 () const
 A virtual implementation for the RealWeakEFTCC class. More...
 
virtual const double getCDMK () const
 The ratio of the real part of the $K$ mixing amplitude over the Standard Model value. More...
 
virtual const double getCepsK () const
 The ratio of the imaginary part of the $K$ mixing amplitude over the Standard Model value. More...
 
const CKMgetCKM () const
 A get method to retrieve the member object of type CKM. More...
 
const double getDAle5Mz () const
 A get method to retrieve the five-flavour hadronic contribution to the electromagnetic coupling, \(\Delta\alpha_{\mathrm{had}}^{(5)}(M_Z^2)\). More...
 
const double getDelGammaZ () const
 A get method to retrieve the theoretical uncertainty in \(\Gamma_Z\), denoted as \(\delta\,\Gamma_Z\). More...
 
const double getDelMw () const
 A get method to retrieve the theoretical uncertainty in \(M_W\), denoted as \(\delta\,M_W\). More...
 
const double getDelR0b () const
 A get method to retrieve the theoretical uncertainty in \(R_b^0\), denoted as \(\delta\,R_b^0\). More...
 
const double getDelR0c () const
 A get method to retrieve the theoretical uncertainty in \(R_c^0\), denoted as \(\delta\,R_c^0\). More...
 
const double getDelR0l () const
 A get method to retrieve the theoretical uncertainty in \(R_l^0\), denoted as \(\delta\,R_l^0\). More...
 
const double getDelSigma0H () const
 A get method to retrieve the theoretical uncertainty in \(\sigma_{Hadron}^0\), denoted as \(\delta\,\sigma_{Hadron}^0\). More...
 
const double getDelSin2th_b () const
 A get method to retrieve the theoretical uncertainty in \(\sin^2\theta_{\rm eff}^{b}\), denoted as \(\delta\sin^2\theta_{\rm eff}^{b}\). More...
 
const double getDelSin2th_l () const
 A get method to retrieve the theoretical uncertainty in \(\sin^2\theta_{\rm eff}^{\rm lept}\), denoted as \(\delta\sin^2\theta_{\rm eff}^{\rm lept}\). More...
 
const double getDelSin2th_q () const
 A get method to retrieve the theoretical uncertainty in \(\sin^2\theta_{\rm eff}^{q\not = b,t}\), denoted as \(\delta\sin^2\theta_{\rm eff}^{q\not = b,t}\). More...
 
const std::string getFlagKappaZ () const
 A method to retrieve the model flag KappaZ. More...
 
const std::string getFlagMw () const
 A method to retrieve the model flag Mw. More...
 
const std::string getFlagRhoZ () const
 A method to retrieve the model flag RhoZ. More...
 
const FlavourgetFlavour () const
 
const double getGF () const
 A get method to retrieve the Fermi constant \(G_\mu\). More...
 
const int getIterationNo () const
 
const ParticlegetLeptons (const QCD::lepton p) const
 A get method to retrieve the member object of a lepton. More...
 
virtual const double getMHl () const
 A get method to retrieve the Higgs mass \(m_h\). More...
 
virtual const double getmq (const QCD::quark q, const double mu) const
 The MSbar running quark mass computed at NLO. More...
 
const double getMuw () const
 A get method to retrieve the matching scale \(\mu_W\) around the weak scale. More...
 
const double getMw () const
 A get method to access the input value of the mass of the \(W\) boson \(M_W\). More...
 
EWSMApproximateFormulaegetMyApproximateFormulae () const
 A get method to retrieve the member pointer of type EWSMApproximateFormulae. More...
 
EWSMcachegetMyEWSMcache () const
 A get method to retrieve the member pointer of type EWSMcache. More...
 
LeptonFlavourgetMyLeptonFlavour () const
 
EWSMOneLoopEWgetMyOneLoopEW () const
 A get method to retrieve the member pointer of type EWSMOneLoopEW,. More...
 
EWSMThreeLoopEWgetMyThreeLoopEW () const
 
EWSMThreeLoopEW2QCDgetMyThreeLoopEW2QCD () const
 
EWSMThreeLoopQCDgetMyThreeLoopQCD () const
 
EWSMTwoFermionsLEP2getMyTwoFermionsLEP2 () const
 A get method to retrieve the member pointer of type EWSMTwoFermionsLEP2. More...
 
EWSMTwoLoopEWgetMyTwoLoopEW () const
 
EWSMTwoLoopQCDgetMyTwoLoopQCD () const
 
const double getMz () const
 A get method to access the mass of the \(Z\) boson \(M_Z\). More...
 
virtual const double getPhiBd () const
 Half the relative phase of the $B_d$ mixing amplitude w.r.t. the Standard Model one. More...
 
virtual const double getPhiBs () const
 Half the relative phase of the $B_s$ mixing amplitude w.r.t. the Standard Model one. More...
 
virtual const StandardModelgetTrueSM () const
 
const gslpp::matrix< gslpp::complex > getUPMNS () const
 A get method to retrieve the object of the PMNS matrix. More...
 
const gslpp::matrix< gslpp::complex > getVCKM () const
 A get method to retrieve the CKM matrix. More...
 
const gslpp::matrix< gslpp::complex > & getYd () const
 A get method to retrieve the Yukawa matrix of the down-type quarks, \(Y_d\). More...
 
const gslpp::matrix< gslpp::complex > & getYe () const
 A get method to retrieve the Yukawa matrix of the charged leptons, \(Y_e\). More...
 
const gslpp::matrix< gslpp::complex > & getYn () const
 A get method to retrieve the Yukawa matrix of the neutrinos, \(Y_\nu\). More...
 
const gslpp::matrix< gslpp::complex > & getYu () const
 A get method to retrieve the Yukawa matrix of the up-type quarks, \(Y_u\). More...
 
virtual const double gLnuN2 () const
 The effective neutrino nucleon LH coupling: gLnuN2. More...
 
virtual const double gRnuN2 () const
 The effective neutrino nucleon RH coupling: gRnuN2. More...
 
virtual const gslpp::complex gV_f (const Particle f) const
 The effective leptonic neutral-current vector coupling \(g_V^l\) in the SM. More...
 
virtual const double gVnue () const
 The effective (muon) neutrino-electron vector coupling: gVnue. More...
 
gslpp::complex I_triangle_1 (const double tau, const double lambda) const
 Loop function entering in the calculation of the effective \(HZ\gamma\) coupling. More...
 
gslpp::complex I_triangle_2 (const double tau, const double lambda) const
 Loop function entering in the calculation of the effective \(HZ\gamma\) coupling. More...
 
virtual bool Init (const std::map< std::string, double > &DPars)
 A method to initialize the model parameters. More...
 
const double intMLL2eeeeus2 (const double s, const double t0, const double t1) const
 
const double intMLR2eeeets2 (const double s, const double t0, const double t1) const
 
const double intMLRtilde2eeeest2 (const double s, const double t0, const double t1) const
 
const double intMRR2eeeeus2 (const double s, const double t0, const double t1) const
 
const bool IsFlagNoApproximateGammaZ () const
 A method to retrieve the model flag NoApproximateGammaZ. More...
 
const bool IsFlagWithoutNonUniversalVC () const
 A method to retrieve the model flag WithoutNonUniversalVC. More...
 
const bool isSMSuccess () const
 A get method to retrieve the success status of the Standard Model update and matching. More...
 
virtual const gslpp::complex kappaZ_f (const Particle f) const
 The effective leptonic neutral-current coupling \(\kappa_Z^l\) in the SM. More...
 
virtual const double LEP2AFBbottom (const double s) const
 
virtual const double LEP2AFBcharm (const double s) const
 
virtual const double LEP2AFBe (const double s) const
 
virtual const double LEP2AFBmu (const double s) const
 
virtual const double LEP2AFBtau (const double s) const
 
virtual const double LEP2dsigmadcosBinE (const double s, const double cos, const double cosmin, const double cosmax) const
 
virtual const double LEP2dsigmadcosBinMu (const double s, const double cos, const double cosmin, const double cosmax) const
 
virtual const double LEP2dsigmadcosBinTau (const double s, const double cos, const double cosmin, const double cosmax) const
 
virtual const double LEP2dsigmadcosE (const double s, const double cos) const
 
virtual const double LEP2dsigmadcosMu (const double s, const double cos) const
 
virtual const double LEP2dsigmadcosTau (const double s, const double cos) const
 
virtual const double LEP2Rbottom (const double s) const
 
virtual const double LEP2Rcharm (const double s) const
 
virtual const double LEP2sigmaBottom (const double s) const
 
virtual const double LEP2sigmaCharm (const double s) const
 
virtual const double LEP2sigmaE (const double s) const
 
virtual const double LEP2sigmaHadron (const double s) const
 
virtual const double LEP2sigmaMu (const double s) const
 
virtual const double LEP2sigmaTau (const double s) const
 
const double MLL2eeff (const Particle f, const double s, const double t) const
 
const double MLR2eeff (const Particle f, const double s) const
 
const double MRL2eeff (const Particle f, const double s) const
 
const double MRR2eeff (const Particle f, const double s, const double t) const
 
virtual const double Mw () const
 The SM prediction for the \(W\)-boson mass in the on-shell scheme, \(M_{W,\mathrm{SM}}\). More...
 
const double Mw_tree () const
 The tree-level mass of the \(W\) boson, \(M_W^{\mathrm{tree}}\). More...
 
const double MwbarFromMw (const double Mw) const
 A method to convert the \(W\)-boson mass in the experimental/running-width scheme to that in the complex-pole/fixed-width scheme. More...
 
const double MwFromMwbar (const double Mwbar) const
 A method to convert the \(W\)-boson mass in the complex-pole/fixed-width scheme to that in the experimental/running-width scheme. More...
 
double Mzbar () const
 The \(Z\)-boson mass \(\overline{M}_Z\) in the complex-pole/fixed-width scheme. More...
 
virtual const double N_nu () const
 The number of neutrinos obtained indirectly from the measurements at the Z pole, \(N_{\nu}\). More...
 
virtual bool PreUpdate ()
 The pre-update method for StandardModel. More...
 
virtual const double Qwemoller (const double q2, const double y) const
 The computation of the electron's weak charge. More...
 
virtual const double Qwn () const
 The computation of the neutron weak charge: Qwn. More...
 
virtual const double Qwp () const
 The computation of the proton weak charge: Qwp. More...
 
virtual const double R0_f (const Particle f) const
 The ratio \(R_\ell^0=\Gamma(Z\to {\rm hadrons})/\Gamma(Z\to \ell^+ \ell^-)\). More...
 
virtual const double R_inv () const
 The ratio of the invisible and leptonic (electron) decay widths of the \(Z\) boson, \(R_{inv}\). More...
 
virtual const double rho_GammaW (const Particle fi, const Particle fj) const
 EW radiative corrections to the width of \(W \to f_i \bar{f}_j\), denoted as \(\rho^W_{ij}\). More...
 
virtual const gslpp::complex rhoZ_f (const Particle f) const
 The effective leptonic neutral-current coupling \(\rho_Z^l\) in the SM. More...
 
virtual const double Ruc () const
 
virtual const double RWc () const
 The ratio \(R_{W,c)=\Gamma(W\to c + X)/\Gamma(W\to had)\). More...
 
virtual const double RWlilj (const Particle li, const Particle lj) const
 The lepton universality ratio \(R_{W,l_i/l_j)=\Gamma(W\to l_i \nu_i)/\Gamma(W\to l_j \nu_j)\). More...
 
virtual const double RZlilj (const Particle li, const Particle lj) const
 The lepton universality ratio \(R_{Z,l_i/l_j)=\Gamma(Z\to l_i^+ l_i^-)/\Gamma(Z\to l_j^+ l_j^-)\). More...
 
const double s02 () const
 The square of the sine of the weak mixing angle \(s_0^2\) defined without weak radiative corrections. More...
 
void setCKM (const CKM &CKMMatrix)
 A set method to change the CKM matrix. More...
 
virtual bool setFlag (const std::string name, const bool value)
 A method to set a flag of StandardModel. More...
 
void setFlagCacheInStandardModel (bool FlagCacheInStandardModel)
 A set method to change the model flag CacheInStandardModel of StandardModel. More...
 
void setFlagNoApproximateGammaZ (bool FlagNoApproximateGammaZ)
 
bool setFlagSigmaForAFB (const bool flagSigmaForAFB_i)
 
bool setFlagSigmaForR (const bool flagSigmaForR_i)
 
virtual bool setFlagStr (const std::string name, const std::string value)
 A method to set a flag of StandardModel. More...
 
void setRequireCKM (bool requireCKM)
 A set method to change the value of requireCKM. More...
 
void setSMSuccess (bool success) const
 A set method to change the success status of the Standard Model update and matching. More...
 
void setYd (const gslpp::matrix< gslpp::complex > &Yd)
 A set method to set the Yukawa matrix of the down-type quarks, \(Y_d\). More...
 
void setYe (const gslpp::matrix< gslpp::complex > &Ye)
 A set method to set the Yukawa matrix of the charged leptons, \(Y_e\). More...
 
void setYu (const gslpp::matrix< gslpp::complex > &Yu)
 A set method to set the Yukawa matrix of the up-type quarks, \(Y_u\). More...
 
virtual const double sigma0_had () const
 The hadronic cross section for \(e^+e^- \to Z \to \mathrm{hadrons}\) at the \(Z\)-pole, \(\sigma_h^0\). More...
 
virtual const double SigmaeeHee (const double sqrt_s, const double Pe, const double Pp) const
 The \(\sigma(e^+ e^- \to e^+ e^- H)\) in the Standard Model. More...
 
virtual const double SigmaeeHvv (const double sqrt_s, const double Pe, const double Pp) const
 The \(\sigma(e^+ e^- \to \nu \bar{\nu} H)\) in the Standard Model. More...
 
virtual const double SigmaeeZH (const double sqrt_s, const double Pe, const double Pp) const
 The \(\sigma(e^+ e^- \to Z H)\) in the Standard Model. More...
 
virtual const double sin2thetaEff (const Particle f) const
 The effective weak mixing angle \(\sin^2\theta_{\rm eff}^{\,\ell}\) for \(Z\ell\bar{\ell}\) at the the \(Z\)-mass scale. More...
 
 StandardModel ()
 The default constructor. More...
 
const double sW2 () const
 
virtual const double sW2 (const double Mw_i) const
 The square of the sine of the weak mixing angle in the on-shell scheme, denoted as \(s_W^2\). More...
 
const double sW2_MSbar_Approx () const
 The (approximated formula for the) square of the sine of the weak mixing angle in the MSbar scheme, denoted as \(\hat{s}_{W}^2\). See: PDG 22, R.L. Workman et al. (Particle Data Group), Prog. Theor. Exp. Phys. 2022, 083C01 (2022) More...
 
const double sW2_ND () const
 The square of the sine of the weak mixing angle in the MSbar-ND scheme (w/o decoupling $\alpha\ln(m_t/M_Z)$ terms), denoted as \(\hat{s}_{ND}^2\). See: PDG 22, R.L. Workman et al. (Particle Data Group), Prog. Theor. Exp. Phys. 2022, 083C01 (2022) (eq. 10.13a/10.13b) More...
 
virtual const double TauLFU_gmuge () const
 The computation of the LFU ratio \(g_\mu/ g_e \). More...
 
virtual const double TauLFU_gtauge () const
 The computation of the LFU ratio \(g_\tau/ g_e \). More...
 
virtual const double TauLFU_gtaugmu () const
 The computation of the LFU ratio \(g_\tau/ g_\mu \). More...
 
virtual const double TauLFU_gtaugmuK () const
 The computation of the LFU ratio \(\left(g_\tau/ g_\mu\right)_K \). More...
 
virtual const double TauLFU_gtaugmuPi () const
 The computation of the LFU ratio \(\left(g_\tau/ g_\mu\right)_\pi \). More...
 
virtual const double ThetaLnuN () const
 The effective neutrino nucleon LH parameter: ThetaLnuN. More...
 
virtual const double ThetaRnuN () const
 The effective neutrino nucleon RH parameter: ThetaRnuN. More...
 
const double tovers2 (const double cosmin, const double cosmax) const
 
const double uovers2 (const double cosmin, const double cosmax) const
 
virtual bool Update (const std::map< std::string, double > &DPars)
 The update method for StandardModel. More...
 
const double v () const
 The Higgs vacuum expectation value. More...
 
virtual ~StandardModel ()
 The default destructor. More...
 
- Public Member Functions inherited from QCD
const double AboveTh (const double mu) const
 The active flavour threshold above the scale \(\mu\) as defined in QCD::Thresholds(). More...
 
void addParameters (std::vector< std::string > params_i)
 A method to add parameters that are specific to only one set of observables. More...
 
const double Als (const double mu, const int Nf_in, const orders order=FULLNLO) const
 Computes the running strong coupling \(\alpha_s(\mu)\) with \(N_f\) active flavours in the \(\overline{\mathrm{MS}}\) scheme. In the cases of LO, NLO and FULLNLO, the coupling is computed with AlsWithInit(). On the other hand, in the cases of NNLO and FULLNNLO, the coupling is computed with AlsWithLambda(). More...
 
const double Als (const double mu, const orders order=FULLNLO, const bool Nf_thr=true) const
 
const double Als4 (const double mu) const
 The value of \(\alpha_s^{\mathrm{FULLNLO}}\) at any scale \(\mu\) with the number of flavours \(n_f = 4\). More...
 
const double AlsByOrder (const double mu, const int Nf_in, const orders order=FULLNLO) const
 
const double AlsByOrder (const double mu, const orders order=FULLNLO, bool Nf_thr=true) const
 
const double AlsOLD (const double mu, const orders order=FULLNLO) const
 Computes the running strong coupling \(\alpha_s(\mu)\) in the \(\overline{\mathrm{MS}}\) scheme. In the cases of LO, NLO and FULLNNLO, the coupling is computed with AlsWithInit(). On the other hand, in the cases of NNLO and FULLNNLO, the coupling is computed with AlsWithLambda(). More...
 
const double AlsWithInit (const double mu, const double alsi, const double mu_i, const int nf, const orders order) const
 Computes the running strong coupling \(\alpha_s(\mu)\) from \(\alpha_s(\mu_i)\) in the \(\overline{\mathrm{MS}}\) scheme, where it is forbidden to across a flavour threshold in the RG running from \(\mu_i\) to \(\mu\). More...
 
const double AlsWithLambda (const double mu, const orders order) const
 Computes the running strong coupling \(\alpha_s(\mu)\) in the \(\overline{\mathrm{MS}}\) scheme with the use of \(\Lambda_{\rm QCD}\). More...
 
const double BelowTh (const double mu) const
 The active flavour threshold below the scale \(\mu\) as defined in QCD::Thresholds(). More...
 
const double Beta0 (const double nf) const
 The \(\beta_0(n_f)\) coefficient for a certain number of flavours \(n_f\). More...
 
const double Beta1 (const double nf) const
 The \(\beta_1(n_f)\) coefficient for a certain number of flavours \(n_f\). More...
 
const double Beta2 (const double nf) const
 The \(\beta_2(n_f)\) coefficient for a certain number of flavours \(n_f\). More...
 
const double Beta3 (const double nf) const
 The \(\beta_3(n_f)\) coefficient for a certain number of flavours \(n_f\). More...
 
void CacheShift (double cache[][5], int n) const
 A member used to manage the caching for this class. More...
 
void CacheShift (int cache[][5], int n) const
 
const orders FullOrder (orders order) const
 Return the FULLORDER enum corresponding to order. More...
 
const double Gamma0 (const double nf) const
 The \(\gamma_0\) coefficient used to compute the running of a mass. More...
 
const double Gamma1 (const double nf) const
 The \(\gamma_1\) coefficient used to compute the running of a mass. More...
 
const double Gamma2 (const double nf) const
 The \(\gamma_2\) coefficient used to compute the running of a mass. More...
 
const double getAlsM () const
 A get method to access the value of \(\alpha_s(M_{\alpha_s})\). More...
 
const BParametergetBBd () const
 For getting the bag parameters corresponding to the operator basis \(O_1 -O_5\) in \(\Delta b = 2\) process in the \(B_d\) meson system. More...
 
const BParametergetBBd_subleading () const
 For getting the subleading bag parameters \(R_2 - R_3\) in \(\Delta b = 2\) process in the \(B_d\) meson system. More...
 
const BParametergetBBs () const
 For getting the bag parameters corresponding to the operator basis \(O_1 -O_5\) in \(\Delta b = 2\) process in the \(B_s\) meson system. More...
 
const BParametergetBBs_subleading () const
 For getting the subleading bag parameters \(R_2 - R_3\) in \(\Delta b = 2\) process in the \(B_s\) meson system. More...
 
const BParametergetBD () const
 For getting the bag parameters corresponding to the operator basis \(O_1 -O_5\) in \(\Delta c = 2\) process in the \(D^0\) meson system. More...
 
const BParametergetBK () const
 For getting the bag parameters corresponding to the operator basis \(O_1 -O_5\) in \(\Delta s = 2\) process in the \(K^0\) meson system. More...
 
const BParametergetBKd1 () const
 
const BParametergetBKd3 () const
 
const double getCF () const
 A get method to access the Casimir factor of QCD. More...
 
const double getMAls () const
 A get method to access the mass scale \(M_{\alpha_s}\) at which the strong coupling constant measurement is provided. More...
 
const MesongetMesons (const QCD::meson m) const
 A get method to access a meson as an object of the type Meson. More...
 
const double getMtpole () const
 A get method to access the pole mass of the top quark. More...
 
const double getMub () const
 A get method to access the threshold between five- and four-flavour theory in GeV. More...
 
const double getMuc () const
 A get method to access the threshold between four- and three-flavour theory in GeV. More...
 
const double getMut () const
 A get method to access the threshold between six- and five-flavour theory in GeV. More...
 
const double getNc () const
 A get method to access the number of colours \(N_c\). More...
 
const double getOptionalParameter (std::string name) const
 A method to get parameters that are specific to only one set of observables. More...
 
const ParticlegetQuarks (const QCD::quark q) const
 A get method to access a quark as an object of the type Particle. More...
 
std::vector< std::string > getUnknownParameters ()
 A method to get the vector of the parameters that have been specified in the configuration file but not being used. More...
 
void initializeBParameter (std::string name_i) const
 A method to initialize B Parameter and the corresponding meson. More...
 
void initializeMeson (QCD::meson meson_i) const
 A method to initialize a meson. More...
 
bool isQCDsuccess () const
 A getter for the QCDsuccess flag. More...
 
const double logLambda (const double nf, orders order) const
 Computes \(\ln\Lambda_\mathrm{QCD}\) with nf flavours in GeV. More...
 
const double Mbar2Mp (const double mbar, const quark q, const orders order=FULLNNLO) const
 Converts the \(\overline{\mathrm{MS}}\) mass \(m(m)\) to the pole mass. More...
 
const double Mofmu2Mbar (const double m, const double mu, const quark q) const
 Converts a quark running mass at an arbitrary scale to the corresponding \(\overline{\mathrm{MS}}\) mass \(m(m)\). More...
 
const double Mp2Mbar (const double mp, const quark q, orders order=FULLNNLO) const
 Converts a quark pole mass to the corresponding \(\overline{\mathrm{MS}}\) mass \(m(m)\). More...
 
const double Mrun (const double mu, const double m, const quark q, const orders order=FULLNNLO) const
 Computes a running quark mass \(m(\mu)\) from \(m(m)\). More...
 
const double Mrun (const double mu_f, const double mu_i, const double m, const quark q, const orders order=FULLNNLO) const
 Runs a quark mass from \(\mu_i\) to \(\mu_f\). More...
 
const double Mrun4 (const double mu_f, const double mu_i, const double m) const
 The running of a mass with the number of flavours \(n_f = 4\). More...
 
const double MS2DRqmass (const double MSbar) const
 Converts a quark mass from the \(\overline{\mathrm{MS}}\) scheme to the \(\overline{\mathrm{DR}}\) scheme. More...
 
const double MS2DRqmass (const double MSscale, const double MSbar) const
 Converts a quark mass from the \(\overline{\mathrm{MS}}\) scheme to the \(\overline{\mathrm{DR}}\) scheme. More...
 
const double Nf (const double mu) const
 The number of active flavour at scale \(\mu\). More...
 
const double NfThresholdCorrections (double mu, double M, double als, int nf, orders order) const
 Threshold corrections in matching \(\alpha_s(n_f+1)\) with \(\alpha_s(n_f)\) from eq. (34) of hep-ph/0512060. More...
 
const std::string orderToString (const orders order) const
 Converts an object of the enum type "orders" to the corresponding string. More...
 
 QCD ()
 Constructor. More...
 
void setComputemt (bool computemt)
 A set method to change the value of computemt. More...
 
void setMtpole (double mtpole_in)
 A method to set the pole mass of the top quark. More...
 
void setNc (double Nc)
 A set method to change the number of colours \(N_c\). More...
 
void setOptionalParameter (std::string name, double value)
 A method to set the parameter value for the parameters that are specific to only one set of observables. More...
 
void setQuarkMass (const quark q, const double mass)
 A set method to change the mass of a quark. More...
 
const double Thresholds (const int i) const
 For accessing the active flavour threshold scales. More...
 
- Public Member Functions inherited from Model
void addMissingModelParameter (const std::string &missingParameterName)
 
std::vector< std::string > getmissingModelParameters ()
 
unsigned int getMissingModelParametersCount ()
 
std::string getModelName () const
 A method to fetch the name of the model. More...
 
const double & getModelParam (std::string name) const
 
bool isModelFWC_DF2 () const
 
bool isModelGeneralTHDM () const
 
bool isModelGeorgiMachacek () const
 
bool IsModelInitialized () const
 A method to check if the model is initialized. More...
 
bool isModelLinearized () const
 
bool isModelNPquadratic () const
 
bool isModelParam (std::string name) const
 
bool isModelSUSY () const
 
bool isModelTHDM () const
 
bool isModelTHDMW () const
 
bool IsUpdateError () const
 A method to check if there was any error in the model update process. More...
 
 Model ()
 The default constructor. More...
 
void raiseMissingModelParameterCount ()
 
void setModelFWC_DF2 ()
 
void setModelGeneralTHDM ()
 
void setModelGeorgiMachacek ()
 
void setModelInitialized (bool ModelInitialized)
 A set method to fix the failure or success of the initialization of the model. More...
 
void setModelLinearized (bool linearized=true)
 
void setModelName (const std::string name)
 A method to set the name of the model. More...
 
void setModelNPquadratic (bool NPquadratic=true)
 
void setModelSUSY ()
 
void setModelTHDM ()
 
void setModelTHDMW ()
 
void setSliced (bool Sliced)
 
void setUpdateError (bool UpdateError)
 A set method to fix the update status as success or failure. More...
 
virtual ~Model ()
 The default destructor. More...
 

Static Public Attributes

static const std::string FlavourWilsonCoefficient_DF2vars [NFlavourWilsonCoefficient_DF2vars]
 
static const int NFlavourWilsonCoefficient_DF2vars = 44
 
- Static Public Attributes inherited from StandardModel
static const double GeVminus2_to_nb = 389379.338
 
static const double Mw_error = 0.00001
 The target accuracy of the iterative calculation of the \(W\)-boson mass in units of GeV. More...
 
static const int NSMvars = 26
 The number of the model parameters in StandardModel. More...
 
static const int NumSMParamsForEWPO = 33
 The number of the SM parameters that are relevant to the EW precision observables. More...
 
static std::string SMvars [NSMvars]
 A string array containing the labels of the model parameters in StandardModel. More...
 
- Static Public Attributes inherited from QCD
static const int NQCDvars = 11
 The number of model parameters in QCD. More...
 
static std::string QCDvars [NQCDvars]
 An array containing the labels under which all QCD parameters are stored in a vector of ModelParameter via InputParser::ReadParameters(). More...
 

Protected Member Functions

virtual void setParameter (const std::string, const double &)
 A method to set the value of a parameter of FlavourWilsonCoefficient_DF2. More...
 
- Protected Member Functions inherited from StandardModel
const double AFB_NoISR_l (const QCD::lepton l_flavor, const double s) const
 
const double AFB_NoISR_q (const QCD::quark q_flavor, const double s) const
 
bool checkEWPOscheme (const std::string scheme) const
 A method to check if a given scheme name in string form is valid. More...
 
virtual void computeCKM ()
 The method to compute the CKM matrix. More...
 
virtual void computeYukawas ()
 The method to compute the Yukawas matrix. More...
 
double Delta_EWQCD (const QCD::quark q) const
 The non-factorizable EW-QCD corrections to the partial widths for \(Z\to q\bar{q}\), denoted as \(\Delta_{\mathrm{EW/QCD}}\). More...
 
const double getIntegrand_AFBnumeratorWithISR_bottom133 (double x) const
 
const double getIntegrand_AFBnumeratorWithISR_bottom167 (double x) const
 
const double getIntegrand_AFBnumeratorWithISR_bottom172 (double x) const
 
const double getIntegrand_AFBnumeratorWithISR_bottom183 (double x) const
 
const double getIntegrand_AFBnumeratorWithISR_bottom189 (double x) const
 
const double getIntegrand_AFBnumeratorWithISR_bottom192 (double x) const
 
const double getIntegrand_AFBnumeratorWithISR_bottom196 (double x) const
 
const double getIntegrand_AFBnumeratorWithISR_bottom200 (double x) const
 
const double getIntegrand_AFBnumeratorWithISR_bottom202 (double x) const
 
const double getIntegrand_AFBnumeratorWithISR_bottom205 (double x) const
 
const double getIntegrand_AFBnumeratorWithISR_bottom207 (double x) const
 
const double getIntegrand_AFBnumeratorWithISR_charm133 (double x) const
 
const double getIntegrand_AFBnumeratorWithISR_charm167 (double x) const
 
const double getIntegrand_AFBnumeratorWithISR_charm172 (double x) const
 
const double getIntegrand_AFBnumeratorWithISR_charm183 (double x) const
 
const double getIntegrand_AFBnumeratorWithISR_charm189 (double x) const
 
const double getIntegrand_AFBnumeratorWithISR_charm192 (double x) const
 
const double getIntegrand_AFBnumeratorWithISR_charm196 (double x) const
 
const double getIntegrand_AFBnumeratorWithISR_charm200 (double x) const
 
const double getIntegrand_AFBnumeratorWithISR_charm202 (double x) const
 
const double getIntegrand_AFBnumeratorWithISR_charm205 (double x) const
 
const double getIntegrand_AFBnumeratorWithISR_charm207 (double x) const
 
const double getIntegrand_AFBnumeratorWithISR_mu130 (double x) const
 
const double getIntegrand_AFBnumeratorWithISR_mu136 (double x) const
 
const double getIntegrand_AFBnumeratorWithISR_mu161 (double x) const
 
const double getIntegrand_AFBnumeratorWithISR_mu172 (double x) const
 
const double getIntegrand_AFBnumeratorWithISR_mu183 (double x) const
 
const double getIntegrand_AFBnumeratorWithISR_mu189 (double x) const
 
const double getIntegrand_AFBnumeratorWithISR_mu192 (double x) const
 
const double getIntegrand_AFBnumeratorWithISR_mu196 (double x) const
 
const double getIntegrand_AFBnumeratorWithISR_mu200 (double x) const
 
const double getIntegrand_AFBnumeratorWithISR_mu202 (double x) const
 
const double getIntegrand_AFBnumeratorWithISR_mu205 (double x) const
 
const double getIntegrand_AFBnumeratorWithISR_mu207 (double x) const
 
const double getIntegrand_AFBnumeratorWithISR_tau130 (double x) const
 
const double getIntegrand_AFBnumeratorWithISR_tau136 (double x) const
 
const double getIntegrand_AFBnumeratorWithISR_tau161 (double x) const
 
const double getIntegrand_AFBnumeratorWithISR_tau172 (double x) const
 
const double getIntegrand_AFBnumeratorWithISR_tau183 (double x) const
 
const double getIntegrand_AFBnumeratorWithISR_tau189 (double x) const
 
const double getIntegrand_AFBnumeratorWithISR_tau192 (double x) const
 
const double getIntegrand_AFBnumeratorWithISR_tau196 (double x) const
 
const double getIntegrand_AFBnumeratorWithISR_tau200 (double x) const
 
const double getIntegrand_AFBnumeratorWithISR_tau202 (double x) const
 
const double getIntegrand_AFBnumeratorWithISR_tau205 (double x) const
 
const double getIntegrand_AFBnumeratorWithISR_tau207 (double x) const
 
const double getIntegrand_dsigmaBox_bottom130 (double x) const
 
const double getIntegrand_dsigmaBox_bottom133 (double x) const
 
const double getIntegrand_dsigmaBox_bottom136 (double x) const
 
const double getIntegrand_dsigmaBox_bottom161 (double x) const
 
const double getIntegrand_dsigmaBox_bottom167 (double x) const
 
const double getIntegrand_dsigmaBox_bottom172 (double x) const
 
const double getIntegrand_dsigmaBox_bottom183 (double x) const
 
const double getIntegrand_dsigmaBox_bottom189 (double x) const
 
const double getIntegrand_dsigmaBox_bottom192 (double x) const
 
const double getIntegrand_dsigmaBox_bottom196 (double x) const
 
const double getIntegrand_dsigmaBox_bottom200 (double x) const
 
const double getIntegrand_dsigmaBox_bottom202 (double x) const
 
const double getIntegrand_dsigmaBox_bottom205 (double x) const
 
const double getIntegrand_dsigmaBox_bottom207 (double x) const
 
const double getIntegrand_dsigmaBox_charm130 (double x) const
 
const double getIntegrand_dsigmaBox_charm133 (double x) const
 
const double getIntegrand_dsigmaBox_charm136 (double x) const
 
const double getIntegrand_dsigmaBox_charm161 (double x) const
 
const double getIntegrand_dsigmaBox_charm167 (double x) const
 
const double getIntegrand_dsigmaBox_charm172 (double x) const
 
const double getIntegrand_dsigmaBox_charm183 (double x) const
 
const double getIntegrand_dsigmaBox_charm189 (double x) const
 
const double getIntegrand_dsigmaBox_charm192 (double x) const
 
const double getIntegrand_dsigmaBox_charm196 (double x) const
 
const double getIntegrand_dsigmaBox_charm200 (double x) const
 
const double getIntegrand_dsigmaBox_charm202 (double x) const
 
const double getIntegrand_dsigmaBox_charm205 (double x) const
 
const double getIntegrand_dsigmaBox_charm207 (double x) const
 
const double getIntegrand_dsigmaBox_down130 (double x) const
 
const double getIntegrand_dsigmaBox_down133 (double x) const
 
const double getIntegrand_dsigmaBox_down136 (double x) const
 
const double getIntegrand_dsigmaBox_down161 (double x) const
 
const double getIntegrand_dsigmaBox_down167 (double x) const
 
const double getIntegrand_dsigmaBox_down172 (double x) const
 
const double getIntegrand_dsigmaBox_down183 (double x) const
 
const double getIntegrand_dsigmaBox_down189 (double x) const
 
const double getIntegrand_dsigmaBox_down192 (double x) const
 
const double getIntegrand_dsigmaBox_down196 (double x) const
 
const double getIntegrand_dsigmaBox_down200 (double x) const
 
const double getIntegrand_dsigmaBox_down202 (double x) const
 
const double getIntegrand_dsigmaBox_down205 (double x) const
 
const double getIntegrand_dsigmaBox_down207 (double x) const
 
const double getIntegrand_dsigmaBox_mu130 (double x) const
 
const double getIntegrand_dsigmaBox_mu133 (double x) const
 
const double getIntegrand_dsigmaBox_mu136 (double x) const
 
const double getIntegrand_dsigmaBox_mu161 (double x) const
 
const double getIntegrand_dsigmaBox_mu167 (double x) const
 
const double getIntegrand_dsigmaBox_mu172 (double x) const
 
const double getIntegrand_dsigmaBox_mu183 (double x) const
 
const double getIntegrand_dsigmaBox_mu189 (double x) const
 
const double getIntegrand_dsigmaBox_mu192 (double x) const
 
const double getIntegrand_dsigmaBox_mu196 (double x) const
 
const double getIntegrand_dsigmaBox_mu200 (double x) const
 
const double getIntegrand_dsigmaBox_mu202 (double x) const
 
const double getIntegrand_dsigmaBox_mu205 (double x) const
 
const double getIntegrand_dsigmaBox_mu207 (double x) const
 
const double getIntegrand_dsigmaBox_strange130 (double x) const
 
const double getIntegrand_dsigmaBox_strange133 (double x) const
 
const double getIntegrand_dsigmaBox_strange136 (double x) const
 
const double getIntegrand_dsigmaBox_strange161 (double x) const
 
const double getIntegrand_dsigmaBox_strange167 (double x) const
 
const double getIntegrand_dsigmaBox_strange172 (double x) const
 
const double getIntegrand_dsigmaBox_strange183 (double x) const
 
const double getIntegrand_dsigmaBox_strange189 (double x) const
 
const double getIntegrand_dsigmaBox_strange192 (double x) const
 
const double getIntegrand_dsigmaBox_strange196 (double x) const
 
const double getIntegrand_dsigmaBox_strange200 (double x) const
 
const double getIntegrand_dsigmaBox_strange202 (double x) const
 
const double getIntegrand_dsigmaBox_strange205 (double x) const
 
const double getIntegrand_dsigmaBox_strange207 (double x) const
 
const double getIntegrand_dsigmaBox_tau130 (double x) const
 
const double getIntegrand_dsigmaBox_tau133 (double x) const
 
const double getIntegrand_dsigmaBox_tau136 (double x) const
 
const double getIntegrand_dsigmaBox_tau161 (double x) const
 
const double getIntegrand_dsigmaBox_tau167 (double x) const
 
const double getIntegrand_dsigmaBox_tau172 (double x) const
 
const double getIntegrand_dsigmaBox_tau183 (double x) const
 
const double getIntegrand_dsigmaBox_tau189 (double x) const
 
const double getIntegrand_dsigmaBox_tau192 (double x) const
 
const double getIntegrand_dsigmaBox_tau196 (double x) const
 
const double getIntegrand_dsigmaBox_tau200 (double x) const
 
const double getIntegrand_dsigmaBox_tau202 (double x) const
 
const double getIntegrand_dsigmaBox_tau205 (double x) const
 
const double getIntegrand_dsigmaBox_tau207 (double x) const
 
const double getIntegrand_dsigmaBox_up130 (double x) const
 
const double getIntegrand_dsigmaBox_up133 (double x) const
 
const double getIntegrand_dsigmaBox_up136 (double x) const
 
const double getIntegrand_dsigmaBox_up161 (double x) const
 
const double getIntegrand_dsigmaBox_up167 (double x) const
 
const double getIntegrand_dsigmaBox_up172 (double x) const
 
const double getIntegrand_dsigmaBox_up183 (double x) const
 
const double getIntegrand_dsigmaBox_up189 (double x) const
 
const double getIntegrand_dsigmaBox_up192 (double x) const
 
const double getIntegrand_dsigmaBox_up196 (double x) const
 
const double getIntegrand_dsigmaBox_up200 (double x) const
 
const double getIntegrand_dsigmaBox_up202 (double x) const
 
const double getIntegrand_dsigmaBox_up205 (double x) const
 
const double getIntegrand_dsigmaBox_up207 (double x) const
 
const double getIntegrand_sigmaWithISR_bottom130 (double x) const
 
const double getIntegrand_sigmaWithISR_bottom133 (double x) const
 
const double getIntegrand_sigmaWithISR_bottom136 (double x) const
 
const double getIntegrand_sigmaWithISR_bottom161 (double x) const
 
const double getIntegrand_sigmaWithISR_bottom167 (double x) const
 
const double getIntegrand_sigmaWithISR_bottom172 (double x) const
 
const double getIntegrand_sigmaWithISR_bottom183 (double x) const
 
const double getIntegrand_sigmaWithISR_bottom189 (double x) const
 
const double getIntegrand_sigmaWithISR_bottom192 (double x) const
 
const double getIntegrand_sigmaWithISR_bottom196 (double x) const
 
const double getIntegrand_sigmaWithISR_bottom200 (double x) const
 
const double getIntegrand_sigmaWithISR_bottom202 (double x) const
 
const double getIntegrand_sigmaWithISR_bottom205 (double x) const
 
const double getIntegrand_sigmaWithISR_bottom207 (double x) const
 
const double getIntegrand_sigmaWithISR_charm130 (double x) const
 
const double getIntegrand_sigmaWithISR_charm133 (double x) const
 
const double getIntegrand_sigmaWithISR_charm136 (double x) const
 
const double getIntegrand_sigmaWithISR_charm161 (double x) const
 
const double getIntegrand_sigmaWithISR_charm167 (double x) const
 
const double getIntegrand_sigmaWithISR_charm172 (double x) const
 
const double getIntegrand_sigmaWithISR_charm183 (double x) const
 
const double getIntegrand_sigmaWithISR_charm189 (double x) const
 
const double getIntegrand_sigmaWithISR_charm192 (double x) const
 
const double getIntegrand_sigmaWithISR_charm196 (double x) const
 
const double getIntegrand_sigmaWithISR_charm200 (double x) const
 
const double getIntegrand_sigmaWithISR_charm202 (double x) const
 
const double getIntegrand_sigmaWithISR_charm205 (double x) const
 
const double getIntegrand_sigmaWithISR_charm207 (double x) const
 
const double getIntegrand_sigmaWithISR_down130 (double x) const
 
const double getIntegrand_sigmaWithISR_down133 (double x) const
 
const double getIntegrand_sigmaWithISR_down136 (double x) const
 
const double getIntegrand_sigmaWithISR_down161 (double x) const
 
const double getIntegrand_sigmaWithISR_down167 (double x) const
 
const double getIntegrand_sigmaWithISR_down172 (double x) const
 
const double getIntegrand_sigmaWithISR_down183 (double x) const
 
const double getIntegrand_sigmaWithISR_down189 (double x) const
 
const double getIntegrand_sigmaWithISR_down192 (double x) const
 
const double getIntegrand_sigmaWithISR_down196 (double x) const
 
const double getIntegrand_sigmaWithISR_down200 (double x) const
 
const double getIntegrand_sigmaWithISR_down202 (double x) const
 
const double getIntegrand_sigmaWithISR_down205 (double x) const
 
const double getIntegrand_sigmaWithISR_down207 (double x) const
 
const double getIntegrand_sigmaWithISR_mu130 (double x) const
 
const double getIntegrand_sigmaWithISR_mu136 (double x) const
 
const double getIntegrand_sigmaWithISR_mu161 (double x) const
 
const double getIntegrand_sigmaWithISR_mu172 (double x) const
 
const double getIntegrand_sigmaWithISR_mu183 (double x) const
 
const double getIntegrand_sigmaWithISR_mu189 (double x) const
 
const double getIntegrand_sigmaWithISR_mu192 (double x) const
 
const double getIntegrand_sigmaWithISR_mu196 (double x) const
 
const double getIntegrand_sigmaWithISR_mu200 (double x) const
 
const double getIntegrand_sigmaWithISR_mu202 (double x) const
 
const double getIntegrand_sigmaWithISR_mu205 (double x) const
 
const double getIntegrand_sigmaWithISR_mu207 (double x) const
 
const double getIntegrand_sigmaWithISR_strange130 (double x) const
 
const double getIntegrand_sigmaWithISR_strange133 (double x) const
 
const double getIntegrand_sigmaWithISR_strange136 (double x) const
 
const double getIntegrand_sigmaWithISR_strange161 (double x) const
 
const double getIntegrand_sigmaWithISR_strange167 (double x) const
 
const double getIntegrand_sigmaWithISR_strange172 (double x) const
 
const double getIntegrand_sigmaWithISR_strange183 (double x) const
 
const double getIntegrand_sigmaWithISR_strange189 (double x) const
 
const double getIntegrand_sigmaWithISR_strange192 (double x) const
 
const double getIntegrand_sigmaWithISR_strange196 (double x) const
 
const double getIntegrand_sigmaWithISR_strange200 (double x) const
 
const double getIntegrand_sigmaWithISR_strange202 (double x) const
 
const double getIntegrand_sigmaWithISR_strange205 (double x) const
 
const double getIntegrand_sigmaWithISR_strange207 (double x) const
 
const double getIntegrand_sigmaWithISR_tau130 (double x) const
 
const double getIntegrand_sigmaWithISR_tau136 (double x) const
 
const double getIntegrand_sigmaWithISR_tau161 (double x) const
 
const double getIntegrand_sigmaWithISR_tau172 (double x) const
 
const double getIntegrand_sigmaWithISR_tau183 (double x) const
 
const double getIntegrand_sigmaWithISR_tau189 (double x) const
 
const double getIntegrand_sigmaWithISR_tau192 (double x) const
 
const double getIntegrand_sigmaWithISR_tau196 (double x) const
 
const double getIntegrand_sigmaWithISR_tau200 (double x) const
 
const double getIntegrand_sigmaWithISR_tau202 (double x) const
 
const double getIntegrand_sigmaWithISR_tau205 (double x) const
 
const double getIntegrand_sigmaWithISR_tau207 (double x) const
 
const double getIntegrand_sigmaWithISR_up130 (double x) const
 
const double getIntegrand_sigmaWithISR_up133 (double x) const
 
const double getIntegrand_sigmaWithISR_up136 (double x) const
 
const double getIntegrand_sigmaWithISR_up161 (double x) const
 
const double getIntegrand_sigmaWithISR_up167 (double x) const
 
const double getIntegrand_sigmaWithISR_up172 (double x) const
 
const double getIntegrand_sigmaWithISR_up183 (double x) const
 
const double getIntegrand_sigmaWithISR_up189 (double x) const
 
const double getIntegrand_sigmaWithISR_up192 (double x) const
 
const double getIntegrand_sigmaWithISR_up196 (double x) const
 
const double getIntegrand_sigmaWithISR_up200 (double x) const
 
const double getIntegrand_sigmaWithISR_up202 (double x) const
 
const double getIntegrand_sigmaWithISR_up205 (double x) const
 
const double getIntegrand_sigmaWithISR_up207 (double x) const
 
const double Integrand_AFBnumeratorWithISR_l (double x, const QCD::lepton l_flavor, const double s) const
 
const double Integrand_AFBnumeratorWithISR_q (double x, const QCD::quark q_flavor, const double s) const
 
const double Integrand_dsigmaBox_l (double cosTheta, const QCD::lepton l_flavor, const double s) const
 
const double Integrand_dsigmaBox_q (double cosTheta, const QCD::quark q_flavor, const double s) const
 
const double Integrand_sigmaWithISR_l (double x, const QCD::lepton l_flavor, const double s) const
 
const double Integrand_sigmaWithISR_q (double x, const QCD::quark q_flavor, const double s) const
 
double m_q (const QCD::quark q, const double mu, const orders order=FULLNLO) const
 
double RAq (const QCD::quark q) const
 The radiator factor associated with the final-state QED and QCD corrections to the the axial-vector-current interactions, \(R_A^q(M_Z^2)\). More...
 
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 \(\kappa_Z^f\) from \(\Delta\rho\), \(\delta\rho_{\rm rem}^{f}\) and \(\Delta r_{\mathrm{rem}}\). More...
 
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 \(W\)-boson mass from \(\Delta\rho\) and \(\Delta r_{\mathrm{rem}}\). More...
 
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 \(\rho_Z^f\) from \(\Delta\rho\), \(\delta\rho_{\rm rem}^{f}\) and \(\Delta r_{\mathrm{rem}}\). More...
 
double RVh () const
 The singlet vector corrections to the hadronic \(Z\)-boson width, denoted as \(R_V^h\). More...
 
double RVq (const QCD::quark q) const
 The radiator factor associated with the final-state QED and QCD corrections to the the vector-current interactions, \(R_V^q(M_Z^2)\). More...
 
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 precision. More...
 
const double sigma_NoISR_l (const QCD::lepton l_flavor, const double s) const
 
const double sigma_NoISR_q (const QCD::quark q_flavor, const double s) const
 
double taub () const
 Top-mass corrections to the \(Zb\bar{b}\) vertex, denoted by \(\tau_b\). More...
 
- Protected Member Functions inherited from QCD
const double MassOfNf (int nf) const
 The Mbar mass of the heaviest quark in the theory with Nf active flavour. More...
 

Protected Attributes

Matching< FlavourWilsonCoefficient_DF2Matching, FlavourWilsonCoefficient_DF2FWCM
 The FlavourWilsonCoefficientMatching_DF2 object. More...
 
- Protected Attributes inherited from StandardModel
double A
 The CKM parameter \(A\) in the Wolfenstein parameterization. More...
 
double ale
 The fine-structure constant \(\alpha\). More...
 
double alpha21
 
double alpha31
 
double AlsMz
 The strong coupling constant at the Z-boson mass, \(\alpha_s(M_Z)\). More...
 
bool bSigmaForAFB
 
bool bSigmaForR
 
double dAl5hMz
 The five-flavour hadronic contribution to the electromagnetic coupling, \(\Delta\alpha_{\mathrm{had}}^{(5)}(M_Z^2)\). (Non-input parameter) More...
 
double dAle5Mz
 The five-flavour hadronic contribution to the electromagnetic coupling, \(\Delta\alpha_{\mathrm{had}}^{(5)}(M_Z^2)\), used as input for FlagMWinput = FALSE. More...
 
double delGammaZ
 The theoretical uncertainty in \(\Gamma_Z\), denoted as \(\delta\,\Gamma_Z\), in GeV. More...
 
double delMw
 The theoretical uncertainty in \(M_W\), denoted as \(\delta\,M_W\), in GeV. More...
 
double delR0b
 The theoretical uncertainty in \(R_b^0\), denoted as \(\delta\,R_b^0\). More...
 
double delR0c
 The theoretical uncertainty in \(R_c^0\), denoted as \(\delta\,R_c^0\). More...
 
double delR0l
 The theoretical uncertainty in \(R_l^0\), denoted as \(\delta\,R_l^0\). More...
 
double delsigma0H
 The theoretical uncertainty in \(\sigma_{Hadron}^0\), denoted as \(\delta\,\sigma_{Hadron}^0\) in nb. More...
 
double delSin2th_b
 The theoretical uncertainty in \(\sin^2\theta_{\rm eff}^{b}\), denoted as \(\delta\sin^2\theta_{\rm eff}^{b}\). More...
 
double delSin2th_l
 The theoretical uncertainty in \(\sin^2\theta_{\rm eff}^{\rm lept}\), denoted as \(\delta\sin^2\theta_{\rm eff}^{\rm lept}\). More...
 
double delSin2th_q
 The theoretical uncertainty in \(\sin^2\theta_{\rm eff}^{q\not = b,t}\), denoted as \(\delta\sin^2\theta_{\rm eff}^{q\not = b,t}\). More...
 
double delta
 
double etab
 The CKM parameter \(\bar{\eta}\) in the Wolfenstein parameterization. More...
 
bool flag_order [orders_EW_size]
 An array of internal flags controlling the inclusions of higher-order corrections. More...
 
bool FlagFixMuwMut
 A boolean for the model flag FixMuwMut. More...
 
bool flagLEP2 [NUMofLEP2RCs]
 
double gamma
 \(\gamma \) used as an input for FlagWolfenstein = FALSE More...
 
double GF
 The Fermi constant \(G_\mu\) in \({\rm GeV}^{-2}\). More...
 
double lambda
 The CKM parameter \(\lambda\) in the Wolfenstein parameterization. More...
 
Particle leptons [6]
 An array of Particle objects for the leptons. More...
 
double mHl
 The Higgs mass \(m_h\) in GeV. More...
 
double muw
 A matching scale \(\mu_W\) around the weak scale in GeV. More...
 
double Mw_inp
 The mass of the \(W\) boson in GeV used as input for FlagMWinput = TRUE. More...
 
CKM myCKM
 An object of type CKM. More...
 
PMNS myPMNS
 
double Mz
 The mass of the \(Z\) boson in GeV. More...
 
bool requireCKM
 An internal flag to control whether the CKM matrix has to be recomputed. More...
 
bool requireYe
 An internal flag to control whether the charged-lepton Yukawa matrix has to be recomputed. More...
 
bool requireYn
 An internal flag to control whether the neutrino Yukawa matrix has to be recomputed. More...
 
double rhob
 The CKM parameter \(\bar{\rho}\) in the Wolfenstein parameterization. More...
 
double s12
 
double s13
 
double s23
 
Flavour SMFlavour
 An object of type Flavour. More...
 
Matching< StandardModelMatching, StandardModelSMM
 An object of type Matching. More...
 
double Vcb
 \(\vert V_{cb} \vert \) used as an input for FlagWolfenstein = FALSE More...
 
double Vub
 \(\vert V_{ub} \vert \) used as an input for FlagWolfenstein = FALSE More...
 
double Vud
 \(\vert V_{ud} \vert \) used as an input for FlagWolfenstein = FALSE and FlagUseVud = TRUE More...
 
double Vus
 \(\vert V_{us} \vert \) used as an input for FlagWolfenstein = FALSE More...
 
gslpp::matrix< gslpp::complex > Yd
 The Yukawa matrix of the down-type quarks. More...
 
gslpp::matrix< gslpp::complex > Ye
 The Yukawa matrix of the charged leptons. More...
 
gslpp::matrix< gslpp::complex > Yn
 The Yukawa matrix of the neutrinos. More...
 
gslpp::matrix< gslpp::complex > Yu
 The Yukawa matrix of the up-type quarks. More...
 
- Protected Attributes inherited from QCD
double AlsM
 The strong coupling constant at the mass scale MAls, \(\alpha_s(M_{\alpha_s})\). More...
 
double CA
 
double CF
 
bool computemt
 Switch for computing the \(\overline{\mathrm{MS}}\) mass of the top quark. More...
 
double dAdA_NA
 
double dFdA_NA
 
double dFdF_NA
 
bool FlagMpole2MbarNumeric
 A flag to determine whether the pole mass to \(\over \mathrm{MS}\) mass conversion is done numerically. More...
 
bool FlagMtPole
 A flag to determine whether the pole mass of the top quark is used as input. More...
 
double MAls
 The mass scale in GeV at which the strong coupling measurement is provided. More...
 
double mtpole
 The pole mass of the top quark. More...
 
double mub
 The threshold between five- and four-flavour theory in GeV. More...
 
double muc
 The threshold between four- and three-flavour theory in GeV. More...
 
double mut
 The threshold between six- and five-flavour theory in GeV. More...
 
double NA
 
double Nc
 The number of colours. More...
 
bool QCDsuccess =true
 
Particle quarks [6]
 The vector of all SM quarks. More...
 
bool requireYd
 Switch for generating the Yukawa couplings to the down-type quarks. More...
 
bool requireYu
 Switch for generating the Yukawa couplings to the up-type quarks. More...
 
double TF
 
- Protected Attributes inherited from Model
bool isSliced = false
 A boolean set to true if the current istance is a slice of an extended object. More...
 
std::map< std::string, std::reference_wrapper< const double > > ModelParamMap
 
bool UpdateError = false
 A boolean set to false if update is successful. More...
 

Private Attributes

gslpp::vector< gslpp::complex > C_bd
 
gslpp::vector< gslpp::complex > C_bs
 The complex Wilson Coefficients. More...
 
gslpp::vector< gslpp::complex > C_c
 
gslpp::vector< gslpp::complex > C_s
 
double imC1_bd
 
double imC1_bs
 
double imC1_c
 
double imC1_s
 
double imC2_bd
 
double imC2_bs
 
double imC2_c
 
double imC2_s
 
double imC3_bd
 
double imC3_bs
 
double imC3_c
 
double imC3_s
 
double imC4_bd
 
double imC4_bs
 
double imC4_c
 
double imC4_s
 
double imC5_bd
 The imaginary parts of the Wilson Coefficients. More...
 
double imC5_bs
 The imaginary parts of the Wilson Coefficients. More...
 
double imC5_c
 The imaginary parts of the Wilson Coefficients. More...
 
double imC5_s
 The imaginary parts of the Wilson Coefficients. More...
 
double reC1_bd
 
double reC1_bs
 
double reC1_c
 
double reC1_s
 
double reC2_bd
 
double reC2_bs
 
double reC2_c
 
double reC2_s
 
double reC3_bd
 
double reC3_bs
 
double reC3_c
 
double reC3_s
 
double reC4_bd
 
double reC4_bs
 
double reC4_c
 
double reC4_s
 
double reC5_bd
 The real parts of the Wilson Coefficients. More...
 
double reC5_bs
 The real parts of the Wilson Coefficients. More...
 
double reC5_c
 The real parts of the Wilson Coefficients. More...
 
double reC5_s
 The real parts of the Wilson Coefficients. More...
 
double WCscale_bd
 
double WCscale_bs
 The scale of the Wilson Coefficients. More...
 
double WCscale_c
 
double WCscale_s
 

Additional Inherited Members

- Public Types inherited from StandardModel
enum  LEP2RCs { Weak = 0 , WeakBox , ISR , QEDFSR , QCDFSR , NUMofLEP2RCs }
 
enum  orders_EW { EW1 = 0 , EW1QCD1 , EW1QCD2 , EW2 , EW2QCD1 , EW3 , orders_EW_size }
 An enumerated type representing perturbative orders of radiative corrections to EW precision observables. More...
 
- Public Types inherited from QCD
enum  lepton { NEUTRINO_1 , ELECTRON , NEUTRINO_2 , MU , NEUTRINO_3 , TAU , NOLEPTON }
 An enum type for leptons. More...
 
enum  meson { P_0 , P_P , K_0 , K_P , D_0 , D_P , D_S , B_D , B_P , B_S , B_C , PHI , K_star , K_star_P , K_S , D_star_P , RHO , RHO_P , OMEGA , MESON_END }
 An enum type for mesons. More...
 
enum  quark { UP , DOWN , CHARM , STRANGE , TOP , BOTTOM }
 An enum type for quarks. More...
 

Constructor & Destructor Documentation

◆ FlavourWilsonCoefficient_DF2()

FlavourWilsonCoefficient_DF2::FlavourWilsonCoefficient_DF2 ( )

FlavourWilsonCoefficient constructor.

Definition at line 17 of file FlavourWilsonCoefficient_DF2.cpp.

17 : StandardModel(),
18 FWCM(*this), C_s(5,0.), C_c(5,0.), C_bd(5,0.), C_bs(5,0.) {
19
20 SMM.setObj((StandardModelMatching&) FWCM.getObj());
21 ModelParamMap.insert(std::make_pair("reC1_s", std::cref(reC1_s)));
22 ModelParamMap.insert(std::make_pair("reC2_s", std::cref(reC2_s)));
23 ModelParamMap.insert(std::make_pair("reC3_s", std::cref(reC3_s)));
24 ModelParamMap.insert(std::make_pair("reC4_s", std::cref(reC4_s)));
25 ModelParamMap.insert(std::make_pair("reC5_s", std::cref(reC5_s)));
26 ModelParamMap.insert(std::make_pair("imC1_s", std::cref(imC1_s)));
27 ModelParamMap.insert(std::make_pair("imC2_s", std::cref(imC2_s)));
28 ModelParamMap.insert(std::make_pair("imC3_s", std::cref(imC3_s)));
29 ModelParamMap.insert(std::make_pair("imC4_s", std::cref(imC4_s)));
30 ModelParamMap.insert(std::make_pair("imC5_s", std::cref(imC5_s)));
31 ModelParamMap.insert(std::make_pair("WCscale_s", std::cref(WCscale_s)));
32 ModelParamMap.insert(std::make_pair("reC1_c", std::cref(reC1_c)));
33 ModelParamMap.insert(std::make_pair("reC2_c", std::cref(reC2_c)));
34 ModelParamMap.insert(std::make_pair("reC3_c", std::cref(reC3_c)));
35 ModelParamMap.insert(std::make_pair("reC4_c", std::cref(reC4_c)));
36 ModelParamMap.insert(std::make_pair("reC5_c", std::cref(reC5_c)));
37 ModelParamMap.insert(std::make_pair("imC1_c", std::cref(imC1_c)));
38 ModelParamMap.insert(std::make_pair("imC2_c", std::cref(imC2_c)));
39 ModelParamMap.insert(std::make_pair("imC3_c", std::cref(imC3_c)));
40 ModelParamMap.insert(std::make_pair("imC4_c", std::cref(imC4_c)));
41 ModelParamMap.insert(std::make_pair("imC5_c", std::cref(imC5_c)));
42 ModelParamMap.insert(std::make_pair("WCscale_c", std::cref(WCscale_c)));
43 ModelParamMap.insert(std::make_pair("reC1_bd", std::cref(reC1_bd)));
44 ModelParamMap.insert(std::make_pair("reC2_bd", std::cref(reC2_bd)));
45 ModelParamMap.insert(std::make_pair("reC3_bd", std::cref(reC3_bd)));
46 ModelParamMap.insert(std::make_pair("reC4_bd", std::cref(reC4_bd)));
47 ModelParamMap.insert(std::make_pair("reC5_bd", std::cref(reC5_bd)));
48 ModelParamMap.insert(std::make_pair("imC1_bd", std::cref(imC1_bd)));
49 ModelParamMap.insert(std::make_pair("imC2_bd", std::cref(imC2_bd)));
50 ModelParamMap.insert(std::make_pair("imC3_bd", std::cref(imC3_bd)));
51 ModelParamMap.insert(std::make_pair("imC4_bd", std::cref(imC4_bd)));
52 ModelParamMap.insert(std::make_pair("imC5_bd", std::cref(imC5_bd)));
53 ModelParamMap.insert(std::make_pair("WCscale_bd", std::cref(WCscale_bd)));
54 ModelParamMap.insert(std::make_pair("reC1_bs", std::cref(reC1_bs)));
55 ModelParamMap.insert(std::make_pair("reC2_bs", std::cref(reC2_bs)));
56 ModelParamMap.insert(std::make_pair("reC3_bs", std::cref(reC3_bs)));
57 ModelParamMap.insert(std::make_pair("reC4_bs", std::cref(reC4_bs)));
58 ModelParamMap.insert(std::make_pair("reC5_bs", std::cref(reC5_bs)));
59 ModelParamMap.insert(std::make_pair("imC1_bs", std::cref(imC1_bs)));
60 ModelParamMap.insert(std::make_pair("imC2_bs", std::cref(imC2_bs)));
61 ModelParamMap.insert(std::make_pair("imC3_bs", std::cref(imC3_bs)));
62 ModelParamMap.insert(std::make_pair("imC4_bs", std::cref(imC4_bs)));
63 ModelParamMap.insert(std::make_pair("imC5_bs", std::cref(imC5_bs)));
64 ModelParamMap.insert(std::make_pair("WCscale_bs", std::cref(WCscale_bs)));
65}
double reC5_bd
The real parts of the Wilson Coefficients.
gslpp::vector< gslpp::complex > C_bd
double imC5_c
The imaginary parts of the Wilson Coefficients.
Matching< FlavourWilsonCoefficient_DF2Matching, FlavourWilsonCoefficient_DF2 > FWCM
The FlavourWilsonCoefficientMatching_DF2 object.
gslpp::vector< gslpp::complex > C_c
double imC5_bd
The imaginary parts of the Wilson Coefficients.
gslpp::vector< gslpp::complex > C_s
double imC5_bs
The imaginary parts of the Wilson Coefficients.
gslpp::vector< gslpp::complex > C_bs
The complex Wilson Coefficients.
double reC5_s
The real parts of the Wilson Coefficients.
double reC5_c
The real parts of the Wilson Coefficients.
double WCscale_bs
The scale of the Wilson Coefficients.
double imC5_s
The imaginary parts of the Wilson Coefficients.
double reC5_bs
The real parts of the Wilson Coefficients.
std::map< std::string, std::reference_wrapper< const double > > ModelParamMap
Definition: Model.h:280
StandardModel()
The default constructor.
Matching< StandardModelMatching, StandardModel > SMM
An object of type Matching.
A class for the matching in the Standard Model.

Member Function Documentation

◆ CheckParameters()

bool FlavourWilsonCoefficient_DF2::CheckParameters ( const std::map< std::string, double > &  DPars)
virtual

A method to check if all the mandatory parameters for FlavourWilsonCoefficient_DF2 have been provided in model initialization.

Parameters
[in]DParsa map of the parameters that are being updated in the Monte Carlo run (including parameters that are varied and those that are held constant)
Returns
a boolean that is true if the execution is successful

Reimplemented from StandardModel.

Definition at line 199 of file FlavourWilsonCoefficient_DF2.cpp.

199 {
200 for (int i = 0; i < NFlavourWilsonCoefficient_DF2vars; i++) {
201 if (DPars.find(FlavourWilsonCoefficient_DF2vars[i]) == DPars.end()) {
202 std::cout << "ERROR: missing mandatory FlavourWilsonCoefficient_DF2 parameter " << FlavourWilsonCoefficient_DF2vars[i] << std::endl;
205 }
206 }
208}
std::map< std::string, double > DPars
Definition: Minimal.cpp:11
static const std::string FlavourWilsonCoefficient_DF2vars[NFlavourWilsonCoefficient_DF2vars]
void addMissingModelParameter(const std::string &missingParameterName)
Definition: Model.h:250
void raiseMissingModelParameterCount()
Definition: Model.h:260
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...

◆ getC_bd()

gslpp::vector< gslpp::complex > FlavourWilsonCoefficient_DF2::getC_bd ( ) const
inline

Definition at line 146 of file FlavourWilsonCoefficient_DF2.h.

147 {
148 return C_bd;
149 }

◆ getC_bs()

gslpp::vector< gslpp::complex > FlavourWilsonCoefficient_DF2::getC_bs ( ) const
inline

Definition at line 151 of file FlavourWilsonCoefficient_DF2.h.

152 {
153 return C_bs;
154 }

◆ getC_c()

gslpp::vector< gslpp::complex > FlavourWilsonCoefficient_DF2::getC_c ( ) const
inline

Definition at line 156 of file FlavourWilsonCoefficient_DF2.h.

157 {
158 return C_c;
159 }

◆ getC_s()

gslpp::vector< gslpp::complex > FlavourWilsonCoefficient_DF2::getC_s ( ) const
inline

Definition at line 161 of file FlavourWilsonCoefficient_DF2.h.

162 {
163 return C_s;
164 }

◆ getMatching()

virtual FlavourWilsonCoefficient_DF2Matching & FlavourWilsonCoefficient_DF2::getMatching ( ) const
inlinevirtual

A get method to access the member reference of type FlavourWilsonCoefficient_DF2Matching.

Returns
a reference to a FlavourWilsonCoefficient_DF2Matching object

Reimplemented from StandardModel.

Definition at line 141 of file FlavourWilsonCoefficient_DF2.h.

142 {
143 return FWCM.getObj();
144 }

◆ getWCscale_bd()

double FlavourWilsonCoefficient_DF2::getWCscale_bd ( ) const
inline
Returns
the scale at which the NP Wilson coefficients in the \(bd\) sector are defined

Definition at line 170 of file FlavourWilsonCoefficient_DF2.h.

171 {
172 return WCscale_bd;
173 }

◆ getWCscale_bs()

double FlavourWilsonCoefficient_DF2::getWCscale_bs ( ) const
inline
Returns
the scale at which the NP Wilson coefficients in the \(bs\) sector are defined

Definition at line 179 of file FlavourWilsonCoefficient_DF2.h.

180 {
181 return WCscale_bs;
182 }

◆ getWCscale_c()

double FlavourWilsonCoefficient_DF2::getWCscale_c ( ) const
inline
Returns
the scale at which the NP Wilson coefficients in the \(cu\) sector are defined

Definition at line 188 of file FlavourWilsonCoefficient_DF2.h.

189 {
190 return WCscale_c;
191 }

◆ getWCscale_s()

double FlavourWilsonCoefficient_DF2::getWCscale_s ( ) const
inline
Returns
the scale at which the NP Wilson coefficients in the \(sd\) sector are defined

Definition at line 197 of file FlavourWilsonCoefficient_DF2.h.

198 {
199 return WCscale_s;
200 }

◆ InitializeModel()

bool FlavourWilsonCoefficient_DF2::InitializeModel ( )
virtual

A method to initialize the model.

This method, called via InputParser::ReadParameters(), allocates memory to the pointers defined in the current class.

Returns
a boolean that is true if model initialization is successful

Reimplemented from StandardModel.

Definition at line 71 of file FlavourWilsonCoefficient_DF2.cpp.

72{
75 return(true);
76}
void setModelFWC_DF2()
Definition: Model.h:186
void setModelInitialized(bool ModelInitialized)
A set method to fix the failure or success of the initialization of the model.
Definition: Model.h:145
virtual bool InitializeModel()
A method to initialize the model.

◆ PostUpdate()

bool FlavourWilsonCoefficient_DF2::PostUpdate ( )
virtual

The post-update method for FlavourWilsonCoefficient_DF2.

This method runs all the procedures that are need to be executed after the model is successfully updated.

Returns
a boolean that is true if the execution is successful

Reimplemented from StandardModel.

Definition at line 78 of file FlavourWilsonCoefficient_DF2.cpp.

79{
80 if(!StandardModel::PostUpdate()) return (false);
81
82 C_s.assign(0,gslpp::complex(reC1_s, imC1_s));
83 C_s.assign(1,gslpp::complex(reC2_s, imC2_s));
84 C_s.assign(2,gslpp::complex(reC3_s, imC3_s));
85 C_s.assign(3,gslpp::complex(reC4_s, imC4_s));
86 C_s.assign(4,gslpp::complex(reC5_s, imC5_s));
87 C_c.assign(0,gslpp::complex(reC1_c, imC1_c));
88 C_c.assign(1,gslpp::complex(reC2_c, imC2_c));
89 C_c.assign(2,gslpp::complex(reC3_c, imC3_c));
90 C_c.assign(3,gslpp::complex(reC4_c, imC4_c));
91 C_c.assign(4,gslpp::complex(reC5_c, imC5_c));
92 C_bd.assign(0,gslpp::complex(reC1_bd, imC1_bd));
93 C_bd.assign(1,gslpp::complex(reC2_bd, imC2_bd));
94 C_bd.assign(2,gslpp::complex(reC3_bd, imC3_bd));
95 C_bd.assign(3,gslpp::complex(reC4_bd, imC4_bd));
96 C_bd.assign(4,gslpp::complex(reC5_bd, imC5_bd));
97 C_bs.assign(0,gslpp::complex(reC1_bs, imC1_bs));
98 C_bs.assign(1,gslpp::complex(reC2_bs, imC2_bs));
99 C_bs.assign(2,gslpp::complex(reC3_bs, imC3_bs));
100 C_bs.assign(3,gslpp::complex(reC4_bs, imC4_bs));
101 C_bs.assign(4,gslpp::complex(reC5_bs, imC5_bs));
102
103 return (true);
104}
virtual bool PostUpdate()
The post-update method for StandardModel.

◆ setParameter()

void FlavourWilsonCoefficient_DF2::setParameter ( const std::string  name,
const double &  value 
)
protectedvirtual

A method to set the value of a parameter of FlavourWilsonCoefficient_DF2.

Parameters
[in]namename of a model parameter
[in]valuethe value to be assigned to the parameter specified by name

Reimplemented from StandardModel.

Definition at line 106 of file FlavourWilsonCoefficient_DF2.cpp.

106 {
107 if(name.compare("reC1_s") == 0)
108 reC1_s = value;
109 else if(name.compare("reC2_s") == 0)
110 reC2_s = value;
111 else if(name.compare("reC3_s") == 0)
112 reC3_s = value;
113 else if(name.compare("reC4_s") == 0)
114 reC4_s = value;
115 else if(name.compare("reC5_s") == 0)
116 reC5_s = value;
117 else if(name.compare("imC1_s") == 0)
118 imC1_s = value;
119 else if(name.compare("imC2_s") == 0)
120 imC2_s = value;
121 else if(name.compare("imC3_s") == 0)
122 imC3_s = value;
123 else if(name.compare("imC4_s") == 0)
124 imC4_s = value;
125 else if(name.compare("imC5_s") == 0)
126 imC5_s = value;
127 else if(name.compare("WCscale_s") == 0)
128 WCscale_s = value;
129 else if(name.compare("reC1_c") == 0)
130 reC1_c = value;
131 else if(name.compare("reC2_c") == 0)
132 reC2_c = value;
133 else if(name.compare("reC3_c") == 0)
134 reC3_c = value;
135 else if(name.compare("reC4_c") == 0)
136 reC4_c = value;
137 else if(name.compare("reC5_c") == 0)
138 reC5_c = value;
139 else if(name.compare("imC1_c") == 0)
140 imC1_c = value;
141 else if(name.compare("imC2_c") == 0)
142 imC2_c = value;
143 else if(name.compare("imC3_c") == 0)
144 imC3_c = value;
145 else if(name.compare("imC4_c") == 0)
146 imC4_c = value;
147 else if(name.compare("imC5_c") == 0)
148 imC5_c = value;
149 else if(name.compare("WCscale_c") == 0)
150 WCscale_c = value;
151 else if(name.compare("reC1_bd") == 0)
152 reC1_bd = value;
153 else if(name.compare("reC2_bd") == 0)
154 reC2_bd = value;
155 else if(name.compare("reC3_bd") == 0)
156 reC3_bd = value;
157 else if(name.compare("reC4_bd") == 0)
158 reC4_bd = value;
159 else if(name.compare("reC5_bd") == 0)
160 reC5_bd = value;
161 else if(name.compare("imC1_bd") == 0)
162 imC1_bd = value;
163 else if(name.compare("imC2_bd") == 0)
164 imC2_bd = value;
165 else if(name.compare("imC3_bd") == 0)
166 imC3_bd = value;
167 else if(name.compare("imC4_bd") == 0)
168 imC4_bd = value;
169 else if(name.compare("imC5_bd") == 0)
170 imC5_bd = value;
171 else if(name.compare("WCscale_bd") == 0)
172 WCscale_bd = value;
173 else if(name.compare("reC1_bs") == 0)
174 reC1_bs = value;
175 else if(name.compare("reC2_bs") == 0)
176 reC2_bs = value;
177 else if(name.compare("reC3_bs") == 0)
178 reC3_bs = value;
179 else if(name.compare("reC4_bs") == 0)
180 reC4_bs = value;
181 else if(name.compare("reC5_bs") == 0)
182 reC5_bs = value;
183 else if(name.compare("imC1_bs") == 0)
184 imC1_bs = value;
185 else if(name.compare("imC2_bs") == 0)
186 imC2_bs = value;
187 else if(name.compare("imC3_bs") == 0)
188 imC3_bs = value;
189 else if(name.compare("imC4_bs") == 0)
190 imC4_bs = value;
191 else if(name.compare("imC5_bs") == 0)
192 imC5_bs = value;
193 else if(name.compare("WCscale_bs") == 0)
194 WCscale_bs = value;
195 else
197}
std::string name
The name of the model.
Definition: Model.h:285
virtual void setParameter(const std::string name, const double &value)
A method to set the value of a parameter of StandardModel.

Member Data Documentation

◆ C_bd

gslpp::vector<gslpp::complex> FlavourWilsonCoefficient_DF2::C_bd
private

Definition at line 221 of file FlavourWilsonCoefficient_DF2.h.

◆ C_bs

gslpp::vector<gslpp::complex> FlavourWilsonCoefficient_DF2::C_bs
private

The complex Wilson Coefficients.

Definition at line 221 of file FlavourWilsonCoefficient_DF2.h.

◆ C_c

gslpp::vector<gslpp::complex> FlavourWilsonCoefficient_DF2::C_c
private

Definition at line 221 of file FlavourWilsonCoefficient_DF2.h.

◆ C_s

gslpp::vector<gslpp::complex> FlavourWilsonCoefficient_DF2::C_s
private

Definition at line 221 of file FlavourWilsonCoefficient_DF2.h.

◆ FlavourWilsonCoefficient_DF2vars

const std::string FlavourWilsonCoefficient_DF2::FlavourWilsonCoefficient_DF2vars
static
Initial value:
= {
"reC1_s","reC2_s","reC3_s","reC4_s","reC5_s",
"reC1_c","reC2_c","reC3_c","reC4_c","reC5_c","reC1_bd","reC2_bd","reC3_bd","reC4_bd","reC5_bd",
"reC1_bs","reC2_bs","reC3_bs","reC4_bs","reC5_bs","imC1_s","imC2_s","imC3_s","imC4_s","imC5_s",
"imC1_c","imC2_c","imC3_c","imC4_c","imC5_c","imC1_bd","imC2_bd","imC3_bd","imC4_bd","imC5_bd",
"imC1_bs","imC2_bs","imC3_bs","imC4_bs","imC5_bs","WCscale_s","WCscale_c","WCscale_bd","WCscale_bs"}

Definition at line 105 of file FlavourWilsonCoefficient_DF2.h.

◆ FWCM

Matching<FlavourWilsonCoefficient_DF2Matching,FlavourWilsonCoefficient_DF2> FlavourWilsonCoefficient_DF2::FWCM
mutableprotected

The FlavourWilsonCoefficientMatching_DF2 object.

Definition at line 210 of file FlavourWilsonCoefficient_DF2.h.

◆ imC1_bd

double FlavourWilsonCoefficient_DF2::imC1_bd
private

Definition at line 219 of file FlavourWilsonCoefficient_DF2.h.

◆ imC1_bs

double FlavourWilsonCoefficient_DF2::imC1_bs
private

Definition at line 220 of file FlavourWilsonCoefficient_DF2.h.

◆ imC1_c

double FlavourWilsonCoefficient_DF2::imC1_c
private

Definition at line 218 of file FlavourWilsonCoefficient_DF2.h.

◆ imC1_s

double FlavourWilsonCoefficient_DF2::imC1_s
private

Definition at line 217 of file FlavourWilsonCoefficient_DF2.h.

◆ imC2_bd

double FlavourWilsonCoefficient_DF2::imC2_bd
private

Definition at line 219 of file FlavourWilsonCoefficient_DF2.h.

◆ imC2_bs

double FlavourWilsonCoefficient_DF2::imC2_bs
private

Definition at line 220 of file FlavourWilsonCoefficient_DF2.h.

◆ imC2_c

double FlavourWilsonCoefficient_DF2::imC2_c
private

Definition at line 218 of file FlavourWilsonCoefficient_DF2.h.

◆ imC2_s

double FlavourWilsonCoefficient_DF2::imC2_s
private

Definition at line 217 of file FlavourWilsonCoefficient_DF2.h.

◆ imC3_bd

double FlavourWilsonCoefficient_DF2::imC3_bd
private

Definition at line 219 of file FlavourWilsonCoefficient_DF2.h.

◆ imC3_bs

double FlavourWilsonCoefficient_DF2::imC3_bs
private

Definition at line 220 of file FlavourWilsonCoefficient_DF2.h.

◆ imC3_c

double FlavourWilsonCoefficient_DF2::imC3_c
private

Definition at line 218 of file FlavourWilsonCoefficient_DF2.h.

◆ imC3_s

double FlavourWilsonCoefficient_DF2::imC3_s
private

Definition at line 217 of file FlavourWilsonCoefficient_DF2.h.

◆ imC4_bd

double FlavourWilsonCoefficient_DF2::imC4_bd
private

Definition at line 219 of file FlavourWilsonCoefficient_DF2.h.

◆ imC4_bs

double FlavourWilsonCoefficient_DF2::imC4_bs
private

Definition at line 220 of file FlavourWilsonCoefficient_DF2.h.

◆ imC4_c

double FlavourWilsonCoefficient_DF2::imC4_c
private

Definition at line 218 of file FlavourWilsonCoefficient_DF2.h.

◆ imC4_s

double FlavourWilsonCoefficient_DF2::imC4_s
private

Definition at line 217 of file FlavourWilsonCoefficient_DF2.h.

◆ imC5_bd

double FlavourWilsonCoefficient_DF2::imC5_bd
private

The imaginary parts of the Wilson Coefficients.

Definition at line 219 of file FlavourWilsonCoefficient_DF2.h.

◆ imC5_bs

double FlavourWilsonCoefficient_DF2::imC5_bs
private

The imaginary parts of the Wilson Coefficients.

Definition at line 220 of file FlavourWilsonCoefficient_DF2.h.

◆ imC5_c

double FlavourWilsonCoefficient_DF2::imC5_c
private

The imaginary parts of the Wilson Coefficients.

Definition at line 218 of file FlavourWilsonCoefficient_DF2.h.

◆ imC5_s

double FlavourWilsonCoefficient_DF2::imC5_s
private

The imaginary parts of the Wilson Coefficients.

Definition at line 217 of file FlavourWilsonCoefficient_DF2.h.

◆ NFlavourWilsonCoefficient_DF2vars

const int FlavourWilsonCoefficient_DF2::NFlavourWilsonCoefficient_DF2vars = 44
static

Definition at line 103 of file FlavourWilsonCoefficient_DF2.h.

◆ reC1_bd

double FlavourWilsonCoefficient_DF2::reC1_bd
private

Definition at line 215 of file FlavourWilsonCoefficient_DF2.h.

◆ reC1_bs

double FlavourWilsonCoefficient_DF2::reC1_bs
private

Definition at line 216 of file FlavourWilsonCoefficient_DF2.h.

◆ reC1_c

double FlavourWilsonCoefficient_DF2::reC1_c
private

Definition at line 214 of file FlavourWilsonCoefficient_DF2.h.

◆ reC1_s

double FlavourWilsonCoefficient_DF2::reC1_s
private

Definition at line 213 of file FlavourWilsonCoefficient_DF2.h.

◆ reC2_bd

double FlavourWilsonCoefficient_DF2::reC2_bd
private

Definition at line 215 of file FlavourWilsonCoefficient_DF2.h.

◆ reC2_bs

double FlavourWilsonCoefficient_DF2::reC2_bs
private

Definition at line 216 of file FlavourWilsonCoefficient_DF2.h.

◆ reC2_c

double FlavourWilsonCoefficient_DF2::reC2_c
private

Definition at line 214 of file FlavourWilsonCoefficient_DF2.h.

◆ reC2_s

double FlavourWilsonCoefficient_DF2::reC2_s
private

Definition at line 213 of file FlavourWilsonCoefficient_DF2.h.

◆ reC3_bd

double FlavourWilsonCoefficient_DF2::reC3_bd
private

Definition at line 215 of file FlavourWilsonCoefficient_DF2.h.

◆ reC3_bs

double FlavourWilsonCoefficient_DF2::reC3_bs
private

Definition at line 216 of file FlavourWilsonCoefficient_DF2.h.

◆ reC3_c

double FlavourWilsonCoefficient_DF2::reC3_c
private

Definition at line 214 of file FlavourWilsonCoefficient_DF2.h.

◆ reC3_s

double FlavourWilsonCoefficient_DF2::reC3_s
private

Definition at line 213 of file FlavourWilsonCoefficient_DF2.h.

◆ reC4_bd

double FlavourWilsonCoefficient_DF2::reC4_bd
private

Definition at line 215 of file FlavourWilsonCoefficient_DF2.h.

◆ reC4_bs

double FlavourWilsonCoefficient_DF2::reC4_bs
private

Definition at line 216 of file FlavourWilsonCoefficient_DF2.h.

◆ reC4_c

double FlavourWilsonCoefficient_DF2::reC4_c
private

Definition at line 214 of file FlavourWilsonCoefficient_DF2.h.

◆ reC4_s

double FlavourWilsonCoefficient_DF2::reC4_s
private

Definition at line 213 of file FlavourWilsonCoefficient_DF2.h.

◆ reC5_bd

double FlavourWilsonCoefficient_DF2::reC5_bd
private

The real parts of the Wilson Coefficients.

Definition at line 215 of file FlavourWilsonCoefficient_DF2.h.

◆ reC5_bs

double FlavourWilsonCoefficient_DF2::reC5_bs
private

The real parts of the Wilson Coefficients.

Definition at line 216 of file FlavourWilsonCoefficient_DF2.h.

◆ reC5_c

double FlavourWilsonCoefficient_DF2::reC5_c
private

The real parts of the Wilson Coefficients.

Definition at line 214 of file FlavourWilsonCoefficient_DF2.h.

◆ reC5_s

double FlavourWilsonCoefficient_DF2::reC5_s
private

The real parts of the Wilson Coefficients.

Definition at line 213 of file FlavourWilsonCoefficient_DF2.h.

◆ WCscale_bd

double FlavourWilsonCoefficient_DF2::WCscale_bd
private

Definition at line 222 of file FlavourWilsonCoefficient_DF2.h.

◆ WCscale_bs

double FlavourWilsonCoefficient_DF2::WCscale_bs
private

The scale of the Wilson Coefficients.

Definition at line 222 of file FlavourWilsonCoefficient_DF2.h.

◆ WCscale_c

double FlavourWilsonCoefficient_DF2::WCscale_c
private

Definition at line 222 of file FlavourWilsonCoefficient_DF2.h.

◆ WCscale_s

double FlavourWilsonCoefficient_DF2::WCscale_s
private

Definition at line 222 of file FlavourWilsonCoefficient_DF2.h.


The documentation for this class was generated from the following files: