a Code for the Combination of Indirect and Direct Constraints on High Energy Physics Models Logo
QCD.h
Go to the documentation of this file.
1/*
2 * Copyright (C) 2012 HEPfit Collaboration
3 *
4 *
5 * For the licensing terms see doc/COPYING.
6 */
7
8#ifndef QCD_H
9#define QCD_H
10
11#include "Model.h"
12#include "Meson.h"
13#include "BParameter.h"
14#include "OrderScheme.h"
15#define MEPS 1.e-10 // mass precision
16
304class QCD : public Model {
305public:
306
310 enum lepton {
317 NOLEPTON
318 };
319
323 enum quark {
329 BOTTOM
330 };
331
332
336 enum meson {
356 MESON_END
357 };
358
359 static const int NQCDvars = 11;
360
365 static std::string QCDvars[NQCDvars];
366
370 QCD();
371
377 const std::string orderToString(const orders order) const;
378
380 // Parameters
381
386 virtual bool Init(const std::map<std::string, double>& DPars);
387
394 virtual bool PreUpdate();
395
403 virtual bool Update(const std::map<std::string, double>& DPars);
404
415 virtual bool PostUpdate();
416
424 virtual bool CheckParameters(const std::map<std::string, double>& DPars);
425
431 void addParameters(std::vector<std::string> params_i);
432
437 void initializeBParameter(std::string name_i) const;
438
443 void initializeMeson(QCD::meson meson_i) const;
444
450 const double getOptionalParameter(std::string name) const
451 {
452 return optionalParameters.at(name);
453 }
454
460 void setOptionalParameter(std::string name, double value)
461 {
462 optionalParameters[name] = value;
463 }
464
469 std::vector<std::string> getUnknownParameters()
470 {
471 return unknownParameters;
472 }
473
475 // Flags
476
483 virtual bool setFlag(const std::string name, const bool value);
484
491 virtual bool setFlagStr(const std::string name, const std::string value);
492
497 virtual bool CheckFlags() const;
498
499
501 // get and set methods for class members
502
507 const double getNc() const
508 {
509 return Nc;
510 }
511
516 void setNc(double Nc)
517 {
518 this->Nc = Nc;
519 }
520
526 const Meson& getMesons(const QCD::meson m) const
527 {
528 return mesonsMap.at(m);
529 }
530
536 const Particle& getQuarks(const QCD::quark q) const
537 {
538 return quarks[q];
539 }
540
546 void setQuarkMass(const quark q, const double mass)
547 {
548 quarks[q].setMass(mass);
549 }
550
555 const double getAlsM() const
556 {
557 return AlsM;
558 }
559
564 const double getMAls() const
565 {
566 return MAls;
567 }
568
573 const double getMut() const
574 {
575 return mut;
576 }
577
582 const double getMub() const
583 {
584 return mub;
585 }
586
591 const double getMuc() const
592 {
593 return muc;
594 }
595
600 const double getMtpole() const
601 {
602 return mtpole;
603 }
604
609 void setMtpole(double mtpole_in)
610 {
611 mtpole = mtpole_in;
612 }
613
618 const double getCF() const
619 {
620 return CF;
621 }
622
629 const BParameter& getBBd() const
630 {
631 return BParameterMap.at("BBd");
632 }
633
640 const BParameter& getBBs() const
641 {
642 return BParameterMap.at("BBs");
643 }
644
652 {
653 return BParameterMap.at("BBd_subleading");
654 }
655
663 {
664 return BParameterMap.at("BBs_subleading");
665 }
666
673 const BParameter& getBD() const
674 {
675 return BParameterMap.at("BD");
676 }
677
684 const BParameter& getBK() const
685 {
686 return BParameterMap.at("BK");
687 }
688
692 const BParameter& getBKd1() const
693 {
694 return BParameterMap.at("BKd1");
695 }
696
700 const BParameter& getBKd3() const
701 {
702 return BParameterMap.at("BKd3");
703 }
704
705
707
714 const double Thresholds(const int i) const;
715
722 const double AboveTh(const double mu) const;
723
730 const double BelowTh(const double mu) const;
731
737 const double Nf(const double mu) const;
738
749 const double NfThresholdCorrections(double mu, double M, double als, int nf, orders order) const;
750
756 const orders FullOrder(orders order) const;
757
759
765 const double Beta0(const double nf) const;
766
772 const double Beta1(const double nf) const;
773
779 const double Beta2(const double nf) const;
780
786 const double Beta3(const double nf) const;
787
800 const double AlsWithInit(const double mu, const double alsi, const double mu_i, const int nf,
801 const orders order) const;
802
811 const double AlsWithLambda(const double mu, const orders order) const;
812
825 const double AlsOLD(const double mu, const orders order = FULLNLO) const;
826 const double Als(const double mu, const orders order = FULLNLO, const bool Nf_thr = true) const;
827
839 const double Als(const double mu, const int Nf_in, const orders order = FULLNLO) const;
840
841 const double AlsByOrder(const double mu, const orders order = FULLNLO, bool Nf_thr = true) const;
842 const double AlsByOrder(const double mu, const int Nf_in, const orders order = FULLNLO) const;
843
850 const double logLambda(const double nf, orders order) const;
851
860 const double Als4(const double mu) const;
861
871 const double Mrun4(const double mu_f, const double mu_i, const double m) const;
872
874
880 const double Gamma0(const double nf) const;
881
887 const double Gamma1(const double nf) const;
888
894 const double Gamma2(const double nf) const;
895
904 const double Mrun(const double mu, const double m, const quark q, const orders order = FULLNNLO) const;
905
915 const double Mrun(const double mu_f, const double mu_i, const double m, const quark q,
916 const orders order = FULLNNLO) const;
917
919
929 const double Mbar2Mp(const double mbar, const quark q, const orders order = FULLNNLO) const;
930
938 const double Mp2Mbar(const double mp, const quark q, orders order = FULLNNLO) const;
939
947 const double Mofmu2Mbar(const double m, const double mu, const quark q) const;
948
956 const double MS2DRqmass(const double MSscale, const double MSbar) const;
957
964 const double MS2DRqmass(const double MSbar) const;
965
967
972 void CacheShift(double cache[][5], int n) const;
973 void CacheShift(int cache[][5], int n) const;
974
979 bool isQCDsuccess() const
980 {
981 return QCDsuccess;
982 }
983
989 {
990 this->computemt = computemt;
991 }
992
993protected:
994
995 mutable bool QCDsuccess=true; // A flag to check if the current calculation is successful
996
1002 virtual void setParameter(const std::string name, const double& value);
1003
1009 const double MassOfNf(int nf) const;
1010
1016
1017 // model parameters
1018 double AlsM;
1019 double MAls;
1020 double mtpole;
1021 double mut;
1022 double mub;
1023 double muc;
1024
1025 double Nc;
1026 double TF,CA,CF,dFdF_NA,dAdA_NA,dFdA_NA,NA; //SU(N)-related quantities
1028
1029private:
1030 mutable std::map<std::string, BParameter> BParameterMap;
1031
1032 double zeta2;
1033 double zeta3;
1034 mutable bool computeFBd;
1035 mutable bool computeFBp;
1036 mutable bool computeBd;
1037 mutable bool computeBs;
1038 static const int CacheSize = 5;
1039 mutable double als_cache[9][CacheSize];
1040 mutable double logLambda5_cache[4][CacheSize];
1041 mutable double logLambdaNLO_cache[9][CacheSize];
1042 mutable double mrun_cache[11][CacheSize];
1043 mutable double mp2mbar_cache[6][CacheSize];
1045 std::map<std::string, double> optionalParameters;
1046 std::vector<std::string> unknownParameters;
1047 mutable std::map<const QCD::meson, Meson> mesonsMap;
1048 bool FlagCsi;
1050
1057 const double Mp2Mbar_pole(const double mp, const orders order = FULLNNLO) const;
1058
1066 const double Mp2Mbar_bar(const double mp, const quark q, const orders order = FULLNNLO) const;
1067
1073 double DeltaMass(double x) const;
1074
1082 const double AlsWithLambda(const double mu, const double logLambda, const orders order) const;
1083
1092 const double ZeroNf6NLO(double *logLambda6, double *logLambda5_in) const;
1093
1101 const double ZeroNf5(double *logLambda5, double *order) const;
1102
1111 const double ZeroNf4NLO(double *logLambda4, double *logLambda5_in) const;
1112
1121 const double ZeroNf3NLO(double *logLambda3, double *logLambda4_in) const;
1122
1128 const double logLambda5(orders order) const;
1129
1137 const double logLambdaNLO(const double nfNEW, const double nfORG, const double logLambdaORG) const;
1138
1149 const double logLambda(const double muMatching, const double mf,
1150 const double nfNEW, const double nfORG,
1151 const double logLambdaORG, orders order) const;
1152
1159 const double threCorrForMass(const double nf_f, const double nf_i) const;
1160
1170 const double MrunTMP(const double mu_f, const double mu_i, const double m, const int nf, const orders order) const;
1171
1179 const double Mp2MbarTMP(double *mu, double *params) const;
1180
1188 const double Mofmu2MbarTMP(double *mu, double *params) const;
1189};
1190
1191#endif /* QCD_H */
std::map< std::string, double > DPars
Definition: Minimal.cpp:11
@ FULLNNLO
Definition: OrderScheme.h:39
@ FULLNLO
Definition: OrderScheme.h:38
A class for the bag parameters.
Definition: BParameter.h:151
A class for mesons.
Definition: Meson.h:315
A class for the template of models.
Definition: Model.h:26
std::string name
The name of the model.
Definition: Model.h:285
A class for particles.
Definition: Particle.h:26
void setMass(double mass)
A set method to fix the particle mass.
Definition: Particle.h:70
A class for parameters related to QCD, hadrons and quarks.
Definition: QCD.h:304
const double getMuc() const
A get method to access the threshold between four- and three-flavour theory in GeV.
Definition: QCD.h:591
double mp2mbar_cache[6][CacheSize]
Cache for pole mass to msbar mass conversion.
Definition: QCD.h:1043
static const int CacheSize
Defines the depth of the cache.
Definition: QCD.h:1038
const double Als4(const double mu) const
The value of at any scale with the number of flavours .
Definition: QCD.cpp:657
const double logLambdaNLO(const double nfNEW, const double nfORG, const double logLambdaORG) const
used for computation of at FULLNLO.
Definition: QCD.cpp:1108
bool isQCDsuccess() const
A getter for the QCDsuccess flag.
Definition: QCD.h:979
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:421
const BParameter & getBK() const
For getting the bag parameters corresponding to the operator basis in process in the meson system.
Definition: QCD.h:684
std::vector< std::string > getUnknownParameters()
A method to get the vector of the parameters that have been specified in the configuration file but n...
Definition: QCD.h:469
meson
An enum type for mesons.
Definition: QCD.h:336
@ OMEGA
Definition: QCD.h:355
@ K_P
Definition: QCD.h:340
@ PHI
Definition: QCD.h:348
@ MESON_END
Definition: QCD.h:356
@ K_star
Definition: QCD.h:349
@ D_P
Definition: QCD.h:342
@ B_C
Definition: QCD.h:347
@ B_P
Definition: QCD.h:345
@ K_star_P
Definition: QCD.h:350
@ P_0
Definition: QCD.h:337
@ RHO_P
Definition: QCD.h:354
@ D_S
Definition: QCD.h:343
@ D_star_P
Definition: QCD.h:352
@ K_0
Definition: QCD.h:339
@ K_S
Definition: QCD.h:351
@ B_D
Definition: QCD.h:344
@ RHO
Definition: QCD.h:353
@ D_0
Definition: QCD.h:341
@ P_P
Definition: QCD.h:338
@ B_S
Definition: QCD.h:346
bool requireYu
Switch for generating the Yukawa couplings to the up-type quarks.
Definition: QCD.h:1012
void addParameters(std::vector< std::string > params_i)
A method to add parameters that are specific to only one set of observables.
Definition: QCD.cpp:199
const 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:1062
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:460
const double Mp2Mbar_pole(const double mp, const orders order=FULLNNLO) const
Converts the top pole mass to the corresponding mass using eq. (16) of hep-ph/0004189.
Definition: QCD.cpp:1683
double mut
The threshold between six- and five-flavour theory in GeV.
Definition: QCD.h:1021
double MAls
The mass scale in GeV at which the strong coupling measurement is provided.
Definition: QCD.h:1019
const BParameter & getBKd1() const
Definition: QCD.h:692
const double Mrun4(const double mu_f, const double mu_i, const double m) const
The running of a mass with the number of flavours .
Definition: QCD.cpp:1531
double CF
Definition: QCD.h:1026
const double MrunTMP(const double mu_f, const double mu_i, const double m, const int nf, const orders order) const
A function to calculate the running of the mass between flavour thresholds.
Definition: QCD.cpp:1489
virtual bool PostUpdate()
The post-update method for QCD.
Definition: QCD.cpp:158
const double getCF() const
A get method to access the Casimir factor of QCD.
Definition: QCD.h:618
std::map< std::string, BParameter > BParameterMap
Definition: QCD.h:1030
const double Beta2(const double nf) const
The coefficient for a certain number of flavours .
Definition: QCD.cpp:611
bool computeBd
Switch for computing from .
Definition: QCD.h:1036
static const int NQCDvars
The number of model parameters in QCD.
Definition: QCD.h:359
double logLambdaNLO_cache[9][CacheSize]
Definition: QCD.h:1041
bool FlagCsi
A flag to determine whether and or (false) and (default, true) are used as inputs.
Definition: QCD.h:1048
void setMtpole(double mtpole_in)
A method to set the pole mass of the top quark.
Definition: QCD.h:609
bool unknownParameterWarning
A flag to stop the unknown parameter warning after the first time.
Definition: QCD.h:1044
void CacheShift(int cache[][5], int n) const
bool FlagMpole2MbarNumeric
A flag to determine whether the pole mass to mass conversion is done numerically.
Definition: QCD.h:1015
double zeta2
computed with the GSL.
Definition: QCD.h:1032
const BParameter & getBKd3() const
Definition: QCD.h:700
void setComputemt(bool computemt)
A set method to change the value of computemt.
Definition: QCD.h:988
const double getOptionalParameter(std::string name) const
A method to get parameters that are specific to only one set of observables.
Definition: QCD.h:450
const BParameter & getBD() const
For getting the bag parameters corresponding to the operator basis in process in the meson system.
Definition: QCD.h:673
double Nc
The number of colours.
Definition: QCD.h:1025
double zeta3
computed with the GSL.
Definition: QCD.h:1033
const Meson & getMesons(const QCD::meson m) const
A get method to access a meson as an object of the type Meson.
Definition: QCD.h:526
double muc
The threshold between four- and three-flavour theory in GeV.
Definition: QCD.h:1023
const 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:1047
const double Gamma0(const double nf) const
The coefficient used to compute the running of a mass.
Definition: QCD.cpp:1271
orders realorder
Definition: QCD.h:1049
const double AlsWithLambda(const double mu, const orders order) const
Computes the running strong coupling in the scheme with the use of .
Definition: QCD.cpp:704
const double getAlsM() const
A get method to access the value of .
Definition: QCD.h:555
double CA
Definition: QCD.h:1026
bool FlagMtPole
A flag to determine whether the pole mass of the top quark is used as input.
Definition: QCD.h:1014
virtual bool PreUpdate()
The pre-update method for QCD.
Definition: QCD.cpp:130
const BParameter & getBBd_subleading() const
For getting the subleading bag parameters in process in the meson system.
Definition: QCD.h:651
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:365
const 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:1617
const double MS2DRqmass(const double MSscale, const double MSbar) const
Converts a quark mass from the scheme to the scheme.
Definition: QCD.cpp:1725
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:1045
std::map< const QCD::meson, Meson > mesonsMap
The map of defined mesons.
Definition: QCD.h:1047
double mrun_cache[11][CacheSize]
Cache for running quark mass.
Definition: QCD.h:1042
virtual bool setFlag(const std::string name, const bool value)
A method to set a flag of QCD.
Definition: QCD.cpp:479
const double logLambda(const double nf, orders order) const
Computes with nf flavours in GeV.
Definition: QCD.cpp:1229
const double Beta1(const double nf) const
The coefficient for a certain number of flavours .
Definition: QCD.cpp:606
const double Gamma1(const double nf) const
The coefficient used to compute the running of a mass.
Definition: QCD.cpp:1276
const double Mrun(const double mu, const double m, const quark q, const orders order=FULLNNLO) const
Computes a running quark mass from .
Definition: QCD.cpp:1353
QCD()
Constructor.
Definition: QCD.cpp:29
void setNc(double Nc)
A set method to change the number of colours .
Definition: QCD.h:516
quark
An enum type for quarks.
Definition: QCD.h:323
@ UP
Definition: QCD.h:324
@ BOTTOM
Definition: QCD.h:329
@ TOP
Definition: QCD.h:328
@ DOWN
Definition: QCD.h:325
@ STRANGE
Definition: QCD.h:327
@ CHARM
Definition: QCD.h:326
const 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:1052
bool computeFBd
Switch for computing from .
Definition: QCD.h:1034
virtual bool setFlagStr(const std::string name, const std::string value)
A method to set a flag of QCD.
Definition: QCD.cpp:511
virtual bool CheckFlags() const
A method to check the sanity of the set of model flags.
Definition: QCD.cpp:517
const BParameter & getBBs_subleading() const
For getting the subleading bag parameters in process in the meson system.
Definition: QCD.h:662
virtual bool Init(const std::map< std::string, double > &DPars)
Initializes the QCD parameters found in the argument.
Definition: QCD.cpp:120
const double Mofmu2Mbar(const double m, const double mu, const quark q) const
Converts a quark running mass at an arbitrary scale to the corresponding mass .
Definition: QCD.cpp:1738
const double Mp2Mbar(const double mp, const quark q, orders order=FULLNNLO) const
Converts a quark pole mass to the corresponding mass .
Definition: QCD.cpp:1672
const double BelowTh(const double mu) const
The active flavour threshold below the scale as defined in QCD::Thresholds().
Definition: QCD.cpp:559
const double logLambda5(orders order) const
for .
Definition: QCD.cpp:1067
double TF
Definition: QCD.h:1026
virtual void setParameter(const std::string name, const double &value)
A method to set the value of a parameter of QCD.
Definition: QCD.cpp:343
const double Beta0(const double nf) const
The coefficient for a certain number of flavours .
Definition: QCD.cpp:601
const double getMut() const
A get method to access the threshold between six- and five-flavour theory in GeV.
Definition: QCD.h:573
double als_cache[9][CacheSize]
Cache for .
Definition: QCD.h:1039
const std::string orderToString(const orders order) const
Converts an object of the enum type "orders" to the corresponding string.
Definition: QCD.cpp:95
const double Mofmu2MbarTMP(double *mu, double *params) const
The member used for finding the numerical solution to the mass from mass.
Definition: QCD.cpp:1730
void setQuarkMass(const quark q, const double mass)
A set method to change the mass of a quark.
Definition: QCD.h:546
const double Beta3(const double nf) const
The coefficient for a certain number of flavours .
Definition: QCD.cpp:618
bool computeFBp
Switch for computing from .
Definition: QCD.h:1035
double dFdF_NA
Definition: QCD.h:1026
bool computemt
Switch for computing the mass of the top quark.
Definition: QCD.h:1011
bool requireYd
Switch for generating the Yukawa couplings to the down-type quarks.
Definition: QCD.h:1013
const double Thresholds(const int i) const
For accessing the active flavour threshold scales.
Definition: QCD.cpp:524
double dAdA_NA
Definition: QCD.h:1026
const double getMtpole() const
A get method to access the pole mass of the top quark.
Definition: QCD.h:600
const double AboveTh(const double mu) const
The active flavour threshold above the scale as defined in QCD::Thresholds().
Definition: QCD.cpp:547
void initializeBParameter(std::string name_i) const
A method to initialize B Parameter and the corresponding meson.
Definition: QCD.cpp:211
bool QCDsuccess
Definition: QCD.h:995
const BParameter & getBBs() const
For getting the bag parameters corresponding to the operator basis in process in the meson system.
Definition: QCD.h:640
const double Mp2Mbar_bar(const double mp, const quark q, const orders order=FULLNNLO) const
Converts the bottom pole mass to the corresponding mass by iteration.
Definition: QCD.cpp:1625
const BParameter & getBBd() const
For getting the bag parameters corresponding to the operator basis in process in the meson system.
Definition: QCD.h:629
const double Nf(const double mu) const
The number of active flavour at scale .
Definition: QCD.cpp:571
const double AlsWithInit(const double mu, const double alsi, const double mu_i, const int nf, const orders order) const
Computes the running strong coupling from in the scheme, where it is forbidden to across a flavour...
Definition: QCD.cpp:627
const 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:709
const orders FullOrder(orders order) const
Return the FULLORDER enum corresponding to order.
Definition: QCD.cpp:728
lepton
An enum type for leptons.
Definition: QCD.h:310
@ NEUTRINO_2
Definition: QCD.h:313
@ NEUTRINO_1
Definition: QCD.h:311
@ MU
Definition: QCD.h:314
@ ELECTRON
Definition: QCD.h:312
@ NOLEPTON
Definition: QCD.h:317
@ NEUTRINO_3
Definition: QCD.h:315
@ TAU
Definition: QCD.h:316
double logLambda5_cache[4][CacheSize]
Definition: QCD.h:1040
const Particle & getQuarks(const QCD::quark q) const
A get method to access a quark as an object of the type Particle.
Definition: QCD.h:536
double dFdA_NA
Definition: QCD.h:1026
const double getNc() const
A get method to access the number of colours .
Definition: QCD.h:507
const double AlsByOrder(const double mu, const orders order=FULLNLO, bool Nf_thr=true) const
Definition: QCD.cpp:804
double DeltaMass(double x) const
A method to compute the correction due to light quarks to the conversion between pole and mass.
Definition: QCD.cpp:1611
Particle quarks[6]
The vector of all SM quarks.
Definition: QCD.h:1027
virtual bool Update(const std::map< std::string, double > &DPars)
The update method for QCD.
Definition: QCD.cpp:136
const double MassOfNf(int nf) const
The Mbar mass of the heaviest quark in the theory with Nf active flavour.
Definition: QCD.cpp:745
const double Als(const double mu, const orders order=FULLNLO, const bool Nf_thr=true) const
Definition: QCD.cpp:762
double mtpole
The pole mass of the top quark.
Definition: QCD.h:1020
const 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:1286
const double Gamma2(const double nf) const
The coefficient used to compute the running of a mass.
Definition: QCD.cpp:1281
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:1046
void initializeMeson(QCD::meson meson_i) const
A method to initialize a meson.
Definition: QCD.cpp:280
double AlsM
The strong coupling constant at the mass scale MAls, .
Definition: QCD.h:1018
const double getMub() const
A get method to access the threshold between five- and four-flavour theory in GeV.
Definition: QCD.h:582
const 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:1057
const double AlsOLD(const double mu, const orders order=FULLNLO) const
Computes the running strong coupling in the scheme. In the cases of LO, NLO and FULLNNLO,...
Definition: QCD.cpp:871
bool computeBs
Switch for computing from .
Definition: QCD.h:1037
double mub
The threshold between five- and four-flavour theory in GeV.
Definition: QCD.h:1022
const double getMAls() const
A get method to access the mass scale at which the strong coupling constant measurement is provided.
Definition: QCD.h:564
const double Mbar2Mp(const double mbar, const quark q, const orders order=FULLNNLO) const
Converts the mass to the pole mass.
Definition: QCD.cpp:1552
double NA
Definition: QCD.h:1026
void CacheShift(double cache[][5], int n) const
A member used to manage the caching for this class.
A class for , the pole mass of the top quark.
Definition: masses.h:164
orders
An enum type for orders in QCD.
Definition: OrderScheme.h:33