LHAPDF  6.5.4
Public Member Functions | Protected Attributes | List of all members
LHAPDF::AlphaS Class Referenceabstract

Calculator interface for computing alpha_s(Q2) in various ways. More...

#include <AlphaS.h>

Inheritance diagram for LHAPDF::AlphaS:
LHAPDF::AlphaS_Analytic LHAPDF::AlphaS_Ipol LHAPDF::AlphaS_ODE

Public Member Functions

 AlphaS ()
 Base class constructor for default param setup.
 
virtual ~AlphaS ()
 Destructor.
 
alpha_s values
double alphasQ (double q) const
 Calculate alphaS(Q)
 
virtual double alphasQ2 (double q2) const =0
 

Protected Member Functions

Calculating beta function values
double _beta (int i, int nf) const
 
std::vector< double > _betas (int nf) const
 

Protected Attributes

int _qcdorder
 Order of QCD evolution (expressed as number of loops)
 
double _mz
 Mass of the Z boson in GeV.
 
double _alphas_mz
 Value of alpha_s(MZ)
 
double _mreference
 Reference mass in GeV.
 
double _alphas_reference
 Value of alpha_s(reference mass)
 
bool _customref
 Decides whether to use custom reference values or fall back on MZ/AlphaS_MZ.
 
std::map< int, double > _quarkmasses
 
std::map< int, double > _flavorthresholds
 
FlavorScheme _flavorscheme
 The flavor scheme in use.
 
int _fixflav
 The allowed numbers of flavours in a fixed scheme.
 

alpha_s metadata

enum  FlavorScheme { FIXED , VARIABLE }
 Enum of flavor schemes.
 
int numFlavorsQ (double q) const
 Calculate the number of active flavours at energy scale Q.
 
virtual int numFlavorsQ2 (double q2) const
 Calculate the number of active flavours at energy scale Q2.
 
double quarkMass (int id) const
 Get a quark mass by PDG code.
 
void setQuarkMass (int id, double value)
 Set quark masses by PDG code. More...
 
double quarkThreshold (int id) const
 Get a flavor scale threshold by PDG code. More...
 
void setQuarkThreshold (int id, double value)
 Set a flavor threshold by PDG code (= quark masses by default) More...
 
int orderQCD ()
 
void setOrderQCD (int order)
 Set the order of QCD (expressed as number of loops) More...
 
void setMZ (double mz)
 Set the Z mass used in this alpha_s. More...
 
void setAlphaSMZ (double alphas)
 Set the alpha_s(MZ) used in this alpha_s. More...
 
void setMassReference (double mref)
 Set the Z mass used in this alpha_s. More...
 
void setAlphaSReference (double alphas)
 Set the alpha_s(MZ) used in this alpha_s. More...
 
virtual void setLambda (unsigned int, double)
 Set the {i}th Lambda constant for i active flavors. More...
 
virtual std::string type () const =0
 Get the implementation type of this AlphaS.
 
void setFlavorScheme (FlavorScheme scheme, int nf=-1)
 Set flavor scheme of alpha_s solver.
 
FlavorScheme flavorScheme () const
 Get flavor scheme.
 

Detailed Description

Calculator interface for computing alpha_s(Q2) in various ways.

The design of the AlphaS classes is that they are substitutible (cf. polymorphism) and are entirely non-dependent on the PDF and Info objects: hence they can be used by external code that actually doesn't want to do anything at all with PDFs, but which just wants to do some alpha_s interpolation.

Member Function Documentation

◆ _beta()

double LHAPDF::AlphaS::_beta ( int  i,
int  nf 
) const
protected

Calculate the i'th beta function given the number of active flavours Currently limited to 0 <= i <= 3 Calculated using the MSbar scheme

◆ _betas()

std::vector<double> LHAPDF::AlphaS::_betas ( int  nf) const
protected

Calculate a vector of beta functions given the number of active flavours Currently returns a 4-element vector of beta0 – beta3

◆ alphasQ2()

virtual double LHAPDF::AlphaS::alphasQ2 ( double  q2) const
pure virtual

Calculate alphaS(Q2)

Todo:
Throw error in this base method if Q < Lambda?

Implemented in LHAPDF::AlphaS_ODE, LHAPDF::AlphaS_Ipol, and LHAPDF::AlphaS_Analytic.

◆ orderQCD()

int LHAPDF::AlphaS::orderQCD ( )
inline

Get the order of QCD (expressed as number of loops)

Used explicitly in the analytic and ODE solvers.

◆ quarkThreshold()

double LHAPDF::AlphaS::quarkThreshold ( int  id) const

Get a flavor scale threshold by PDG code.

Used in the analytic and ODE solvers.

◆ setAlphaSMZ()

void LHAPDF::AlphaS::setAlphaSMZ ( double  alphas)
inline

Set the alpha_s(MZ) used in this alpha_s.

Used in the ODE solver.

◆ setAlphaSReference()

void LHAPDF::AlphaS::setAlphaSReference ( double  alphas)
inline

Set the alpha_s(MZ) used in this alpha_s.

Used in the ODE solver.

◆ setLambda()

virtual void LHAPDF::AlphaS::setLambda ( unsigned int  ,
double   
)
inlinevirtual

Set the {i}th Lambda constant for i active flavors.

Used in the analytic solver.

Reimplemented in LHAPDF::AlphaS_Analytic.

◆ setMassReference()

void LHAPDF::AlphaS::setMassReference ( double  mref)
inline

Set the Z mass used in this alpha_s.

Used in the ODE solver.

◆ setMZ()

void LHAPDF::AlphaS::setMZ ( double  mz)
inline

Set the Z mass used in this alpha_s.

Used in the ODE solver.

◆ setOrderQCD()

void LHAPDF::AlphaS::setOrderQCD ( int  order)
inline

Set the order of QCD (expressed as number of loops)

Used in the analytic and ODE solvers.

◆ setQuarkMass()

void LHAPDF::AlphaS::setQuarkMass ( int  id,
double  value 
)

Set quark masses by PDG code.

Used in the analytic and ODE solvers.

◆ setQuarkThreshold()

void LHAPDF::AlphaS::setQuarkThreshold ( int  id,
double  value 
)

Set a flavor threshold by PDG code (= quark masses by default)

Used in the analytic and ODE solvers.

Member Data Documentation

◆ _quarkmasses

std::map<int, double> LHAPDF::AlphaS::_quarkmasses
protected

Masses of quarks in GeV Used for working out flavor thresholds and the number of quarks that are active at energy scale Q.


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