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

A class for parameters related to QCD, hadrons and quarks. More...

#include <QCD.h>

+ Inheritance diagram for QCD:

Detailed Description

A class for parameters related to QCD, hadrons and quarks.

Author
HEPfit Collaboration

This class is a Model class that assigns and updates parameters related to and derived from QCD. A complete list of parameters in the QCD class can be found below. This class includes, but is not limited to, the running of the strong coupling constant (Full NNLO), running of the quark masses and conversions between pole mass and \(\overline{\mathrm{MS}}\) mass. All hadronization parameters like the bag parameters for the mesons and their decay constants are assigned and updated by this class.

Initialization

The constructor QCD() sets the charge and isospin of the quarks. It also sets the mass scale of the light quarks UP, DOWN and STRANGE to 2 GeV. The cache is initialized too along with the computation of \(\zeta(2)\) and \(\zeta(3)\).

The initializations and updates of the model parameters and flags are explained below.

Model parameters

The model parameters of QCD are summarized below:

Label LaTeX symbol Description
AlsM \(\alpha_s(M_{\alpha_s})\) The strong coupling constant at the scale \(M_{\alpha_s}\).
MAls \(M_{\alpha_s}\) The mass scale in GeV at which the strong coupling constant measurement is provided. Must be in the 5 flavour energy range.
mup \(m_{u}(2\,\mathrm{GeV})\) The \(\overline{\mathrm{MS}}\) mass of the up quark at 2 GeV, \(m_u(2\,\mathrm{GeV})\), in GeV.
mdown \(m_{d}(2\,\mathrm{GeV})\) The \(\overline{\mathrm{MS}}\) mass of the down quark at 2 GeV, \(m_d(2\,\mathrm{GeV})\), in GeV.
mcharm \(m_{c}(m_c)\) The \(\overline{\mathrm{MS}}\) scale-invariant mass of the charm quark, \(m_c(m_c)\), in GeV.
mstrange \(m_{s}(2\,\mathrm{GeV})\) The \(\overline{\mathrm{MS}}\) mass of the strange quark at 2 GeV , \(m_s(2\,\mathrm{GeV})\), in GeV.
mtop \(m_{t}\) The pole mass of the top quark in GeV.
mbottom \(m_{b}(m_b)\) The \(\overline{\mathrm{MS}}\) scale-invariant mass of the bottom quark, \(m_b(m_b)\), in GeV.
muc \(\mu_c\) The threshold between three- and four-flavour theory in GeV.
mub \(\mu_b\) The threshold between four- and five-flavour theory in GeV.
mut \(\mu_t\) The threshold between five- and six-flavour theory in GeV.
Attention
The parameters AlsM and MAls are not used in StandardModel and the model classes inherited from it.

The set of the model parameters are initialized and updated with the methods Init() and Update(), respectively, where the former calls the latter. In Update(), the methods PreUpdate() and PostUpdate() are called to run all the procedures that are need to be executed before and after the model parameters are updated. The \(\overline{\mathrm{MS}}\) mass for the top quark is computed and the scale set in PostUpdate() with the updated parameters. Inside the Update() method, the individual model parameter is assigned with the protected member function setParameter().

Computation of the strong coupling constant \(\alpha_s\):

The strong coupling constant \(\alpha_s\) at an arbitrary scale can be computed with the member functions:

  • AlsWithInit(const double mu, const double alsi, const double mu_i, const orders order),
  • AlsWithLambda(const double mu, const orders order),

where another function

  • Als(const double mu, const orders order = FULLNLO)

calls AlsWithInit() for order=LO/FULLNLO(NLO), and AlsWithLambda() for order=FULLNNLO(NNLO).

The function AlsWithInit() computes \(\alpha_s(\mu)\) with a given initial value \(\alpha_s(\mu_i)\):

\[ \alpha_s(\mu)=\frac{\alpha_s(\mu_i)}{v(\mu)},\qquad \alpha_s(\mu)=\frac{\alpha_s(\mu_i)}{v(\mu)} \left[1-\frac{\beta_1}{\beta_0}\frac{\alpha_s(\mu_i)}{4\pi} \frac{\ln v(\mu)}{v(\mu)}\right], \]

at LO and FULLNLO (NLO) respectively, where

\[ v(\mu) = 1- \beta_0\frac{\alpha_s(\mu_i)}{2\pi}\ln\frac{\mu_i}{\mu}, \]

and the one-loop and two-loop beta functions are given by

\[ \beta_0 = \frac{11N_c-2N_f}{3},\qquad \beta_1 = \frac{34}{3}N_c^2-\frac{10}{3}N_cN_f-2C_FN_f. \]

The function AlsWithLambda() computes \(\alpha_s(\mu)\) with the use of \(\Lambda_{\rm QCD}\), where the value of \(\Lambda_{\rm QCD}\) for \(N_f=5\) is derived from \(\alpha_s(M_{\alpha_s})\) by solving the equation (see e.g., [63], which follows the convention in [49] for \(\ln(\mu^2/\Lambda^2)\)):

\[ \frac{\alpha_s(\mu)}{4\pi}=\frac{1}{\beta_0L}-\frac{\beta_1\ln L}{\beta_0^3L^2}+ \frac{1}{\left(\beta_0L\right)^3}\left[\frac{\beta_1^2}{\beta_0^2}\left(\ln^2L-\ln L-1\right)+ \frac{\beta_2}{\beta_0}\right]+O\left(\frac{\ln^3L}{L^4}\right), \]

where \(L\equiv\ln(\mu^2/\Lambda^2)\), and the three-loop beta function is given by

\[ \beta_2 = \frac{2857}{54}N_c^3+C_F^2N_f-\frac{205}{18}C_FN_cN_f -\frac{1415}{54}N_c^2N_f+\frac{11}{9}C_FN_f^2 + \frac{79}{54}N_cN_f^2. \]

For \(N_f < 5\), \(\Lambda_{QCD}\) can be obtained by solving the following matching condition at \(\mu\) [63] [59] :

\begin{eqnarray} \beta_0^{'}\ln\frac{\Lambda^{'2}}{\Lambda^2} &=& (\beta_0^{'}-\beta_0)L+\left(\frac{\beta_1^{'}}{\beta_0^{'}}- \frac{\beta_1}{\beta_0}\right)\ln L-\frac{\beta_1^{'}}{\beta_0^{'}}\ln\frac{\beta_0^{'}}{\beta_0} - C_1 \\ &&+ \frac{1}{\beta_0L}\left[\frac{\beta_1}{\beta_0}\left(\frac{\beta_1^{'}}{\beta_0^{'}}-\frac{\beta_1}{\beta_0} \right)\ln L + \frac{\beta_1^{'2}}{\beta_0^{'2}}-\frac{\beta_1^2}{\beta_0^2}-\frac{\beta_2^{'2}}{\beta_0^{'2}} +\frac{\beta_2^2}{\beta_0^2}+\frac{\beta_1^{'}}{\beta_0^{'}}C_1-C_1^2-C_2\right]+ O\left(\frac{\ln^2L}{L^2}\right), \end{eqnarray}

where the primed (unprimed) quantities refer to those pertaining to \(N_f-1\) \((N_f)\). The terms \(C_1\) and \(C_2\) are given by

\[ C_1 = \frac{2}{3}\ln\frac{\mu^2}{\mu_f^2},\qquad C_2 = -16\left(\frac{1}{36}\ln^2\frac{\mu^2}{\mu_f^2}-\frac{19}{24}\ln\frac{\mu^2}{\mu_f^2}+\frac{11}{72}\right). \]

where \(\mu_f=m_f(\mu_f)\) denotes the \(\overline{\mathrm{MS}}\) invariant mass of the \(N_f\)-th flavour. Moreover, the matching condition at a flavour threshold is given by

\[ \alpha_s^{(N_f-1)}(\mu) = (\zeta_g^f)^2\alpha_s^{(N_f)}(\mu), \]

where

\[ (\zeta_g^f)^2 = 1+\frac{\alpha_s^{(N_f)}(\mu)}{\pi}\left(-\frac{1}{6}\ln\frac{\mu^2}{\mu_f^2}\right)+ \left(\frac{\alpha_s^{(N_f)}(\mu)}{\pi}\right)^2\left(\frac{1}{36}\ln^2\frac{\mu^2}{\mu_f^2} - \frac{19}{24}\ln\frac{\mu^2}{\mu_f^2} + \frac{11}{72}\right) + O\left(\left(\frac{\alpha_s}{\pi} \right)^3\right). \]

For the top quark mass, the pole mass \(m_t\) is used as an input instead of the \(\overline{\mathrm{MS}}\) invariant mass. Then the \(\overline{\mathrm{MS}}\) invariant mass is computed from the pole mass with the computed value of \(\alpha_s^{(6)}(m_t)\), which in turn is computed from \(\alpha_s^{(5)}(M_Z)\). Hence, the matching condition from \(N_f = 6\) to \(N_f = 5\) has to be free from the \(\overline{\mathrm{MS}}\) invariant mass of the top quark. In this case we use the following matching condition

\[ \alpha_s^{(6)}(\mu_t) = \frac{1}{(\zeta_g^{OS,t})^2}\alpha_s^{(5)}(\mu_t), \]

where

\[ \frac{1}{(\zeta_g^{OS,t})^2} = 1+\frac{\alpha_s^{(5)}(\mu)}{\pi}\left(\frac{1}{6}\ln\frac{\mu^2}{m_t^2}\right)+ \left(\frac{\alpha_s^{(5)}(\mu)}{\pi}\right)^2\left(\frac{1}{36}\ln^2\frac{\mu^2}{m_t^2} + \frac{19}{24}\ln\frac{\mu^2}{m_t^2} + \frac{7}{24}\right) + O\left(\left(\frac{\alpha_s}{\pi}\right)^3\right). \]

Besides, \(\Lambda_{QCD}\) for \(N_f=6\) is derived from that for \(N_f=5\) with the relation:

\begin{eqnarray} \beta_0\ln\frac{\Lambda^{2}}{\Lambda^{'2}} &=& (\beta_0-\beta_0^{'})L^{'}+\left(\frac{\beta_1}{\beta_0}- \frac{\beta_1^{'}}{\beta_0^{'}}\right)\ln L^{'}-\frac{\beta_1}{\beta_0}\ln\frac{\beta_0}{\beta_0^{'}}-C_1^{'} \\ && + \frac{1}{\beta_0^{'}L^{'}}\left[\frac{\beta_1^{'}}{\beta_0^{'}}\left(\frac{\beta_1}{\beta_0}-\frac{\beta_1^{'}} {\beta_0^{'}}\right)\ln L^{'}+\frac{\beta_1^{2}}{\beta_0^{2}}-\frac{\beta_1^{'2}}{\beta_0^{'2}}- \frac{\beta_2^{2}}{\beta_0^{2}}+\frac{\beta_2^{'2}}{\beta_0^{'2}}+\frac{\beta_1}{\beta_0}C_1^{'}-C_1^{'2}- C_2^{'}\right]+O\left(\frac{\ln^2L}{L^2}\right), \end{eqnarray}

where

\[ C_1^{'} = -\frac{2}{3}\ln\frac{\mu^2}{m_t^2},\qquad C_2^{'} = -16\left(\frac{1}{36}\ln^2\frac{\mu^2}{m_t^2}+\frac{19}{24}\ln\frac{\mu^2}{m_t^2}+\frac{7}{24}\right). \]

Computation of the running of the quark masses:

In the \(\overline{\mathrm{MS}}\) scheme the quark mass at a scale \(\mu\) is given by (see e.g., [64])

\[ m_q(\mu) = m_q(\mu_0)\left[\frac{\alpha_s(\mu)}{\alpha_s(\mu_0)}\right]^{\frac{\gamma^{(0)}_m}{2\beta_0}} \left\{1+A_1\frac{\alpha_s(\mu)-\alpha_s(\mu_0)}{4\pi}+\frac{A_1^2}{2}\left(\frac{\alpha_s(\mu)-\alpha_s(\mu_0)}{4\pi} \right)^2+\frac{A_2}{2}\left[\left(\frac{\alpha_s(\mu)}{4\pi}\right)^2-\left(\frac{\alpha_s(\mu_0)}{4\pi} \right)^2\right]\right\}, \]

where

\[ A_1 = \frac{\gamma^{(1)}_m}{2\beta_0} - \frac{\beta_1\gamma^{(0)}_m}{2\beta_0^2},\qquad A_2 = \frac{\beta_1^2\gamma^{(0)}_m}{2\beta_0^3}-\frac{\beta_2\gamma^{(0)}_m}{2\beta_0^2}- \frac{\beta_1\gamma^{(1)}_m}{2\beta_0^2}+\frac{\gamma^{(2)}_m}{2\beta_0}, \]

and

\[ \gamma_m^{(0)} = 6C_F,\qquad \gamma_m^{(1)} = C_F\left(3C_F+\frac{97}{3}N_c-\frac{10}{3}N_f\right),\\ \gamma_m^{(2)} = 129C_F^3-\frac{129}{2}C_F^2N_c+\frac{11413}{54}C_FN_c^2+C_F^2N_f(-46+48\zeta(3))+ C_FN_CN_f\left(-\frac{556}{27}-48\zeta(3)\right)-\frac{70}{27}C_FN_f^2, \]

in the \(\overline{\mathrm{MS}}\) scheme. The threshold conditions are given by [59]

\[ m_q^{(N_f-1)}(\mu)=\zeta_m^f\, m_q^{(N_f)}(\mu),\qquad m_q^{(N_f)}(\mu)=\frac{1}{\zeta_m^f}\, m_q^{(N_f-1)}(\mu), \]

where

\[ \zeta^f_m=1+\left(\frac{\alpha_s^{(N_f)}(\mu)}{\pi}\right)^2\left(\frac{1}{12}\ln^2\frac{\mu^2}{\mu_f^2} -\frac{5}{36}\ln\frac{\mu^2}{\mu_f^2}+\frac{89}{432}\right),\\ \frac{1}{\zeta^f_m}=1+\left(\frac{\alpha_s^{(N_f-1)}(\mu)}{\pi}\right)^2\left(-\frac{1}{12}\ln^2\frac{\mu^2} {\mu_f^2}+\frac{5}{36}\ln\frac{\mu^2}{\mu_f^2}-\frac{89}{432}\right), \]

with \(\mu_f=m_f(\mu_f)\).

Pole mass vs. \( \overline{\mathrm{MS}} \) scale invariant mass:

The pole mass \(M_q\) of a heavy quark is related to the \(\overline{\mathrm{MS}}\) scale invariant mass \(m_q(m_q)\) as [59]

\[ m_q(m_q)=M_q\left\{1-\frac{4}{3}\frac{\alpha_s(M_q)}{\pi}+\left[-\frac{2251}{268}-2\zeta(2) -\frac{2}{3}\zeta(2)\ln2+\frac{\zeta(3)}{6}+\frac{n_l}{3}\left(\zeta(2)+\frac{71}{48}\right) -\frac{4}{3}\sum_{1\le i \le n_l}\Delta\left(\frac{M_i}{M_q}\right)\right]\left(\frac{\alpha_s(M_q)}{\pi} \right)^2\right\}, \]

where \(n_l\) is the number of the light quarks, and \(\Delta(x)\) is given by

\[ \Delta(x)=\frac{\pi^2}{8}x-0.597x^2+0.230x^3+O\left(x^4\right). \]

It should be noted that the above formula requires the light quark pole masses which are not well defined. The pole mass in terms of the \(\overline{\mathrm{MS}}\) mass is given by

\[ M_q=m_q(m_q)\left\{1+\frac{4}{3}\frac{\alpha_s(m_q)}{\pi}+\left[\frac{307}{32}+2\zeta(2) +\frac{2}{3}\zeta(2)\ln2-\frac{\zeta(3)}{6}-\frac{n_l}{3}\left(\zeta(2)+\frac{71}{48}\right) +\frac{4}{3}\sum_{1\le i \le n_l}\Delta\left(\frac{m_i}{m_q}\right)\right]\left(\frac{\alpha_s(m_q)}{\pi} \right)^2\right\} \]

To get the pole mass of the light quarks we solve the above equation numerically and derive it from the corresponding \(\overline{\mathrm{MS}}\) mass.

Definition at line 304 of file QCD.h.

Public Types

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...
 

Public Member Functions

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
 
virtual bool CheckFlags () const
 A method to check the sanity of the set of model flags. More...
 
virtual bool CheckParameters (const std::map< std::string, double > &DPars)
 A method to check if all the mandatory parameters for QCD have been provided in model initialization. More...
 
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...
 
virtual bool Init (const std::map< std::string, double > &DPars)
 Initializes the QCD parameters found in the argument. 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...
 
virtual bool PostUpdate ()
 The post-update method for QCD. More...
 
virtual bool PreUpdate ()
 The pre-update method for QCD. More...
 
 QCD ()
 Constructor. More...
 
virtual bool setFlag (const std::string name, const bool value)
 A method to set a flag of QCD. More...
 
virtual bool setFlagStr (const std::string name, const std::string value)
 A method to set a flag of QCD. 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...
 
virtual bool Update (const std::map< std::string, double > &DPars)
 The update method for QCD. 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 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

double MassOfNf (int nf) const
 The Mbar mass of the heaviest quark in the theory with Nf active flavour. More...
 
virtual void setParameter (const std::string name, const double &value)
 A method to set the value of a parameter of QCD. More...
 

Protected Attributes

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...
 

Private Member Functions

double AlsWithLambda (const double mu, const double logLambda, const orders order) const
 The strong coupling constant computed with using \(\Lambda_{\rm QCD}\). More...
 
double logLambda (const double muMatching, const double mf, const double nfNEW, const double nfORG, const double logLambdaORG, orders order) const
 \(\log(\Lambda_{\rm QCD})\) used for computation of \(\alpha_s\) at FULLNNLO. More...
 
double logLambda5 (orders order) const
 \(\log(\Lambda_{\rm QCD})\) for \(n_f = 5\). More...
 
double logLambdaNLO (const double nfNEW, const double nfORG, const double logLambdaORG) const
 \(\log(\Lambda_{\rm QCD})\) used for computation of \(\alpha_s\) at FULLNLO. More...
 
double Mp2MbarTMP (double *mu, double *params) const
 The member used for finding the numerical solution to the pole mass from the \(\overline{\rm MS}\) mass. More...
 
double MrunTMP (const double mu_f, const double mu_i, const double m, const orders order) const
 A function to calculate the running of the mass between flavour thresholds. More...
 
double threCorrForMass (const double nf_f, const double nf_i) const
 The threshold correction for running of a mass when crossing a flavour threshold. More...
 
double ZeroNf3NLO (double *logLambda3, double *logLambda4_in) const
 A member for calculating the difference in \(\alpha_s^{\mathrm{FULLNLO}}\) across the three-four flavour threshold using AlsWithLambda(). More...
 
double ZeroNf4NLO (double *logLambda4, double *logLambda5_in) const
 A member for calculating the difference in \(\alpha_s^{\mathrm{FULLNLO}}\) across the four-five flavour threshold using AlsWithLambda(). More...
 
double ZeroNf5 (double *logLambda5, double *order) const
 A member for calculating the difference in \(\alpha_s\) using AlsWithLambda() and the input value of \(\alpha_s(M_{\alpha_s})\) given as a model parameter. More...
 
double ZeroNf6NLO (double *logLambda6, double *logLambda5_in) const
 A member for calculating the difference in \(\alpha_s^{\mathrm{FULLNLO}}\) across the six-five flavour threshold using AlsWithLambda(). More...
 

Private Attributes

double als_cache [9][CacheSize]
 Cache for \(\alpha_s\). More...
 
std::map< std::string, BParameterBParameterMap
 
bool computeBd
 Switch for computing \(B_{B_d}\) from \(B_{B_s}\). More...
 
bool computeBs
 Switch for computing \(B_{B_s}\) from \(F_{B_s}\sqrt{B_{B_s}}\). More...
 
bool computeFBd
 Switch for computing \(F_{B_d}\) from \(F_{B_s}\). More...
 
bool computeFBp
 Switch for computing \(F_{B^+}\) from \(F_{B_s}\). More...
 
bool FlagCsi
 A flag to determine whether \(B_{B_s}\) and \(B_{B_s}/B_{B_d}\) or \(F_{B_s}\sqrt{B_{B_s}}\) (false) and \(\xi \equiv F_{B_s}\sqrt{B_{B_s}}/(F_{B_d}\sqrt{B_{B_d}})\) (default, true) are used as inputs. More...
 
double logLambda5_cache [4][CacheSize]
 
double logLambdaNLO_cache [9][CacheSize]
 
std::map< const QCD::meson, MesonmesonsMap
 The map of defined mesons. More...
 
double mp2mbar_cache [5][CacheSize]
 Cache for pole mass to msbar mass conversion. More...
 
double mrun_cache [10][CacheSize]
 Cache for running quark mass. More...
 
std::map< std::string, double > optionalParameters
 A map for containing the list and values of the parameters that are used only by a specific set of observables. More...
 
orders realorder
 
std::vector< std::string > unknownParameters
 A vector for containing the names of the parameters that are not being used but specified in the configuration file. More...
 
bool unknownParameterWarning
 A flag to stop the unknown parameter warning after the first time. More...
 
double zeta2
 \(\zeta(2)\) computed with the GSL. More...
 
double zeta3
 \(\zeta(3)\) computed with the GSL. More...
 

Static Private Attributes

static const int CacheSize = 5
 Defines the depth of the cache. More...
 

Member Enumeration Documentation

◆ lepton

An enum type for leptons.

Enumerator
NEUTRINO_1 

The 1st-generation neutrino

ELECTRON 

Electron

NEUTRINO_2 

The 2nd-generation neutrino

MU 

Muon

NEUTRINO_3 

The 3rd-generation neutrino

TAU 

Tau

NOLEPTON 

a lepton when none is needed

Definition at line 310 of file QCD.h.

310  {
311  NEUTRINO_1,
312  ELECTRON,
313  NEUTRINO_2,
314  MU,
315  NEUTRINO_3,
316  TAU,
317  NOLEPTON
318  };

◆ meson

enum QCD::meson

An enum type for mesons.

Enumerator
P_0 

\(\pi^0\) meson

P_P 

\(\pi^\pm\) meson

K_0 

\(K^0\) meson

K_P 

\(K^\pm\) meson

D_0 

\(D^0\) meson

D_P 

\(D^\pm\) meson

B_D 

\(B_d\) meson

B_P 

\(B^\pm\) meson

B_S 

\(B_s\) meson

B_C 

\(B_c\) meson

PHI 

\(\phi\) meson

K_star 

\(K^*\) meson

K_star_P 

\(K^{*,\pm}\) meson

D_star_P 

\(D^{*,\pm}\) meson

RHO 

\(\rho\) meson

RHO_P 

\(\rho\) meson

OMEGA 

\(\omega\) meson

MESON_END 

The size of this enum.

Definition at line 336 of file QCD.h.

336  {
337  P_0,
338  P_P,
339  K_0,
340  K_P,
341  D_0,
342  D_P,
343  B_D,
344  B_P,
345  B_S,
346  B_C,
347  PHI,
348  K_star,
349  K_star_P,
350  D_star_P,
351  RHO,
352  RHO_P,
353  OMEGA,
354  MESON_END
355  };

◆ quark

enum QCD::quark

An enum type for quarks.

Enumerator
UP 

Up quark

DOWN 

Down quark

CHARM 

Charm quark

STRANGE 

Strange quark

TOP 

Top quark

BOTTOM 

Bottom quark

Definition at line 323 of file QCD.h.

323  {
324  UP,
325  DOWN,
326  CHARM,
327  STRANGE,
328  TOP,
329  BOTTOM
330  };

Constructor & Destructor Documentation

◆ QCD()

QCD::QCD ( )

Constructor.

Definition at line 30 of file QCD.cpp.

31 {
32  FlagCsi = true;
33  computeFBd = false;
34  computeFBp = false;
35  computeBd = false;
36  computeBs = false;
37  Nc = 3.;
38  TF = 0.5;
39  CF = Nc / 2. - 1. / (2. * Nc);
40  CA = Nc;
41  dFdA_NA = Nc*(Nc*Nc+6.)/48.;
42  dAdA_NA = Nc*Nc*(Nc*Nc+36.)/24.;
43  dFdF_NA = (Nc*Nc-6.+18./Nc/Nc)/96.;
44  NA = Nc*Nc-1.;
45 
46  // Particle(std::string name, double mass, double mass_scale = 0., double width = 0., double charge = 0.,double isospin = 0.);
47  quarks[UP] = Particle("UP", 0., 2., 0., 2. / 3., .5);
48  quarks[CHARM] = Particle("CHARM", 0., 0., 0., 2. / 3., .5);
49  quarks[TOP] = Particle("TOP", 0., 0., 0., 2. / 3., .5);
50  quarks[DOWN] = Particle("DOWN", 0., 2., 0., -1. / 3., -.5);
51  quarks[STRANGE] = Particle("STRANGE", 0., 2., 0., -1. / 3., -.5);
52  quarks[BOTTOM] = Particle("BOTTOM", 0., 0., 0., -1. / 3., -.5);
53 
56  for (int i = 0; i < CacheSize; i++) {
57  for (int j = 0; j < 8; j++)
58  als_cache[j][i] = 0.;
59  for (int j = 0; j < 4; j++)
60  logLambda5_cache[j][i] = 0.;
61  for (int j = 0; j < 10; j++)
62  mrun_cache[j][i] = 0.;
63  for (int j = 0; j < 5; j++)
64  mp2mbar_cache[j][i] = 0.;
65  }
66 
67  ModelParamMap.insert(std::make_pair("AlsM", std::cref(AlsM)));
68  ModelParamMap.insert(std::make_pair("MAls", std::cref(MAls)));
69  ModelParamMap.insert(std::make_pair("mup", std::cref(quarks[UP].getMass())));
70  ModelParamMap.insert(std::make_pair("mdown", std::cref(quarks[DOWN].getMass())));
71  ModelParamMap.insert(std::make_pair("mcharm", std::cref(quarks[CHARM].getMass())));
72  ModelParamMap.insert(std::make_pair("mstrange", std::cref(quarks[STRANGE].getMass())));
73  ModelParamMap.insert(std::make_pair("mtop", std::cref(mtpole)));
74  ModelParamMap.insert(std::make_pair("mbottom", std::cref(quarks[BOTTOM].getMass())));
75  ModelParamMap.insert(std::make_pair("muc", std::cref(muc)));
76  ModelParamMap.insert(std::make_pair("mub", std::cref(mub)));
77  ModelParamMap.insert(std::make_pair("mut", std::cref(mut)));
78 
80  realorder = LO;
81 }

Member Function Documentation

◆ AboveTh()

double QCD::AboveTh ( const double  mu) const

The active flavour threshold above the scale \(\mu\) as defined in QCD::Thresholds().

Parameters
[in]mua scale \(\mu\) in GeV
Returns
the higher active flavour threshold

Definition at line 420 of file QCD.cpp.

421 {
422  int i;
423  for (i = 4; i >= 0; i--)
424  if (mu < Thresholds(i)) return (Thresholds(i));
425 
426  throw std::runtime_error("Error in QCD::AboveTh()");
427 }

◆ addParameters()

void QCD::addParameters ( std::vector< std::string >  params_i)

A method to add parameters that are specific to only one set of observables.

Parameters
[in]params_ia vector of parameters to be added (including parameters that are varied and those that are held constant)

Definition at line 182 of file QCD.cpp.

183 {
184  for (std::vector<std::string>::iterator it = params_i.begin(); it < params_i.end(); it++) {
185  if (optionalParameters.find(*it) == optionalParameters.end()){
186  optionalParameters[*it] = 0.;
187  ModelParamMap.insert(std::make_pair(*it, std::cref(optionalParameters[*it])));
188  }
189  }
190 }

◆ Als()

double QCD::Als ( const double  mu,
const orders  order = FULLNLO,
bool  Nf_thr = true 
) const
virtual

Definition at line 637 of file QCD.cpp.

638 {
639  switch (order)
640  {
641  case LO:
642  realorder = order;
643  return AlsByOrder(mu, LO, Nf_thr);
644  case FULLNLO:
645  realorder = order;
646  return (AlsByOrder(mu, LO, Nf_thr) + AlsByOrder(mu, NLO, Nf_thr));
647  case FULLNNLO:
648  realorder = order;
649  return (AlsByOrder(mu, LO, Nf_thr) + AlsByOrder(mu, NLO, Nf_thr) + AlsByOrder(mu, NNLO, Nf_thr));
650  case FULLNNNLO:
651  realorder = order;
652  return (AlsByOrder(mu, LO, Nf_thr) + AlsByOrder(mu, NLO, Nf_thr) + AlsByOrder(mu, NNLO, Nf_thr) + AlsByOrder(mu, NNNLO, Nf_thr));
653  default:
654  throw std::runtime_error("QCD::Als(): " + orderToString(order) + " is not implemented.");
655  }
656 }

◆ Als4()

double QCD::Als4 ( const double  mu) const

The value of \(\alpha_s^{\mathrm{FULLNLO}}\) at any scale \(\mu\) with the number of flavours \(n_f = 4\).

Parameters
[in]muthe scale at which \(\alpha_s\) has to be computed
Returns
\(\alpha_s^{\mathrm{FULLNLO}}(\mu)\) with \(n_f = 4\)
Attention
Temporary function waiting for the implementation of NNLO exact.

Definition at line 536 of file QCD.cpp.

537 {
538  double v = 1. - Beta0(4.) * AlsM / 2. / M_PI * log(MAls / mu);
539  return (AlsM / v * (1. - Beta1(4.) / Beta0(4.) * AlsM / 4. / M_PI * log(v) / v));
540 }

◆ AlsByOrder()

double QCD::AlsByOrder ( const double  mu,
const orders  order = FULLNLO,
bool  Nf_thr = true 
) const
virtual

Definition at line 658 of file QCD.cpp.

659 {
660  int i, nfAls = (int) Nf(MAls), nfmu = Nf_thr ? (int) Nf(mu) : nfAls;
661  double als, alstmp, mutmp;
662  orders fullord;
663 
664  for (i = 0; i < CacheSize; ++i)
665  if ((mu == als_cache[0][i]) && ((double) order == als_cache[1][i]) &&
666  (AlsM == als_cache[2][i]) && (MAls == als_cache[3][i]) &&
667  (mut == als_cache[4][i]) && (mub == als_cache[5][i]) &&
668  (muc == als_cache[6][i]) && (Nf_thr == (bool) als_cache[7][i]))
669  return als_cache[8][i];
670 
671  switch (order)
672  {
673  case LO:
674  case NLO:
675  case NNLO:
676  case NNNLO:
677  if (nfAls == nfmu)
678  als = AlsWithInit(mu, AlsM, MAls, order);
679  fullord = FullOrder(order);
680  if (nfAls > nfmu) {
681  mutmp = BelowTh(MAls);
682  alstmp = AlsWithInit(mutmp, AlsM, MAls, realorder);
683  alstmp *= (1. - NfThresholdCorrections(mutmp, MassOfNf(nfAls), alstmp, nfAls, fullord));
684  for (i = nfAls - 1; i > nfmu; i--) {
685  mutmp = BelowTh(mutmp - MEPS);
686  alstmp = AlsWithInit(mutmp, alstmp, AboveTh(mutmp) - MEPS, realorder);
687  alstmp *= (1. - NfThresholdCorrections(mutmp, MassOfNf(i), alstmp, i, fullord));
688  }
689  als = AlsWithInit(mu, alstmp, AboveTh(mu) - MEPS, order);
690  }
691 
692  if (nfAls < nfmu) {
693  mutmp = AboveTh(MAls) - MEPS;
694  alstmp = AlsWithInit(mutmp, AlsM, MAls, realorder);
695  alstmp *= (1. + NfThresholdCorrections(mutmp, MassOfNf(nfAls + 1), alstmp, nfAls + 1, fullord));
696  for (i = nfAls + 1; i < nfmu; i++) {
697  mutmp = AboveTh(mutmp) - MEPS;
698  alstmp = AlsWithInit(mutmp, alstmp, BelowTh(mutmp) + MEPS, realorder);
699  alstmp *= (1. + NfThresholdCorrections(mutmp, MassOfNf(i + 1), alstmp, i + 1, fullord));
700  }
701  als = AlsWithInit(mu, alstmp, BelowTh(mu) + MEPS, order);
702  }
703 
704  CacheShift(als_cache, 9);
705  als_cache[0][0] = mu;
706  als_cache[1][0] = (double) order;
707  als_cache[2][0] = AlsM;
708  als_cache[3][0] = MAls;
709  als_cache[4][0] = mut;
710  als_cache[5][0] = mub;
711  als_cache[6][0] = muc;
712  als_cache[7][0] = (int) Nf_thr;
713  als_cache[8][0] = als;
714 
715  return als;
716  default:
717  throw std::runtime_error("QCD::Als(): " + orderToString(order) + " is not implemented.");
718  }
719 }

◆ AlsOLD()

double QCD::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().

Parameters
[in]muthe scale \(\mu\) in GeV
[in]orderorder in the \(\alpha_s\) expansion as defined in OrderScheme
[in]order_qedorder in the \(\alpha_e\) expansion as defined in OrderScheme. Default to NO_QED.
[in]Nf_thrtrue (default): \(n_f\) = Nf(mu), false: \(n_f\) = Nf(AlsM)
Returns
the strong coupling constant \(\alpha_s(\mu)\) in the \(\overline{\mathrm{MS}}\) scheme

Definition at line 721 of file QCD.cpp.

722 {
723  int i;
724  for (i = 0; i < CacheSize; ++i)
725  if ((mu == als_cache[0][i]) && ((double) order == als_cache[1][i]) &&
726  (AlsM == als_cache[2][i]) && (MAls == als_cache[3][i]) &&
727  (mut == als_cache[4][i]) && (mub == als_cache[5][i]) &&
728  (muc == als_cache[6][i]))
729  return als_cache[7][i];
730 
731  double nfmu = Nf(mu), nfz = Nf(MAls), mu_thre1, mu_thre2, Als_tmp, mf;
732  double als;
733 
734  switch (order) {
735  case LO:
736  case FULLNLO:
737  case NLO:
738  if (nfmu == nfz)
739  als = AlsWithInit(mu, AlsM, MAls, order);
740  else if (nfmu > nfz) {
741  if (order == NLO)
742  throw std::runtime_error("NLO is not implemented in QCD::Als(mu,order).");
743  if (nfmu == nfz + 1.) {
744  mu_thre1 = AboveTh(MAls); // mut
745  Als_tmp = AlsWithInit(mu_thre1 - MEPS, AlsM, MAls, order);
746  if (order == FULLNLO) {
747  mf = getQuarks(TOP).getMass(); //mf = mtpole;
748  Als_tmp = (1. + Als_tmp / M_PI * log(mu_thre1 / mf) / 3.) * Als_tmp;
749  }
750  als = AlsWithInit(mu, Als_tmp, mu_thre1 + MEPS, order);
751  } else
752  throw std::runtime_error("Error in QCD::Als(mu,order)");
753  } else {
754  if (order == NLO)
755  throw std::runtime_error("NLO is not implemented in QCD::Als(mu,order).");
756  if (nfmu == nfz - 1.) {
757  mu_thre1 = BelowTh(MAls); // mub
758  Als_tmp = AlsWithInit(mu_thre1 + MEPS, AlsM, MAls, order);
759  if (order == FULLNLO) {
760  mf = getQuarks(BOTTOM).getMass();
761  Als_tmp = (1. - Als_tmp / M_PI * log(mu_thre1 / mf) / 3.) * Als_tmp;
762  }
763  als = AlsWithInit(mu, Als_tmp, mu_thre1 - MEPS, order);
764  } else if (nfmu == nfz - 2.) {
765  mu_thre1 = BelowTh(MAls); // mub
766  mu_thre2 = AboveTh(mu); // muc
767  Als_tmp = Als(mu_thre1 + MEPS, order);
768  if (order == FULLNLO) {
769  mf = getQuarks(BOTTOM).getMass();
770  Als_tmp = (1. - Als_tmp / M_PI * log(mu_thre1 / mf) / 3.) * Als_tmp;
771  }
772  Als_tmp = AlsWithInit(mu_thre2, Als_tmp, mu_thre1 - MEPS, order);
773  if (order == FULLNLO) {
774  mf = getQuarks(CHARM).getMass();
775  Als_tmp = (1. - Als_tmp / M_PI * log(mu_thre2 / mf) / 3.) * Als_tmp;
776  }
777  als = AlsWithInit(mu, Als_tmp, mu_thre2 - MEPS, order);
778  } else
779  throw std::runtime_error("Error in QCD::Als(mu,order)");
780  }
781  break;
782  case NNLO:
783 // case NNNLO:
784  case FULLNNLO:
785 // case FULLNNNLO:
786  /* alpha_s(mu) computed with Lambda_QCD for Nf=nfmu */
787  als = AlsWithLambda(mu, order);
788  break;
789  default:
790  throw std::runtime_error(orderToString(order) + " is not implemented in QCD::Als(mu,order).");
791  }
792 
793  CacheShift(als_cache, 8);
794  als_cache[0][0] = mu;
795  als_cache[1][0] = (double) order;
796  als_cache[2][0] = AlsM;
797  als_cache[3][0] = MAls;
798  als_cache[4][0] = mut;
799  als_cache[5][0] = mub;
800  als_cache[6][0] = muc;
801  als_cache[7][0] = als;
802 
803  return als;
804 }

◆ AlsWithInit()

double QCD::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\).

Parameters
[in]mua scale \(\mu\) in GeV
[in]alsithe initial value for the coupling at the scale given below
[in]mu_ithe initial scale \(\mu_i\) in GeV
[in]orderLO, NLO or FULLNLO in the \(\alpha_s\) expansion defined in OrderScheme
Returns
the strong coupling constant \(\alpha_s(\mu)\) in the \(\overline{\mathrm{MS}}\) scheme

Definition at line 498 of file QCD.cpp.

500 {
501  double nf = Nf(mu_i);
502 // double nf = Nf(mu);
503 // if (nf != Nf(mu_i))
504 // throw std::runtime_error("Error in QCD::AlsWithInit().");
505 
506  double b1_b0 = Beta1(nf)/Beta0(nf);
507  double v = 1. - Beta0(nf) * alsi / 2. / M_PI * log(mu_i / mu);
508  double logv = log(v);
509 
510  switch (order) {
511  case LO:
512  return (alsi / v);
513  case NLO:
514  return (- alsi * alsi / 4. / M_PI / v / v * b1_b0 * logv );
515  case NNLO:
516  return (alsi * alsi * alsi / 4. / 4. / M_PI /M_PI / v / v / v * (
517  Beta2(nf) / Beta0(nf) * (1. - v) + b1_b0 * b1_b0 * (logv * logv -
518  logv + v - 1.)));
519  case NNNLO:
520  return (alsi * alsi * alsi * alsi / 4. / 4. / 4. / M_PI /M_PI / M_PI /
521  v / v / v / v * ( Beta3(nf) / Beta0(nf) * (1. - v * v) / 2. +
522  b1_b0 * Beta2(nf) / Beta0(nf) * ((2. * v - 3.) * logv + v * v -
523  v) + b1_b0 * b1_b0 * b1_b0 * (-logv * logv * logv + 2.5 *
524  logv * logv + 2. * (1. - v) * logv - 0.5 * (v - 1.) * (v - 1.))));
525  case FULLNLO:
526  return (AlsWithInit(mu,alsi,mu_i,LO) + AlsWithInit(mu,alsi,mu_i,NLO));
527  case FULLNNLO:
528  return(AlsWithInit(mu,alsi,mu_i,LO)+AlsWithInit(mu,alsi,mu_i,NLO)+AlsWithInit(mu,alsi,mu_i,NNLO));
529  case FULLNNNLO:
530  return(AlsWithInit(mu,alsi,mu_i,LO)+AlsWithInit(mu,alsi,mu_i,NLO)+AlsWithInit(mu,alsi,mu_i,NNLO)+AlsWithInit(mu,alsi,mu_i,NNNLO));
531  default:
532  throw std::runtime_error("QCD::AlsWithInit(): " + orderToString(order) + " is not implemented.");
533  }
534 }

◆ AlsWithLambda() [1/2]

double QCD::AlsWithLambda ( const double  mu,
const double  logLambda,
const orders  order 
) const
private

The strong coupling constant computed with using \(\Lambda_{\rm QCD}\).

Parameters
[in]muthe scale of the strong coupling constant
[in]logLambda\(\log(\Lambda_{\rm QCD})\)
[in]orderthe QCD order at which \(\alpha_s\) is required
Returns
\(\alpha_s(\mu)\) for the specified order

Definition at line 542 of file QCD.cpp.

544 {
545  double nf = Nf(mu);
546  double L = 2. * (log(mu) - logLambda);
547 
548  // LO contribution
549  double b0 = Beta0(nf);
550  double b0L = b0*L;
551  double alsLO = 4. * M_PI / b0L;
552  if (order == LO) return alsLO;
553 
554  // NLO contribution
555  double b1 = Beta1(nf);
556  double log_L = log(L);
557  double alsNLO = 4. * M_PI / b0L * (-b1 * log_L / b0 / b0L);
558  if (order == NLO) return alsNLO;
559  if (order == FULLNLO) return (alsLO + alsNLO);
560 
561  // NNLO contribution
562  double b2 = Beta2(nf);
563  double alsNNLO = 4. * M_PI / b0L * (1. / b0L / b0L
564  * (b1 * b1 / b0 / b0 * (log_L * log_L - log_L - 1.) + b2 / b0));
565  if (order == NNLO) return alsNNLO;
566  if (order == FULLNNLO) return (alsLO + alsNLO + alsNNLO);
567 
568  // NNNLO contribution
569  double b3 = Beta3(nf);
570  double alsNNNLO = 4. * M_PI / b0L * (-1. / b0L / b0L / b0L
571  * (b1 * b1 * b1 / b0 / b0 / b0 * (log_L * log_L * log_L - 5./2. * log_L * log_L -2. * log_L - 0.5) + 3. * b1 * b2 * log_L / b0 / b0 - 0.5 * b3 / b0));
572  if (order == NNNLO) return alsNNNLO;
573  if (order == FULLNNNLO) return (alsLO + alsNLO + alsNNLO + alsNNNLO);
574 
575  throw std::runtime_error(orderToString(order) + " is not implemented in QCD::AlsWithLambda().");
576 }

◆ AlsWithLambda() [2/2]

double QCD::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}\).

Parameters
[in]muA scale \(\mu\) in GeV
[in]orderLO, NLO, FULLNLO, NNLO or FULLNNLO in the \(\alpha_s\) expansion defined in OrderScheme
Returns
the strong coupling constant \(\alpha_s(\mu)\) in the \(\overline{\mathrm{MS}}\) scheme

Definition at line 578 of file QCD.cpp.

579 {
580  return AlsWithLambda(mu, logLambda(Nf(mu), order), order);
581 }

◆ BelowTh()

double QCD::BelowTh ( const double  mu) const

The active flavour threshold below the scale \(\mu\) as defined in QCD::Thresholds().

Parameters
[in]mua scale \(\mu\) in GeV
Returns
the lower active flavour threshold

Definition at line 429 of file QCD.cpp.

430 {
431  int i;
432  for (i = 0; i < 5; i++)
433  if (mu >= Thresholds(i)) return (Thresholds(i));
434 
435  throw std::runtime_error("Error in QCD::BelowTh()");
436 }

◆ Beta0()

double QCD::Beta0 ( const double  nf) const

The \(\beta_0(n_f)\) coefficient for a certain number of flavours \(n_f\).

Parameters
[in]nfthe number of active flavours \(n_f\)
Returns
\(\beta_0(n_f)\)

Definition at line 466 of file QCD.cpp.

467 {
468  return ( (11. * CA - 4. * TF * nf) / 3. );
469 }

◆ Beta1()

double QCD::Beta1 ( const double  nf) const

The \(\beta_1(n_f)\) coefficient for a certain number of flavours \(n_f\).

Parameters
[in]nfthe number of active flavours \(n_f\)
Returns
\(\beta_1(n_f)\)

Definition at line 471 of file QCD.cpp.

472 {
473  return ( 34./3. * CA * CA - ( 20./3. * CA + 4. * CF )* TF * nf);
474 }

◆ Beta2()

double QCD::Beta2 ( const double  nf) const

The \(\beta_2(n_f)\) coefficient for a certain number of flavours \(n_f\).

Parameters
[in]nfthe number of active flavours \(n_f\)
Returns
\(\beta_2(n_f)\)

Definition at line 476 of file QCD.cpp.

477 {
478  return ( 2857./54. * CA * CA * CA - (1415./27. * CA * CA + 205./9. * CF * CA -
479  2. * CF * CF) * TF * nf +
480  (158./27. * CA + 44./9. * CF ) * TF * TF * nf * nf );
481 }

◆ Beta3()

double QCD::Beta3 ( const double  nf) const

The \(\beta_3(n_f)\) coefficient for a certain number of flavours \(n_f\).

Parameters
[in]nfthe number of active flavours \(n_f\)
Returns
\(\beta_3(n_f)\)

Definition at line 484 of file QCD.cpp.

485 {
486  return ( CA * CF * TF * TF * nf * nf * (17152./243. + 448./9. * zeta3) +
487  CA * CF * CF * TF * nf * (-4204./27. + 352./9. * zeta3) +
488  424./243. * CA * TF * TF * TF * nf * nf * nf + CA * CA * CF * TF * nf *
489  (7073./243. - 656./9. * zeta3) + CA * CA * TF * TF * nf * nf * (7930./81.+
490  224./9. * zeta3) + 1232./243. * CF * TF * TF * TF * nf * nf * nf +
491  CA * CA * CA * TF * nf * (-39143./81. + 136./3. * zeta3) + CA * CA * CA *
492  CA * (150653./486. - 44./9. * zeta3) + CF * CF * TF * TF * nf * nf *
493  (1352./27. - 704./9. * zeta3) + 46. * CF * CF * CF * TF * nf +
494  nf * dFdA_NA * (512./9. - 1664./3. * zeta3) + nf * nf * dFdF_NA * (
495  -704./9. + 512./3. * zeta3) + dAdA_NA * (-80./9. + 704./3. * zeta3) );
496 }

◆ CacheShift() [1/2]

void QCD::CacheShift ( double  cache[][5],
int  n 
) const

A member used to manage the caching for this class.

Parameters
[in]cachethe cache to be moved
[in]nthe dimension of the cache to be shifted

◆ CacheShift() [2/2]

void QCD::CacheShift ( int  cache[][5],
int  n 
) const

◆ CheckFlags()

bool QCD::CheckFlags ( ) const
virtual

A method to check the sanity of the set of model flags.

Returns
a boolean that is true if the set of model flags is sane

Implements Model.

Reimplemented in StandardModel.

Definition at line 399 of file QCD.cpp.

400 {
401  return (true);
402 }

◆ CheckParameters()

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

A method to check if all the mandatory parameters for QCD 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

Implements Model.

Reimplemented in NPSMEFTd6, StandardModel, HiggsKigen, GeneralTHDM, HiggsChiral, THDMW, HiggsKvgenKfgen, HiggsKvKfgen, THDM, GeorgiMachacek, NPZbbbarLinearized, NPZbbbar, RealWeakEFTLFV, FlavourWilsonCoefficient, NPEpsilons, FlavourWilsonCoefficient_DF2, HiggsKvKf, NPSTUZbbbarLR, NPEpsilons_pureNP, RealWeakEFTCC, NPSTU, CMFV, NPDF2, myModel, and myModel.

Definition at line 335 of file QCD.cpp.

336 {
337  for (int i = 0; i < NQCDvars; i++)
338  if (DPars.find(QCDvars[i]) == DPars.end()) {
339  std::cout << "ERROR: missing mandatory QCD parameter " << QCDvars[i] << std::endl;
342  }
343  for (std::map<std::string, double>::iterator it = optionalParameters.begin(); it != optionalParameters.end(); it++) {
344  if (DPars.find(it->first) == DPars.end()) {
345  std::cout << "ERROR: missing optional parameter " << it->first << std::endl;
347  addMissingModelParameter(it->first);
348  }
349  }
350  if (!BParameterMap.empty()) {
351  for (std::map<std::string, BParameter>::iterator it = BParameterMap.begin(); it != BParameterMap.end(); it++) {
352  std::vector<std::string> parameters = it->second.parameterList(it->first);
353  for (std::vector<std::string>::iterator it1 = parameters.begin(); it1 != parameters.end(); it1++) {
354  if (DPars.find(*it1) == DPars.end()) {
355  std::cout << "ERROR: missing parameter for " << it->first << ": " << *it1 << std::endl;
358  }
359  }
360  }
361  }
362  if (!mesonsMap.empty()) {
363  for (std::map<const QCD::meson, Meson>::iterator it = mesonsMap.begin(); it != mesonsMap.end(); it++) {
364  std::vector<std::string> parameters = it->second.parameterList(it->second.getName());
365  for (std::vector<std::string>::iterator it1 = parameters.begin(); it1 != parameters.end(); it1++) {
366  if (DPars.find(*it1) == DPars.end()) {
367  std::cout << "ERROR: missing parameter for " << mesonsMap.at(it->first).getName() << ": " << *it1 << std::endl;
370  }
371  }
372  }
373  }
374  if (getMissingModelParametersCount() > 0) return false;
375  else return true;
376 }

◆ FullOrder()

orders QCD::FullOrder ( orders  order) const

Return the FULLORDER enum corresponding to order.

Parameters
[in]orderof the expansion in \(\alpha_s\)
Returns
the FULLORDER enum corresponding to order

Definition at line 603 of file QCD.cpp.

604 {
605  switch(order)
606  {
607  case LO:
608  return(LO);
609  case NLO:
610  return(FULLNLO);
611  case NNLO:
612  return(FULLNNLO);
613  case NNNLO:
614  return(FULLNNNLO);
615  default:
616  throw std::runtime_error("QCD::FullOrder(): " + orderToString(order) + " is not implemented.");
617  }
618 }

◆ Gamma0()

double QCD::Gamma0 ( const double  nf) const

The \(\gamma_0\) coefficient used to compute the running of a mass.

Parameters
[in]nfthe number of active flavours \(n_f\)
Returns
the \(\gamma_0\) coefficient.

Definition at line 1008 of file QCD.cpp.

1009 {
1010  return ( 6. * CF);
1011 }

◆ Gamma1()

double QCD::Gamma1 ( const double  nf) const

The \(\gamma_1\) coefficient used to compute the running of a mass.

Parameters
[in]nfthe number of active flavours \(n_f\)
Returns
the \(\gamma_1\) coefficient

Definition at line 1013 of file QCD.cpp.

1014 {
1015  return ( CF * (3. * CF + 97. / 3. * Nc - 10. / 3. * nf));
1016 }

◆ Gamma2()

double QCD::Gamma2 ( const double  nf) const

The \(\gamma_2\) coefficient used to compute the running of a mass.

Parameters
[in]nfThe number of active flavours \(n_f\)
Returns
the \(\gamma_2\) coefficient

Definition at line 1018 of file QCD.cpp.

1019 {
1020  return ( 129. * CF * CF * CF - 129. / 2. * CF * CF * Nc + 11413. / 54. * CF * Nc * Nc
1021  + CF * CF * nf * (-46. + 48. * zeta3) + CF * Nc * nf * (-556. / 27. - 48. * zeta3)
1022  - 70. / 27. * CF * nf * nf);
1023 }

◆ getAlsM()

double QCD::getAlsM ( ) const
inline

A get method to access the value of \(\alpha_s(M_{\alpha_s})\).

Returns
the strong coupling constant at \(M_{\alpha_s}\), \(\alpha_s(M_{\alpha_s})\)

Definition at line 543 of file QCD.h.

544  {
545  return AlsM;
546  }

◆ getBBd()

BParameter QCD::getBBd ( ) const
inline

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.

Returns
the vector of bag parameters

Definition at line 608 of file QCD.h.

609  {
610  return BParameterMap.at("BBd");
611  }

◆ getBBs()

BParameter QCD::getBBs ( ) const
inline

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.

Returns
the vector of bag parameters

Definition at line 619 of file QCD.h.

620  {
621  return BParameterMap.at("BBs");
622  }

◆ getBD()

BParameter QCD::getBD ( ) const
inline

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.

Returns
the vector of bag parameters

Definition at line 630 of file QCD.h.

631  {
632  return BParameterMap.at("BD");
633  }

◆ getBK()

BParameter QCD::getBK ( ) const
inline

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.

Returns
the vector of bag parameters

Definition at line 641 of file QCD.h.

642  {
643  return BParameterMap.at("BK");
644  }

◆ getBKd1()

BParameter QCD::getBKd1 ( ) const
inline
Returns

Definition at line 649 of file QCD.h.

650  {
651  return BParameterMap.at("BKd1");
652  }

◆ getBKd3()

BParameter QCD::getBKd3 ( ) const
inline
Returns

Definition at line 657 of file QCD.h.

658  {
659  return BParameterMap.at("BKd3");
660  }

◆ getCF()

double QCD::getCF ( ) const
inline

A get method to access the Casimir factor of QCD.

Returns
the Casimir factor

Definition at line 597 of file QCD.h.

598  {
599  return CF;
600  }

◆ getMAls()

double QCD::getMAls ( ) const
inline

A get method to access the mass scale \(M_{\alpha_s}\) at which the strong coupling constant measurement is provided.

Returns
the mass scale in GeV \(M_{\alpha_s}\) at which the strong coupling constant measurement is provided

Definition at line 552 of file QCD.h.

553  {
554  return MAls;
555  }

◆ getMesons()

Meson QCD::getMesons ( const QCD::meson  m) const
inline

A get method to access a meson as an object of the type Meson.

Parameters
[in]mthe name of a meson
Returns
the object of the meson specified in the argument

Definition at line 524 of file QCD.h.

525  {
526  return mesonsMap.at(m);
527  }

◆ getMtpole()

double QCD::getMtpole ( ) const
inline

A get method to access the pole mass of the top quark.

Returns
the pole mass of the top quark \(m_t^{pole}\)

Definition at line 588 of file QCD.h.

589  {
590  return mtpole;
591  }

◆ getMub()

double QCD::getMub ( ) const
inline

A get method to access the threshold between five- and four-flavour theory in GeV.

Returns
the threshold \(\mu_b\)

Definition at line 570 of file QCD.h.

571  {
572  return mub;
573  }

◆ getMuc()

double QCD::getMuc ( ) const
inline

A get method to access the threshold between four- and three-flavour theory in GeV.

Returns
the threshold \(\mu_c\)

Definition at line 579 of file QCD.h.

580  {
581  return muc;
582  }

◆ getMut()

double QCD::getMut ( ) const
inline

A get method to access the threshold between six- and five-flavour theory in GeV.

Returns
the threshold \(\mu_t\)

Definition at line 561 of file QCD.h.

562  {
563  return mut;
564  }

◆ getNc()

double QCD::getNc ( ) const
inline

A get method to access the number of colours \(N_c\).

Returns
the number of colours

Definition at line 505 of file QCD.h.

506  {
507  return Nc;
508  }

◆ getOptionalParameter()

double QCD::getOptionalParameter ( std::string  name) const
inline

A method to get parameters that are specific to only one set of observables.

Parameters
[in]namethe name of the parameter
Returns
a double that is the value of the parameter

Definition at line 448 of file QCD.h.

449  {
450  return optionalParameters.at(name);
451  }

◆ getQuarks()

Particle QCD::getQuarks ( const QCD::quark  q) const
inline

A get method to access a quark as an object of the type Particle.

Parameters
[in]qthe name of a quark
Returns
the object of the quark found in the argument

Definition at line 534 of file QCD.h.

535  {
536  return quarks[q];
537  }

◆ getUnknownParameters()

std::vector<std::string> QCD::getUnknownParameters ( )
inline

A method to get the vector of the parameters that have been specified in the configuration file but not being used.

Returns
a vector of strings that contain the names of the parameters

Definition at line 467 of file QCD.h.

468  {
469  return unknownParameters;
470  }

◆ Init()

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

Initializes the QCD parameters found in the argument.

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

Implements Model.

Reimplemented in StandardModel, GeneralTHDM, THDMW, THDM, GeorgiMachacek, RealWeakEFTLFV, FlavourWilsonCoefficient, RealWeakEFTCC, myModel, and myModel.

Definition at line 107 of file QCD.cpp.

108 {
109  bool check = CheckParameters(DPars);
110  if (!check) return (check);
111  check *= Update(DPars);
112  unknownParameterWarning = false;
113  return (check);
114 }

◆ initializeBParameter()

void QCD::initializeBParameter ( std::string  name_i) const

A method to initialize B Parameter and the corresponding meson.

Parameters
[in]name_iname of the B parameters set

Definition at line 192 of file QCD.cpp.

193 {
194  if (BParameterMap.find(name_i) != BParameterMap.end()) return;
195 
196  if (name_i.compare("BBs") == 0 || name_i.compare("BBd") == 0) {
197  BParameterMap.insert(std::pair<std::string, BParameter >("BBs", BParameter(5, "BBs")));
198  BParameterMap.at("BBs").setFlagCsi(FlagCsi);
199  BParameterMap.at("BBs").ModelParameterMapInsert(ModelParamMap);
200  computeBs = true;
202 
203  BParameterMap.insert(std::pair<std::string, BParameter >("BBd", BParameter(5, "BBd")));
204  BParameterMap.at("BBd").setFlagCsi(FlagCsi);
205  BParameterMap.at("BBd").ModelParameterMapInsert(ModelParamMap);
206  computeBd = true;
208  return;
209  }
210  if (name_i.compare("BD") == 0) {
211  BParameterMap.insert(std::pair<std::string, BParameter >(name_i, BParameter(5, name_i)));
212  BParameterMap.at(name_i).ModelParameterMapInsert(ModelParamMap);
214  return;
215  }
216  if (name_i.compare("BK") == 0) {
217  BParameterMap.insert(std::pair<std::string, BParameter >(name_i, BParameter(5, name_i)));
218  BParameterMap.at(name_i).ModelParameterMapInsert(ModelParamMap);
220  return;
221  }
222  if (name_i.compare("BKd1") == 0) {
223  BParameterMap.insert(std::pair<std::string, BParameter >(name_i, BParameter(10, name_i)));
224  BParameterMap.at(name_i).ModelParameterMapInsert(ModelParamMap);
226  return;
227  }
228  if (name_i.compare("BKd3") == 0) {
229  BParameterMap.insert(std::pair<std::string, BParameter >(name_i, BParameter(10, name_i)));
230  BParameterMap.at(name_i).ModelParameterMapInsert(ModelParamMap);
232  return;
233  }
234 }

◆ initializeMeson()

void QCD::initializeMeson ( QCD::meson  meson_i) const

A method to initialize a meson.

Parameters
[in]meson_ithe enumerator corresponding to the meson

Definition at line 236 of file QCD.cpp.

237 {
238  if (mesonsMap.find(meson_i) != mesonsMap.end()) return;
239 
240  mesonsMap.insert(std::pair<const QCD::meson, Meson>(meson_i, Meson()));
241 
242  if (meson_i == QCD::P_0) mesonsMap.at(meson_i).setName("P_0");
243  else if (meson_i == QCD::P_P) mesonsMap.at(meson_i).setName("P_P");
244  else if (meson_i == QCD::K_0) mesonsMap.at(meson_i).setName("K_0");
245  else if (meson_i == QCD::K_P) mesonsMap.at(meson_i).setName("K_P");
246  else if (meson_i == QCD::D_0) mesonsMap.at(meson_i).setName("D_0");
247  else if (meson_i == QCD::D_P) mesonsMap.at(meson_i).setName("D_P");
248  else if (meson_i == QCD::B_D) mesonsMap.at(meson_i).setName("B_D");
249  else if (meson_i == QCD::B_P) mesonsMap.at(meson_i).setName("B_P");
250  else if (meson_i == QCD::B_S) mesonsMap.at(meson_i).setName("B_S");
251  else if (meson_i == QCD::B_C) mesonsMap.at(meson_i).setName("B_C");
252  else if (meson_i == QCD::PHI) mesonsMap.at(meson_i).setName("PHI");
253  else if (meson_i == QCD::K_star) mesonsMap.at(meson_i).setName("K_star");
254  else if (meson_i == QCD::K_star_P) mesonsMap.at(meson_i).setName("K_star_P");
255  else if (meson_i == QCD::D_star_P) mesonsMap.at(meson_i).setName("D_star_P");
256  else if (meson_i == QCD::RHO) mesonsMap.at(meson_i).setName("RHO");
257  else if (meson_i == QCD::RHO_P) mesonsMap.at(meson_i).setName("RHO_P");
258  else if (meson_i == QCD::OMEGA) mesonsMap.at(meson_i).setName("OMEGA");
259  else {
260  std::stringstream out;
261  out << meson_i;
262  throw std::runtime_error("QCD::initializeMeson() meson " + out.str() + " not implemented");
263  }
264 
265  if (meson_i == QCD::B_D) computeFBd = true;
266  if (meson_i == QCD::B_P) computeFBp = true;
267 
269 
270  mesonsMap.at(meson_i).ModelParameterMapInsert(ModelParamMap);
271 }

◆ logLambda() [1/2]

double QCD::logLambda ( const double  muMatching,
const double  mf,
const double  nfNEW,
const double  nfORG,
const double  logLambdaORG,
orders  order 
) const
private

\(\log(\Lambda_{\rm QCD})\) used for computation of \(\alpha_s\) at FULLNNLO.

Parameters
[in]muMatchingthe scale at which the matching is done during crossing a flavour threshold
[in]mfthe mass of the quark sitting at the flavour threshold being crossed
[in]nfNEWthe number of flavours after crossing the flavour threshold, \(n_{f}^{\mathrm{NEW}}\)
[in]nfORGthe number of flavours before crossing the flavour threshold, \(n_{f}^{\mathrm{ORG}}\)
[in]logLambdaORGthe value of \(\log(\Lambda_{\rm QCD})\) with \(n_f = n_{f}^{\mathrm{ORG}}\)
[in]orderthe QCD order of the calculation
Returns
\(\log(\Lambda_{\rm QCD})\) for \(n_f = n_{f}^{\mathrm{NEW}}\)

Definition at line 914 of file QCD.cpp.

917 {
918  if (fabs(nfNEW - nfORG) != 1.)
919  throw std::runtime_error("Error in QCD::logLambda()");
920  if (order == NLO) order = FULLNLO;
921  if (order == NNLO) order = FULLNNLO;
922 
923  /* We do not use the codes below for FULLNLO, since threshold corrections
924  * can be regarded as an NNLO effect as long as setting the matching scale
925  * to be close to the mass scale of the decoupling quark. In order to use
926  * the relation als^{nf+1} = als^{nf} exactly, we use logLambdaNLO method.
927  */
928  if (order == FULLNLO)
929  return logLambdaNLO(nfNEW, nfORG, logLambdaORG);
930 
931  double logMuMatching = log(muMatching);
932  double L = 2. * (logMuMatching - logLambdaORG);
933  double rNEW = 0.0, rORG = 0.0, log_mu2_mf2 = 0.0, log_L = 0.0;
934  double C1 = 0.0, C2 = 0.0; // threshold corrections
935  double logLambdaNEW;
936 
937  // LO contribution
938  logLambdaNEW = 1. / 2. / Beta0(nfNEW)
939  *(Beta0(nfNEW) - Beta0(nfORG)) * L + logLambdaORG;
940 
941  // NLO contribution
942  if (order == FULLNLO || order == FULLNNLO) {
943  rNEW = Beta1(nfNEW) / Beta0(nfNEW);
944  rORG = Beta1(nfORG) / Beta0(nfORG);
945  log_mu2_mf2 = 2. * (logMuMatching - log(mf));
946  log_L = log(L);
947  if (nfNEW < nfORG)
948  C1 = 2. / 3. * log_mu2_mf2;
949  else
950  C1 = -2. / 3. * log_mu2_mf2;
951  logLambdaNEW += 1. / 2. / Beta0(nfNEW)
952  *((rNEW - rORG) * log_L
953  - rNEW * log(Beta0(nfNEW) / Beta0(nfORG)) - C1);
954  }
955 
956  // NNLO contribution
957  if (order == FULLNNLO) {
958  if (nfNEW == 5. && nfORG == 6.)
959  C2 = -16. * (log_mu2_mf2 * log_mu2_mf2 / 36. - 19. / 24. * log_mu2_mf2 - 7. / 24.);
960  else if (nfNEW == 6. && nfORG == 5.)
961  C2 = -16. * (log_mu2_mf2 * log_mu2_mf2 / 36. + 19. / 24. * log_mu2_mf2 + 7. / 24.);
962  else {
963  if (nfNEW < nfORG)
964  C2 = -16. * (log_mu2_mf2 * log_mu2_mf2 / 36. - 19. / 24. * log_mu2_mf2 + 11. / 72.);
965  else
966  C2 = -16. * (log_mu2_mf2 * log_mu2_mf2 / 36. + 19. / 24. * log_mu2_mf2 - 11. / 72.);
967  }
968  logLambdaNEW += 1. / 2. / Beta0(nfNEW) / Beta0(nfORG) / L
969  * (rORG * (rNEW - rORG) * log_L + rNEW * rNEW - rORG * rORG
970  - Beta2(nfNEW) / Beta0(nfNEW) + Beta2(nfORG) / Beta0(nfORG)
971  + rNEW * C1 - C1 * C1 - C2);
972  }
973 
974  return logLambdaNEW;
975 }

◆ logLambda() [2/2]

double QCD::logLambda ( const double  nf,
orders  order 
) const

Computes \(\ln\Lambda_\mathrm{QCD}\) with nf flavours in GeV.

Parameters
[in]nfthe number of active flavours \(n_f\)
[in]orderLO, NLO, FULLNLO, NNLO or FULLNNLO in the \(\alpha_s\) expansion defined in OrderScheme
Returns
\(\ln\Lambda_\mathrm{QCD}\) with nf flavours in GeV

Definition at line 977 of file QCD.cpp.

978 {
979  if (order == NLO) order = FULLNLO;
980  if (order == NNLO) order = FULLNNLO;
981 
982  double muMatching, mf, logLambdaORG, logLambdaNEW;
983  if (nf == 5.)
984  return logLambda5(order);
985  else if (nf == 6.) {
986  muMatching = Thresholds(1); // mut
987  /* matching condition from Nf=5 to Nf=6 is given in terms of the top pole mass. */
988  mf = mtpole; // top pole mass
989  return logLambda(muMatching, mf, 6., 5., logLambda5(order), order);
990  } else if (nf == 4. || nf == 3.) {
991  muMatching = Thresholds(2); // mub
992  mf = getQuarks(BOTTOM).getMass(); // m_b(m_b)
993  logLambdaORG = logLambda5(order);
994  logLambdaNEW = logLambda(muMatching, mf, 4., 5., logLambdaORG, order);
995  if (nf == 3.) {
996  muMatching = Thresholds(3); // muc
997  mf = getQuarks(CHARM).getMass(); // m_c(m_c)
998  logLambdaORG = logLambdaNEW;
999  logLambdaNEW = logLambda(muMatching, mf, 3., 4., logLambdaORG, order);
1000  }
1001  return logLambdaNEW;
1002  } else
1003  throw std::runtime_error("Error in QCD::logLambda()");
1004 }

◆ logLambda5()

double QCD::logLambda5 ( orders  order) const
private

\(\log(\Lambda_{\rm QCD})\) for \(n_f = 5\).

Parameters
[in]orderthe QCD order of the computation
Returns
\(\log(\Lambda_{\rm QCD}^{(5)})\)

Definition at line 829 of file QCD.cpp.

830 {
831  if (order == NLO) order = FULLNLO;
832  if (order == NNLO) order = FULLNNLO;
833 
834  for (int i = 0; i < CacheSize; ++i)
835  if ((AlsM == logLambda5_cache[0][i])
836  && (MAls == logLambda5_cache[1][i])
837  && ((double) order == logLambda5_cache[2][i]))
838  return logLambda5_cache[3][i];
839 
841  logLambda5_cache[0][0] = AlsM;
842  logLambda5_cache[1][0] = MAls;
843  logLambda5_cache[2][0] = (double) order;
844 
845  if (order == LO)
846  logLambda5_cache[3][0] = log(MAls) - 2. * M_PI / Beta0(5.) / AlsM;
847  else {
848  double xmin = -4., xmax = -0.2;
849  TF1 f = TF1("f", this, &QCD::ZeroNf5, xmin, xmax, 1, "QCD", "zeroNf5");
850 
851  ROOT::Math::WrappedTF1 wf1(f);
852  double ledouble = (double) order;
853  wf1.SetParameters(&ledouble);
854 
855  ROOT::Math::BrentRootFinder brf;
856  brf.SetFunction(wf1, xmin, xmax);
857 
858  if (brf.Solve()) logLambda5_cache[3][0] = brf.Root();
859  else
860  throw std::runtime_error("Error in QCD::logLambda5()");
861  }
862  return ( logLambda5_cache[3][0]);
863 }

◆ logLambdaNLO()

double QCD::logLambdaNLO ( const double  nfNEW,
const double  nfORG,
const double  logLambdaORG 
) const
private

\(\log(\Lambda_{\rm QCD})\) used for computation of \(\alpha_s\) at FULLNLO.

Parameters
[in]nfNEWthe number of flavours after crossing the flavour threshold, \(n_{f}^{\mathrm{NEW}}\)
[in]nfORGthe number of flavours before crossing the flavour threshold, \(n_{f}^{\mathrm{ORG}}\)
[in]logLambdaORGthe value of \(\log(\Lambda_{\rm QCD})\) with \(n_f = n_{f}^{\mathrm{ORG}}\)
Returns
\(\log(\Lambda_{\rm QCD})\) for \(n_f = n_{f}^{\mathrm{NEW}}\)

Definition at line 865 of file QCD.cpp.

867 {
868  for (int i = 0; i < CacheSize; ++i)
869  if ((AlsM == logLambdaNLO_cache[0][i])
870  && (MAls == logLambdaNLO_cache[1][i])
871  && (mut == logLambdaNLO_cache[2][i])
872  && (mub == logLambdaNLO_cache[3][i])
873  && (muc == logLambdaNLO_cache[4][i])
874  && (nfNEW == logLambdaNLO_cache[5][i])
875  && (nfORG == logLambdaNLO_cache[6][i])
876  && (logLambdaORG == logLambdaNLO_cache[7][i]))
877  return logLambdaNLO_cache[8][i];
878 
880  logLambdaNLO_cache[0][0] = AlsM;
881  logLambdaNLO_cache[1][0] = MAls;
882  logLambdaNLO_cache[2][0] = mut;
883  logLambdaNLO_cache[3][0] = mub;
884  logLambdaNLO_cache[4][0] = muc;
885  logLambdaNLO_cache[5][0] = nfNEW;
886  logLambdaNLO_cache[6][0] = nfORG;
887  logLambdaNLO_cache[7][0] = logLambdaORG;
888 
889  double xmin = -4., xmax = -0.2;
890 
891  TF1 f;
892  if (nfNEW == 6. && nfORG == 5.) {
893  f = TF1("f", this, &QCD::ZeroNf6NLO, xmin, xmax, 1, "QCD", "zeroNf6NLO");
894  } else if (nfNEW == 4. && nfORG == 5.) {
895  f = TF1("f", this, &QCD::ZeroNf4NLO, xmin, xmax, 1, "QCD", "zeroNf4NLO");
896  } else if (nfNEW == 3. && nfORG == 4.) {
897  f = TF1("f", this, &QCD::ZeroNf3NLO, xmin, xmax, 1, "QCD", "zeroNf3NLO");
898  } else
899  throw std::runtime_error("Error in QCD::logLambdaNLO()");
900 
901  ROOT::Math::WrappedTF1 wf1(f);
902  wf1.SetParameters(&logLambdaORG);
903 
904  ROOT::Math::BrentRootFinder brf;
905  brf.SetFunction(wf1, xmin, xmax);
906 
907  if (brf.Solve()) logLambdaNLO_cache[8][0] = brf.Root();
908  else
909  throw std::runtime_error("Error in QCD::logLambdaNLO()");
910 
911  return ( logLambdaNLO_cache[8][0]);
912 }

◆ MassOfNf()

double QCD::MassOfNf ( int  nf) const
protected

The Mbar mass of the heaviest quark in the theory with Nf active flavour.

Parameters
[in]Nfthe number of active flavour
Returns
MSbar \(m_q(m_q)\)

Definition at line 620 of file QCD.cpp.

621 {
622  switch(nf)
623  {
624  case 6:
625  return(quarks[TOP].getMass());
626  case 5:
627  return(quarks[BOTTOM].getMass());
628  case 4:
629  return(quarks[CHARM].getMass());
630  case 3:
631  return(quarks[STRANGE].getMass());
632  default:
633  throw std::runtime_error("QCD::MassOfNf(): no running masses for light quarks");
634  }
635 }

◆ Mbar2Mp()

double QCD::Mbar2Mp ( const double  mbar,
const orders  order = FULLNNLO 
) const

Converts the \(\overline{\mathrm{MS}}\) mass \(m(m)\) to the pole mass.

Parameters
[in]mbarthe \(\overline{\mathrm{MS}}\) mass \(m(m)\) in GeV
[in]orderLO, NLO, FULLNLO, NNLO or FULLNNLO in the \(\alpha_s\) expansion defined in OrderScheme
Returns
the pole mass in GeV
Attention
Can only be used for conversion of mass of the top and bottom quarks.

Definition at line 1221 of file QCD.cpp.

1222 {
1223  // LO contribution
1224  double MpLO = mbar;
1225  if (order == LO) return MpLO;
1226 
1227  // alpha_s(mbar)/pi
1228  orders orderForAls;
1229  if (order == NLO || order == FULLNLO) orderForAls = FULLNLO;
1230  if (order == NNLO || order == FULLNNLO) orderForAls = FULLNNLO;
1231  double a = Als(mbar + MEPS, orderForAls) / M_PI;
1232 
1233  // NLO contribution
1234  double MpNLO = mbar * 4. / 3. * a;
1235  if (order == NLO) return MpNLO;
1236  if (order == FULLNLO) return (MpLO + MpNLO);
1237 
1238  // NNLO contribution
1239  double nl, x;
1240  if (mbar < 3.)
1241  throw std::runtime_error("QCD::Mbar2Mp() can convert only top and bottom masses");
1242  else if (mbar < 50.) {
1243  // for the b quark
1244  nl = 4.;
1245  /* simply adding m_s(2 GeV) and m_c(m_c) */
1246  x = (quarks[STRANGE].getMass() + quarks[CHARM].getMass()) / mbar;
1247  } else {
1248  // for the top quark
1249  nl = 5.;
1250  /* simply adding m_s(2 GeV), m_c(m_c) and m_b(m_b) */
1251  x = (quarks[STRANGE].getMass() + quarks[CHARM].getMass()
1252  + quarks[BOTTOM].getMass()) / mbar;
1253  }
1254  double Delta = M_PI * M_PI / 8. * x - 0.597 * x * x + 0.230 * x * x*x;
1255  double MpNNLO = mbar * (307. / 32. + 2. * zeta2 + 2. / 3. * zeta2 * log(2.0) - zeta3 / 6.
1256  - nl / 3. * (zeta2 + 71. / 48.) + 4. / 3. * Delta) * a*a;
1257  if (order == NNLO) return MpNNLO;
1258  if (order == FULLNNLO) return (MpLO + MpNLO + MpNNLO);
1259 
1260  throw std::runtime_error(orderToString(order) + " is not implemented in QCD::Mbar2Mp().");
1261 }

◆ Mp2Mbar()

double QCD::Mp2Mbar ( const double  mp,
const orders  order = FULLNNLO 
) const

Converts a quark pole mass to the corresponding \(\overline{\mathrm{MS}}\) mass \(m(m)\).

Parameters
[in]mpthe pole mass of the bottom or top quark in GeV
[in]orderLO, NLO, FULLNLO, NNLO or FULLNNLO in the \(\alpha_s\) expansion defined in OrderScheme
Returns
the \(\overline{\mathrm{MS}}\) mass \(m(m)\) in GeV

Definition at line 1270 of file QCD.cpp.

1271 {
1272  if (order == NLO || order == NNLO)
1273  throw std::runtime_error(orderToString(order) + " is not implemented in QCD::Mp2Mbar().");
1274 
1275  int i;
1276  double ms = quarks[STRANGE].getMass(), mc = quarks[CHARM].getMass();
1277  double alsmp = Als(mp, order);
1278  for (i = 0; i < CacheSize; ++i)
1279  if (alsmp == mp2mbar_cache[0][i] && ms == mp2mbar_cache[1][i] &&
1280  mc == mp2mbar_cache[2][i] && (double) order == mp2mbar_cache[3][i])
1281  return mp2mbar_cache[4][i];
1282 
1284  mp2mbar_cache[0][0] = alsmp;
1285  mp2mbar_cache[1][0] = ms;
1286  mp2mbar_cache[2][0] = mc;
1287  mp2mbar_cache[3][0] = (double) order;
1288 
1289  TF1 f("f", this, &QCD::Mp2MbarTMP, mp / 2., 2. * mp, 2, "QCD", "mp2mbara");
1290 
1291  ROOT::Math::WrappedTF1 wf1(f);
1292  double params[2];
1293  params[0] = mp;
1294  params[1] = (double) order;
1295  wf1.SetParameters(params);
1296 
1297  ROOT::Math::BrentRootFinder brf;
1298 
1299  brf.SetFunction(wf1, .7 * mp, 1.3 * mp);
1300  if (brf.Solve())
1301  mp2mbar_cache[4][0] = brf.Root();
1302  else
1303  throw std::runtime_error("error in QCD::mp2mbar");
1304 
1305  return (mp2mbar_cache[4][0]);
1306 }

◆ Mp2MbarTMP()

double QCD::Mp2MbarTMP ( double *  mu,
double *  params 
) const
private

The member used for finding the numerical solution to the pole mass from the \(\overline{\rm MS}\) mass.

Parameters
[in]mua pointer to the \(\overline{\rm MS}\) mass
[in]paramsa pointer to a vector containing the pole mass and the QCD order of the computation
Returns
the difference in the pole mass and the pole mass as computed from the \(\overline{\rm MS}\) mass

Definition at line 1263 of file QCD.cpp.

1264 {
1265  double mp = params[0];
1266  orders order = (orders) params[1];
1267  return (mp - Mbar2Mp(*mu, order));
1268 }

◆ Mrun() [1/2]

double QCD::Mrun ( const double  mu,
const double  m,
const orders  order = FULLNNLO 
) const

Computes a running quark mass \(m(\mu)\) from \(m(m)\).

Parameters
[in]mua scale \(\mu\) in GeV
[in]mthe \(\overline{\mathrm{MS}}\) mass \(m(m)\) in GeV
[in]orderLO, NLO, FULLNLO, NNLO or FULLNNLO in the \(\alpha_s\) expansion defined in OrderScheme
Returns
the running quark mass \(m(\mu)\) in GeV

Definition at line 1064 of file QCD.cpp.

1065 {
1066  return Mrun(mu, m, m, order);
1067 }

◆ Mrun() [2/2]

double QCD::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\).

Parameters
[in]mu_fa scale \(\mu_f\) in GeV
[in]mu_ia scale \(\mu_i\) in GeV
[in]mthe \(\overline{\mathrm{MS}}\) mass \(m(\mu_i)\) in GeV
[in]orderLO, NLO, FULLNLO, NNLO or FULLNNLO in the \(\alpha_s\) expansion defined in OrderScheme
Returns
the running quark mass \(m(\mu_f)\) in GeV

Definition at line 1069 of file QCD.cpp.

1071 {
1072  // Note: When the scale evolves across a flavour threshold, the definitions
1073  // of the outputs for "NLO" and "NNLO" become complicated.
1074 
1075  int i;
1076  for (i = 0; i < CacheSize; ++i) {
1077  if ((mu_f == mrun_cache[0][i]) && (mu_i == mrun_cache[1][i]) &&
1078  (m == mrun_cache[2][i]) && ((double) order == mrun_cache[3][i]) &&
1079  (AlsM == mrun_cache[4][i]) && (MAls == mrun_cache[5][i]) &&
1080  (mut == mrun_cache[6][i]) && (mub == mrun_cache[7][i]) &&
1081  (muc == mrun_cache[8][i]))
1082  return mrun_cache[9][i];
1083  }
1084 
1085  double nfi = Nf(mu_i), nff = Nf(mu_f);
1086  double mu_threshold, mu_threshold2, mu_threshold3, mrun;
1087  if (nff == nfi)
1088  mrun = MrunTMP(mu_f, mu_i, m, order);
1089  else if (nff > nfi) {
1090  if (order == NLO || order == NNLO)
1091  throw std::runtime_error(orderToString(order) + " is not implemented in QCD::Mrun(mu_f,mu_i,m,order)");
1092  mu_threshold = AboveTh(mu_i);
1093  mrun = MrunTMP(mu_threshold - MEPS, mu_i, m, order);
1094  if (order == FULLNNLO)
1095  mrun *= threCorrForMass(nfi + 1., nfi); // threshold corrections
1096  if (nff == nfi + 1.) {
1097  mrun = MrunTMP(mu_f, mu_threshold + MEPS, mrun, order);
1098  } else if (nff == nfi + 2.) {
1099  mu_threshold2 = BelowTh(mu_f);
1100  mrun = MrunTMP(mu_threshold2 - MEPS, mu_threshold + MEPS, mrun, order);
1101  if (order == FULLNNLO)
1102  mrun *= threCorrForMass(nff, nfi + 1.); // threshold corrections
1103  mrun = MrunTMP(mu_f, mu_threshold2 + MEPS, mrun, order);
1104  } else if (nff == nfi + 3.) {
1105  mu_threshold2 = AboveTh(mu_threshold);
1106  mrun = MrunTMP(mu_threshold2 - MEPS, mu_threshold + MEPS, mrun, order);
1107  if (order == FULLNNLO)
1108  mrun *= threCorrForMass(nfi + 2., nfi + 1.); // threshold corrections
1109  mu_threshold3 = BelowTh(mu_f);
1110  mrun = MrunTMP(mu_threshold3 - MEPS, mu_threshold2 + MEPS, mrun, order);
1111  if (order == FULLNNLO)
1112  mrun *= threCorrForMass(nff, nfi + 2.); // threshold corrections
1113  mrun = MrunTMP(mu_f, mu_threshold3 + MEPS, mrun, order);
1114  } else
1115  throw std::runtime_error("Error in QCD::Mrun(mu_f,mu_i,m,order)");
1116  } else {
1117  if (order == NLO || order == NNLO)
1118  throw std::runtime_error(orderToString(order) + " is not implemented in QCD::Mrun(mu_f,mu_i,m,order)");
1119  mu_threshold = BelowTh(mu_i);
1120  mrun = MrunTMP(mu_threshold + MEPS, mu_i, m, order);
1121  if (order == FULLNNLO)
1122  mrun *= threCorrForMass(nfi - 1., nfi); // threshold corrections
1123  if (nff == nfi - 1.)
1124  mrun = MrunTMP(mu_f, mu_threshold - MEPS, mrun, order);
1125  else if (nff == nfi - 2.) {
1126  mu_threshold2 = AboveTh(mu_f);
1127  mrun = MrunTMP(mu_threshold2 + MEPS, mu_threshold - MEPS, mrun, order);
1128  if (order == FULLNNLO)
1129  mrun *= threCorrForMass(nff, nfi - 1.); // threshold corrections
1130  mrun = MrunTMP(mu_f, mu_threshold2 - MEPS, mrun, order);
1131  } else
1132  throw std::runtime_error("Error in QCD::Mrun(mu_f,mu_i,m,order)");
1133  }
1134 
1135  if (mrun < 0.0) {
1136  std::stringstream out;
1137  out << "QCD::Mrun(): A quark mass becomes tachyonic in QCD::Mrun("
1138  << mu_f << ", " << mu_i << ", " << m << ", " << orderToString(order) << ")"
1139  << std::endl
1140  << " Als(" << mu_i << ", " << orderToString(order) << ")/(4pi)="
1141  << Als(mu_i, order) / (4. * M_PI) << std::endl
1142  << " Als(" << mu_f << ", " << orderToString(order) << ")/(4pi)="
1143  << Als(mu_f, order) / (4. * M_PI);
1144  throw std::runtime_error(out.str());
1145  }
1146 
1147  CacheShift(mrun_cache, 10);
1148  mrun_cache[0][0] = mu_f;
1149  mrun_cache[1][0] = mu_i;
1150  mrun_cache[2][0] = m;
1151  mrun_cache[3][0] = (double) order;
1152  mrun_cache[4][0] = AlsM;
1153  mrun_cache[5][0] = MAls;
1154  mrun_cache[6][0] = mut;
1155  mrun_cache[7][0] = mub;
1156  mrun_cache[8][0] = muc;
1157  mrun_cache[9][0] = mrun;
1158 
1159  return mrun;
1160 }

◆ Mrun4()

double QCD::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\).

Parameters
[in]mu_fthe final scale \(\mu_f\) to which the mass is run
[in]mu_ithe initial scale \(\mu_i\) from which the mass is run
[in]mthe mass at the scale \(\mu_i\)
Returns
the mass at the scale \(\mu_f\) with \(n_f = 4\)
Attention
Temporary function waiting for the implementation of NNLO exact.

Definition at line 1199 of file QCD.cpp.

1200 {
1201  double nf = 4.;
1202 
1203  // alpha_s/(4pi)
1204  double ai = Als4(mu_i) / (4. * M_PI);
1205  double af = Als4(mu_f) / (4. * M_PI);
1206 
1207  // LO contribution
1208  double b0 = Beta0(nf), g0 = Gamma0(nf);
1209  double mLO = m * pow(af / ai, g0 / (2. * b0));
1210 
1211  // NLO contribution
1212  double b1 = Beta1(nf), g1 = Gamma1(nf);
1213  double A1 = g1 / (2. * b0) - b1 * g0 / (2. * b0 * b0);
1214  double mNLO = mLO * A1 * (af - ai);
1215  return (mLO + mNLO);
1216 
1217 }

◆ MrunTMP()

double QCD::MrunTMP ( const double  mu_f,
const double  mu_i,
const double  m,
const orders  order 
) const
private

A function to calculate the running of the mass between flavour thresholds.

Parameters
[in]mu_fthe final scale \(\mu_f\) to which the mass if run
[in]mu_ithe initial scale \(\mu_i\) from which the mass if run
[in]mthe mass at the scale \(\mu_i\)
[in]orderthe QCD order at which the running is being calculated
Returns
the mass run from \(\mu_i\) to \(\mu_f\)

Definition at line 1162 of file QCD.cpp.

1164 {
1165  double nf = Nf(mu_f);
1166  if (nf != Nf(mu_i))
1167  throw std::runtime_error("Error in QCD::MrunTMP().");
1168 
1169  // alpha_s/(4pi)
1170  orders orderForAls;
1171  if (order == LO) orderForAls = LO;
1172  if (order == NLO || order == FULLNLO) orderForAls = FULLNLO;
1173  if (order == NNLO || order == FULLNNLO) orderForAls = FULLNNLO;
1174  double ai = Als(mu_i, orderForAls) / (4. * M_PI);
1175  double af = Als(mu_f, orderForAls) / (4. * M_PI);
1176 
1177  // LO contribution
1178  double b0 = Beta0(nf), g0 = Gamma0(nf);
1179  double mLO = m * pow(af / ai, g0 / (2. * b0));
1180  if (order == LO) return mLO;
1181 
1182  // NLO contribution
1183  double b1 = Beta1(nf), g1 = Gamma1(nf);
1184  double A1 = g1 / (2. * b0) - b1 * g0 / (2. * b0 * b0);
1185  double mNLO = mLO * A1 * (af - ai);
1186  if (order == NLO) return mNLO;
1187  if (order == FULLNLO) return (mLO + mNLO);
1188 
1189  // NNLO contribution
1190  double b2 = Beta2(nf), g2 = Gamma2(nf);
1191  double A2 = b1 * b1 * g0 / (2. * b0 * b0 * b0) - b2 * g0 / (2. * b0 * b0) - b1 * g1 / (2. * b0 * b0) + g2 / (2. * b0);
1192  double mNNLO = mLO * (A1 * A1 / 2. * (af - ai)*(af - ai) + A2 / 2. * (af * af - ai * ai));
1193  if (order == NNLO) return mNNLO;
1194  if (order == FULLNNLO) return (mLO + mNLO + mNNLO);
1195 
1196  throw std::runtime_error(orderToString(order) + " is not implemented in QCD::MrunTMP()");
1197 }

◆ MS2DRqmass() [1/2]

double QCD::MS2DRqmass ( const double  MSbar) const

Converts a quark mass from the \(\overline{\mathrm{MS}}\) scheme to the \(\overline{\mathrm{DR}}\) scheme.

Parameters
[in]MSbarThe \(\overline{\mathrm{MS}}\) mass \(m(m)\) in GeV
Returns
The \(\overline{\mathrm{DR}}\) mass \(m(m)\) in GeV

Definition at line 1308 of file QCD.cpp.

1309 {
1310  return (MSbar / (1. + Als(MSbar, FULLNLO) / 4. / M_PI * CF));
1311 }

◆ MS2DRqmass() [2/2]

double QCD::MS2DRqmass ( const double  MSscale,
const double  MSbar 
) const

Converts a quark mass from the \(\overline{\mathrm{MS}}\) scheme to the \(\overline{\mathrm{DR}}\) scheme.

Parameters
[in]MSscalethe scale at which the \(\overline{\mathrm{MS}}\) mass is defined
[in]MSbarthe \(\overline{\mathrm{MS}}\) mass \(m(m)\) in GeV
Returns
the \(\overline{\mathrm{DR}}\) mass \(m(m)\) in GeV

Definition at line 1313 of file QCD.cpp.

1314 {
1315  return (MSbar / (1. + Als(MSscale, FULLNLO) / 4. / M_PI * CF));
1316 }

◆ Nf()

double QCD::Nf ( const double  mu) const

The number of active flavour at scale \(\mu\).

Parameters
[in]mua scale \(\mu\) in GeV
Returns
active N_f

Definition at line 438 of file QCD.cpp.

439 {
440  int i;
441  for (i = 1; i < 5; i++)
442  if (mu >= Thresholds(i))
443  return (7. - (double) i);
444 
445  throw std::runtime_error("Error in QCD::Nf()");
446 }

◆ NfThresholdCorrections()

double QCD::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.

Parameters
[in]muthe matching scale
[in]Mthe running quark mass
[in]alsvalue of \(\alpha_s(\mu)\) in the \(n_f\) flavour theory
[in]nfnumber of active flavour
[in]orderorder of the expansion in \(\alpha_s\)
Returns
Threshold correction (without the leading term equal to 1)

Definition at line 583 of file QCD.cpp.

584 {
585  double lmM = 2.*log(mu/M), res = 0.;
586  switch(order)
587  {
588  case FULLNNNLO:
589  res += als*als*als/M_PI/M_PI/M_PI*(-564731./124416. + 82043./27648.*gslpp_special_functions::zeta(3) +
590  2191./576.*lmM + 511./576.*lmM*lmM + lmM*lmM*lmM/216. + ((double) nf - 1.) * (2633./31104. - 67./576.*lmM + lmM*lmM/36.));
591  case FULLNNLO:
592  res += als*als/M_PI/M_PI*(-11./72. + 19./24.*lmM + lmM*lmM/36.);
593  case FULLNLO:
594  res += als/M_PI*lmM/6.;
595  case LO:
596  break;
597  default:
598  throw std::runtime_error("QCD::ThresholdCorrections(): order not implemented.");
599  }
600  return(res);
601 }

◆ orderToString()

std::string QCD::orderToString ( const orders  order) const

Converts an object of the enum type "orders" to the corresponding string.

Parameters
[in]orderan object of the enum type "orders"
Returns
the string of the given "order"

Definition at line 83 of file QCD.cpp.

84 {
85  switch (order) {
86  case LO:
87  return "LO";
88  case NLO:
89  return "NLO";
90  case FULLNLO:
91  return "FULLNLO";
92  case NNLO:
93  return "NNLO";
94  case FULLNNLO:
95  return "FULLNNLO";
96  case NNNLO:
97  return "NNNLO";
98  case FULLNNNLO:
99  return "FULLNNNLO";
100  default:
101  throw std::runtime_error("QCD::orderToString(): order not implemented.");
102  }
103 }

◆ PostUpdate()

bool QCD::PostUpdate ( )
virtual

The post-update method for QCD.

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

  • computing the decay constant \(F_{B_D}\) from \(F_{B_s}\)
  • computing the bag parameters \(B_{B_d}\) from \(B_{B_s}\)
  • computing the \(\overline{\rm MS}\) mass of the top quark at the \(\overline{\rm MS}\) mass, \(m_t^{\overline{\rm MS}}(m_t^{\overline{\rm MS}})\) and setting the scale at the same value.
    Returns
    a boolean that is true if the execution is successful

Implements Model.

Reimplemented in NPSMEFTd6, StandardModel, GeneralTHDM, THDMW, HiggsKigen, THDM, GeorgiMachacek, NPZbbbarLinearized, NPZbbbar, RealWeakEFTLFV, FlavourWilsonCoefficient, NPEpsilons, FlavourWilsonCoefficient_DF2, HiggsChiral, RealWeakEFTCC, NPbase, myModel, and myModel.

Definition at line 143 of file QCD.cpp.

144 {
145  if (computeFBd) mesonsMap.at(QCD::B_D).setDecayconst(mesonsMap.at(QCD::B_S).getDecayconst() / mesonsMap.at(QCD::B_D).getFBsoFBd());
146  if (computeFBp) mesonsMap.at(QCD::B_P).setDecayconst(mesonsMap.at(QCD::B_S).getDecayconst() / mesonsMap.at(QCD::B_P).getFBsoFBd()); /**** FOR NOW FB+ = FBd ****/
147  if (computeBs && FlagCsi)
148  {
149  BParameterMap.at("BBs").setBpars(0, BParameterMap.at("BBs").getFBsSqrtBBs1() * BParameterMap.at("BBs").getFBsSqrtBBs1() / mesonsMap.at(QCD::B_S).getDecayconst() / mesonsMap.at(QCD::B_S).getDecayconst());
150  BParameterMap.at("BBs").setBpars(1, BParameterMap.at("BBs").getFBsSqrtBBs2() * BParameterMap.at("BBs").getFBsSqrtBBs2() / mesonsMap.at(QCD::B_S).getDecayconst() / mesonsMap.at(QCD::B_S).getDecayconst());
151  BParameterMap.at("BBs").setBpars(2, BParameterMap.at("BBs").getFBsSqrtBBs3() * BParameterMap.at("BBs").getFBsSqrtBBs3() / mesonsMap.at(QCD::B_S).getDecayconst() / mesonsMap.at(QCD::B_S).getDecayconst());
152  BParameterMap.at("BBs").setBpars(3, BParameterMap.at("BBs").getFBsSqrtBBs4() * BParameterMap.at("BBs").getFBsSqrtBBs4() / mesonsMap.at(QCD::B_S).getDecayconst() / mesonsMap.at(QCD::B_S).getDecayconst());
153  BParameterMap.at("BBs").setBpars(4, BParameterMap.at("BBs").getFBsSqrtBBs5() * BParameterMap.at("BBs").getFBsSqrtBBs5() / mesonsMap.at(QCD::B_S).getDecayconst() / mesonsMap.at(QCD::B_S).getDecayconst());
154  }
155  if (computeBd) {
156  if (FlagCsi) {
157  BParameterMap.at("BBd").setBpars(0, mesonsMap.at(QCD::B_D).getFBsoFBd() * mesonsMap.at(QCD::B_D).getFBsoFBd() * BParameterMap.at("BBs").getBpars()(0) / BParameterMap.at("BBd").getcsi() / BParameterMap.at("BBd").getcsi());
158  BParameterMap.at("BBd").setBBsoBBd(BParameterMap.at("BBs").getBpars()(0) / BParameterMap.at("BBd").getBpars()(0));
159  BParameterMap.at("BBd").setBpars(1, BParameterMap.at("BBd").getFBdSqrtBBd2() * BParameterMap.at("BBd").getFBdSqrtBBd2() / mesonsMap.at(QCD::B_D).getDecayconst() / mesonsMap.at(QCD::B_D).getDecayconst());
160  BParameterMap.at("BBd").setBpars(2, BParameterMap.at("BBd").getFBdSqrtBBd3() * BParameterMap.at("BBd").getFBdSqrtBBd3() / mesonsMap.at(QCD::B_D).getDecayconst() / mesonsMap.at(QCD::B_D).getDecayconst());
161  BParameterMap.at("BBd").setBpars(3, BParameterMap.at("BBd").getFBdSqrtBBd4() * BParameterMap.at("BBd").getFBdSqrtBBd4() / mesonsMap.at(QCD::B_D).getDecayconst() / mesonsMap.at(QCD::B_D).getDecayconst());
162  BParameterMap.at("BBd").setBpars(4, BParameterMap.at("BBd").getFBdSqrtBBd5() * BParameterMap.at("BBd").getFBdSqrtBBd5() / mesonsMap.at(QCD::B_D).getDecayconst() / mesonsMap.at(QCD::B_D).getDecayconst());
163  } else
164  BParameterMap.at("BBd").setBpars(0, BParameterMap.at("BBs").getBpars()(0) / BParameterMap.at("BBd").getBBsoBBd());
165  }
166  if (computemt) {
168 #if SUSYFIT_DEBUG & 2
169 // std::cout << "WARNING: using NLO mt for debugging purpose!"<< std::endl;
170 // quarks[TOP].setMass(Mp2Mbar(mtpole, FULLNLO));
171  mut = quarks[TOP].getMass();
172  std::cout << "postupdate: " << mtpole << std::endl;
173  std::cout << "postupdate: " << Mp2Mbar(mtpole, FULLNNLO) << std::endl;
174  std::cout << "mut set to " << mut << std::endl;
175 #endif
176  quarks[TOP].setMass_scale(quarks[TOP].getMass());
177  }
178 
179  return (true);
180 }

◆ PreUpdate()

bool QCD::PreUpdate ( )
virtual

The pre-update method for QCD.

This method resets the internal flags requireYu, requireYd, computeBd, computeFBd and computemt before updating the model parameters with the method Update().

Returns
a boolean that is true if the execution is successful

Implements Model.

Reimplemented in StandardModel, GeneralTHDM, THDMW, THDM, GeorgiMachacek, RealWeakEFTLFV, FlavourWilsonCoefficient, RealWeakEFTCC, myModel, and myModel.

Definition at line 116 of file QCD.cpp.

117 {
118  requireYu = false;
119  requireYd = false;
120  computemt = false;
121 
122  return (true);
123 }

◆ setFlag()

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

A method to set a flag of QCD.

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

Implements Model.

Reimplemented in GeneralTHDM, GeorgiMachacek, NPSMEFT6dtopquark, NPSMEFTd6, StandardModel, HiggsKigen, THDM, RealWeakEFTLFV, FlavourWilsonCoefficient, NPEpsilons, HiggsChiral, RealWeakEFTCC, myModel, and myModel.

Definition at line 380 of file QCD.cpp.

381 {
382  bool res = false;
383  if (name.compare("FlagCsi") == 0) {
384  FlagCsi = value;
385  if (computeBs) BParameterMap.at("BBs").setFlagCsi(FlagCsi);
386  if (computeBd) BParameterMap.at("BBd").setFlagCsi(FlagCsi);
387  res = true;
388  }
389 
390  return res;
391 }

◆ setFlagStr()

bool QCD::setFlagStr ( const std::string  name,
const std::string  value 
)
virtual

A method to set a flag of QCD.

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

Implements Model.

Reimplemented in StandardModel, GeneralTHDM, THDMW, and THDM.

Definition at line 393 of file QCD.cpp.

394 {
395  std::cout << "WARNING: wrong name or value for ModelFlag " << name << std::endl;
396  return (false);
397 }

◆ setNc()

void QCD::setNc ( double  Nc)
inline

A set method to change the number of colours \(N_c\).

Parameters
[in]Ncthe number of colours

Definition at line 514 of file QCD.h.

515  {
516  this->Nc = Nc;
517  }

◆ setOptionalParameter()

void QCD::setOptionalParameter ( std::string  name,
double  value 
)
inline

A method to set the parameter value for the parameters that are specific to only one set of observables.

Parameters
[in]namethe name of the parameter
[in]valuethe value of the parameter

Definition at line 458 of file QCD.h.

459  {
460  optionalParameters[name] = value;
461  }

◆ setParameter()

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

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

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

Implements Model.

Reimplemented in GeneralTHDM, THDMW, GeorgiMachacek, FlavourWilsonCoefficient_DF2, FlavourWilsonCoefficient, myModel, myModel, CMFV, NPSMEFTd6, StandardModel, HiggsKigen, HiggsChiral, THDM, NPEpsilons, HiggsKvgenKfgen, HiggsKvKfgen, HiggsKvKf, NPZbbbar, NPEpsilons_pureNP, NPZbbbarLinearized, NPSTU, NPSTUZbbbarLR, NPSMEFT6dtopquark, and NPDF2.

Definition at line 273 of file QCD.cpp.

274 {
275  int notABparameter = 0;
276  int notAMesonParameter = 0;
277 
278  if (name.compare("AlsM") == 0) {
279  AlsM = value;
280  computemt = true;
281  requireYu = true;
282  requireYd = true;
283  } else if (name.compare("MAls") == 0) {
284  MAls = value;
285  computemt = true;
286  requireYu = true;
287  requireYd = true;
288  } else if (name.compare("mup") == 0) {
289  if (value < MEPS) UpdateError = true;
290  quarks[UP].setMass(value);
291  requireYu = true;
292  } else if (name.compare("mdown") == 0) {
293  if (value < MEPS) UpdateError = true;
294  quarks[DOWN].setMass(value);
295  requireYd = true;
296  } else if (name.compare("mcharm") == 0) {
297  quarks[CHARM].setMass(value);
298  quarks[CHARM].setMass_scale(value);
299  requireYu = true;
300  } else if (name.compare("mstrange") == 0) {
301  if (value < MEPS) UpdateError = true;
302  quarks[STRANGE].setMass(value);
303  requireYd = true;
304  } else if (name.compare("mtop") == 0) {
305  mtpole = value;
306  requireYu = true;
307  computemt = true;
308  } else if (name.compare("mbottom") == 0) {
309  quarks[BOTTOM].setMass(value);
310  quarks[BOTTOM].setMass_scale(value);
311  requireYd = true;
312  } else if (name.compare("muc") == 0)
313  muc = value;
314  else if (name.compare("mub") == 0)
315  mub = value;
316  else if (name.compare("mut") == 0)
317  mut = value;
318  else if (optionalParameters.find(name) != optionalParameters.end())
319  setOptionalParameter(name, value);
320  else {
321  if (!BParameterMap.empty())
322  for (std::map<std::string, BParameter>::iterator it = BParameterMap.begin(); it != BParameterMap.end(); it++)
323  if(it->second.setParameter(name, value))
324  notABparameter += 1;
325  if (!mesonsMap.empty())
326  for (std::map<const QCD::meson, Meson>::iterator it = mesonsMap.begin(); it != mesonsMap.end(); it++)
327  if(it->second.setParameter(name, value))
328  notAMesonParameter += 1;
329  if (unknownParameterWarning && !isSliced && notABparameter == 0 && notAMesonParameter == 0)
330  if (std::find(unknownParameters.begin(), unknownParameters.end(), name) == unknownParameters.end()) unknownParameters.push_back(name);
331  }
332 
333 }

◆ threCorrForMass()

double QCD::threCorrForMass ( const double  nf_f,
const double  nf_i 
) const
private

The threshold correction for running of a mass when crossing a flavour threshold.

Parameters
[in]nf_fthe number of flavours \(n_f\) after crossing the threshold
[in]nf_ithe number of flavours \(n_i\) before crossing the threshold
Returns
the threshold correction factor

Definition at line 1025 of file QCD.cpp.

1026 {
1027  if (fabs(nf_f - nf_i) != 1.)
1028  throw std::runtime_error("Error in QCD::threCorrForMass()");
1029 
1030  double mu_threshold, mf, log_mu2_mf2;
1031  if (nf_f > nf_i) {
1032  if (nf_f == 6.) {
1033  mu_threshold = mut;
1034  mf = quarks[TOP].getMass(); // m_t(m_t)
1035  } else if (nf_f == 5.) {
1036  mu_threshold = mub;
1037  mf = quarks[BOTTOM].getMass(); // m_b(m_b)
1038  } else if (nf_f == 4.) {
1039  mu_threshold = muc;
1040  mf = quarks[CHARM].getMass(); // m_c(m_c)
1041  } else
1042  throw std::runtime_error("Error in QCD::threCorrForMass()");
1043  log_mu2_mf2 = 2. * log(mu_threshold / mf);
1044  return (1. + pow(Als(mu_threshold - MEPS, FULLNNLO) / M_PI, 2.)
1045  *(-log_mu2_mf2 * log_mu2_mf2 / 12. + 5. / 36. * log_mu2_mf2 - 89. / 432.));
1046  } else {
1047  if (nf_i == 6.) {
1048  mu_threshold = mut;
1049  mf = quarks[TOP].getMass(); // m_t(m_t)
1050  } else if (nf_i == 5.) {
1051  mu_threshold = mub;
1052  mf = quarks[BOTTOM].getMass(); // m_b(m_b)
1053  } else if (nf_i == 4.) {
1054  mu_threshold = muc;
1055  mf = quarks[CHARM].getMass(); // m_c(m_c)
1056  } else
1057  throw std::runtime_error("Error in QCD::threCorrForMass()");
1058  log_mu2_mf2 = 2. * log(mu_threshold / mf);
1059  return (1. + pow(Als(mu_threshold + MEPS, FULLNNLO) / M_PI, 2.)
1060  *(log_mu2_mf2 * log_mu2_mf2 / 12. - 5. / 36. * log_mu2_mf2 + 89. / 432.));
1061  }
1062 }

◆ Thresholds()

double QCD::Thresholds ( const int  i) const

For accessing the active flavour threshold scales.

Parameters
[in]ithe index referring to active flavour thresholds.
Returns
the threshold scale: 1.0E10 (i = 0), \(\mu_t\) (i = 1), \(\mu_b\) (i = 2), \(\mu_c\) (i = 3) and 0. (default)

Definition at line 406 of file QCD.cpp.

407 {
408  if (!(mut > mub && mub > muc))
409  throw std::runtime_error("inverted thresholds in QCD::Thresholds()!");
410 
411  switch (i) {
412  case 0: return (1.0E10);
413  case 1: return (mut);
414  case 2: return (mub);
415  case 3: return (muc);
416  default: return (0.);
417  }
418 }

◆ Update()

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

The update method for QCD.

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 (including parameters that are varied and those that are held constant)
Returns
a boolean that is true if the execution is successful

Implements Model.

Reimplemented in StandardModel, GeneralTHDM, THDMW, THDM, GeorgiMachacek, RealWeakEFTLFV, FlavourWilsonCoefficient, RealWeakEFTCC, NPbase, myModel, and myModel.

Definition at line 125 of file QCD.cpp.

126 {
127  if (!PreUpdate()) return (false);
128 
129  UpdateError = false;
130 
131  for (std::map<std::string, double>::const_iterator it = DPars.begin(); it != DPars.end(); it++)
132  {
133  setParameter(it->first, it->second);
134  }
135 
136  if (UpdateError) return (false);
137 
138  if (!PostUpdate()) return (false);
139 
140  return (true);
141 }

◆ ZeroNf3NLO()

double QCD::ZeroNf3NLO ( double *  logLambda3,
double *  logLambda4_in 
) const
private

A member for calculating the difference in \(\alpha_s^{\mathrm{FULLNLO}}\) across the three-four flavour threshold using AlsWithLambda().

Parameters
[in]logLambda3\(\log(\Lambda_{\rm QCD}^{(3)})\)
[in]logLambda4_in\(\log(\Lambda_{\rm QCD}^{(4)})\)
Returns
the difference \(\alpha_s^{\mathrm{FULLNLO}}(\mu_c-\epsilon)\) - \(\alpha_s^{\mathrm{FULLNLO}}(\mu_c+\epsilon)\) with \(\epsilon = 10^{-10}\)

Definition at line 823 of file QCD.cpp.

824 {
825  return ( AlsWithLambda(muc - 1.e-10, *logLambda3, FULLNLO)
826  - AlsWithLambda(muc + 1.e-10, *logLambda4_in, FULLNLO));
827 }

◆ ZeroNf4NLO()

double QCD::ZeroNf4NLO ( double *  logLambda4,
double *  logLambda5_in 
) const
private

A member for calculating the difference in \(\alpha_s^{\mathrm{FULLNLO}}\) across the four-five flavour threshold using AlsWithLambda().

Parameters
[in]logLambda4\(\log(\Lambda_{\rm QCD}^{(4)})\)
[in]logLambda5_in\(\log(\Lambda_{\rm QCD}^{(5)})\)
Returns
the difference \(\alpha_s^{\mathrm{FULLNLO}}(\mu_b-\epsilon)\) - \(\alpha_s^{\mathrm{FULLNLO}}(\mu_b+\epsilon)\) with \(\epsilon = 10^{-10}\)

Definition at line 817 of file QCD.cpp.

818 {
819  return ( AlsWithLambda(mub - 1.e-10, *logLambda4, FULLNLO)
820  - AlsWithLambda(mub + 1.e-10, *logLambda5_in, FULLNLO));
821 }

◆ ZeroNf5()

double QCD::ZeroNf5 ( double *  logLambda5,
double *  order 
) const
private

A member for calculating the difference in \(\alpha_s\) using AlsWithLambda() and the input value of \(\alpha_s(M_{\alpha_s})\) given as a model parameter.

Parameters
[in]logLambda5\(\log(\Lambda_{\rm QCD}^{(5)})\)
[in]orderthe QCD order of the calculation
Returns
AlsWithLambda( \(M_{\alpha_s}\), *logLambda5, *order) - \(\alpha_s(M_{\alpha_s})\)

Definition at line 812 of file QCD.cpp.

813 {
814  return ( AlsWithLambda(MAls, *logLambda5, (orders) * order) - AlsM);
815 }

◆ ZeroNf6NLO()

double QCD::ZeroNf6NLO ( double *  logLambda6,
double *  logLambda5_in 
) const
private

A member for calculating the difference in \(\alpha_s^{\mathrm{FULLNLO}}\) across the six-five flavour threshold using AlsWithLambda().

Parameters
[in]logLambda6\(\log(\Lambda_{\rm QCD}^{(6)})\)
[in]logLambda5_in\(\log(\Lambda_{\rm QCD}^{(5)})\)
Returns
the difference \(\alpha_s^{\mathrm{FULLNLO}}(\mu_t+\epsilon)\) - \(\alpha_s^{\mathrm{FULLNLO}}(\mu_t-\epsilon)\) with \(\epsilon = 10^{-10}\)

Definition at line 806 of file QCD.cpp.

807 {
808  return ( AlsWithLambda(mut + 1.e-10, *logLambda6, FULLNLO)
809  - AlsWithLambda(mut - 1.e-10, *logLambda5_in, FULLNLO));
810 }

Member Data Documentation

◆ als_cache

double QCD::als_cache[9][CacheSize]
mutableprivate

Cache for \(\alpha_s\).

Definition at line 946 of file QCD.h.

◆ AlsM

double QCD::AlsM
protected

The strong coupling constant at the mass scale MAls, \(\alpha_s(M_{\alpha_s})\).

Definition at line 925 of file QCD.h.

◆ BParameterMap

std::map<std::string, BParameter> QCD::BParameterMap
mutableprivate

Definition at line 937 of file QCD.h.

◆ CA

double QCD::CA
protected

Definition at line 933 of file QCD.h.

◆ CacheSize

const int QCD::CacheSize = 5
staticprivate

Defines the depth of the cache.

Definition at line 945 of file QCD.h.

◆ CF

double QCD::CF
protected

Definition at line 933 of file QCD.h.

◆ computeBd

bool QCD::computeBd
mutableprivate

Switch for computing \(B_{B_d}\) from \(B_{B_s}\).

Definition at line 943 of file QCD.h.

◆ computeBs

bool QCD::computeBs
mutableprivate

Switch for computing \(B_{B_s}\) from \(F_{B_s}\sqrt{B_{B_s}}\).

Definition at line 944 of file QCD.h.

◆ computeFBd

bool QCD::computeFBd
mutableprivate

Switch for computing \(F_{B_d}\) from \(F_{B_s}\).

Definition at line 941 of file QCD.h.

◆ computeFBp

bool QCD::computeFBp
mutableprivate

Switch for computing \(F_{B^+}\) from \(F_{B_s}\).

Definition at line 942 of file QCD.h.

◆ computemt

bool QCD::computemt
protected

Switch for computing the \(\overline{\mathrm{MS}}\) mass of the top quark.

Definition at line 920 of file QCD.h.

◆ dAdA_NA

double QCD::dAdA_NA
protected

Definition at line 933 of file QCD.h.

◆ dFdA_NA

double QCD::dFdA_NA
protected

Definition at line 933 of file QCD.h.

◆ dFdF_NA

double QCD::dFdF_NA
protected

Definition at line 933 of file QCD.h.

◆ FlagCsi

bool QCD::FlagCsi
private

A flag to determine whether \(B_{B_s}\) and \(B_{B_s}/B_{B_d}\) or \(F_{B_s}\sqrt{B_{B_s}}\) (false) and \(\xi \equiv F_{B_s}\sqrt{B_{B_s}}/(F_{B_d}\sqrt{B_{B_d}})\) (default, true) are used as inputs.

Definition at line 955 of file QCD.h.

◆ logLambda5_cache

double QCD::logLambda5_cache[4][CacheSize]
mutableprivate

Definition at line 947 of file QCD.h.

◆ logLambdaNLO_cache

double QCD::logLambdaNLO_cache[9][CacheSize]
mutableprivate

Definition at line 948 of file QCD.h.

◆ MAls

double QCD::MAls
protected

The mass scale in GeV at which the strong coupling measurement is provided.

Definition at line 926 of file QCD.h.

◆ mesonsMap

std::map<const QCD::meson, Meson> QCD::mesonsMap
mutableprivate

The map of defined mesons.

Definition at line 954 of file QCD.h.

◆ mp2mbar_cache

double QCD::mp2mbar_cache[5][CacheSize]
mutableprivate

Cache for pole mass to msbar mass conversion.

Definition at line 950 of file QCD.h.

◆ mrun_cache

double QCD::mrun_cache[10][CacheSize]
mutableprivate

Cache for running quark mass.

Definition at line 949 of file QCD.h.

◆ mtpole

double QCD::mtpole
protected

The pole mass of the top quark.

Definition at line 927 of file QCD.h.

◆ mub

double QCD::mub
protected

The threshold between five- and four-flavour theory in GeV.

Definition at line 929 of file QCD.h.

◆ muc

double QCD::muc
protected

The threshold between four- and three-flavour theory in GeV.

Definition at line 930 of file QCD.h.

◆ mut

double QCD::mut
protected

The threshold between six- and five-flavour theory in GeV.

Definition at line 928 of file QCD.h.

◆ NA

double QCD::NA
protected

Definition at line 933 of file QCD.h.

◆ Nc

double QCD::Nc
protected

The number of colours.

Definition at line 932 of file QCD.h.

◆ NQCDvars

const int QCD::NQCDvars = 11
static

The number of model parameters in QCD.

Definition at line 357 of file QCD.h.

◆ optionalParameters

std::map<std::string, double> QCD::optionalParameters
private

A map for containing the list and values of the parameters that are used only by a specific set of observables.

Definition at line 952 of file QCD.h.

◆ QCDvars

std::string QCD::QCDvars
static
Initial value:
= {
"AlsM", "MAls",
"mup", "mdown", "mcharm", "mstrange", "mtop", "mbottom",
"muc", "mub", "mut"
}

An array containing the labels under which all QCD parameters are stored in a vector of ModelParameter via InputParser::ReadParameters().

Definition at line 363 of file QCD.h.

◆ quarks

Particle QCD::quarks[6]
protected

The vector of all SM quarks.

Definition at line 934 of file QCD.h.

◆ realorder

orders QCD::realorder
mutableprivate

Definition at line 956 of file QCD.h.

◆ requireYd

bool QCD::requireYd
protected

Switch for generating the Yukawa couplings to the down-type quarks.

Definition at line 922 of file QCD.h.

◆ requireYu

bool QCD::requireYu
protected

Switch for generating the Yukawa couplings to the up-type quarks.

Definition at line 921 of file QCD.h.

◆ TF

double QCD::TF
protected

Definition at line 933 of file QCD.h.

◆ unknownParameters

std::vector<std::string> QCD::unknownParameters
private

A vector for containing the names of the parameters that are not being used but specified in the configuration file.

Definition at line 953 of file QCD.h.

◆ unknownParameterWarning

bool QCD::unknownParameterWarning
private

A flag to stop the unknown parameter warning after the first time.

Definition at line 951 of file QCD.h.

◆ zeta2

double QCD::zeta2
private

\(\zeta(2)\) computed with the GSL.

Definition at line 939 of file QCD.h.

◆ zeta3

double QCD::zeta3
private

\(\zeta(3)\) computed with the GSL.

Definition at line 940 of file QCD.h.


The documentation for this class was generated from the following files:
QCD::TAU
Definition: QCD.h:316
gslpp_special_functions::zeta
double zeta(int i)
Definition: gslpp_special_functions.cpp:20
BParameter
A class for the bag parameters.
Definition: BParameter.h:151
QCD::NEUTRINO_3
Definition: QCD.h:315
QCD::ZeroNf6NLO
double ZeroNf6NLO(double *logLambda6, double *logLambda5_in) const
A member for calculating the difference in across the six-five flavour threshold using AlsWithLambda...
Definition: QCD.cpp:806
QCD::computeFBd
bool computeFBd
Switch for computing from .
Definition: QCD.h:941
QCD::logLambdaNLO
double logLambdaNLO(const double nfNEW, const double nfORG, const double logLambdaORG) const
used for computation of at FULLNLO.
Definition: QCD.cpp:865
QCD::BOTTOM
Definition: QCD.h:329
QCD::NOLEPTON
Definition: QCD.h:317
QCD::ZeroNf5
double ZeroNf5(double *logLambda5, double *order) const
A member for calculating the difference in using AlsWithLambda() and the input value of given as a ...
Definition: QCD.cpp:812
Meson
A class for mesons.
Definition: Meson.h:315
QCD::CacheSize
static const int CacheSize
Defines the depth of the cache.
Definition: QCD.h:945
QCD::Als
virtual double Als(const double mu, const orders order=FULLNLO, bool Nf_thr=true) const
Definition: QCD.cpp:637
QCD::mub
double mub
The threshold between five- and four-flavour theory in GeV.
Definition: QCD.h:929
QCD::AlsByOrder
virtual double AlsByOrder(const double mu, const orders order=FULLNLO, bool Nf_thr=true) const
Definition: QCD.cpp:658
Particle
A class for particles.
Definition: Particle.h:26
QCD::unknownParameterWarning
bool unknownParameterWarning
A flag to stop the unknown parameter warning after the first time.
Definition: QCD.h:951
QCD::RHO_P
Definition: QCD.h:352
QCD::MAls
double MAls
The mass scale in GeV at which the strong coupling measurement is provided.
Definition: QCD.h:926
QCD::setOptionalParameter
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 observabl...
Definition: QCD.h:458
QCD::B_S
Definition: QCD.h:345
QCD::Beta1
double Beta1(const double nf) const
The coefficient for a certain number of flavours .
Definition: QCD.cpp:471
QCD::mesonsMap
std::map< const QCD::meson, Meson > mesonsMap
The map of defined mesons.
Definition: QCD.h:954
QCD::Nf
double Nf(const double mu) const
The number of active flavour at scale .
Definition: QCD.cpp:438
FULLNNNLO
Definition: OrderScheme.h:39
QCD::initializeMeson
void initializeMeson(QCD::meson meson_i) const
A method to initialize a meson.
Definition: QCD.cpp:236
NNNLO
Definition: OrderScheme.h:36
QCD::computeFBp
bool computeFBp
Switch for computing from .
Definition: QCD.h:942
QCD::computemt
bool computemt
Switch for computing the mass of the top quark.
Definition: QCD.h:920
QCD::dFdA_NA
double dFdA_NA
Definition: QCD.h:933
LO
Definition: OrderScheme.h:33
QCD::UP
Definition: QCD.h:324
QCD::D_P
Definition: QCD.h:342
QCD::Mp2Mbar
double Mp2Mbar(const double mp, const orders order=FULLNNLO) const
Converts a quark pole mass to the corresponding mass .
Definition: QCD.cpp:1270
Model::addMissingModelParameter
void addMissingModelParameter(const std::string &missingParameterName)
Definition: Model.h:232
QCD::CHARM
Definition: QCD.h:326
QCD::TF
double TF
Definition: QCD.h:933
QCD::als_cache
double als_cache[9][CacheSize]
Cache for .
Definition: QCD.h:946
QCD::dFdF_NA
double dFdF_NA
Definition: QCD.h:933
QCD::logLambda
double logLambda(const double nf, orders order) const
Computes with nf flavours in GeV.
Definition: QCD.cpp:977
QCD::FlagCsi
bool FlagCsi
A flag to determine whether and or (false) and (default, true) are used as inputs.
Definition: QCD.h:955
QCD::NEUTRINO_2
Definition: QCD.h:313
QCD::setParameter
virtual void setParameter(const std::string name, const double &value)
A method to set the value of a parameter of QCD.
Definition: QCD.cpp:273
gslpp::log
complex log(const complex &z)
Definition: gslpp_complex.cpp:342
QCD::P_P
Definition: QCD.h:338
QCD::ELECTRON
Definition: QCD.h:312
QCD::orderToString
std::string orderToString(const orders order) const
Converts an object of the enum type "orders" to the corresponding string.
Definition: QCD.cpp:83
Model::UpdateError
bool UpdateError
A boolean set to false if update is successful.
Definition: Model.h:254
QCD::MassOfNf
double MassOfNf(int nf) const
The Mbar mass of the heaviest quark in the theory with Nf active flavour.
Definition: QCD.cpp:620
QCD::mp2mbar_cache
double mp2mbar_cache[5][CacheSize]
Cache for pole mass to msbar mass conversion.
Definition: QCD.h:950
QCD::mtpole
double mtpole
The pole mass of the top quark.
Definition: QCD.h:927
QCD::Mbar2Mp
double Mbar2Mp(const double mbar, const orders order=FULLNNLO) const
Converts the mass to the pole mass.
Definition: QCD.cpp:1221
QCD::mrun_cache
double mrun_cache[10][CacheSize]
Cache for running quark mass.
Definition: QCD.h:949
Model::ModelParamMap
std::map< std::string, std::reference_wrapper< const double > > ModelParamMap
Definition: Model.h:262
QCD::CheckParameters
virtual bool CheckParameters(const std::map< std::string, double > &DPars)
A method to check if all the mandatory parameters for QCD have been provided in model initialization.
Definition: QCD.cpp:335
QCD::Beta0
double Beta0(const double nf) const
The coefficient for a certain number of flavours .
Definition: QCD.cpp:466
QCD::Beta2
double Beta2(const double nf) const
The coefficient for a certain number of flavours .
Definition: QCD.cpp:476
QCD::PostUpdate
virtual bool PostUpdate()
The post-update method for QCD.
Definition: QCD.cpp:143
QCD::muc
double muc
The threshold between four- and three-flavour theory in GeV.
Definition: QCD.h:930
QCD::Als4
double Als4(const double mu) const
The value of at any scale with the number of flavours .
Definition: QCD.cpp:536
Particle::getMass
const double & getMass() const
A get method to access the particle mass.
Definition: Particle.h:61
QCD::BParameterMap
std::map< std::string, BParameter > BParameterMap
Definition: QCD.h:937
QCD::Update
virtual bool Update(const std::map< std::string, double > &DPars)
The update method for QCD.
Definition: QCD.cpp:125
QCD::RHO
Definition: QCD.h:351
QCD::dAdA_NA
double dAdA_NA
Definition: QCD.h:933
QCD::K_star
Definition: QCD.h:348
QCD::B_P
Definition: QCD.h:344
QCD::TOP
Definition: QCD.h:328
QCD::OMEGA
Definition: QCD.h:353
Particle::setMass_scale
void setMass_scale(double mass_scale)
A set method to fix the scale at which the particle mass is defined.
Definition: Particle.h:142
QCD::zeta3
double zeta3
computed with the GSL.
Definition: QCD.h:940
gslpp::pow
complex pow(const complex &z1, const complex &z2)
Definition: gslpp_complex.cpp:395
QCD::Gamma1
double Gamma1(const double nf) const
The coefficient used to compute the running of a mass.
Definition: QCD.cpp:1013
Model::getMissingModelParametersCount
unsigned int getMissingModelParametersCount()
Definition: Model.h:247
Model::raiseMissingModelParameterCount
void raiseMissingModelParameterCount()
Definition: Model.h:242
Model::isSliced
bool isSliced
A boolean set to true if the current istance is a slice of an extended object.
Definition: Model.h:264
QCD::Beta3
double Beta3(const double nf) const
The coefficient for a certain number of flavours .
Definition: QCD.cpp:484
QCD::Gamma0
double Gamma0(const double nf) const
The coefficient used to compute the running of a mass.
Definition: QCD.cpp:1008
QCD::zeta2
double zeta2
computed with the GSL.
Definition: QCD.h:939
QCD::threCorrForMass
double threCorrForMass(const double nf_f, const double nf_i) const
The threshold correction for running of a mass when crossing a flavour threshold.
Definition: QCD.cpp:1025
NNLO
Definition: OrderScheme.h:35
QCD::B_C
Definition: QCD.h:346
QCD::CacheShift
void CacheShift(double cache[][5], int n) const
A member used to manage the caching for this class.
QCD::getQuarks
Particle getQuarks(const QCD::quark q) const
A get method to access a quark as an object of the type Particle.
Definition: QCD.h:534
QCD::optionalParameters
std::map< std::string, double > optionalParameters
A map for containing the list and values of the parameters that are used only by a specific set of ob...
Definition: QCD.h:952
QCD::D_0
Definition: QCD.h:341
QCD::K_P
Definition: QCD.h:340
QCD::BelowTh
double BelowTh(const double mu) const
The active flavour threshold below the scale as defined in QCD::Thresholds().
Definition: QCD.cpp:429
QCD::requireYd
bool requireYd
Switch for generating the Yukawa couplings to the down-type quarks.
Definition: QCD.h:922
QCD::Mp2MbarTMP
double Mp2MbarTMP(double *mu, double *params) const
The member used for finding the numerical solution to the pole mass from the mass.
Definition: QCD.cpp:1263
Particle::setMass
void setMass(double mass)
A set method to fix the particle mass.
Definition: Particle.h:70
QCD::D_star_P
Definition: QCD.h:350
QCD::QCDvars
static std::string QCDvars[NQCDvars]
An array containing the labels under which all QCD parameters are stored in a vector of ModelParamete...
Definition: QCD.h:363
orders
orders
An enum type for orders in QCD.
Definition: OrderScheme.h:31
QCD::computeBs
bool computeBs
Switch for computing from .
Definition: QCD.h:944
QCD::PreUpdate
virtual bool PreUpdate()
The pre-update method for QCD.
Definition: QCD.cpp:116
QCD::Thresholds
double Thresholds(const int i) const
For accessing the active flavour threshold scales.
Definition: QCD.cpp:406
QCD::logLambdaNLO_cache
double logLambdaNLO_cache[9][CacheSize]
Definition: QCD.h:948
QCD::Gamma2
double Gamma2(const double nf) const
The coefficient used to compute the running of a mass.
Definition: QCD.cpp:1018
QCD::P_0
Definition: QCD.h:337
QCD::K_star_P
Definition: QCD.h:349
QCD::NA
double NA
Definition: QCD.h:933
QCD::K_0
Definition: QCD.h:339
QCD::ZeroNf4NLO
double ZeroNf4NLO(double *logLambda4, double *logLambda5_in) const
A member for calculating the difference in across the four-five flavour threshold using AlsWithLambd...
Definition: QCD.cpp:817
QCD::Mrun
double Mrun(const double mu, const double m, const orders order=FULLNNLO) const
Computes a running quark mass from .
Definition: QCD.cpp:1064
QCD::STRANGE
Definition: QCD.h:327
QCD::requireYu
bool requireYu
Switch for generating the Yukawa couplings to the up-type quarks.
Definition: QCD.h:921
QCD::FullOrder
orders FullOrder(orders order) const
Return the FULLORDER enum corresponding to order.
Definition: QCD.cpp:603
QCD::PHI
Definition: QCD.h:347
QCD::computeBd
bool computeBd
Switch for computing from .
Definition: QCD.h:943
QCD::B_D
Definition: QCD.h:343
Model::name
std::string name
The name of the model.
Definition: Model.h:267
QCD::AboveTh
double AboveTh(const double mu) const
The active flavour threshold above the scale as defined in QCD::Thresholds().
Definition: QCD.cpp:420
QCD::NfThresholdCorrections
double NfThresholdCorrections(double mu, double M, double als, int nf, orders order) const
Threshold corrections in matching with from eq. (34) of hep-ph/0512060.
Definition: QCD.cpp:583
QCD::Nc
double Nc
The number of colours.
Definition: QCD.h:932
QCD::ZeroNf3NLO
double ZeroNf3NLO(double *logLambda3, double *logLambda4_in) const
A member for calculating the difference in across the three-four flavour threshold using AlsWithLamb...
Definition: QCD.cpp:823
NLO
Definition: OrderScheme.h:34
QCD::AlsWithLambda
double AlsWithLambda(const double mu, const orders order) const
Computes the running strong coupling in the scheme with the use of .
Definition: QCD.cpp:578
QCD::realorder
orders realorder
Definition: QCD.h:956
QCD::MESON_END
Definition: QCD.h:354
QCD::CA
double CA
Definition: QCD.h:933
QCD::unknownParameters
std::vector< std::string > unknownParameters
A vector for containing the names of the parameters that are not being used but specified in the conf...
Definition: QCD.h:953
QCD::AlsWithInit
double AlsWithInit(const double mu, const double alsi, const double mu_i, const orders order) const
Computes the running strong coupling from in the scheme, where it is forbidden to across a flavour...
Definition: QCD.cpp:498
QCD::CF
double CF
Definition: QCD.h:933
FULLNNLO
Definition: OrderScheme.h:38
QCD::DOWN
Definition: QCD.h:325
QCD::logLambda5_cache
double logLambda5_cache[4][CacheSize]
Definition: QCD.h:947
FULLNLO
Definition: OrderScheme.h:37
QCD::MrunTMP
double MrunTMP(const double mu_f, const double mu_i, const double m, const orders order) const
A function to calculate the running of the mass between flavour thresholds.
Definition: QCD.cpp:1162
QCD::NQCDvars
static const int NQCDvars
The number of model parameters in QCD.
Definition: QCD.h:357
QCD::NEUTRINO_1
Definition: QCD.h:311
QCD::quarks
Particle quarks[6]
The vector of all SM quarks.
Definition: QCD.h:934
QCD::logLambda5
double logLambda5(orders order) const
for .
Definition: QCD.cpp:829
QCD::MU
Definition: QCD.h:314
QCD::mut
double mut
The threshold between six- and five-flavour theory in GeV.
Definition: QCD.h:928
QCD::AlsM
double AlsM
The strong coupling constant at the mass scale MAls, .
Definition: QCD.h:925