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

Model for WEFT LFV contributions to \(\Delta F=1\) processes like \( b\to s\) decays. More...

#include <RealWeakEFTLFV.h>

+ Inheritance diagram for RealWeakEFTLFV:

Detailed Description

Model for WEFT LFV contributions to \(\Delta F=1\) processes like \( b\to s\) decays.

Author
HEPfit Collaboration

Model parameters

The model parameters of RealWeakEFTLFV model are summarized below:

Label LaTeX symbol Description
C7, C7p \(C_7\) and \(C_7^\prime\) The Wilson coefficient of the dipole operator \(O_7\) and the chirality flipped \(O_7^\prime\)
C8, C8p \(C_8\) and \(C_8^\prime\) The Wilson coefficient of the chromomagnetic operator \(O_8\) and the chirality flipped \(O_8^\prime\)
C9_11, C9p_11 \(C_{9,e}\) and \(C_{9,e}^\prime\) The Wilson coefficient of the semileptonic operator \(O_{9,e}\) and the chirality flipped \(O_{9,e}^\prime\) coupling to electrons
C10_11, C10p_11 \(C_{10,e}\) and \(C_{10,e}^\prime\) The Wilson coefficient of the semileptonic operator \(O_{10,e}\) and the chirality flipped \(O_{10,e}^\prime\) coupling to electrons
CS_11, CSp_11 \(C_{S,e}\) and \(C_{S,e}^\prime\) The Wilson coefficient of the scalar operator \(O_{S,e}\) and the chirality flipped \(O_{S,e}^\prime\) coupling to electrons
CP_11, CPp_11 \(C_{P,e}\) and \(C_{P,e}^\prime\) The Wilson coefficient of the pseudo-scalar operator \(O_{P,e}\) and the chirality flipped \(O_{P,e}^\prime\) coupling to electrons
C9_22, C9p_22 \(C_{9,\mu}\) and \(C_{9,\mu}^\prime\) The Wilson coefficient of the semileptonic operator \(O_{9,\mu}\) and the chirality flipped \(O_{9,\mu}^\prime\) coupling to muons
C10_22, C10p_22 \(C_{10,\mu}\) and \(C_{10,\mu}^\prime\) The Wilson coefficient of the semileptonic operator \(O_{10,\mu}\) and the chirality flipped \(O_{10,\mu}^\prime\) coupling to muons
CS_22, CSp_22 \(C_{S,\mu}\) and \(C_{S,\mu}^\prime\) The Wilson coefficient of the scalar operator \(O_{S,\mu}\) and the chirality flipped \(O_{S,\mu}^\prime\) coupling to muons
CP_22, CPp_22 \(C_{P,\mu}\) and \(C_{P,\mu}^\prime\) The Wilson coefficient of the pseudo-scalar operator \(O_{P,\mu}\) and the chirality flipped \(O_{P,\mu}^\prime\) coupling to muons

Definition at line 97 of file RealWeakEFTLFV.h.

Public Member Functions

virtual bool CheckParameters (const std::map< std::string, double > &DPars)
 A method to check if all the mandatory parameters for RealWeakEFTLFV have been provided in model initialization. More...
 
double getC7 () const
 
double getC7p () const
 
double getC8 () const
 
virtual RealWeakEFTLFVMatchinggetMatching () const
 A get method to access the member reference of type RealWeakEFTLFVMatching. More...
 
virtual bool Init (const std::map< std::string, double > &DPars)
 Initializes the RealWeakEFTLFV parameters found in the argument. More...
 
virtual bool InitializeModel ()
 The post-update method for RealWeakEFTLFV. More...
 
virtual bool PostUpdate ()
 The post-update method for RealWeakEFTLFV. More...
 
virtual bool PreUpdate ()
 The pre-update method for RealWeakEFTLFV. More...
 
 RealWeakEFTLFV ()
 RealWeakEFTLFV constructor. More...
 
virtual bool setFlag (const std::string name, const bool value)
 A method to set a flag of RealWeakEFTLFV. More...
 
virtual bool Update (const std::map< std::string, double > &DPars)
 The update method for RealWeakEFTLFV. More...
 
 ~RealWeakEFTLFV ()
 RealWeakEFTLFV destructor. More...
 
- Public Member Functions inherited from StandardModel
virtual 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 double AFB (const Particle f) 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...
 
double ale_OS (const double mu, orders order=FULLNLO) const
 The running electromagnetic coupling \(\alpha(\mu)\) in the on-shell scheme. More...
 
double alphaMz () const
 The electromagnetic coupling at the \(Z\)-mass scale, \(\alpha(M_Z^2)=\alpha/(1-\Delta\alpha(M_Z^2))\). More...
 
double Als (double mu, orders order=FULLNLO, bool qed_flag=false, bool Nf_thr=true) const
 The running QCD coupling \(\alpha(\mu)\) in the \(\overline{MS}\) scheme including QED corrections. More...
 
double AlsByOrder (double mu, orders order=FULLNLO, bool qed_flag=false, bool Nf_thr=true) 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...
 
double Beta_e (int nm, unsigned int nf) const
 QED beta function coefficients - eq. (36) hep-ph/0512066. More...
 
double Beta_s (int nm, unsigned int nf) const
 QCD beta function coefficients including QED corrections - eq. (36) hep-ph/0512066. More...
 
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...
 
double computeBrHtobb () const
 The Br \((H\to bb)\) in the Standard Model. More...
 
double computeBrHtocc () const
 The Br \((H\to cc)\) in the Standard Model. More...
 
double computeBrHtogaga () const
 The Br \((H\to\gamma\gamma)\) in the Standard Model. More...
 
double computeBrHtogg () const
 The Br \((H\to gg)\) in the Standard Model. More...
 
double computeBrHtomumu () const
 The Br \((H\to \mu\mu)\) in the Standard Model. More...
 
double computeBrHtoss () const
 The Br \((H\to ss)\) in the Standard Model. More...
 
double computeBrHtotautau () const
 The Br \((H\to \tau\tau)\) in the Standard Model. More...
 
double computeBrHtoWW () const
 The Br \((H\to WW)\) in the Standard Model. More...
 
double computeBrHtoZga () const
 The Br \((H\to Z\gamma)\) in the Standard Model. More...
 
double computeBrHtoZZ () const
 The Br \((H\to ZZ)\) in the Standard Model. More...
 
double computeBrHtoZZinv () const
 The Br \((H\to ZZ \to inv)\) 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...
 
double computeGammaHgaga_tt () const
 The top loop contribution to \(H\to\gamma\gamma\) in the Standard Model. More...
 
double computeGammaHgaga_tW () const
 The mixed \(t-W\) loop contribution to \(H\to\gamma\gamma\) in the Standard Model. More...
 
double computeGammaHgaga_WW () const
 The \(W\) loop contribution to \(H\to\gamma\gamma\) in the Standard Model. More...
 
double computeGammaHgg_bb () const
 The bottom loop contribution to \(H\to gg\) in the Standard Model. More...
 
double computeGammaHgg_tb () const
 The top-bottom interference contribution to \(H\to gg\) in the Standard Model. More...
 
double computeGammaHgg_tt () const
 The top loop contribution to \(H\to gg\) in the Standard Model. More...
 
double computeGammaHTotal () const
 The Higgs total width in the Standard Model. More...
 
double computeGammaHZga_tt () const
 The top loop contribution to \(H\to Z\gamma\) in the Standard Model. More...
 
double computeGammaHZga_tW () const
 The mixed \(t-W\) loop contribution to \(H\to Z\gamma\) in the Standard Model. More...
 
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. [137]. More...
 
double computeSigmaggH (const double sqrt_s) const
 The ggH cross section in the Standard Model. More...
 
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...
 
double computeSigmaggH_tb (const double sqrt_s) const
 The top-bottom interference contribution to the ggH cross section in the Standard Model. More...
 
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...
 
double computeSigmattH (const double sqrt_s) const
 The ttH production cross section in the Standard Model. More...
 
double computeSigmaVBF (const double sqrt_s) const
 The VBF cross section in the Standard Model. More...
 
double computeSigmaWF (const double sqrt_s) const
 The W fusion contribution \(\sigma_{WF}\) to higgs-production cross section in the Standard Model. More...
 
double computeSigmaWH (const double sqrt_s) const
 The WH production cross section in the Standard Model. More...
 
double computeSigmaZF (const double sqrt_s) const
 The Z fusion contribution \(\sigma_{ZF}\) to higgs-production cross section in the Standard Model. More...
 
double computeSigmaZH (const double sqrt_s) const
 The ZH production cross section in the Standard Model. More...
 
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 double cW2 () const
 
virtual 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...
 
double DeltaAlpha () const
 The total corrections to the electromagnetic coupling \(\alpha\) at the \(Z\)-mass scale, denoted as \(\Delta\alpha(M_Z^2)\). More...
 
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...
 
double DeltaAlphaLepton (const double s) const
 Leptonic contribution to the electromagnetic coupling \(\alpha\), denoted as \(\Delta\alpha_{\mathrm{lept}}(s)\). More...
 
double DeltaAlphaTop (const double s) const
 Top-quark contribution to the electromagnetic coupling \(\alpha\), denoted as \(\Delta\alpha_{\mathrm{top}}(s)\). More...
 
virtual 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 double DeltaR () const
 The SM prediction for \(\Delta r\) derived from that for the \(W\) boson mass. More...
 
virtual double DeltaRbar () const
 The SM prediction for \(\Delta \overline{r}\) derived from that for the \(W\)-boson mass. More...
 
virtual 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 double epsilon1 () const
 The SM contribution to the epsilon parameter \(\varepsilon_1\). More...
 
virtual double epsilon2 () const
 The SM contribution to the epsilon parameter \(\varepsilon_2\). More...
 
virtual double epsilon3 () const
 The SM contribution to the epsilon parameter \(\varepsilon_3\). More...
 
virtual double epsilonb () const
 The SM contribution to the epsilon parameter \(\varepsilon_b\). More...
 
virtual gslpp::complex gA_f (const Particle f) const
 The effective leptonic neutral-current axial-vector coupling \(g_A^l\) in the SM. More...
 
virtual double Gamma_had () const
 The hadronic decay width of the \(Z\) boson, \(\Gamma_{h}\). More...
 
virtual double Gamma_inv () const
 The invisible partial decay width of the \(Z\) boson, \(\Gamma_{\mathrm{inv}}\). More...
 
virtual double Gamma_Z () const
 The total decay width of the \(Z\) boson, \(\Gamma_Z\). More...
 
virtual double GammaW () const
 The total width of the \(W\) boson, \(\Gamma_W\). More...
 
virtual 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 double GammaZ (const Particle f) const
 The \(Z\to \ell\bar{\ell}\) partial decay width, \(\Gamma_\ell\). More...
 
double getAle () const
 A get method to retrieve the fine-structure constant \(\alpha\). More...
 
double getAlsMz () const
 A get method to access the value of \(\alpha_s(M_Z)\). More...
 
virtual double getCBd () const
 The ratio of the absolute value of the $B_d$ mixing amplitude over the Standard Model value. More...
 
virtual double getCBs () const
 The ratio of the absolute value of the $B_s$ mixing amplitude over the Standard Model value. More...
 
virtual double getCCC1 () const
 A virtual implementation for the RealWeakEFTCC class. More...
 
virtual double getCCC2 () const
 A virtual implementation for the RealWeakEFTCC class. More...
 
virtual double getCCC3 () const
 A virtual implementation for the RealWeakEFTCC class. More...
 
virtual double getCCC4 () const
 A virtual implementation for the RealWeakEFTCC class. More...
 
virtual double getCCC5 () const
 A virtual implementation for the RealWeakEFTCC class. More...
 
virtual double getCDMK () const
 The ratio of the real part of the $K$ mixing amplitude over the Standard Model value. More...
 
virtual double getCepsK () const
 The ratio of the imaginary part of the $K$ mixing amplitude over the Standard Model value. More...
 
CKM getCKM () const
 A get method to retrieve the member object of type CKM. More...
 
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...
 
double getDelGammaZ () const
 A get method to retrieve the theoretical uncertainty in \(\Gamma_Z\), denoted as \(\delta\,\Gamma_Z\). More...
 
double getDelMw () const
 A get method to retrieve the theoretical uncertainty in \(M_W\), denoted as \(\delta\,M_W\). More...
 
double getDelR0b () const
 A get method to retrieve the theoretical uncertainty in \(R_b^0\), denoted as \(\delta\,R_b^0\). More...
 
double getDelR0c () const
 A get method to retrieve the theoretical uncertainty in \(R_c^0\), denoted as \(\delta\,R_c^0\). More...
 
double getDelR0l () const
 A get method to retrieve the theoretical uncertainty in \(R_l^0\), denoted as \(\delta\,R_l^0\). More...
 
double getDelSigma0H () const
 A get method to retrieve the theoretical uncertainty in \(\sigma_{Hadron}^0\), denoted as \(\delta\,\sigma_{Hadron}^0\). More...
 
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...
 
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...
 
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...
 
std::string getFlagKappaZ () const
 A method to retrieve the model flag KappaZ. More...
 
std::string getFlagMw () const
 A method to retrieve the model flag Mw. More...
 
std::string getFlagRhoZ () const
 A method to retrieve the model flag RhoZ. More...
 
const FlavourgetFlavour () const
 
double getGF () const
 A get method to retrieve the Fermi constant \(G_\mu\). More...
 
int getIterationNo () const
 
Particle getLeptons (const QCD::lepton p) const
 A get method to retrieve the member object of a lepton. More...
 
virtual double getMHl () const
 A get method to retrieve the Higgs mass \(m_h\). More...
 
virtual double getmq (const QCD::quark q, const double mu) const
 
double getMuw () const
 A get method to retrieve the matching scale \(\mu_W\) around the weak scale. 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
 
EWSMTwoLoopEWgetMyTwoLoopEW () const
 
EWSMTwoLoopQCDgetMyTwoLoopQCD () const
 
double getMz () const
 A get method to access the mass of the \(Z\) boson \(M_Z\). More...
 
virtual double getPhiBd () const
 Half the relative phase of the $B_d$ mixing amplitude w.r.t. the Standard Model one. More...
 
virtual double getPhiBs () const
 Half the relative phase of the $B_s$ mixing amplitude w.r.t. the Standard Model one. More...
 
virtual StandardModel getTrueSM () const
 
gslpp::matrix< gslpp::complexgetUPMNS () const
 A get method to retrieve the object of the PMNS matrix. More...
 
gslpp::matrix< gslpp::complexgetVCKM () const
 A get method to retrieve the CKM matrix. More...
 
gslpp::matrix< gslpp::complexgetYd () const
 A get method to retrieve the Yukawa matrix of the down-type quarks, \(Y_d\). More...
 
gslpp::matrix< gslpp::complexgetYe () const
 A get method to retrieve the Yukawa matrix of the charged leptons, \(Y_e\). More...
 
gslpp::matrix< gslpp::complexgetYn () const
 A get method to retrieve the Yukawa matrix of the neutrinos, \(Y_\nu\). More...
 
gslpp::matrix< gslpp::complexgetYu () const
 A get method to retrieve the Yukawa matrix of the up-type quarks, \(Y_u\). More...
 
virtual gslpp::complex gV_f (const Particle f) const
 The effective leptonic neutral-current vector coupling \(g_V^l\) in the SM. More...
 
bool IsFlagNoApproximateGammaZ () const
 A method to retrieve the model flag NoApproximateGammaZ. More...
 
bool IsFlagWithoutNonUniversalVC () const
 A method to retrieve the model flag WithoutNonUniversalVC. More...
 
virtual gslpp::complex kappaZ_f (const Particle f) const
 The effective leptonic neutral-current coupling \(\kappa_Z^l\) in the SM. More...
 
virtual double Mw () const
 The SM prediction for the \(W\)-boson mass in the on-shell scheme, \(M_{W,\mathrm{SM}}\). More...
 
virtual double Mw_tree () const
 The tree-level mass of the \(W\) boson, \(M_W^{\mathrm{tree}}\). More...
 
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...
 
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 double N_nu () const
 The number of neutrinos obtained indirectly from the measurements at the Z pole, \(N_{\nu}\). More...
 
virtual double R0_f (const Particle f) const
 The ratio \(R_\ell^0=\Gamma(Z\to {\rm hadrons})/\Gamma(Z\to \ell^+ \ell^-)\). More...
 
virtual double R_inv () const
 The ratio of the invisible and leptonic (electron) decay widths of the \(Z\) boson, \(R_{inv}\). More...
 
virtual 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 gslpp::complex rhoZ_f (const Particle f) const
 The effective leptonic neutral-current coupling \(\rho_Z^l\) in the SM. More...
 
double s02 () const
 The square of the sine of the weak mixing angle \(s_0^2\) defined without weak radiative corrections. More...
 
void setFlagCacheInStandardModel (bool FlagCacheInStandardModel)
 A set method to change the model flag CacheInStandardModel of StandardModel. More...
 
void setFlagNoApproximateGammaZ (bool FlagNoApproximateGammaZ)
 
virtual bool setFlagStr (const std::string name, const std::string value)
 A method to set a flag of StandardModel. More...
 
virtual 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 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...
 
double sW2 () const
 
virtual 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...
 
virtual double v () const
 The Higgs vacuum expectation value. More...
 
virtual ~StandardModel ()
 The default destructor. More...
 
- Public Member Functions inherited from QCD
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...
 
virtual double Als (const double mu, const orders order=FULLNLO, bool Nf_thr=true) 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...
 
virtual double AlsByOrder (const double mu, const orders order=FULLNLO, bool Nf_thr=true) 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...
 
double AlsWithInit (const double mu, const double alsi, const double mu_i, 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...
 
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...
 
double BelowTh (const double mu) const
 The active flavour threshold below the scale \(\mu\) as defined in QCD::Thresholds(). More...
 
double Beta0 (const double nf) const
 The \(\beta_0(n_f)\) coefficient for a certain number of flavours \(n_f\). More...
 
double Beta1 (const double nf) const
 The \(\beta_1(n_f)\) coefficient for a certain number of flavours \(n_f\). More...
 
double Beta2 (const double nf) const
 The \(\beta_2(n_f)\) coefficient for a certain number of flavours \(n_f\). More...
 
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
 
orders FullOrder (orders order) const
 Return the FULLORDER enum corresponding to order. More...
 
double Gamma0 (const double nf) const
 The \(\gamma_0\) coefficient used to compute the running of a mass. More...
 
double Gamma1 (const double nf) const
 The \(\gamma_1\) coefficient used to compute the running of a mass. More...
 
double Gamma2 (const double nf) const
 The \(\gamma_2\) coefficient used to compute the running of a mass. More...
 
double getAlsM () const
 A get method to access the value of \(\alpha_s(M_{\alpha_s})\). More...
 
BParameter getBBd () 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...
 
BParameter getBBs () 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...
 
BParameter getBD () 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...
 
BParameter getBK () 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...
 
BParameter getBKd1 () const
 
BParameter getBKd3 () const
 
double getCF () const
 A get method to access the Casimir factor of QCD. More...
 
double getMAls () const
 A get method to access the mass scale \(M_{\alpha_s}\) at which the strong coupling constant measurement is provided. More...
 
Meson getMesons (const QCD::meson m) const
 A get method to access a meson as an object of the type Meson. More...
 
double getMtpole () const
 A get method to access the pole mass of the top quark. More...
 
double getMub () const
 A get method to access the threshold between five- and four-flavour theory in GeV. More...
 
double getMuc () const
 A get method to access the threshold between four- and three-flavour theory in GeV. More...
 
double getMut () const
 A get method to access the threshold between six- and five-flavour theory in GeV. More...
 
double getNc () const
 A get method to access the number of colours \(N_c\). More...
 
double getOptionalParameter (std::string name) const
 A method to get parameters that are specific to only one set of observables. More...
 
Particle getQuarks (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...
 
double logLambda (const double nf, orders order) const
 Computes \(\ln\Lambda_\mathrm{QCD}\) with nf flavours in GeV. More...
 
double Mbar2Mp (const double mbar, const orders order=FULLNNLO) const
 Converts the \(\overline{\mathrm{MS}}\) mass \(m(m)\) to the pole mass. More...
 
double Mp2Mbar (const double mp, const orders order=FULLNNLO) const
 Converts a quark pole mass to the corresponding \(\overline{\mathrm{MS}}\) mass \(m(m)\). More...
 
double Mrun (const double mu, const double m, const orders order=FULLNNLO) const
 Computes a running quark mass \(m(\mu)\) from \(m(m)\). More...
 
double Mrun (const double mu_f, const double mu_i, const double m, const orders order=FULLNNLO) const
 Runs a quark mass from \(\mu_i\) to \(\mu_f\). More...
 
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...
 
double MS2DRqmass (const double MSbar) const
 Converts a quark mass from the \(\overline{\mathrm{MS}}\) scheme to the \(\overline{\mathrm{DR}}\) scheme. More...
 
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...
 
double Nf (const double mu) const
 The number of active flavour at scale \(\mu\). More...
 
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...
 
std::string orderToString (const orders order) const
 Converts an object of the enum type "orders" to the corresponding string. More...
 
 QCD ()
 Constructor. 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...
 
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 isModelGeneralTHDM () const
 
bool isModelGeorgiMachacek () const
 
bool IsModelInitialized () const
 A method to check if the model is initialized. More...
 
bool isModelLinearized () 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 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 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 int NRealWeakEFTLFVvars = 21
 
static const std::string RealWeakEFTLFVvars [NRealWeakEFTLFVvars]
 
- 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...
 

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, B_D, B_P, B_S, B_C, PHI, K_star, K_star_P, 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...
 
- Protected Member Functions inherited from StandardModel
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 Yukawa matrices. 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...
 
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...
 
virtual void setParameter (const std::string name, const double &value)
 A method to set the value of a parameter of StandardModel. More...
 
double taub () const
 Top-mass corrections to the \(Zb\bar{b}\) vertex, denoted by \(\tau_b\). More...
 
- Protected Member Functions inherited from QCD
double MassOfNf (int nf) const
 The Mbar mass of the heaviest quark in the theory with Nf active flavour. 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...
 
double dAle5Mz
 The five-flavour hadronic contribution to the electromagnetic coupling, \(\Delta\alpha_{\mathrm{had}}^{(5)}(M_Z^2)\). 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...
 
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...
 
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 Vus
 \(\vert V_{us} \vert \) used as an input for FlagWolfenstein = FALSE More...
 
gslpp::matrix< gslpp::complexYd
 The Yukawa matrix of the down-type quarks. More...
 
gslpp::matrix< gslpp::complexYe
 The Yukawa matrix of the charged leptons. More...
 
gslpp::matrix< gslpp::complexYn
 The Yukawa matrix of the neutrinos. More...
 
gslpp::matrix< gslpp::complexYu
 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
 
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...
 
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
 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
 A boolean set to false if update is successful. More...
 

Constructor & Destructor Documentation

◆ RealWeakEFTLFV()

RealWeakEFTLFV::RealWeakEFTLFV ( )

RealWeakEFTLFV constructor.

Definition at line 15 of file RealWeakEFTLFV.cpp.

15  : StandardModel(), ReWEFTM(*this) {
16 
17  SMM.setObj((StandardModelMatching&) ReWEFTM.getObj());
18  ModelParamMap.insert(std::make_pair("C7", std::cref(C7)));
19  ModelParamMap.insert(std::make_pair("C7p", std::cref(C7p)));
20  ModelParamMap.insert(std::make_pair("C8", std::cref(C8)));
21  ModelParamMap.insert(std::make_pair("C8p", std::cref(C8p)));
22 
23  ModelParamMap.insert(std::make_pair("C9_11", std::cref(C9_11)));
24  ModelParamMap.insert(std::make_pair("C9p_11", std::cref(C9p_11)));
25  ModelParamMap.insert(std::make_pair("C10_11", std::cref(C10_11)));
26  ModelParamMap.insert(std::make_pair("C10p_11", std::cref(C10p_11)));
27  ModelParamMap.insert(std::make_pair("CS_11", std::cref(CS_11)));
28  ModelParamMap.insert(std::make_pair("CSp_11", std::cref(CSp_11)));
29  ModelParamMap.insert(std::make_pair("CP_11", std::cref(CP_11)));
30  ModelParamMap.insert(std::make_pair("CPp_11", std::cref(CPp_11)));
31 
32  ModelParamMap.insert(std::make_pair("C9_22", std::cref(C9_22)));
33  ModelParamMap.insert(std::make_pair("C9p_22", std::cref(C9p_22)));
34  ModelParamMap.insert(std::make_pair("C10_22", std::cref(C10_22)));
35  ModelParamMap.insert(std::make_pair("C10p_22", std::cref(C10p_22)));
36  ModelParamMap.insert(std::make_pair("CS_22", std::cref(CS_22)));
37  ModelParamMap.insert(std::make_pair("CSp_22", std::cref(CSp_22)));
38  ModelParamMap.insert(std::make_pair("CP_22", std::cref(CP_22)));
39  ModelParamMap.insert(std::make_pair("CPp_22", std::cref(CPp_22)));
40 
41  ModelParamMap.insert(std::make_pair("WCscale", std::cref(WCscale)));
42 }

◆ ~RealWeakEFTLFV()

RealWeakEFTLFV::~RealWeakEFTLFV ( )

RealWeakEFTLFV destructor.

Definition at line 44 of file RealWeakEFTLFV.cpp.

44  {
45  if (IsModelInitialized()) {
46  }
47 }

Member Function Documentation

◆ CheckParameters()

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

A method to check if all the mandatory parameters for RealWeakEFTLFV 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 143 of file RealWeakEFTLFV.cpp.

143  {
144  for (int i = 0; i < NRealWeakEFTLFVvars; i++) {
145  if (DPars.find(RealWeakEFTLFVvars[i]) == DPars.end()) {
146  std::cout << "ERROR: missing mandatory RealWeakEFTLFV parameter " << RealWeakEFTLFVvars[i] << std::endl;
149  }
150  }
151  return(StandardModel::CheckParameters(DPars));
152 }

◆ getC7()

double RealWeakEFTLFV::getC7 ( ) const
inline
Returns
\( C_7\)

Definition at line 180 of file RealWeakEFTLFV.h.

180  {
181  return C7;
182  }

◆ getC7p()

double RealWeakEFTLFV::getC7p ( ) const
inline
Returns
\( C_7'\)

Definition at line 188 of file RealWeakEFTLFV.h.

188  {
189  return C7p;
190  }

◆ getC8()

double RealWeakEFTLFV::getC8 ( ) const
inline
Returns
\( C_8\)

Definition at line 196 of file RealWeakEFTLFV.h.

196  {
197  return C8;
198  }

◆ getMatching()

virtual RealWeakEFTLFVMatching& RealWeakEFTLFV::getMatching ( ) const
inlinevirtual

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

Returns
a reference to a RealWeakEFTLFVMatching object

Reimplemented from StandardModel.

Definition at line 171 of file RealWeakEFTLFV.h.

172  {
173  return ReWEFTM.getObj();
174  }

◆ Init()

bool RealWeakEFTLFV::Init ( const std::map< std::string, double > &  DPars)
virtual

Initializes the RealWeakEFTLFV parameters found in the argument.

Parameters
[in]DParsa map containing the parameters (all as double) to be used in Monte Carlo

Reimplemented from StandardModel.

Definition at line 58 of file RealWeakEFTLFV.cpp.

58  {
59  return(StandardModel::Init(DPars));
60 }

◆ InitializeModel()

bool RealWeakEFTLFV::InitializeModel ( )
virtual

The post-update method for RealWeakEFTLFV.

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 52 of file RealWeakEFTLFV.cpp.

53 {
55  return(true);
56 }

◆ PostUpdate()

bool RealWeakEFTLFV::PostUpdate ( )
virtual

The post-update method for RealWeakEFTLFV.

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 85 of file RealWeakEFTLFV.cpp.

86 {
87  if(!StandardModel::PostUpdate()) return (false);
88 
89  /* Necessary for updating StandardModel parameters in StandardModelMatching,
90  * and FlavourWC and FlavourWC-derived parameters in FlavourWCMatching */
91  ReWEFTM.getObj().updateRealWeakEFTLFVParameters();
92 
93  return (true);
94 }

◆ PreUpdate()

bool RealWeakEFTLFV::PreUpdate ( )
virtual

The pre-update method for RealWeakEFTLFV.

Returns
a boolean that is true if the execution is successful

Reimplemented from StandardModel.

Definition at line 62 of file RealWeakEFTLFV.cpp.

63 {
64  if(!StandardModel::PreUpdate()) return (false);
65 
66  return (true);
67 }

◆ setFlag()

bool RealWeakEFTLFV::setFlag ( const std::string  name,
const bool  value 
)
virtual

A method to set a flag of RealWeakEFTLFV.

Parameters
[in]namename of a model flag
[in]valuethe boolean to be assigned to the flag specified by name
Returns
a boolean that is true if the execution is successful

Reimplemented from StandardModel.

Definition at line 157 of file RealWeakEFTLFV.cpp.

158 {
159  return StandardModel::setFlag(name,value);
160 }

◆ Update()

bool RealWeakEFTLFV::Update ( const std::map< std::string, double > &  DPars)
virtual

The update method for RealWeakEFTLFV.

This method updates all the model parameters with given DPars.

Parameters
[in]DParsa map of the parameters that are being updated in the Monte Carlo run
Returns
a boolean that is true if the execution is successful

Reimplemented from StandardModel.

Definition at line 69 of file RealWeakEFTLFV.cpp.

69  {
70 
71  if(!PreUpdate()) return (false);
72 
73  UpdateError = false;
74 
75  for (std::map<std::string, double>::const_iterator it = DPars.begin(); it != DPars.end(); it++)
76  setParameter(it->first, it->second);
77 
78  if (UpdateError) return (false);
79 
80  if(!PostUpdate()) return (false);
81 
82  return (true);
83 }

Member Data Documentation

◆ NRealWeakEFTLFVvars

const int RealWeakEFTLFV::NRealWeakEFTLFVvars = 21
static

Definition at line 100 of file RealWeakEFTLFV.h.

◆ RealWeakEFTLFVvars

const std::string RealWeakEFTLFV::RealWeakEFTLFVvars
static
Initial value:
= {"C7", "C7p", "C8", "C8p",
"C9_11", "C9p_11", "C10_11", "C10p_11", "CS_11", "CSp_11", "CP_11", "CPp_11",
"C9_22", "C9p_22", "C10_22", "C10p_22", "CS_22", "CSp_22", "CP_22", "CPp_22",
"WCscale"}

Definition at line 102 of file RealWeakEFTLFV.h.


The documentation for this class was generated from the following files:
StandardModel::setParameter
virtual void setParameter(const std::string name, const double &value)
A method to set the value of a parameter of StandardModel.
Definition: StandardModel.cpp:231
RealWeakEFTLFV::PreUpdate
virtual bool PreUpdate()
The pre-update method for RealWeakEFTLFV.
Definition: RealWeakEFTLFV.cpp:62
Model::addMissingModelParameter
void addMissingModelParameter(const std::string &missingParameterName)
Definition: Model.h:232
RealWeakEFTLFV::RealWeakEFTLFVvars
static const std::string RealWeakEFTLFVvars[NRealWeakEFTLFVvars]
Definition: RealWeakEFTLFV.h:102
StandardModel::CheckParameters
virtual bool CheckParameters(const std::map< std::string, double > &DPars)
A method to check if all the mandatory parameters for StandardModel have been provided in model initi...
Definition: StandardModel.cpp:313
RealWeakEFTLFV::PostUpdate
virtual bool PostUpdate()
The post-update method for RealWeakEFTLFV.
Definition: RealWeakEFTLFV.cpp:85
Matching::setObj
void setObj(T &obji)
Definition: Matching.h:15
Model::UpdateError
bool UpdateError
A boolean set to false if update is successful.
Definition: Model.h:254
StandardModel::SMM
Matching< StandardModelMatching, StandardModel > SMM
An object of type Matching.
Definition: StandardModel.h:2506
RealWeakEFTLFV::NRealWeakEFTLFVvars
static const int NRealWeakEFTLFVvars
Definition: RealWeakEFTLFV.h:100
StandardModel::setFlag
virtual bool setFlag(const std::string name, const bool value)
A method to set a flag of StandardModel.
Definition: StandardModel.cpp:378
Model::ModelParamMap
std::map< std::string, std::reference_wrapper< const double > > ModelParamMap
Definition: Model.h:262
StandardModel::StandardModel
StandardModel()
The default constructor.
Definition: StandardModel.cpp:35
StandardModel::Init
virtual bool Init(const std::map< std::string, double > &DPars)
A method to initialize the model parameters.
Definition: StandardModel.cpp:159
StandardModelMatching
A class for the matching in the Standard Model.
Definition: StandardModelMatching.h:26
StandardModel::PreUpdate
virtual bool PreUpdate()
The pre-update method for StandardModel.
Definition: StandardModel.cpp:172
Model::raiseMissingModelParameterCount
void raiseMissingModelParameterCount()
Definition: Model.h:242
StandardModel::PostUpdate
virtual bool PostUpdate()
The post-update method for StandardModel.
Definition: StandardModel.cpp:199
Model::IsModelInitialized
bool IsModelInitialized() const
A method to check if the model is initialized.
Definition: Model.h:136
StandardModel::InitializeModel
virtual bool InitializeModel()
A method to initialize the model.
Definition: StandardModel.cpp:140
Model::name
std::string name
The name of the model.
Definition: Model.h:267
Model::setModelInitialized
void setModelInitialized(bool ModelInitialized)
A set method to fix the failure or success of the initialization of the model.
Definition: Model.h:145