Python API Reference
A Higher Order Perturbative Parton Evolution Toolkit
HOPPET is a Fortran package for carrying out DGLAP evolution and other common manipulations of parton distribution functions (PDFs).
Citation: G.P. Salam, J. Rojo, ‘A Higher Order Perturbative Parton Evolution Toolkit (HOPPET)’, Comput. Phys. Commun. 180 (2009) 120-156, arXiv:0804.3755
and
A. Karlberg, P. Nason, G.P. Salam, G. Zanderighi & F. Dreyer, arXiv:2510.09310.
Example:
>>> import hoppet as hp
>>> import numpy as np
>>>
>>> def main():
>>> dy = 0.1
>>> nloop = 3
>>> # Start hoppet
>>> hp.Start(dy, nloop)
>>>
>>> asQ0 = 0.35
>>> Q0 = np.sqrt(2.0)
>>> # Do the evolution.
>>> hp.Evolve(asQ0, Q0, nloop, 1.0, hp.BenchmarkPDFunpol, Q0)
>>>
>>> # Evaluate the PDFs at some x values and print them
>>> xvals = [1e-5,1e-4,1e-3,1e-2,0.1,0.3,0.5,0.7,0.9]
>>> Q = 100.0
>>>
>>> print('')
>>> print(' Evaluating PDFs at Q =',Q, ' GeV')
>>> print(' x u-ubar d-dbar 2(ubr+dbr) c+cbar gluon')
>>> for ix in range(9):
>>> pdf_array = hp.Eval(xvals[ix], Q)
>>> print('{:7.1E} {:11.4E} {:11.4E} {:11.4E} {:11.4E} {:11.4E}'.format(
>>> xvals[ix],
>>> pdf_array[6 + 2] - pdf_array[6 - 2],
>>> pdf_array[6 + 1] - pdf_array[6 - 1],
>>> 2 * (pdf_array[6 - 1] + pdf_array[6 - 2]),
>>> pdf_array[6 - 4] + pdf_array[6 + 4],
>>> pdf_array[6 + 0]
>>> ))
>>>
>>> hp.DeleteAll()
For more examples, see https://github.com/hoppet-code/hoppet/tree/master/example_py
- hoppet.AlphaQED(Q)[source]
Return the QED coupling at scale Q
- Parameters:
Q (float) – Scale in GeV
- Return type:
float
- Returns:
The QED coupling
- hoppet.AlphaS(Q)[source]
Return the strong coupling at scale Q
- Parameters:
Q (float) – Scale in GeV
- Return type:
float
- Returns:
The strong coupling
- hoppet.Assign(callback)[source]
Assign a PDF function to HOPPET without evolution.
- Parameters:
callback (callable) – Python function with signature callback(x, Q) returning array of 13 PDF values
The callback should return PDFs in HOPPET order:
[tbar, bbar, cbar, sbar, ubar, dbar, g, d, u, s, c, b, t]
or if QED evolution is included
[tbar, bbar, cbar, sbar, ubar, dbar, g, d, u, s, c, b, t, EMPTY, photon, electron, muon, tau]
- hoppet.BenchmarkPDFunpol(x, Q)[source]
Evaluate the unpolarized benchmark PDF set.
- Parameters:
x (float) – Longitudinal momentum fraction (0 < x < 1)
Q (float) – Momentum scale [GeV]
- Returns:
Array of 13 benchmark PDF values
- Return type:
list
Useful for testing and validation against known results. For a benchmark PDF with photon and leptons look at example_py/tabulation_example_qed.py
- hoppet.CachedEvolve(callback)[source]
Perform a cached evolution.
- Parameters:
callback (callable) – PDF function with signature callback(x, Q) returning array of 13/18 PDF values
More efficient than
Evolve()when doing multiple evolutions. Needs a call toPreEvolve()first.The callback should return PDFs in HOPPET order:
[tbar, bbar, cbar, sbar, ubar, dbar, g, d, u, s, c, b, t]
or if QED evolution is included
[tbar, bbar, cbar, sbar, ubar, dbar, g, d, u, s, c, b, t, EMPTY, photon, electron, muon, tau]
- hoppet.DeleteAll()[source]
Clean up and free all HOPPET internal arrays and memory.
Call this at the end of your program to free resources.
- hoppet.Eval(x, Q)[source]
Evaluate evolved PDFs at given x and Q.
- Parameters:
x (float) – Longitudinal momentum fraction (0 < x < 1)
Q (float) – Momentum scale [GeV]
- Returns:
Array of 13/18 PDF values in HOPPET order: [tbar, bbar, cbar, sbar, ubar, dbar, g, d, u, s, c, b, t] or if QED evolution is included [tbar, bbar, cbar, sbar, ubar, dbar, g, d, u, s, c, b, t, EMPTY, photon, electron, muon, tau]
- Return type:
list
- hoppet.EvalIFlv(x, Q, iflv)[source]
Return xf(x,Q) for the flavour indicated by iflv, which should be one of the hoppet iflv_* constants (iflv_g, iflv_d, iflv_ubar, etc.)
- Parameters:
x (float) – Longitudinal momentum fraction (0 < x < 1)
Q (float) – Scale in GeV
iflv (int) – One of the hoppet.iflv_g, hoppet.iflv_d, etc.
- Return type:
float
- Returns:
xf(x,Q) for the flavour indicated by iflv
- hoppet.EvalSplit(x, Q, iloop, nf)[source]
Return the value of
[P(iloop,nf) otimes pdf] (x,Q)
where P(iloop,nf) is the iloop-splitting function for the given value of nf, and pdf is our internally stored pdf.
The normalisation is such that the nloop dglap evolution equation is
- dpdf/dlnQ^2 = sum_{iloop=1}^nloop
(alphas/(2*pi))^iloop * P(iloop,nf) otimes pdf
Note that each time nf changes relative to a previous call for the same iloop, the convolution has to be repeated for the whole table. So for efficient results when requiring multiple nf values, calls with the same nf value should be grouped together.
In particular, for repeated calls with the same value of nf, the convolutions are carried out only on the first call (i.e. once for each value of iloop). Multiple calls with different values for iloop can be carried out without problems.
Note that iloop can also be of the form ij or ijk, which means P(i)*P(j)*pdf or P(i)*P(j)*P(k)*pdf. The sum of i+j+k is currently bounded to be <= 4.
The number of loops must be consistent with iloop
- Parameters:
x (float) – Longitudinal momentum fraction (0 < x < 1)
Q (float) – Momentum scale [GeV]
iloop (int) – Perturbative order (1=LO, 2=NLO, etc.)
nf (int) – Number of active flavours
- Returns:
Array of 13/18 PDF values in HOPPET order: [tbar, bbar, cbar, sbar, ubar, dbar, g, d, u, s, c, b, t] or if QED evolution is included [tbar, bbar, cbar, sbar, ubar, dbar, g, d, u, s, c, b, t, EMPTY, photon, electron, muon, tau]
- Return type:
list
- hoppet.Evolve(asQ0, Q0alphas, nloop, muR_Q, callback, Q0pdf)[source]
Evolve PDFs from initial scale Q0 to all scales in the grid.
- Parameters:
asQ0 (float) – Strong coupling at initial scale
Q0alphas (float) – Reference scale for alpha_s [GeV]
nloop (int) – Number of loops for evolution
muR_Q (float) – Renormalization scale ratio (muR/Q)
callback (callable) – PDF function with signature callback(x, Q)
Q0pdf (float) – Initial scale for PDF [GeV]
The callback should return PDFs in HOPPET order:
[tbar, bbar, cbar, sbar, ubar, dbar, g, d, u, s, c, b, t]
or if QED evolution is included
[tbar, bbar, cbar, sbar, ubar, dbar, g, d, u, s, c, b, t, EMPTY, photon, electron, muon, tau]
- hoppet.InitStrFct(order_max, separate_orders, xR, xF)[source]
Initialize structure function calculations.
- Parameters:
order_max (int) – Maximum perturbative order <=4
separate_orders (int) – Whether to separate by order
xR (float) – Renormalization scale factor
xF (float) – Factorization scale factor
- hoppet.InitStrFctFlav(order_max, separate_orders, xR, xF)[source]
Initialize flavour-decomposed structure function calculations.
- Parameters:
order_max (int) – Maximum perturbative order <=4
separate_orders (int) – Whether to separate by order
xR (float) – Renormalization scale factor
xF (float) – Factorization scale factor
- hoppet.PreEvolve(asQ0, Q0alphas, nloop, muR_Q, Q0pdf)[source]
Prepare a cached evolution. Once this has been called, one can use hoppetCachedEvolve (C++) or
CachedEvolve()(Python) to carry out the cached evolution of a specific initial condition.- Parameters:
asQ0 (float) – Strong coupling at initial scale Q0alphas
Q0alphas (float) – Initial scale for alpha_s [GeV]
nloop (int) – Number of loops for evolution (1=LO, 2=NLO, 3=NNLO, 4=N3LO)
muR_Q (float) – Ratio of renormalisation scheme to factorisation scale during evolution
Q0pdf (float) – Initial scale for PDF [GeV]
- hoppet.SetApproximateDGLAPN3LO(splitting_variant)[source]
Arrange for the use of various approximate N3LO splitting functions.
- Parameters:
splitting_variant (int) – One of the hoppet.n3lo_splitting_approximation_* options.
To be called before hoppetStart
- hoppet.SetCoupling(asQ0, Q0alphas, nloop)[source]
Set up the strong coupling such that alphas(Q)=alphas_Q, with the given number of loops (nloop).
The user should have set the quark masses or requested a FFN scheme prior to calling this function.
This function is provided mainly for use in conjunction with hoppetAssign (C++)
Assign()(Python). In particular, it has the side effect of modifying the structure of the PDF tables to make sure they know about the mass thresholds.If QED has been requested, a QED coupling will also be set up (its value is not currently configurable from this interface).
If you call hoppetEvolve (C++)
Evolve()(Python), there is no need to separately call this routine.- Parameters:
asQ0 (float) – alphas at the scale Q0
Q0alphas (float) – the scale Q0
nloop (int) – Perturbative order (1: LO, 2: NLO, 3: NNLO, 4: N3LO)
- hoppet.SetExactDGLAP(exact_nfthreshold, exact_splitting)[source]
Arrange for the use of exact NNLO splitting and mass-threshold functions.
- Parameters:
exact_nfthreshold (boolean) – If True use the exact NNLO mass threshold functions. If False use the faster parametrisations
exact_splitting (boolean) – If True use the exact NNLO splitting functions. If False use the faster parametrisations
To be called before hoppetStart
- hoppet.SetFFN(fixed_nf)[source]
Set things up to be a fixed-flavour number scheme with the given fixed_nf number of flavours
- Parameters:
fixed_nf (int) – Value of fixed number of flavours (e.g. 5)
- hoppet.SetMSbarMassVFN(mc, mb, mt)[source]
Set things up to be a variable-flavour number scheme with the given quark (MSbar) masses. Thresholds are crossed at the MSbar masses, both for the coupling and the PDF evolution.
- Parameters:
mc (float) – Charm mass
mb (float) – Bottom mass
mt (float) – Top mass
- hoppet.SetN3LOnfthresholds(variant)[source]
Arrange for the use of N3LO mass thresholds or not.
- Parameters:
variant (int) – One of the hoppet.n3lo_nfthresholds_* options.
To be called before hoppetStart
- hoppet.SetPoleMassVFN(mc, mb, mt)[source]
Set things up to be a variable-flavour number scheme with the given quark (pole) masses. Thresholds are crossed at the pole masses, both for the coupling and the PDF evolution.
- Parameters:
mc (float) – Charm mass
mb (float) – Bottom mass
mt (float) – Top mass
- hoppet.SetQED(withqed, qcdqed, plq)[source]
Enable or disable QED evolution with photon and lepton PDFs.
- Parameters:
withqed (int) – 1 to enable QED evolution, 0 to disable
qcdqed (int) – Treatment of QCD-QED coupling (implementation-specific)
plq (int) – 1 to enable lepton splitting function
Note – Must be called before
Start()as it modifies global PDF settings.
- hoppet.SetSplittingN3LO(splitting_variant)[source]
Change the variant for the N3LO splitting functions.
- Parameters:
splitting_variant (int) – One of the hoppet.n3lo_splitting_* options.
To be called before hoppetStart
- hoppet.SetSplittingNNLO(splitting_variant)[source]
Arrange for the use of various NNLO splitting functions.
- Parameters:
splitting_variant (int) – One of the hoppet.nnlo_splitting_* options.
To be called before hoppetStart
- hoppet.SetVFN(mc, mb, mt)[source]
Set things up to be a variable-flavour number scheme with the given quark (pole) masses. Now deprecated; use hoppetSetPoleMassVFN instead, which is what is being called undder the hood.
- Parameters:
mc (float) – Charm mass
mb (float) – Bottom mass
mt (float) – Top mass
- hoppet.SetYLnlnQInterpOrders(yorder, lnlnQorder)[source]
Override the default interpolation order in y and lnlnQ.
- Parameters:
yorder (int) – The interpolation order in y (default: 5)
lnlnQorder (int) – The interpolation order in lnlnQ (default: 4)
2 corresponds to quadratic, 3 to cubic etc.
- hoppet.Start(dy, nloop)[source]
Initialize HOPPET evolution tables with basic parameters.
- Parameters:
dy (float) – Step size in ln(1/x). Typical range: 0.05-0.2
nloop (int) – Maximum number of loops (1=LO, 2=NLO, 3=NNLO, 4=N3LO)
Note – This is the basic initialization. For more control, use
:param
StartExtended().:
- hoppet.StartExtended(ymax, dy, Qmin, Qmax, dlnlnQ, nloop, order, factscheme)[source]
Initialize HOPPET evolution tables with extended parameters.
- Parameters:
ymax (float) – Highest value of ln(1/x) user wants to access
dy (float) – Internal ln(1/x) grid spacing (0.05-0.2 is sensible)
Qmin (float) – Lower limit of Q range [GeV]
Qmax (float) – Upper limit of Q range [GeV]
dlnlnQ (float) – Internal table spacing in ln(ln(Q)) (e.g. dy/4)
nloop (int) – Maximum number of loops (1=LO, 2=NLO, 3=NNLO, 4=N3LO)
order (int) – Order of numerical interpolation (e.g. -6) factscheme (int): One of the hoppet.factscheme_* factorization scheme identifier (e.g. hoppet.factscheme_MSbar)
This provides full control over the evolution grid setup.
- hoppet.StartStrFct(order_max)[source]
Minimal setup of structure functions
- Parameters:
order_max (int) – highest order in QCD to compute (1: LO, 2: NLO, 3: NNLO, 4: N3LO)
- hoppet.StartStrFctExtended(order_max, nflav, scale_choice, constant_mu, param_coefs, wmass, zmass)[source]
Setup of constants and parameters needed for structure functions
- Parameters:
order_max (int) – highest order in QCD to compute (1: LO, 2: NLO, 3: NNLO, 4: N3LO)
nflav (int) – integer number of flavours (if negative use variable flavour)
scale_choice (int) – (0: fixed scale, 1: use Q, 2: use arbitrary scale)
constant_mu (float) – if scale_choice = scale_choice_fixed (= 0) then this is the fixed scale
param_coefs (boolean) – if .true. use parametrised coefficients functions
wmass (float) – Mass of the W boson
zmass (float) – Mass of the z boson
- hoppet.StrFct(x, Q, muR_in, muF_in)[source]
Calculate structure functions with specified scales.
- Parameters:
x (float) – Longitudinal momentum fraction
Q (float) – Hard scale [GeV]
muR_in (float) – Renormalization scale [GeV]
muF_in (float) – Factorization scale [GeV]
- Returns:
Array of structure function values where the indices correspond to the hoppet.iF1Wp etc. constants
- Return type:
list
- hoppet.StrFctFlav(x, Q, muR_in, muF_in, flav)[source]
Calculate flavour-decomposed structure functions.
- Parameters:
x (float) – Longitudinal momentum fraction
Q (float) – Hard scale [GeV]
muR_in (float) – Renormalization scale [GeV]
muF_in (float) – Factorization scale [GeV]
flav (int) – Flavour index (one of hoppet.iflv* constants)
- Returns:
Array of structure function values for specified flavour where indices correspond to FL, F2, F3
- Return type:
list
- hoppet.StrFctLO(x, Q, muR_in, muF_in)[source]
Calculate leading-order structure functions.
- Parameters:
x (float) – Longitudinal momentum fraction
Q (float) – Hard scale [GeV]
muR_in (float) – Renormalization scale [GeV]
muF_in (float) – Factorization scale [GeV]
- Returns:
Array of LO structure function values where the indices correspond to the hoppet.iF1Wp etc. constants
- Return type:
list
- hoppet.StrFctLOFlav(x, Q, muR_in, muF_in, flav)[source]
Calculate LO flavour-decomposed structure functions.
- Parameters:
x (float) – Longitudinal momentum fraction
Q (float) – Hard scale [GeV]
muR_in (float) – Renormalization scale [GeV]
muF_in (float) – Factorization scale [GeV]
flav (int) – Flavour index
- Returns:
Array of LO structure function values for specified flavour where indices correspond to FL, F2, F3
- Return type:
list
- hoppet.StrFctN3LO(x, Q, muR_in, muF_in)[source]
Calculate next-to-next-to-next-to-leading-order structure functions.
- Parameters:
x (float) – Longitudinal momentum fraction
Q (float) – Hard scale [GeV]
muR_in (float) – Renormalization scale [GeV]
muF_in (float) – Factorization scale [GeV]
- Returns:
Array of N3LO structure function values where the indices correspond to the hoppet.iF1Wp etc. constants
- Return type:
list
- hoppet.StrFctNLO(x, Q, muR_in, muF_in)[source]
Calculate next-to-leading-order structure functions.
- Parameters:
x (float) – Longitudinal momentum fraction
Q (float) – Hard scale [GeV]
muR_in (float) – Renormalization scale [GeV]
muF_in (float) – Factorization scale [GeV]
- Returns:
Array of NLO structure function values where the indices correspond to the hoppet.iF1Wp etc. constants
- Return type:
list
- hoppet.StrFctNLOFlav(x, Q, muR_in, muF_in, flav)[source]
Calculate NLO flavour-decomposed structure functions.
- Parameters:
x (float) – Longitudinal momentum fraction
Q (float) – Hard scale [GeV]
muR_in (float) – Renormalization scale [GeV]
muF_in (float) – Factorization scale [GeV]
flav (int) – Flavour index
- Returns:
Array of NLO structure function values for specified flavour where indices correspond to FL, F2, F3
- Return type:
list
- hoppet.StrFctNNLO(x, Q, muR_in, muF_in)[source]
Calculate next-to-next-to-leading-order structure functions.
- Parameters:
x (float) – Longitudinal momentum fraction
Q (float) – Hard scale [GeV]
muR_in (float) – Renormalization scale [GeV]
muF_in (float) – Factorization scale [GeV]
- Returns:
Array of NNLO structure function values where the indices correspond to the hoppet.iF1Wp etc. constants
- Return type:
list
- hoppet.StrFctNoMu(x, Q)[source]
Calculate structure functions with default scale choices.
- Parameters:
x (float) – Longitudinal momentum fraction
Q (float) – Hard scale [GeV] (used for both muR and muF)
- Returns:
Array of structure function values where the indices correspond to the hoppet.iF1Wp etc. constants
- Return type:
list
- hoppet.StrFctNoMuFlav(x, Q, flav)[source]
Calculate flavour-decomposed structure functions with default scales.
- Parameters:
x (float) – Longitudinal momentum fraction
Q (float) – Hard scale [GeV]
flav (int) – Flavour index
- Returns:
Array of structure function values for specified flavour where indices correspond to FL, F2, F3
- Return type:
list
- hoppet.WriteLHAPDFGrid(basename, pdf_index)[source]
Write out the contents of tables(0) (assumed to be the PDF) in the LHAPDF format
- Parameters:
basename (string) – The basename of the output
pdf_index (int) – The intended index in the LHAPDF meaning. If index is 0 both the PDF and the .info file is saved to disk