eko package

eko.run_dglap(theory_card, operators_card)[source]

This function takes a DGLAP theory configuration dictionary and performs the solution of the DGLAP equations.

The EKO \(\mathbf E_{k,j}(a_s^1\leftarrow a_s^0)\) is determined in order to fullfill the following evolution

\[\mathbf f(x_k,a_s^1) = \mathbf E_{k,j}(a_s^1\leftarrow a_s^0) \mathbf f(x_j,a_s^0)\]
Parameters

setup (dict) – input card - see Input & Output

Returns

output – output dictionary - see Input & Output

Return type

dict

Subpackages

Submodules

eko.basis_rotation module

This module contains the definitions of the Flavor Basis and Evolution Basis.

eko.basis_rotation.ad_projector(ad_lab, nf)[source]

Build a projector (as a numpy array) for the given anomalous dimension sector.

Parameters
  • ad_lab (str) – name of anomalous dimension sector

  • nf (int) – number of light flavors

Returns

proj – projector over the specified sector

Return type

np.ndarray

eko.basis_rotation.ad_projectors(nf)[source]

Build projectors tensor (as a numpy array), collecting all the individual sector projectors.

Parameters

nf (int) – number of light flavors

Returns

projs – projectors tensor

Return type

np.ndarray

eko.basis_rotation.anomalous_dimensions_basis = ('S_qq', 'S_qg', 'S_gq', 'S_gg', 'NS_m', 'NS_p', 'NS_v')

Sorted elements in Anomalous Dimensions Basis as str.

eko.basis_rotation.evol_basis = ('ph', 'S', 'g', 'V', 'V3', 'V8', 'V15', 'V24', 'V35', 'T3', 'T8', 'T15', 'T24', 'T35')

Sorted elements in Evolution Basis as str.

Definition: here.

corresponding PDF : \(\gamma, \Sigma, g, V, V_{3}, V_{8}, V_{15}, V_{24}, V_{35}, T_{3}, T_{8}, T_{15}, T_{24}, T_{35}\)

eko.basis_rotation.evol_basis_pids = (22, 100, 21, 200, 203, 208, 215, 224, 235, 103, 108, 115, 124, 135)

PID representation of evol_basis.

eko.basis_rotation.flavor_basis_names = ('ph', 'tbar', 'bbar', 'cbar', 'sbar', 'ubar', 'dbar', 'g', 'd', 'u', 's', 'c', 'b', 't')

String representation of flavor_basis_pids.

eko.basis_rotation.flavor_basis_pids = (22, -6, -5, -4, -3, -2, -1, 21, 1, 2, 3, 4, 5, 6)

Sorted elements in Flavor Basis as PID.

Definition: here

corresponding PDF : \(\gamma, \bar t, \bar b, \bar c, \bar s, \bar u, \bar d, g, d, u, s, c, b, t\)

eko.basis_rotation.map_ad_to_evolution = {'NS_m': ['V3.V3', 'V8.V8', 'V15.V15', 'V24.V24', 'V35.V35'], 'NS_p': ['T3.T3', 'T8.T8', 'T15.T15', 'T24.T24', 'T35.T35'], 'NS_v': ['V.V'], 'S_gg': ['g.g'], 'S_gq': ['g.S'], 'S_qg': ['S.g'], 'S_qq': ['S.S']}

Map anomalous dimension sectors’ names to their members

eko.basis_rotation.rotate_flavor_to_evolution = array([[ 1,  0,  0,  0,  0,  0,  0,  0,  0,  0,  0,  0,  0,  0],        [ 0,  1,  1,  1,  1,  1,  1,  0,  1,  1,  1,  1,  1,  1],        [ 0,  0,  0,  0,  0,  0,  0,  1,  0,  0,  0,  0,  0,  0],        [ 0, -1, -1, -1, -1, -1, -1,  0,  1,  1,  1,  1,  1,  1],        [ 0,  0,  0,  0,  0, -1,  1,  0, -1,  1,  0,  0,  0,  0],        [ 0,  0,  0,  0,  2, -1, -1,  0,  1,  1, -2,  0,  0,  0],        [ 0,  0,  0,  3, -1, -1, -1,  0,  1,  1,  1, -3,  0,  0],        [ 0,  0,  4, -1, -1, -1, -1,  0,  1,  1,  1,  1, -4,  0],        [ 0,  5, -1, -1, -1, -1, -1,  0,  1,  1,  1,  1,  1, -5],        [ 0,  0,  0,  0,  0,  1, -1,  0, -1,  1,  0,  0,  0,  0],        [ 0,  0,  0,  0, -2,  1,  1,  0,  1,  1, -2,  0,  0,  0],        [ 0,  0,  0, -3,  1,  1,  1,  0,  1,  1,  1, -3,  0,  0],        [ 0,  0, -4,  1,  1,  1,  1,  0,  1,  1,  1,  1, -4,  0],        [ 0, -5,  1,  1,  1,  1,  1,  0,  1,  1,  1,  1,  1, -5]])

Basis rotation matrix between Flavor Basis and Evolution Basis.

eko.beta module

This module contains the QCD beta function coefficients.

See pQCD ingredients.

eko.beta.b(k, nf)[source]

Compute b coefficient.

Parameters
  • k (int) – perturbative order

  • nf (int) – number of active flavors

Returns

b – b_k(nf)

Return type

float

eko.beta.beta(k, nf)[source]

Compute value of a beta coefficients

Parameters
  • k (int) – perturbative order

  • nf (int) – number of active flavors

Returns

beta – beta_k(nf)

Return type

float

eko.beta.beta_0(nf: int)[source]

Computes the first coefficient of the QCD beta function.

Implements Eq. (3.1) of [HRU+17].

Parameters

nf (int) – number of active flavors

Returns

beta_0 – first coefficient of the QCD beta function \(\beta_0^{n_f}\)

Return type

float

eko.beta.beta_1(nf: int)[source]

Computes the second coefficient of the QCD beta function.

Implements Eq. (3.2) of [HRU+17].

Parameters

nf (int) – number of active flavors

Returns

beta_1 – second coefficient of the QCD beta function \(\beta_1^{n_f}\)

Return type

float

eko.beta.beta_2(nf: int)[source]

Computes the third coefficient of the QCD beta function

Implements Eq. (3.3) of [HRU+17].

Parameters

nf (int) – number of active flavors

Returns

beta_2 – third coefficient of the QCD beta function \(\beta_2^{n_f}\)

Return type

float

eko.constants module

This files sets the physical constants.

The constants are:

  • \(N_C\) the number of colors - defaults to \(3\)

  • \(T_R\) the normalization of fundamental generators - defaults to \(1/2\)

  • \(C_A\) second Casimir constant in the adjoint representation - defaults to \(N_C = 3\)

  • \(C_F\) second Casimir constant in the fundamental representation - defaults to \(\frac{N_C^2-1}{2N_C} = 4/3\)

eko.constants.CA = 3.0

second Casimir constant in the adjoint representation

eko.constants.CF = 1.3333333333333333

second Casimir constant in the fundamental representation

eko.constants.NC = 3

the number of colors

eko.constants.TR = 0.5

the normalization of fundamental generators

eko.gamma module

This module contains the QCD gamma function coefficients.

See pQCD ingredients.

eko.gamma.gamma(order, nf)[source]

Compute the value of a gamma coefficient

Parameters
  • order (int) – perturbative order

  • nf (int) – number of active flavors

Returns

gamma – QCD gamma function coefficient

Return type

float

eko.gamma.gamma_0()[source]
eko.gamma.gamma_1(nf)[source]
eko.gamma.gamma_2(nf)[source]
eko.gamma.gamma_3(nf)[source]

eko.interpolation module

Library providing all necessary tools for PDF interpolation.

This library also provides a class to generate the interpolator InterpolatorDispatcher. Upon construction the dispatcher generates a number of functions to evaluate the interpolator.

class eko.interpolation.Area(lower_index, poly_number, block, xgrid)[source]

Bases: object

Class that define each of the area of each of the subgrid interpolators.

Upon construction an array of coefficients is generated.

Parameters
  • lower_index (int) – lower index of the area

  • poly_number (int) – number of polynomial

  • block (tuple(int, int)) – kmin and kmax

  • xgrid (numpy.ndarray) – Grid in x-space from which the interpolators are constructed

class eko.interpolation.BasisFunction(xgrid, poly_number, list_of_blocks, mode_log=True, mode_N=True)[source]

Bases: object

Object containing a list of areas for a given polynomial number defined by (xmin-xmax) and containing a list of coefficients.

Upon construction will generate all areas and generate and compile a function to evaluate in N (or x) the interpolator

Parameters
  • xgrid (numpy.ndarray) – Grid in x-space from which the interpolators are constructed

  • poly_number (int) – number of polynomial

  • list_of_blocks (list(tuple(int, int))) – list of tuples with the (kmin, kmax) values for each area

  • mode_log (bool) – use logarithmic interpolation?

  • mode_N (bool) – if true compiles the function on N, otherwise compiles x

areas_to_const()[source]

Returns an array containing all areas with (xmin, xmax, numpy.array of coefficients)

Returns

area config

Return type

numpy.ndarray

compile_n()[source]

Compiles the function to evaluate the interpolator in N space.

Generates a function evaluate_Nx() with a (N, logx) signature.

\[\tilde p(N)*\exp(- N * \ln(x))\]

The polynomials contain naturally factors of \(\exp(N * j * \ln(x_{min/max}))\) which can be joined with the Mellin inversion factor.

compile_x()[source]

Compiles the function to evaluate the interpolator in x space

evaluate_x(x)[source]

Evaluate basis function in x-space (regardless of the true space).

Parameters

x (float) – evaluated point

Returns

res – p(x)

Return type

float

is_below_x(x)[source]

Are all areas below x?

Parameters

x (float) – reference value

Returns

is_below_x – xmax of highest area <= x?

Return type

bool

class eko.interpolation.InterpolatorDispatcher(xgrid, polynomial_degree, log=True, mode_N=True)[source]

Bases: object

Setups the interpolator.

Upon construction will generate a list of BasisFunction objects. Each of these BasisFunction objects expose a callable method (also accessible as the __call__ method of the class) which will be numba-compiled.

Parameters
  • xgrid_in (numpy.ndarray) – Grid in x-space from which the interpolators are constructed

  • polynomial_degree (int) – degree of the interpolation polynomial

  • log (bool) – Whether it is a log or linear interpolator

  • mode_N (bool) – if true compiles the function on N, otherwise compiles x

classmethod from_dict(operators_card, mode_N=True)[source]

Create object from dictionary.

operators runcard parameters

Name

Type

default

description

interpolation_xgrid

list(float)

[required]

the interpolation grid

interpolation_polynomial_degree

int

4

polynomial degree of the interpolating function

interpolation_is_log

bool

True

use logarithmic interpolation?

Parameters

operators_card (dict) – input configurations

get_interpolation(targetgrid)[source]

Computes interpolation matrix between targetgrid and xgrid.

\[f(targetgrid) = R \cdot f(xgrid)\]
Parameters

targetgrid (array) – grid to interpolate to

Returns

R – interpolation matrix, do be multiplied from the left(!)

Return type

array

to_dict()[source]

Returns the configuration for the underlying xgrid (from which the instance can be created again).

The output dictionary contains a numpy array.

Returns

ret – full grid configuration

Return type

dict

eko.interpolation.evaluate_Nx(N, logx, area_list)[source]

Evaluates a single linear Lagrange interpolator in N-space multiplied by the Mellin-inversion factor.

\[\tilde p(N)*\exp(- N * \ln(x))\]
Parameters
  • N (complex) – Mellin variable

  • logx (float) – logarithm of inversion point

  • area_list (list) – area configuration of basis function

Returns

res – basis function * inversion factor

Return type

float

eko.interpolation.evaluate_x(x, area_list)[source]

Get a single linear Lagrange interpolator in x-space

\[p(x)\]
Parameters
  • x (float) – interpolation point

  • area_list (list) – area configuration of basis function

Returns

res – basis function(x)

Return type

float

eko.interpolation.log_evaluate_Nx(N, logx, area_list)[source]

Evaluates a single logarithmic Lagrange interpolator in N-space multiplied by the Mellin-inversion factor.

\[\tilde p(N)*\exp(- N * \ln(x))\]
Parameters
  • N (complex) – Mellin variable

  • logx (float) – logarithm of inversion point

  • area_list (list) – area configuration of basis function

Returns

res – kernel * inversion factor

Return type

float

eko.interpolation.log_evaluate_x(x, area_list)[source]

Get a single logarithmic Lagrange interpolator in x-space

\[p(x)\]
Parameters
  • x (float) – interpolation point

  • area_list (list) – area configuration of basis function

Returns

res – basis function(x)

Return type

float

eko.interpolation.make_grid(n_low, n_mid, n_high=0, x_min=1e-07, x_low=0.1, x_high=0.9, x_high_max=0.9999)[source]

Creates a log-lin-log-spaced grid.

1.0 is always part of the grid and the final grid is unified, i.e. esp. points that might appear twice at the borders of the regions are unified.

Parameters
  • n_low (int) – points in the small-x region

  • n_mid (int) – points in the medium-x region

  • n_high (int) – points in the large-x region

  • x_min (float) – minimum x (included)

  • x_low (float) – seperation point between small and medium

  • x_high_max (float) – closest point before 1

Returns

xgrid – generated grid

Return type

numpy.ndarray

eko.interpolation.make_lambert_grid(n_pts, x_min=1e-07, x_max=1.0)[source]

Creates a smoothly spaced grid that is linear near 1 and logarithmic near 0.

It is generated by the relation:

\[x(y) = \frac{1}{5} W_{0}(5 \exp(5-y))\]

where \(W_{0}\) is the principle branch of the Lambert W function and \(y\) is a variable which extremes are given as function of \(x\) by the direct relation:

\[y(x) = 5(1-x)-\log(x)\]

This method is implemted in PineAPPL, [CNSZ20] eq 2.11 and relative paragraph.

Parameters
  • n_pts (int) – number of points

  • x_min (float) – minimum x (included)

  • x_max (float) – maximum x (included)

Returns

xgrid – generated grid

Return type

numpy.ndarray

eko.mellin module

This module contains the implementation of the inverse Mellin transformation.

Although this module provides three different path implementations in practice only the Talbot path [AV04]

\[p_{\text{Talbot}}(t) = o + r \cdot ( \theta \cot(\theta) + i\theta)\quad \text{with}~\theta = \pi(2t-1)\]

is used, as it results in the most efficient convergence. The default values for the parameters \(r,o\) are given by \(r = 1/2, o = 0\) for the non-singlet integrals and by \(r = \frac{2}{5} \frac{16}{1 - \ln(x)}, o = 1\) for the singlet sector. Note that the non-singlet kernels evolve poles only up to \(N=0\) whereas the singlet kernels have poles up to \(N=1\).

eko.mellin.Talbot_jac(t, r, o)[source]

Derivative of Talbot path.

\[\frac{dp_{\text{Talbot}}(t)}{dt}\]
Parameters
  • t (float) – way parameter

  • r (float) – scaling parameter - effectively corresponds to the intersection of the path with the real axis

  • o (float) – offset on real axis

Returns

jac – derivative of Talbot path

Return type

complex

eko.mellin.Talbot_path(t, r, o)[source]

Talbot path.

\[p_{\text{Talbot}}(t) = o + r \cdot ( \theta \cot(\theta) + i\theta ), \theta = \pi(2t-1)\]
Parameters
  • t (float) – way parameter

  • r (float) – scaling parameter - effectively corresponds to the intersection of the path with the real axis

  • o (float) – offset on real axis

Returns

path – Talbot path

Return type

complex

eko.mellin.edge_jac(t, m, _c, phi)[source]

Derivative of edged path

\[\frac{dp_{\text{edge}}(t)}{dt}\]
Parameters
  • t (float) – way parameter

  • m (float) – length of the path

  • c (float, optional) – intersection of path with real axis

  • phi (complex, optional) – bended angle

Returns

path – Derivative of edged path

Return type

complex

eko.mellin.edge_path(t, m, c, phi)[source]

Edged path with a given angle.

\[p_{\text{edge}}(t) = c + m\left|t - \frac 1 2\right|\exp(i\phi)\]
Parameters
  • t (float) – way parameter

  • m (float) – length of the path

  • c (float, optional) – intersection of path with real axis

  • phi (complex, optional) – bended angle

Returns

path – Edged path

Return type

complex

eko.mellin.line_jac(_t, m, _c)[source]

Derivative of Textbook path.

\[\frac{dp_{\text{line}}(t)}{dt}\]
Parameters
  • t (float) – way parameter

  • m (float) – scaling parameter

  • o (float) – offset on real axis

Returns

jac – derivative of Textbook path

Return type

complex

eko.mellin.line_path(t, m, c)[source]

Textbook path, i.e. a straight line parallel to the imaginary axis.

\[p_{\text{line}}(t) = c + m \cdot (2t - 1)\]
Parameters
  • t (float) – way parameter

  • m (float) – scaling parameter

  • c (float) – offset on real axis

Returns

path – Textbook path

Return type

complex

eko.member module

class eko.member.MemberName(name)[source]

Bases: object

Operator member name in operator evolution space

Parameters

name (str) – operator name

property input

Returns input flavor name (given by the second part of the name)

property target

Returns target flavor name (given by the first part of the name)

class eko.member.OpMember(value, error)[source]

Bases: object

A single operator for a specific element in evolution basis.

This class provide some basic mathematical operations such as products. It can also be applied to a pdf vector via apply_pdf(). This class will never be exposed to the outside, but will be an internal member of the Operator and PhysicalOperator instances.

Parameters
  • value (np.array) – operator matrix

  • error (np.array) – operator error matrix

copy()[source]
classmethod id_like(other)[source]

Creates an identity operator.

Parameters

other (OpMember) – reference member

Returns

1

Return type

cls

class eko.member.OperatorBase(op_members, q2_final)[source]

Bases: object

Abstract base class to hold a dictionary of interpolation matrices.

Parameters
operation(other)[source]

Choose mathematical operation by rank

Parameters

other (OperatorBase) – operand

Returns

operation to perform (np.matmul or np.multiply)

Return type

callable

static operator_multiply(left, right, operation)[source]

Multiply two operators.

Parameters
  • left (OperatorBase) – left operand

  • right (OperatorBase) – right operand

  • operation (callable) – operation to perform (np.matmul or np.multiply)

Returns

new operator members dictionary

Return type

dict

classmethod promote_names(op_members, q2_final)[source]

Promote string keys to MemberName.

Parameters
  • op_members (dict) – mapping of str onto OpMember

  • q2_final (float) – final scale

class eko.member.ScalarOperator(op_members, q2_final)[source]

Bases: eko.member.OperatorBase

eko.msbar_masses module

This module contains the RGE for the \(\overline{MS}\) masses

eko.msbar_masses.compute_msbar_mass(theory_card)[source]

Compute the \(\overline{MS}\) masses solving the equation \(m_{\bar{MS}}(\mu) = \mu\)

Parameters

theory_card (dict) – theory run card

Returns

masses – list of \(\overline{MS}\) masses squared

Return type

list

eko.msbar_masses.evolve_msbar_mass(m2_ref, q2m_ref, strong_coupling, nf_ref=None, fact_to_ren=1.0, q2_to=None)[source]

Compute the \(\overline{MS}\) mass. If the final scale is not given it solves the equation \(m_{\overline{MS}}(m) = m\) for a fixed number of nf Else perform the mass evolution.

Parameters
  • m2_ref (float) – squared initial mass reference

  • q2m_ref (float) – squared initial scale

  • nf_ref (int, optional (not used when q2_to is given)) – number of active flavours at the scale q2m_ref, where the solution is searched

  • fact_to_ren (float) – \(\mu_F^2/\muR^2\)

  • strong_coupling (eko.strong_coupling.StrongCoupling) – Instance of StrongCoupling able to generate a_s for any q

  • q2_to (float, optional) – scale at which the mass is computed

Returns

m2\(m_{\overline{MS}}(\mu_2)^2\)

Return type

float

eko.msbar_masses.msbar_ker_dispatcher(q2_to, q2m_ref, strong_coupling, fact_to_ren, nf)[source]

Select the \(\overline{MS}\) kernel and compute the strong coupling values

Parameters
Returns

Expanded or exact \(\overline{MS}\) kernel

Return type

ker

eko.msbar_masses.msbar_ker_exact(a0, a1, order, nf)[source]

Exact \(\overline{MS}\) RGE kernel

Parameters
  • a0 (float) – strong coupling at the initial scale

  • a1 (float) – strong coupling at the final scale

  • oreder (int) – perturbative order

  • nf (int) – number of active flavours

Returns

ker – Exact \(\overline{MS}\) kernel:

..math:

k_{exact} = e^{int_{a_s(mu_{h,0}^2)}^{a_s(mu^2)} gamma(a_s) / beta(a_s) da_s}

Return type

float

eko.msbar_masses.msbar_ker_expanded(a0, a1, order, nf)[source]

Expanded \(\overline{MS}\) RGE kernel

Parameters
  • a0 (float) – strong coupling at the initial scale

  • a1 (float) – strong coupling at the final scale

  • oreder (int) – perturbative order

  • nf (int) – number of active flavours

Returns

ker – Expanded \(\overline{MS}\) kernel:

..math:

k_{expanded} &= left (frac{a_s(mu^2)}{a_s(mu_{h,0}^2)} right )^{c_0} frac{j_{exp}(a_s(mu^2))}{j_{exp}(a_s(mu_{h,0}^2))} \ j_{exp}(a_s) &= 1 + a_s left [ c_1 - b_1 c_0 right ] + frac{a_s^2}{2} left [c_2 - c_1 b_1 - b_2 c_0 + b_1^2 c_0 + (c_1 - b_1 c_0)^2 right]

Return type

float

eko.output module

This file contains the output management

class eko.output.Output[source]

Bases: dict

Wrapper for the output to help with application to PDFs and dumping to file.

apply_pdf(lhapdf_like, targetgrid=None, rotate_to_evolution_basis=False)[source]

Apply all available operators to the input PDFs.

Parameters
  • lhapdf_like (object) – object that provides an xfxQ2 callable (as lhapdf and ekomark.toyLH.toyPDF do) (and thus is in flavor basis)

  • targetgrid (list) – if given, interpolates to the PDFs given at targetgrid (instead of xgrid)

  • rotate_to_evolution_basis (bool) – if True rotate to evolution basis

Returns

out_grid – output PDFs and their associated errors for the computed Q2grid

Return type

dict

apply_pdf_flavor(lhapdf_like, targetgrid=None, flavor_rotation=None)[source]

Apply all available operators to the input PDFs.

Parameters
  • lhapdf_like (object) –

    object that provides an xfxQ2 callable (as lhapdf and ekomark.toyLH.toyPDF do) (and thus is in flavor basis)

  • targetgrid (list) – if given, interpolates to the PDFs given at targetgrid (instead of xgrid)

  • flavor_rotation (np.ndarray) – Rotation matrix in flavor space

Returns

out_grid – output PDFs and their associated errors for the computed Q2grid

Return type

dict

dump_tar(tarname)[source]

Writes representation into a tar archive containing:

  • metadata (in YAML)

  • operator (in numpy .npy format)

Parameters

tarname (str) – target file name

dump_yaml(stream=None, binarize=True, skip_q2_grid=False)[source]

Serialize result as YAML.

Parameters
  • stream (None or stream) – if given, dump is written on it

  • binarize (bool) – dump in binary format (instead of list format)

  • skip_q2_grid (bool) – avoid dumping Q2grid (i.e. the actual operators) into the yaml file (default: False)

Returns

dump – result of dump(output, stream), i.e. a string, if no stream is given or Null, if written successfully to stream

Return type

any

dump_yaml_to_file(filename, binarize=True, skip_q2_grid=False)[source]

Writes YAML representation to a file.

Parameters
  • filename (str) – target file name

  • binarize (bool) – dump in binary format (instead of list format)

  • skip_q2_grid (bool) – avoid dumping Q2grid (i.e. the actual operators) into the yaml file (default: False)

Returns

ret – result of dump(output, stream), i.e. Null if written successfully

Return type

any

flavor_reshape(targetbasis=None, inputbasis=None)[source]

Changes the operators to have in the output targetbasis and/or in the input inputbasis.

The operation is in-place.

Parameters
  • targetbasis (numpy.ndarray) – target rotation specified in the flavor basis

  • inputbasis (None or list) – input rotation specified in the flavor basis

get_raw(binarize=True, skip_q2_grid=False)[source]

Serialize result as dict/YAML.

This maps the original numpy matrices to lists.

Parameters

binarize (bool) – dump in binary format (instead of list format)

Returns

out – dictionary which will be written on output

Return type

dict

classmethod load_tar(tarname)[source]

Load tar representation from file (compliant with dump_tar() output).

Parameters

tarname (str) – source tar name

Returns

obj – loaded object

Return type

output

classmethod load_yaml(stream, skip_q2_grid=False)[source]

Load YAML representation from stream

Parameters
  • stream (any) – source stream

  • skip_q2_grid (bool) – avoid loading Q2grid (i.e. the actual operators) from the yaml file (default: False)

Returns

obj – loaded object

Return type

output

classmethod load_yaml_from_file(filename, skip_q2_grid=False)[source]

Load YAML representation from file

Parameters
  • filename (str) – source file name

  • skip_q2_grid (bool) – avoid loading Q2grid (i.e. the actual operators) from the yaml file (default: False)

Returns

obj – loaded object

Return type

output

to_evol(source=True, target=False)[source]

Rotate the operator into evolution basis.

This also assigns also the pids. The operation is in-place.

Parameters
  • source (bool) – rotate on the input tensor

  • target (bool) – rotate on the output tensor

xgrid_reshape(targetgrid=None, inputgrid=None)[source]

Changes the operators to have in the output targetgrid and/or in the input inputgrid.

The operation is in-place.

Parameters
  • targetgrid (None or list) – xgrid for the target

  • inputgrid (None or list) – xgrid for the input

eko.runner module

This file contains the main application class of eko

class eko.runner.Runner(theory_card, operators_card)[source]

Bases: object

Represents a single input configuration.

For details about the configuration, see here

Parameters

setup (dict) – input configurations

banner = '\noooooooooooo oooo    oooo  \\\\ .oooooo.\n`888\'     `8 `888   .8P\'  //////    `Y8b\n 888          888  d8\'   \\\\o\\/////    888\n 888oooo8     88888     \\\\\\\\/8/////   888\n 888    "     888`88b.      888 ///   888\n 888       o  888  `88b.    `88b //  d88\'\no888ooooood8 o888o  o888o     `Y8bood8P\'\n'
get_output()[source]

Collects all data for output (to run the evolution)

Returns

ret – output instance

Return type

eko.output.Output

eko.strong_coupling module

This file contains the QCD beta function coefficients and the handling of the running coupling \(\alpha_s\).

See pQCD ingredients.

class eko.strong_coupling.StrongCoupling(alpha_s_ref, scale_ref, masses, thresholds_ratios, order=0, method='exact', nf_ref=None, max_nf=None, hqm_scheme='POLE')[source]

Bases: object

Computes the strong coupling constant \(a_s\).

Note that

  • all scale parameters (scale_ref and scale_to), have to be given as squared values, i.e. in units of \(\text{GeV}^2\)

  • although, we only provide methods for \(a_s = \frac{\alpha_s(\mu^2)}{4\pi}\) the reference value has to be given in terms of \(\alpha_s(\mu_0^2)\) due to legacy reasons

  • the order refers to the perturbative order of the beta function, thus order=0 means leading order beta function, means evolution with \(\beta_0\), means running at 1-loop - so there is a natural mismatch between order and the number of loops by one unit

Normalization is given by [HRU+17]:

\[\frac{da_s(\mu^2)}{d\ln\mu^2} = \beta(a_s) \ = - \sum\limits_{n=0} \beta_n a_s^{n+2}(\mu^2) \quad \text{with}~ a_s = \frac{\alpha_s(\mu^2)}{4\pi}\]

See pQCD ingredients.

Parameters
  • alpha_s_ref (float) – alpha_s(!) at the reference scale \(\alpha_s(\mu_0^2)\)

  • scale_ref (float) – reference scale \(\mu_0^2\)

  • masses (list(float)) – list with quark masses

  • thresholds_ratios (list(float)) – list with ratios between the mass and the thresholds

  • order (int) – Evaluated order of the beta function: 0 = LO, …

  • method (["expanded", "exact"]) – Applied method to solve the beta function

  • nf_ref (int) – if given, the number of flavors at the reference scale

  • max_nf (int) – if given, the maximum number of flavors

a_s(scale_to, fact_scale=None, nf_to=None)[source]

Computes strong coupling \(a_s(\mu_R^2) = \frac{\alpha_s(\mu_R^2)}{4\pi}\).

Parameters
  • scale_to (float) – final scale to evolve to \(\mu_R^2\)

  • fact_scale (float) – factorization scale (if different from final scale)

Returns

a_s – strong coupling \(a_s(\mu_R^2) = \frac{\alpha_s(\mu_R^2)}{4\pi}\)

Return type

float

compute(as_ref, nf, scale_from, scale_to)[source]

Wrapper in order to pass the computation to the corresponding method (depending on the calculation method).

Parameters
  • as_ref (float) – reference alpha_s

  • nf (int) – value of nf for computing alpha_s

  • scale_from (float) – reference scale

  • scale_to (float) – target scale

Returns

a_s – strong coupling at target scale \(a_s(Q^2)\)

Return type

float

compute_exact(as_ref, nf, scale_from, scale_to)[source]

Compute via RGE.

Parameters
  • as_ref (float) – reference alpha_s

  • nf (int) – value of nf for computing alpha_s

  • scale_from (float) – reference scale

  • scale_to (float) – target scale

Returns

a_s – strong coupling at target scale \(a_s(Q^2)\)

Return type

float

classmethod from_dict(theory_card, masses=None)[source]

Create object from theory dictionary.

Parameters
  • theory_card (dict) – theory dictionary

  • masses (list) – list of \(\overline{MS}\) masses squared or None if POLE masses are used

Returns

cls – created object

Return type

StrongCoupling

property q2_ref

reference scale

eko.strong_coupling.as_expanded(order, as_ref, nf, scale_from, scale_to)[source]

Compute expanded expression.

Parameters
  • order (int) – perturbation order

  • as_ref (float) – reference alpha_s

  • nf (int) – value of nf for computing alpha_s

  • scale_from (float) – reference scale

  • scale_to (float) – target scale

Returns

a_s – coupling at target scale \(a_s(Q^2)\)

Return type

float

eko.strong_coupling.compute_matching_coeffs_down(mass_scheme)[source]

Matching coefficients [SS06] [CKS06] at threshold when moving to a regime with less flavors.

This is the perturbative inverse of matching_coeffs_up and has been obtained via

Module[{f, g, l, sol},
    f[a_] := a + Sum[d[n, k]*L^k*a^(1 + n), {n, 3}, {k, 0, n}];
    g[a_] := a + Sum[c[n, k]*L^k*a^(1 + n), {n, 3}, {k, 0, n}] /. {c[1, 0] -> 0};
    l = CoefficientList[Normal@Series[f[g[a]], {a, 0, 5}], {a, L}];
    sol = First@
        Solve[{l[[3]] == 0, l[[4]] == 0, l[[5]] == 0},
        Flatten@Table[d[n, k], {n, 3}, {k, 0, n}]];
    Do[Print@r, {r, sol}];
    Print@Series[f[g[a]] /. sol, {a, 0, 5}];
    Print@Series[g[f[a]] /. sol, {a, 0, 5}];
]
Parameters

mass_scheme – Heavy quark mass scheme: “POLE” or “MSBAR”

Returns

downward matching coefficient matrix

Return type

matching_coeffs_down

eko.strong_coupling.compute_matching_coeffs_up(mass_scheme)[source]

Matching coefficients [CKS06, SS06, Vog05] at threshold when moving to a regime with more flavors.

\[a_s^{(n_l+1)} = a_s^{(n_l)} + \sum\limits_{n=1} (a_s^{(n_l)})^n \sum\limits_{k=0}^n c_{nl} \log(\mu_R^2/\mu_F^2)\]
Parameters

mass_scheme – Heavy quark mass scheme: “POLE” or “MSBAR”

Returns

forward matching coefficient matrix

Return type

matching_coeffs_down

eko.strong_coupling.strong_coupling_mod_ev(mod_ev)[source]

Map ModEv key to the available strong coupling evolution methods

eko.thresholds module

This module holds the classes that define the FNS.

class eko.thresholds.PathSegment(q2_from, q2_to, nf)[source]

Bases: object

Oriented path in the threshold landscape.

Parameters
  • q2_from (float) – starting point

  • q2_to (float) – final point

  • nf (int) – number of active flavors

property is_backward

True if q2_from bigger than q2_to

property tuple

Tuple representation suitable for hashing.

class eko.thresholds.ThresholdsAtlas(masses, q2_ref=None, nf_ref=None, thresholds_ratios=None, max_nf=None)[source]

Bases: object

Holds information about the thresholds any Q2 has to pass in order to get there from a given q2_ref.

Parameters
  • masses (list(float)) – list of quark masses squared

  • q2_ref (float) – reference scale

  • nf_ref (int) – number of active flavors at the reference scale

  • thresholds_ratios (list(float)) – list of ratios between masses and matching thresholds squared

  • max_nf (int) – maximum number of active flavors

static build_area_walls(masses, thresholds_ratios=None, max_nf=None)[source]

Create the object from the run card.

The thresholds are computed by \((m_q \cdot k_q^{Thr})\).

Returns

threshold list

Return type

list

classmethod ffns(nf, q2_ref=None)[source]

Create a FFNS setup.

The function creates simply sufficient thresholds at 0 (in the beginning), since the number of flavors is determined by counting from below.

Parameters
  • nf (int) – number of light flavors

  • q2_ref (float) – reference scale

classmethod from_dict(theory_card, prefix='k', max_nf_name='MaxNfPdf', masses=None)[source]

Create the object from the run card.

The thresholds are computed by \((m_q \cdot k_q^{Thr})\).

Parameters
  • theory_card (dict) – run card with the keys given at the head of the module

  • prefix (str) – prefix for the ratio parameters

  • masses (list) – list of \(\overline{MS}\) masses squared or None if POLE masses are used

Returns

created object

Return type

ThresholdsAtlas

nf(q2)[source]

Finds the number of flavor active at the given scale.

Parameters

q2 (float) – reference scale

Returns

nf – number of active flavors

Return type

int

path(q2_to, nf_to=None, q2_from=None, nf_from=None)[source]

Get path from q2_from to q2_to.

Parameters
  • q2_to (float) – target value of q2

  • q2_from (float) – starting value of q2

Returns

path – List of PathSegment to go through in order to get from q2_from to q2_to.

Return type

list(PathSegment)

eko.version module