Available functions
The dcoumentation of all available functions is listed below in alphabetical order. Documentation can also be found using the search bar on the left, or via the Julia REPL, e.g:
julia> using QCDNUM
julia> ?QCDNUM.evsgns
Module documentation
QCDNUM.EvolutionParams
— TypeEvolutionParams
Struct for holding all QCDNUM Parameters.
QCDNUM.GridParams
— TypeGridParams
Struct for holding the QCDNUM grid parameters.
QCDNUM.InputPDF
— TypeInputPDF(; func::Function, map::AbstractArray{<:Real})
Struct containing all necessary info to pass a PDF (parton distribution function) into QCDNUM.
QCDNUM.SPLINTParams
— TypeSPLINTParams
Struct for storage of parameters used with SPLINT package of QCDNUM.
QCDNUM.SplineAddresses
— TypeSplineAddresses
Lookup table for addresses of different structure function splines.
QCDNUM.WrappedPDF
— Typeconst WrappedPDF = FunctionWrappers.FunctionWrapper{Float64, Tuple{Int32, Float64}}
QCDNUM.WrappedSpFun
— Typeconst WrappedSpFun = FunctionWrappers.FunctionWrapper{Float64, Tuple{Int32, Float64}}
QCDNUM.allfij
— Methodallfij(iset, ix, iq, n, ichk)
Get all flavour-pdf values for given ix
and iq
grid points.
Arguments
iset::Integer
: pdf set id (1-24).ix::Integer
: x grid point.iq::Integer
: q2 grid point.n::Integer
: number of additional pdfs to be returned.ichk::Integer
: flag to steer error checking. See QCDNUM
docs. ichk = -1
makes code faster, at the risk of not checking certain things.
Returns
pdf::Array{Float64,2}
: pdf values over flavours.
QCDNUM.allfxq
— Methodallfxq(iset, x, qmu2, n, ichk)
Get all flavour-pdf values for a given x
and mu^2
.
Arguments
iset::Integer
: pdf set id (1-24).x::Float64
: input value ofx
.qmu2::Float64
: input value ofqmu2
.n::Integer
: number of additional pdfs to be returned.ichk::Integer
: flag to steer error checking. See QCDNUM
docs. ichk = -1
makes code faster, at the risk of not checking certain things.
Returns
pdf::Array{Float64,2}
: pdf values over flavours.
QCDNUM.altabn
— Methodaltabn(iset, iq, n)
Returns the value of (alphaS / 2pi)^n, properly truncated, at the factorisation scale muF^2.
Arguments
iset::Integer
: identifier of tha active alpha_S table (0) or pdf set (1-24).iq::Integer
: index of q2 grid pointn::Integer
: power of alpha_S for the different perturbative series.
Returns
asn::Float64
: value of (alpha_S / 2pi)^n, 0 if error.ierr::Integer
: set, on exit, to 1 ifiq
is close to or below the value of
Lambda^2, and to 2 if iq
is outside the grid boundaries.
QCDNUM.asfunc
— Methodasfunc(r2)
Evolve alpha_S(mu_R^2)
. Does not use mu^2
grid or weight tables.
Returns
alphas::FLoat64
:alpha_S
value.nf::Integer
: number of flavours at scale r2.ierr::Integer
: error code
QCDNUM.bvalij
— Methodbvalij(iset, id, ix, iq, ichk)
Get the value of a basis pdf at a given point (ix
, iq
).
Arguments
iset::Integer
: pdf set id (1-24)id::Integer
: basis pdf identifier from 0 to 12+n,
where n is the number of additional pdfs in iset.
ix::Integer
: x index.iq::Integer
: qq index.ichk::Integer
: flag to steer error checking.
See allfxq().
Returns
pdf::Float64
: pdf values.
QCDNUM.bvalxq
— Methodbvalxq(iset, id, x, qq, ichk)
Get the value of a basis pdf at a given point (x
, qq
).
Arguments
iset::Integer
: pdf set id (1-24)id::Integer
: basis pdf identifier from 0 to 12+n,
where n is the number of additional pdfs in iset.
x::Float64
: x value.qq::Float64
: qq value.ichk::Integer
: flag to steer error checking.
See allfxq().
Returns
pdf::Float64
: pdf values.
QCDNUM.cpypar
— Methodcpypar(iset)
Copy the evolution parameters of a pdf set to a local array. The array has the following param values:
- 1:
iord
- 2:
alfas
- 3:
r2alf
- 4:
nfix
- 5:
q2c
- 6:
q2b
- 7:
q2t
- 8:
ar
- 9:
br
- 10:
xmin
- 11:
qmin
- 12:
qmax
In addition to the evolution parameters is given the pdf type in array(13): 1 = unpolarised, 2 = polarised, 3 = time-like, 4 = external, 5 = user.
Arguments
iset::Integer
: pdf set identifier.
Returns
array::Vector{Float64}
QCDNUM.dmpwgt
— Methoddmpwgt(itype, lun, filename)
Dump weight tables of a given itype
to filename
.
QCDNUM.dsp_funs1
— Methoddsp_funs1(ia, u, ichk)
Evaluate function for 1D spline.
Possible values of ichk
:
-1
: extrapolate the spline0
: return 01
: throw an error message
Arguments
ia::Integer
: address of splineu::Float64
: x or qqichk::Integer
: defines behaviour when outside
spline range
QCDNUM.dsp_funs2
— Methoddsp_funs2(ia, x, q, ichk)
Evaluate function for 2D spline.
Possible values of ichk
:
-1
: extrapolate the spline0
: return 01
: throw an error message
Arguments
ia::Integer
: address of splinex::Float64
: x valueq::Float64
: qq valueichk::Integer
: defines behaviour when outside
spline range
QCDNUM.dsp_ints1
— Methoddsp_ints1(ia, u1, u2)
Evaluate integral of spline between u1 and u2.
The integration limits must lie inside the spline range.
QCDNUM.dsp_ints2
— Methoddsp_ints2(ia, x1, x2, q1, q2, rs, np)
Evaluate integral of spline between x1, x2, q1 and q2. Also takes care of rscut and integrations makes use of N-point Gauss quadrature, as defined by the choice of np.
The integration limits must lie inside the spline range.
QCDNUM.dsp_rscut
— Methoddsp_rscut(ia)
Get the root(s) cut for the spline at ia
.
QCDNUM.dsp_rsmax
— Methoddsp_rsmax(ia, rsc)
Get the root(s) cut limit for the spline at ia
and cut rsc
.
QCDNUM.dsp_spgetval
— Methoddsp_spsgtval(ia, i)
Get some extra info, val
, from the spline ia
.
Arguments
ia::Integer
: spline address.i::Integer
: storage index, runs from 1-100.
Returns
val::Float64
: value to store.
QCDNUM.dsp_uread
— Methoddsp_uread(i)
Read something from the reserved user space.
Arguments
i::Integer
: where to read, from 1 to nuser
Returns
val::Float64
: what is read
QCDNUM.evolfg
— Functionevolfg(itype::Integer, func::Union{Base.CFunction, Ptr{Nothing}}, def::Array{Float64}, iq0::Integer)
evolfg(itype::Integer, func::WrappedPDF, def::Array{Float64}, iq0::Integer)
evolfg(itype::Integer, func::InputPDF, iq0::Integer)
Evolve the flavour pdf set.
Arguments
itype::Integer
: select un-polarised (1), polarised (2) or
time-like (3) evolution.
func::Union{Base.CFunction, Ptr{Nothing}}
: User-defined function
that returns input x * f_j(x)
at iq0
. j
is from 0
to 2 * nf
.
def::Array{Float64}
: input array containing the contribution of
quark species i
to the input distribution j
.
iq0::Integer
: grid index of the starting scalemu_0^2
.
Returns
epsi::Float64
: max deviation of the quadratic spline interpolation
from linear interpolation mid-between grid points.
QCDNUM.evolve
— Methodevolve(input_pdf, qcdnum_params)
High-level interface to QCD evolution with QCDNUM.
QCDNUM.evsgns
— Functionevsgns(itype::Integer, func::Union{Base.CFunction, Ptr{Nothing}}, isns::Array{Int32,1}, n::Integer, iq0::Integer)
evsgns(itype::Integer, func::WrappedPDF, isns::Array{Int32,1}, n::Integer, iq0::Integer)
Evolve an arbitrary set of single/non-singlet pdfs. The evolution can only run in FFNS or MFNS mode, as it is not possible to correctly match at the thresholds as in evolfg.
Arguments
The arguments are as for evolfg, expect def::Array{Float64} is replaced with:
isns::Array{Int32,1}
: Input int array specifing the
evolution type. Entries can be (+1, -1, +-2) corresponding to singlet, valence non-singlet and +/- q_ns singlets respectively.
n::Integer
: Number of singlet/non-singlet pdfs to evolve
Returns
epsi::Float64
: Maximum deviation of the quadratic spline from
linear interpolation mid-between the grid points.
QCDNUM.extpdf
— Methodextpdf(fun, iset, n, offset)
Import a pdfset from an external source.
Arguments
fun::Union{Base.CFunction, Ptr{Nothing}}
: User-defined function with the signature
fun(ipdf::Integer, x::Float64, qq::Float64, first::UInt8)::Float64 specifying the values at x and qq of pdfset ipdf.
iset::Integer
: Pdfset identifier, between 1 and 24.n::Integer
: Number of pdf tables in addition to gluon tables.offset::Float64
: Relative offset at the thresholds mu_h^2, used
to catch matching discontinuities.
Returns
epsi::Float64
: Maximum deviation of the quadratic spline from
linear interpolation mid-between the grid points.
QCDNUM.fflist
— Methodfflist(iset, c, isel, x, q, ichk)
A fast routine to generate a list of interpolated pdfs.
Arguments
iset::Integer
: pdf set identifier [1-24]c::Array{Float64}
: coefficients of quarks/anti-quarksisel::Integer
: selection flagx::Array{Float64}
: list of x values.q::Array{Float64}
: list of q2 values.ichk::Integer
: flag to steer error checking. Seeallfxq()
.
Returns
f::Array{Float64}
: output list of pdf values.
QCDNUM.ffromr
— Methodrfromf(fscale2)
Convert the renormalisation scale, muR^2, to the factorisation scale, muF^2.
QCDNUM.fillwt
— Methodfillwt(itype)
Fill weight tables for all order and number of flavours. itype
is used to select un-polarised pdfs (1), polarised pdfs (2) or fragmentation functions (3).
Returns
nwds::Integer
: number of words used in memory.
QCDNUM.fsplne
— Methodfsplne(iset, id, x, iq)
Spline interpolation of a basis PDF in x, at the grid point iq
. Provided as a diagnostic tool to investigate possible quadratic spline oscillations.
Arguments
iset::Integer
: PDF set identifier [1-24] id::Integer
: Identifier of a basis pdf |e^±| x::Float64
: x value iq::Integer
: Index of a Q^2 grid point
Returns
pdf::Float64
: pdf value
QCDNUM.ftable
— Methodfftabl(ist, c, isel, x, q, n, ichk)
A fast routine to generate a table of interpolated pdfs.
iset::Integer
: pdf set identifier [1-24]c::Array{Float64}
: coefficients of quarks/anti-quarksisel::Integer
: selection flagx::Array{Float64}
: list of x values.q::Array{Float64}
: list of q2 values.ichk::Integer
: flag to steer error checking. Seeallfxq()
.
Returns
table::Matrix{Float64}
: output table of pdf values.
QCDNUM.fvalij
— Methodfvalij(iset, id, ix, iq, ichk)
Get the value of a flavour momentum density at a given point (ix
, iq
).
Arguments
iset::Integer
: pdf set id (1-24)id::Integer
: basis pdf identifier from -6 to 6+n,
where n is the number of additional pdfs in iset.
ix::Integer
: x index.iq::Intger
: qq index.ichk::Integer
: flag to steer error checking.
See allfxq().
Returns
pdf::Float64
: pdf values.
QCDNUM.fvalxq
— Methodfvalxq(iset, id, x, qq, ichk)
Get the value of a flavour momentum density at a given point (x
, qq
).
Arguments
iset::Integer
: pdf set id (1-24)id::Integer
: basis pdf identifier from -6 to 6+n,
where n is the number of additional pdfs in iset.
x::Float64
: x value.qq::Float64
: qq value.ichk::Integer
: flag to steer error checking.
See allfxq().
Returns
pdf::Float64
: pdf values.
QCDNUM.getabr
— Methodgetabr()
Get the relation between the factorisation scale muF^2 and the renormalisation scale muR^2.
muR^2 = ar
muF^2 + br
.
QCDNUM.getalf
— Methodgetalf()
Set the starting value of alpha_S
and the starting renormalisation scale r2
. By default alpha_S(m_Z^2) = 0.118
.
QCDNUM.getcbt
— Methodgetcbt()
Return the current threshold settings for the FFNS or VFNS.
Returns
nfix::Integer
: number of flavours in the FFNS. For VNFS set
nfix = 0
.
qc/b/t::Float64
: q2 values of the heavy flavour thresholds
in the VFNS. Ignored if FFNS.
QCDNUM.getint
— Methodgetint(param)
Get QCDNUM integer parameters.
Arguments
param::String
: Name of parameter. Can be "iter" (number of
evolutions in backwards iteration), "tlmc" (time-like matching conditions), "nopt" (number of perturbative terms) or "edbg" (evolution loop debug printout).
Returns
ival::Int32
: Integer value that is read
QCDNUM.getlim
— Methodgetlim(iset)
Read current grid boundary values for a pdf set iset
.
Returns
xmin::Float64
: min x boundaryqmin::Float64
: min q2 boundaryqmax::Float64
: max q2 boundary
QCDNUM.getord
— Methodgetord()
Get order of perturbative QCD calculations. iord
= 1, 2, 3 for LO, NLO and NNLO respectively. By default, iord = 2
.
QCDNUM.getval
— Methodgetval(param)
Get QCDNUM parameters.
Arguments
param::String
: Name of parameter. Can be "null" (result of
calc that cannot be performed), "epsi" (tolerance level in float comparison |x-y| < epsi), "epsg" (numerical accuracy of Gauss integration in weight table calc), "elim" (allowed diff between quadratic and linear spline interpolation mid-between x grid points - to disable, set elim<=0), "alim" (Max allowed value of alpha_s(mu^2)), "qmin" (smallest possible boundary of mu^2 grid) "qmax" (largest possible boundary of mu^2 grid).
Returns
val::Float64
: Value.
QCDNUM.gqcopy
— Methodgqcopy(n)
Copy the current mu^2
grid into an array of length n
.
QCDNUM.gqmake
— Methodgqmake(qarr, wgt, n, nqin)
Define a logarithmically-spaced mu_F^2
grid on which the parton densities are evolved.
Arguments
qarr::Array{Float64,1}
: input array containingn
values ofmu^2
in ascending order. The lower edge of the grid should be above 0.1 GeV^2.
wgt::Array{Float64,1}
: relative grip point density in each region
defined by qarr
.
n::Integer
: number of values inqarr
andwgt
(n>=2).nqin::Integer
: requested number of grid points.
Returns
nqout::Integer
: number of generated grid points.
QCDNUM.grpars
— Methodgrpars()
Get the current grid definitions.
Returns
nx::Integer
: number of points in x grid. xmi::Float64
: lower boundary of x grid. xma::Float64
: upper boundary of x grid. nq::Integer
: number of points in qq grid. qmi::Float64
: lower boundary of qq grid. qma::Float64
: upper boundary of qq grid. iord::Integer
: order of spline interpolation.
QCDNUM.gxcopy
— Methodgxcopy(n)
Copy the current x
grid into an array of length n
.
QCDNUM.gxmake
— Methodgxmake(xmin, iwt, n, nxin, iord)
Define a logarithmically-spaced x
grid.
Arguments
xmin::Array{Float64,1}
: an input array containitnn
values
of x
in ascending order. xmin[1]
defines the lower end of the grid and other values define approx positions where the point density will change according to the values set in iwt
.
iwt::Array{Int32,1}
: input integer weights given in ascending order
and must always be an integer multiple of the previous weight.
n::Integer
: number of values specified inxmin
andiwt
.nxin::Integer
: Requested number of grid points.iord::Integer
: iord = 2(3) for linear(quadratic) spline interpolation.
Returns
nxout::Integer
: the number of generated grid points.
QCDNUM.ievtyp
— Methodievtype(iset)
Get the pdf evolution type for set iset::Integer
.
QCDNUM.init
— Methodinit()
High-level default initialisation for QCDNUM.
QCDNUM.iqfrmq
— Methodiqfrmq(q2)
Get grid point index of the closest grid point at or below q2
value.
QCDNUM.isp_s2make
— Methodisp_s2make(istepx, istepq)
Create a spline object in memory, and return the address. Every istep-th grid point is taken as a node point of the spline and the grid boundaries are always included as node points.
Arguments
istepx::Integer
: steps taken in sampling the QCDNUM x gridistepq::Integer
: steps taken in sampling the QCDNUM qq grid
QCDNUM.isp_s2user
— Methodisp_s2user(xarr, nx, qarr, nq)
Set your own node points in the spline in case the automatic sampling fails.
The routine will discard points outside the x-qq evolution grid, round the remaining nodes down to the nearest grid-point and then sort them in ascending order, discarding equal values. Thus you are allowed to enter un-sorted scattered arrays.
Arguments
xarr::Array{Float64}
: array of x valuesnx::Integer
: length of xarrqarr::Array{Float64}
: array of qq valuesnq::Integer
: length of qarr
Returns
iasp::Integer
: address of the spline object
QCDNUM.isp_splinetype
— Methodisp_splinetype(ia)
Get the type of spline at address ia
.
Possible types are:
-1
: x spline0
: not a spline1
: qq spline2
: x-qq spline
QCDNUM.isp_spread
— Methodssp_spread(filename)
Read spline from filename
and return address ia
.
QCDNUM.isp_spsize
— Methodisp_spsize(ia)
Get used space for spline ia
, or total memory size when ia = 0
.
QCDNUM.isp_spvers
— Methodisp_spvers()
SPLINT version as a date.
QCDNUM.isp_sqmake
— Methodisp_sqmake(istepq)
Create a 1D qq spline object in memory, and return the address. Every istep-th grid point is taken as a node point of the spline and the grid boundaries are always included as node points.
Arguments
istepq::Integer
: steps taken in sampling the QCDNUM qq grid
QCDNUM.isp_squser
— Methodisp_squser(qarr, nq)
Set your own node points in the 1D qq spline in case the automatic sampling fails.
Arguments
qarr::Array{Float64}
: array of qq valuesnq::Integer
: length of qarr
Returns
iasp::Integer
: address of the spline object
QCDNUM.isp_sxmake
— Methodisp_sxmake(istepx)
Create a 1D x spline object in memory, and return the address. Every istep-th grid point is taken as a node point of the spline and the grid boundaries are always included as node points.
Arguments
istepx::Integer
: steps taken in sampling the QCDNUM x grid
QCDNUM.isp_sxuser
— Methodisp_sxuser(xarr, nx)
Set your own node points in the 1D x spline in case the automatic sampling fails.
Arguments
xarr::Array{Float64}
: array of x valuesnx::Integer
: length of xarr
Returns
iasp::Integer
: address of the spline object
QCDNUM.ixfrmx
— Methodixfrmx(x)
Get grid point index of closest grid point at or below x
value.
QCDNUM.keygrp
— Methodkeygrp(iset, igroup)
Returns the parameter key of the pdf group in iset
. Useful to check if parameters match. Like a more specific keypar.
igroup
can be:
- 1: order
- 2: alpha_S
- 3: fnsandthresholds
- 4: scale
- 5: cutes
- 6: all
QCDNUM.keypar
— Methodkeypar(iset)
Returns the parameter key of the pdf in iset
. Useful to check if parameters match with.
QCDNUM.load_params
— FunctionQCDNUM.load_params(file_name::String)
QCDNUM.load_params(::Type{HDF5.File}, file_name::String)
QCDNUM.load_params(src::HDF5.H5DataStore)
Load stored QCDNUM or SPLINT parameters.
QCDNUM.make_grid
— Methodmake_grid(grid_params)
High-level interface to build QCDNUM grid from GridParams.
QCDNUM.mixfns
— Methodmixfns(nfix, r2c, r2b, r2t)
Select the MFNS mode and set thresholds on mu_R^2.
Arguments
nfix::Integer
: Fixed number of flavours for MFNS. Can be
in the range [3, 6].
r2c/b/t::Float64
: Thresholds defined on the renormalisation scale
mu_R^2 for c, b and t.
QCDNUM.nfrmiq
— Methodnfrmiq(iset, iq)
Returns the number of active flavours nf
at a q2 grid point iq
.
Arguments
iset::Integer
: pdf set identifier.iq::Integer
: q2 grid point.
Returns
nf::Integer
: number of active flavours. 0 ifiq
is
outside the q2 grid.
ithresh::Integer
: threshold indicator that is set to
+1(-1) if iq
is at a threshold with the larger (smaller) number of flavours, 0 otherwise.
QCDNUM.nptabs
— Methodnptabs(iset)
Get the number of pdf tables in set iset::Integer
.
QCDNUM.nwused
— Methodnwused()
Returns the szie nwtot
of the QCDNUM store (the parameter nwf0 in qcdnum.inc) and the number of words used, nwuse
.
The output argument ndummy
is not used at present.
QCDNUM.nxtlun
— Methodnxtlun(lmin)
Get next free logical unit number above max(lmin, 10)
. Returns 0 if there is no free logical unit. Can be called before or after qcinit
. Handy if you want to open a file on a unit that is guaranteed to be free.
QCDNUM.pdfcpy
— Methodpdfcpy(iset1, iset2)
Copy pdf set.
QCDNUM.pullcp
— Methodpushcp()
Pull the current parameters from LIFO stack (load stashed parameters).
QCDNUM.pushcp
— Methodpushcp()
Push the current parameters to LIFO stack (temporarily stash them).
QCDNUM.qcdnum_lock
— Methodqcdnum_lock()::AbstractLock
Get the global lock for QCDNUM operations.
The QCDNUM Fortran package itself is not thread-safe, so QCDNUM.jl uses a global lock do guard against parallel calls to the Fortran library.
QCDNUM.qcinit
— Methodqcinit(lun, filename)
Initialise QCDNUM - should be called before anything else.
Arguments
lun::Integer
: the output logical unit number. When set to 6,
the QCDNUM messages appear on the standard output. When set to -6, the QCDNUM banner printout is suppressed.
filename::String
: the output filename to store log. Irrelevant
when lun
is set to 6/-6
QCDNUM.qfrmiq
— Methodqfrmiq(iq)
Get q2
value at grid point iq
.
QCDNUM.qqatiq
— Methodqqatiq(q, iq)
Check if q2
coincides with a grid point iq
.
Returns
out::Bool
: true
or false
QCDNUM.qstore
— Methodqstore(action, i, val)
QCDNUM reserves 500 words of memory for user use.
Arguments
action::String
: Can be "write", "read", "lock" or "unlock".i::Integer
: Where to read/write in the store.val::Float64
: What to write in the store.
QCDNUM.qstore
— Methodqstore(action, i)
QCDNUM reserves 500 words of memory for user use.
Arguments
action::String
: Can be "write", "read", "lock" or "unlock".i::Integer
: Where to read/write in the store.
Returns
val::Float64
: What is read from the store.
QCDNUM.readwt
— Methodreadwt(lun, filename)
Read weight tables from filename.
QCDNUM.rfromf
— Methodrfromf(fscale2)
Convert the factorisation scale, muF^2, to the renormalisation scale, muR^2.
QCDNUM.save_params
— FunctionQCDNUM.save_params(file_name::String, params::Union{EvolutionParams,SPLINTParams})
QCDNUM.save_params(::Type{HDF5.File}, file_name::String, params::Union{EvolutionParams,SPLINTParams})
QCDNUM.save_params(trg::HDF5.H5DataStore, params::Union{EvolutionParams,SPLINTParams})
Store the QCDNUM or SPLINT parameters for reproducibility.
QCDNUM.setabr
— Methodsetabr(ar, br)
Define the relation between the factorisation scale muF^2 and the renormalisation scale muR^2.
muR^2 = ar
muF^2 + br
.
QCDNUM.setalf
— Methodsetalf(alfs, r2)
Set the starting value of alpha_S
and the starting renormalisation scale r2
. By default alpha_S(m_Z^2) = 0.118
.
QCDNUM.setcbt
— Methodsetcbt(nfix, iqc, iqb, iqt)
Select FFNS or VFNS, and thresholds on mu_F^2
if necessary.
Arguments
nfix::Integer
: number of flavours in the FFNS. For VNFS set
nfix = 0
.
iqc/b/t::Integer
: grid indices of the heavy flavour thresholds
in the VFNS. Ignored if FFNS.
QCDNUM.setint
— Methodsetint(param, ival)
Set QCDNUM integer parameters.
Arguments
param::String
: Name of parameter. Can be "iter" (number of
evolutions in backwards iteration), "tlmc" (time-like matching conditions), "nopt" (number of perturbative terms) or "edbg" (evolution loop debug printout).
ival::Integer
: Value to set.
QCDNUM.setlim
— Methodsetlim(ixmin, iqmin, iqmax)
Restrict the range of a pdf evolution or import to only part of the x-q2 grid.
ixmin
, iqmin
and iqmax
re-define the range of the grid. To release a cut, eneter a value of 0. Fatal error if the cuts result in a kinematic domain that is too small or empty.
QCDNUM.setlun
— Methodsetlun(lun, filename)
Redirect the QCDNUM messages. Arguments are the same as for qcinit
, but it can be called at any point in the program after qcinit
.
Arguments
lun::Integer
: the output logical unit number. When set to 6,
the QCDNUM messages appear on the standard output. When set to -6, the QCDNUM banner printout is suppressed.
filename::String
: the output filename to store log. Irrelevant
when lun
is set to 6/-6
QCDNUM.setord
— Methodsetord(iord)
Set order of perturbative QCD calculations. iord
= 1, 2, 3 for LO, NLO and NNLO respectively. By default, iord = 2
.
QCDNUM.setval
— Methodsetval(param, val)
Set QCDNUM parameters.
Arguments
param::String
: Name of parameter. Can be "null" (result of
calc that cannot be performed), "epsi" (tolerance level in float comparison |x-y| < epsi), "epsg" (numerical accuracy of Gauss integration in weight table calc), "elim" (allowed diff between quadratic and linear spline interpolation mid-between x grid points - to disable, set elim<=0), "alim" (Max allowed value of alpha_s(mu^2)), "qmin" (smallest possible boundary of mu^2 grid) "qmax" (largest possible boundary of mu^2 grid).
val::Float64
: Value to set.
QCDNUM.splchk
— Methodsplchk(iset, id, iq)
Return for a basis pdf at a given Q^2 grid point the maximum deviation between a linear interpolation and the spline interpolation used by QCDNUM.
To be used if having issues with evolfg or extpdf.
Arguments
iset::Integer
: PDF set identifier [1-24] id::Integer
: Identifier of a basis PDF |e^±| iq::Integer
: Index of a Q^2 grid point
Returns
epsi::Float64
: Value of the max deviation
QCDNUM.splint_init
— Methodsplint_init()
High-level interface to splint initialisation.
QCDNUM.ssp_erase
— Methodssp_erase(ia)
Clear the memory from ia
onwards. ia
=0 can also be used to erase all spline objects in memory.
QCDNUM.ssp_extrapu
— Methodssp_extrapu(ia, n)
Define the extrapolation at the kinematic limit for a spline at address ia
.
The extrapolation index n
can be:
0
: constant1
: linear2
: quadratic3
: cubic
QCDNUM.ssp_extrapv
— Methodssp_extrapv(ia, n)
Define the extrapolation at the kinematic limit for a spline at address ia
.
The extrapolation index n
can be:
0
: constant1
: linear2
: quadratic3
: cubic
QCDNUM.ssp_nprint
— Methodssp_nprint(ia)
Print a list of nodes and grids for the spline at address ia
.
QCDNUM.ssp_s2f123
— Methodssp_s2f123(ia, iset, def, istf, rs)
Fast structure function input for 2D splines over x and qq.
Arguments
ia::Integer
: address of the spline objectiset::Integer
: QCDNUM pdf-set indexdef::Array{Float64}
: Array of (anti-)quark
coefficients
istf::Integer
: structure function index
1<=>FL, 2<=>F2, 3<=>xF3, 4<=>FL'
rs::Float64
: sqrt(s) cut - 0 for no kinematic cut
QCDNUM.ssp_s2fill
— Methodssp_s2fill(iasp::Integer, fun::Union{Base.CFunction, Ptr{Nothing}}, rs::Float64)
ssp_s2fill(iasp::Integer, fun::WrappedSpFun, rs::Float64)
Fill the spline object by passing a function. The function must have the signature fun(ix::Integer, iq::Integer, first::Bool).
Arguments
iasp::Integer
: address of the spline objectfun::Union{Base.CFunction, Ptr{Nothing}}
: function to be splinedrs::Float64
: set a sqrt(s) cut - 0 for no kinematic cut
QCDNUM.ssp_spdump
— Methodssp_spdump(ia, filename)
Dump spline at address ia
to filename
.
QCDNUM.ssp_spinit
— Methodssp_spinit(nuser)
Initialise SPLINT - should be called before other SPLINT functions.
Arguments
nuser::Integer
: the number of words reserved for user storage
QCDNUM.ssp_splims
— Methodssp_splims(ia)
Get node limits of spline at address ia
.
Here, u and v refer to the x and qq dimensions respectively.
Returns tuple containing
nu::Integer
: number of nodes in u directionu1::Float64
: lower u limitu2::Float64
: upper u limitnv::Integer
: number of nodes in v directionv1::Float64
: lower v limitv2::Float64
: upper v limitn::Integer
: number of active nodes below kinematic limit
QCDNUM.ssp_spsetval
— Methodssp_spsetval(ia, i, val)
Store some extra info, val
, along with the spline at ia
.
Arguments
ia::Integer
: spline address.i::Integer
: storage index, runs from 1-100.val::Float64
: value to store.
QCDNUM.ssp_sqf123
— Methodssp_sq123(ia, iset, def, istf, ix)
Fast structure function input for splines over qq.
Arguments
ia::Integer
: address of the spline objectiset::Integer
: QCDNUM pdf-set indexdef::Array{Float64}
: Array of (anti-)quark
coefficients
istf::Integer
: structure function index
1<=>FL, 2<=>F2, 3<=>xF3, 4<=>FL'
ix::Integer
: index of x value
QCDNUM.ssp_sqfill
— Methodssp_sqfill(iasp, fun, ix)
Fill the 1D qq spline object by passing a function. The function must have the signature fun(ix::Integer, iq::Integer, first::Bool).
Arguments
iasp::Integer
: address of the spline objectfun::Union{Base.CFunction, Ptr{Nothing}}
: function to be splinedix::Integer
: fixed ix value to pass to fun
QCDNUM.ssp_sxf123
— Methodssp_sxf123(ia, iset, def, istf, iq)
Fast structure function input for splines over x.
Arguments
ia::Integer
: address of the spline objectiset::Integer
: QCDNUM pdf-set indexdef::Array{Float64}
: Array of (anti-)quark
coefficients
istf::Integer
: structure function index
1<=>FL, 2<=>F2, 3<=>xF3, 4<=>FL'
iq::Integer
: index of qq value
QCDNUM.ssp_sxfill
— Methodssp_sxfill(iasp, fun, iq)
Fill the 1D x spline object by passing a function. The function must have the signature fun(ix::Integer, iq::Integer, first::Bool).
Arguments
iasp::Integer
: address of the spline objectfun::Union{Base.CFunction, Ptr{Nothing}}
: function to be splinediq::Integer
: fixed iq value to pass to fun
QCDNUM.ssp_unodes
— Methodssp_unodes(ia, n, nu)
Copy u-nodes from spline at address ia
to local array.
Arguments
ia::Integer
: address of splinen::Integer
: dimension of array to copy tonu::Integer
: number of u-nodes copied
Returns
array::Array{Float64}
: array of u-nodes
QCDNUM.ssp_uwrite
— Methodssp_uwrite(i, val)
Write something to the reserved user space.
Arguments
i::Integer
: where to write, from 1 to nuserval::Float
: what to write
QCDNUM.ssp_vnodes
— Methodssp_vnodes(ia, n, nv)
Copy v-nodes from spline at address ia
to local array.
Arguments
ia::Integer
: address of splinen::Integer
: dimension of array to copy tonv::Integer
: number of v-nodes copied
Returns
array::Array{Float64}
: array of v-nodes
QCDNUM.sumfij
— Methodsumfij(iset, c, isel, ix, iq, ichk)
Return the gluon or a weighted sum of quark densities, depending on the selection flag, isel
.
isel
values:
0: Gluon density |xg> 1: Linear combination c
summed over active flavours 2-8: Specific singlet/non-singlet quark component 9: Intrinsic heavy flavours 12+i: Additional pdf |xf_i> in iset
See QCDNUM manual for more information.
Arguments
iset::Integer
: pdf set id (1-24).c::Array{Float64}
: Coefficients of quarks/anti-quarksisel::Integer
: Selection flagix::Integer
: x grid point.iq::Integer
: q2 grid point.n::Integer
: number of additional pdfs to be returned.ichk::Integer
: flag to steer error checking. See QCDNUM
docs. ichk = -1
makes code faster, at the risk of not checking certain things.
Returns
pdf::Float64
: pdf value.
QCDNUM.sumfxq
— Methodsumfxq(iset, c, isel, x, qmu2, ichk)
Return the gluon or a weighted sum of quark densities, depending on the selection flag, isel
.
isel
values:
0: Gluon density |xg> 1: Linear combination c
summed over active flavours 2-8: Specific singlet/non-singlet quark component 9: Intrinsic heavy flavours 12+i: Additional pdf |xf_i> in iset
See QCDNUM manual for more information.
Arguments
iset::Integer
: pdf set id (1-24).c::Array{Float64}
: Coefficients of quarks/anti-quarksisel::Integer
: Selection flagx::Float64
: input value ofx
.qmu2::Float64
: input value ofqmu2
.n::Integer
: number of additional pdfs to be returned.ichk::Integer
: flag to steer error checking. See QCDNUM
docs. ichk = -1
makes code faster, at the risk of not checking certain things.
Returns
pdf::Float64
: pdf value.
QCDNUM.usepar
— Methodusepar(iset)
Activate the parameters of iset
, ie. copy them to iset = 0
and re-initialise the active look-up tables.
QCDNUM.usrpdf
— Methodusrpdf(fun, iset, n, offset)
Create a user-defined type-5 pdfset (same type as the output of evsgns).
Arguments
fun::Union{Base.CFunction, Ptr{Nothing}}
: User-defined function with the signature
fun(ipdf::Integer, x::Float64, qq::Float64, first::UInt8)::Float64 specifying the values at x and qq of pdfset ipdf.
iset::Integer
: Pdfset identifier, between 1 and 24.n::Integer
: Number of pdf tables in addition to gluon tables.offset::Float64
: Relative offset at the thresholds mu_h^2, used
to catch matching discontinuities.
Returns
epsi::Float64
: Maximum deviation of the quadratic spline from
linear interpolation mid-between the grid points.
QCDNUM.wtfile
— Methodwtfile(itype, filename)
Maintains an up-to-date weight table in filename.
QCDNUM.xfrmix
— Methodxfrmix(ix)
Get x
value at grid point ix
.
QCDNUM.xxatix
— Methodxxatix(x, ix)
Check if x
coincides with a grid point ix
.
Returns
out::Bool
: true
or false
QCDNUM.zmfillw
— Methodzmfillw()
Fill weight tables for zero-mass structure function calculations.
QCDNUM.zmstfun
— Methodzmstfun()
Calculate a structure function from a linear combination of parton densities.
Arguments
istf::Integer
: structure function index
where (1,2,3,4) = (FL, F2, xF3, fL^').
def::Array{Float64}
: coeffs of the quark linear combination
for which the structure function is to be calculated.
x::Array{Float64}
: list ofx
values.Q2::Array{Float64}
: list ofQ2
values.n::Integer
: number of items inx
,Q2
andf
.ichk::Integer
: flag for grid boundary checks. See QCDNUM docs.
Returns
f::Array{Float64}
: list of structure functions.
QCDNUM.zmwords
— Methodzmwords()
Check the number of words available in the ZMSTF workspace and the number of words used.
QCDNUM.@qlccall
— MacroQCDNUM.@qlccall expr
Equivalent to @lock qcdnum_lock() @ccall ...
.