eko package

Evolution Kernel Operators.

Subpackages

Submodules

eko.basis_rotation module

Contains the definitions of the Flavor Basis and Evolution Basis.

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.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.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.unified_evol_basis = ('g', 'ph', 'S', 'Sdelta', 'V', 'Vdelta', 'Td3', 'Vd3', 'Tu3', 'Vu3', 'Td8', 'Vd8', 'Tu8', 'Vu8')

Sorted elements in Unified Evolution Basis as str.

Definition: here.

corresponding PDF : \(g, \gamma, \Sigma, \Sigma_\Delta, V, V_\Delta, T^d_3, V^d_3, T^u_3, V^u_3, T^d_8, V^d_8, T^u_8, V^u_8\)

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.unified_evol_basis_pids = (21, 22, 100, 101, 200, 201, 105, 205, 104, 204, 110, 210, 109, 209)

PID representation of unified_evol_basis.

The notation used for the non singlet compunents is the following: pid_ns(u) = pid_ns + 1, pid_ns(d) = pid_ns + 2.

eko.basis_rotation.anomalous_dimensions_basis = ((100, 100), (100, 21), (21, 100), (21, 21), (10201, 0), (10101, 0), (10200, 0))

Sorted elements in Anomalous Dimensions Basis as str.

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.basis_rotation.rotate_flavor_to_unified_evolution = array([[ 0,  0,  0,  0,  0,  0,  0,  1,  0,  0,  0,  0,  0,  0],        [ 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,  1, -1,  1, -1,  1, -1,  0, -1,  1, -1,  1, -1,  1],        [ 0, -1, -1, -1, -1, -1, -1,  0,  1,  1,  1,  1,  1,  1],        [ 0, -1,  1, -1,  1, -1,  1,  0, -1,  1, -1,  1, -1,  1],        [ 0,  0,  0,  0, -1,  0,  1,  0,  1,  0, -1,  0,  0,  0],        [ 0,  0,  0,  0,  1,  0, -1,  0,  1,  0, -1,  0,  0,  0],        [ 0,  0,  0, -1,  0,  1,  0,  0,  0,  1,  0, -1,  0,  0],        [ 0,  0,  0,  1,  0, -1,  0,  0,  0,  1,  0, -1,  0,  0],        [ 0,  0, -2,  0,  1,  0,  1,  0,  1,  0,  1,  0, -2,  0],        [ 0,  0,  2,  0, -1,  0, -1,  0,  1,  0,  1,  0, -2,  0],        [ 0, -2,  0,  1,  0,  1,  0,  0,  0,  1,  0,  1,  0, -2],        [ 0,  2,  0, -1,  0, -1,  0,  0,  0,  1,  0,  1,  0, -2]])

Basis rotation matrix between Flavor Basis and Unified Evolution Basis.

eko.basis_rotation.map_ad_to_evolution = {(21, 21): ['g.g'], (21, 100): ['g.S'], (100, 21): ['S.g'], (100, 100): ['S.S'], (10101, 0): ['T3.T3', 'T8.T8', 'T15.T15', 'T24.T24', 'T35.T35'], (10200, 0): ['V.V'], (10201, 0): ['V3.V3', 'V8.V8', 'V15.V15', 'V24.V24', 'V35.V35']}

Map anomalous dimension sectors’ names to their members

eko.basis_rotation.ad_projector(ad_lab, nf, qed)[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

  • qed (bool) – activate qed

Returns:

proj – projector over the specified sector

Return type:

np.ndarray

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

Select light flavors for a given ad_lab.

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

  • nf (int) – number of light flavors

Returns:

sector – list of sector for a given ad_lab

Return type:

list

eko.basis_rotation.ad_projectors(nf, qed)[source]

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

Parameters:
  • nf (int) – number of light flavors

  • qed (bool) – activate qed

Returns:

projs – projectors tensor

Return type:

np.ndarray

eko.basis_rotation.intrinsic_unified_evol_labels(nf)[source]

Collect all labels in the intrinsic unified evolution basis.

Parameters:

nf (int) – number of light flavors

Returns:

labels – active distributions

Return type:

list(str)

eko.beta module

Compute the QCD beta function coefficients.

See pQCD ingredients.

eko.beta.beta_qcd_as2(nf)[source]

Compute the first coefficient of the QCD beta function.

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

Parameters:

nf (int) – number of active flavors

Returns:

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

Return type:

float

eko.beta.beta_qed_aem2(nf, nl)[source]

Compute the first coefficient of the QED beta function.

Implements Eq. (7) of [Sur96].

Parameters:
  • nf (int) – number of active flavors

  • nl (int) – number of leptons partecipating to alphaem running

Returns:

beta_qed_aem2 – first coefficient of the QED beta function \(\\beta_qed_aem2^{n_f}\)

Return type:

float

eko.beta.beta_qcd_as3(nf)[source]

Compute the second coefficient of the QCD beta function.

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

Parameters:

nf (int) – number of active flavors

Returns:

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

Return type:

float

eko.beta.beta_qed_aem3(nf, nl)[source]

Compute the second coefficient of the QED beta function.

Implements Eq. (7) of [Sur96].

Parameters:
  • nf (int) – number of active flavors

  • nl (int) – number of leptons partecipating to alphaem running

Returns:

beta_qed_aem3 – second coefficient of the QED beta function \(\\beta_qed_aem3^{n_f}\)

Return type:

float

eko.beta.beta_qcd_as2aem1(nf)[source]

Compute the first QED correction of the QCD beta function.

Implements Eq. (7) of [Sur96].

Parameters:

nf (int) – number of active flavors

Returns:

beta_as2aem1 – first QED correction of the QCD beta function \(\\beta_as2aem1^{n_f}\)

Return type:

float

eko.beta.beta_qed_aem2as1(nf)[source]

Compute the first QCD correction of the QED beta function.

Implements Eq. (7) of [Sur96].

Parameters:

nf (int) – number of active flavors

Returns:

beta_aem2as1 – first QCD correction of the QED beta function \(\\beta_aem2as1^{n_f}\)

Return type:

float

eko.beta.beta_qcd_as4(nf)[source]

Compute the third coefficient of the QCD beta function.

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

Parameters:

nf (int) – number of active flavors

Returns:

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

Return type:

float

eko.beta.beta_qcd_as5(nf)[source]

Compute the fourth coefficient of the QCD beta function.

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

Parameters:

nf (int) – number of active flavors

Returns:

beta_qcd_as5 – fourth coefficient of the QCD beta function \(\\beta_qcd_as5^{n_f}\)

Return type:

float

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

Compute value of a beta_qcd coefficients.

Parameters:
  • k (tuple(int, int)) – perturbative order

  • nf (int) – number of active flavors

Returns:

beta_qcd – beta_qcd_k(nf)

Return type:

float

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

Compute value of a beta_qed coefficients.

Parameters:
  • k (tuple(int, int)) – perturbative order

  • nf (int) – number of active flavors

  • nl (int) – number of leptons partecipating to alphaem running

Returns:

beta_qed – beta_qed_k(nf)

Return type:

float

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

Compute b_qcd coefficient.

Parameters:
  • k (tuple(int, int)) – perturbative order

  • nf (int) – number of active flavors

Returns:

b_qcd – b_qcd_k(nf)

Return type:

float

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

Compute b_qed coefficient.

Parameters:
  • k (tuple(int, int)) – perturbative order

  • nf (int) – number of active flavors

  • nl (int) – number of leptons partecipating to alphaem running

Returns:

b_qed – b_qed_k(nf)

Return type:

float

eko.constants module

Sets the physical constants.

eko.constants.NC = 3

The number of colors.

eko.constants.TR = 0.5

The normalization of fundamental generators.

Defaults to \(T_R = 1/2\).

eko.constants.CA = 3.0

Second Casimir constant in the adjoint representation.

Defaults to \(N_C = 3\).

eko.constants.CF = 1.3333333333333333

Second Casimir constant in the fundamental representation.

Defaults to \(\frac{N_C^2-1}{2N_C} = 4/3\).

eko.constants.MTAU = 1.777

Mass of the tau.

eko.constants.eu2 = 0.4444444444444444

Up quarks charge squared.

eko.constants.ed2 = 0.1111111111111111

Down quarks charge squared.

eko.constants.zeta2 = 1.6449340668482264

\(\zeta(2)\)

eko.constants.zeta3 = 1.2020569031595942

\(\zeta(3)\)

eko.constants.zeta4 = 1.0823232337111381

\(\zeta(4)\)

eko.constants.zeta5 = 1.03692775514337

\(\zeta(5)\)

eko.constants.log2 = 0.6931471805599453

\(\ln(2)\)

eko.constants.li4half = 0.517479

\(Li_{4}(1/2)\)

eko.constants.update_colors(nc)[source]

Update the number of colors to \(NC = nc\).

The Casimirs for a generic value of \(NC\) are consistenly updated as well.

Parameters:

nc (int) – Number of colors

eko.constants.uplike_flavors(nf)[source]

Compute the number of up flavors.

Parameters:

nf (int) – Number of active flavors

Returns:

nu

Return type:

int

eko.constants.charge_combinations(nf)[source]

Compute the combination of charges.

Parameters:

nf (int) – Number of active flavors

Returns:

  • e2avg (float)

  • vue2m (float)

  • vde2m (float)

eko.couplings module

Manage running (and fixed) couplings.

Manage QCD coupling \(\alpha_s\) and QED coupling \(\alpha\). We provide an interface to access them simultaneously and provide several strategies to solve the associated RGE.

See pQCD ingredients.

eko.couplings.couplings_mod_ev(mod_ev: EvolutionMethod) CouplingEvolutionMethod[source]

Map ModEv key to the available strong coupling evolution methods.

Parameters:

mod_ev (str) – evolution mode

Returns:

coupling mode

Return type:

str

eko.couplings.exact_lo(ref, beta0, lmu)[source]

Compute expanded solution at LO.

Parameters:
  • ref (float) – reference value of the coupling

  • beta0 (float) – first coefficient of the beta function

  • lmu (float) – logarithm of the ratio between target and reference scales

Returns:

coupling at target scale \(a(\mu_R^2)\)

Return type:

float

eko.couplings.expanded_nlo(ref, beta0, b1, lmu)[source]

Compute expanded solution at NLO.

Implement the default expression for NLO expanded solution, e.g. the one implemented in APFEL.

Parameters:
  • ref (float) – reference value of the coupling

  • beta0 (float) – first coefficient of the beta function

  • b1 (float) – second coefficient of the b function

  • lmu (float) – logarithm of the ratio between target and reference scales

Returns:

coupling at target scale \(a(\mu_R^2)\)

Return type:

float

eko.couplings.expanded_nnlo(ref, beta0, b1, b2, lmu)[source]

Compute expanded solution at NNLO.

Implement the NNLO expanded solution from APFEL (not the default expression)

Parameters:
  • ref (float) – reference value of the coupling

  • beta0 (float) – first coefficient of the beta function

  • b1 (float) – second coefficient of the b function

  • b2 (float) – third coefficient of the b function

  • lmu (float) – logarithm of the ratio between target and reference scales

Returns:

coupling at target scale \(a(\mu_R^2)\)

Return type:

float

eko.couplings.expanded_n3lo(ref, beta0, b1, b2, b3, lmu)[source]

Compute expanded solution at N3LO.

Implement the N3LO expanded solution obtained via iterated solution of the RGE [Rot]

Parameters:
  • ref (float) – reference value of the coupling

  • beta0 (float) – first coefficient of the beta function

  • b1 (float) – second coefficient of the b function

  • b2 (float) – third coefficient of the b function

  • b3 (float) – fourth coefficient of the b function

  • lmu (float) – logarithm of the ratio between target and reference scales

Returns:

coupling at target scale \(a(\mu_R^2)\)

Return type:

float

eko.couplings.expanded_qcd(ref, order, beta0, b_vec, lmu)[source]

Compute QCD expanded solution at a given order.

Parameters:
  • ref (float) – reference value of the strong coupling

  • order (int) – QCD order

  • beta0 (float) – first coefficient of the beta function

  • b_vec (list) – list of b function coefficients (including b0)

  • lmu (float) – logarithm of the ratio between target and reference scales

Returns:

strong coupling at target scale \(a_s(\mu_R^2)\)

Return type:

float

eko.couplings.expanded_qed(ref, order, beta0, b_vec, lmu)[source]

Compute QED expanded solution at a given order.

Parameters:
  • ref (float) – reference value of the QED coupling

  • order (int) – QED order

  • beta0 (float) – first coefficient of the beta function

  • b_vec (list) – list of b function coefficients (including b0)

  • lmu (float) – logarithm of the ratio between target and reference scales

Returns:

QED coupling at target scale \(a_em(\mu_R^2)\)

Return type:

float

eko.couplings.couplings_expanded_alphaem_running(order, couplings_ref, nf, nl, scale_from, scale_to, decoupled_running)[source]

Compute coupled expanded expression of the couplings for running alphaem.

Implement Eqs. (17-18) from [Sur96]

Parameters:
  • order (tuple(int, int)) – perturbation order

  • couplings_ref (numpy.ndarray) – reference alpha_s and alpha

  • nf (int) – value of nf for computing the couplings

  • nl (int) – number of leptons partecipating to alphaem running

  • scale_from (float) – reference scale

  • scale_to (float) – target scale

  • decoupled_running (bool) – whether the running of the couplings is decoupled or not

Returns:

couplings at target scale \(a(\mu_R^2)\)

Return type:

numpy.ndarray

eko.couplings.couplings_expanded_fixed_alphaem(order, couplings_ref, nf, scale_from, scale_to)[source]

Compute coupled expanded expression of the couplings for fixed alphaem.

Parameters:
  • order (tuple(int, int)) – perturbation order

  • couplings_ref (numpy.ndarray) – reference alpha_s and alpha

  • nf (int) – value of nf for computing the couplings

  • scale_from (float) – reference scale

  • scale_to (float) – target scale

Returns:

couplings at target scale \(a(\mu_R^2)\)

Return type:

numpy.ndarray

class eko.couplings.Couplings(couplings: CouplingsInfo, order: Tuple[int, int], method: CouplingEvolutionMethod, masses: List[float], hqm_scheme: QuarkMassScheme, thresholds_ratios: Iterable[float])[source]

Bases: object

Compute the strong and electromagnetic coupling constants \(a_s, a_{em}\).

Note that

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

  • the order refers to the perturbative order of the beta functions, i.e. the number of loops for QCD and QED respectively. QCD is always running with at least 1 loop, while QED might not run at all (so 0 loop).

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:
  • couplings – reference configuration

  • order – Number of loops in beta functions (QCD, QED)

  • method – Applied method to solve the beta functions

  • masses – list with quark masses squared

  • hqm_scheme – heavy quark mass scheme

  • thresholds_ratios – list with ratios between the matching scales and the mass squared

property mu2_ref

Return reference scale.

unidimensional_exact(beta0, b_vec, u, a_ref, method, rtol)[source]

Compute single coupling via decoupled RGE.

Parameters:
  • beta0 (float) – first coefficient of the beta function

  • b_vec (list) – list of b function coefficients (including b0)

  • u (float) – \(log(scale_to / scale_from)\)

  • a_ref (float) – reference alpha_s or alpha

  • method (string) – method for solving the RGE

  • rtol (float) – relative acuracy of the solution

Returns:

coupling at target scale \(a(Q^2)\)

Return type:

float

compute_exact_alphaem_running(a_ref, nf, nl, scale_from, scale_to)[source]

Compute couplings via RGE with running alphaem.

Parameters:
  • a_ref (numpy.ndarray) – reference alpha_s and alpha

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

  • nl (int) – number of leptons partecipating to alphaem running

  • scale_from (float) – reference scale

  • scale_to (float) – target scale

Returns:

couplings at target scale \(a(Q^2)\)

Return type:

numpy.ndarray

compute_exact_fixed_alphaem(a_ref, nf, scale_from, scale_to)[source]

Compute couplings via RGE with fixed alphaem.

Parameters:
  • a_ref (numpy.ndarray) – reference alpha_s and alpha

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

  • scale_from (float) – reference scale

  • scale_to (float) – target scale

Returns:

couplings at target scale \(a(Q^2)\)

Return type:

numpy.ndarray

compute(a_ref, nf, nl, scale_from, scale_to)[source]

Compute actual couplings.

This is a wrapper in order to pass the computation to the corresponding method (depending on the calculation method).

Parameters:
  • a_ref (numpy.ndarray) – reference a

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

  • nl (int) – number of leptons partecipating to alphaem running

  • scale_from (float) – reference scale

  • scale_to (float) – target scale

Returns:

couplings at target scale \(a(Q^2)\)

Return type:

numpy.ndarray

a(scale_to, nf_to=None)[source]

Compute couplings \(a_i(\mu_R^2) = \frac{\alpha_i(\mu_R^2)}{4\pi}\).

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

  • nf_to (int) – final nf value

Returns:

couplings \(a_i(\mu_R^2) = \frac{\alpha_i(\mu_R^2)}{4\pi}\)

Return type:

numpy.ndarray

a_s(scale_to, nf_to=None)[source]

Compute strong coupling.

The strong oupling uses the normalization \(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\)

  • nf_to (int) – final nf value

Returns:

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

Return type:

float

a_em(scale_to, nf_to=None)[source]

Compute electromagnetic coupling.

The electromagnetic oupling uses the normalization \(a_em(\mu_R^2) = \frac{\alpha_em(\mu_R^2)}{4\pi}\).

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

  • nf_to (int) – final nf value

Returns:

a_em – couplings \(a_em(\mu_R^2) = \frac{\alpha_em(\mu_R^2)}{4\pi}\)

Return type:

float

eko.couplings.compute_matching_coeffs_up(mass_scheme, nf)[source]

Compute the upward matching coefficients.

The matching coefficients [CKS06, SS06, Vog05] are normalized at threshold when moving to a regime with more flavors.

We follow notation of [Vog05] (eq 2.43) for POLE masses.

The inverse \(\overline{MS}\) values instead are given in [SS06] (eq 3.1) multiplied by a factor of 4 (and 4^2 …).

\[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 (str) – Heavy quark mass scheme: “POLE” or “MSBAR”

  • nf (int) – number of active flavors in the lower patch

Returns:

forward matching coefficient matrix

Return type:

numpy.ndarray

eko.couplings.compute_matching_coeffs_down(mass_scheme, nf)[source]

Compute the downward matching coefficients.

This is the inverse function to compute_matching_coeffs_up().

Parameters:
  • mass_scheme (str) – Heavy quark mass scheme: “POLE” or “MSBAR”

  • nf (int) – number of active flavors in the lower patch

Returns:

downward matching coefficient matrix

Return type:

numpy.ndarray

eko.couplings.invert_matching_coeffs(c_up)[source]

Compute the perturbative inverse of the matching conditions.

They can be obtained e.g. via Mathematica by:

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:

c_up (numpy.ndarray) – forward matching coefficient matrix

Returns:

downward matching coefficient matrix

Return type:

numpy.ndarray

eko.gamma module

The QCD gamma function coefficients.

See pQCD ingredients.

eko.gamma.gamma_qcd_as1()[source]

Compute the first coefficient of the QCD gamma function.

Implements 15 of [VLvR97].

Returns:

gamma_0 – first coefficient of the QCD gamma function \(\gamma_{m,0}^{n_f}\)

Return type:

float

eko.gamma.gamma_qcd_as2(nf)[source]

Compute the second coefficient of the QCD gamma function.

Implements 15 of [VLvR97].

Parameters:

nf (int) – number of active flavors

Returns:

gamma_1 – second coefficient of the QCD gamma function \(\gamma_{m,1}^{n_f}\)

Return type:

float

eko.gamma.gamma_qcd_as3(nf)[source]

Compute the third coefficient of the QCD gamma function.

Implements 15 of [VLvR97].

Parameters:

nf (int) – number of active flavors

Returns:

gamma_2 – third coefficient of the QCD gamma function \(\gamma_{m,2}^{n_f}\)

Return type:

float

eko.gamma.gamma_qcd_as4(nf)[source]

Compute the fourth coefficient of the QCD gamma function.

Implements 15 of [VLvR97].

Parameters:

nf (int) – number of active flavors

Returns:

gamma_3 – fourth coefficient of the QCD gamma function \(\gamma_{m,3}^{n_f}\)

Return type:

float

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

Compute the value of a QCD gamma coefficient.

Parameters:
  • order (int) – perturbative order

  • nf (int) – number of active flavors

Returns:

gammaQCD gamma function coefficient

Return type:

float

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

Define area of 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

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

Evaluate logarithmic interpolator in N-space.

A single logarithmic Lagrange interpolator is evaluated 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.evaluate_Nx(N, logx, area_list)[source]

Evaluate linear interpolator in N-space.

A single linear Lagrange interpolator is evaluated 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_grid(N, is_log, logx, area_list)[source]

Evaluate interpolator in N-space.

Parameters:
  • N (complex) – Mellin variable

  • is_log (boolean) – is a logarithmic interpolation

  • logx (float) – logarithm of inversion point

  • area_list (list) – area configuration of basis function

Returns:

pj – basis function * inversion factor

Return type:

float

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

Evaluate 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_x(x, area_list)[source]

Evaluate 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

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

Bases: object

Represent an element of the polynomial basis.

It contains a list of areas for a given polynomial number defined by (xmin-xmax) which in turn contain a list of coefficients.

Upon construction it 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 on x

is_below_x(x)[source]

Check all areas to be below specified threshold.

Parameters:

x (float) – reference value

Returns:

is_below_x – xmax of highest area <= x?

Return type:

bool

areas_to_const()[source]

Return an array containing all areas.

The area format is (xmin, xmax, numpy.array of coefficients).

Returns:

area config

Return type:

numpy.ndarray

compile_x()[source]

Compile 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

compile_n()[source]

Compile 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.

class eko.interpolation.XGrid(xgrid: Sequence | ndarray[Any, dtype[_ScalarType_co]], log: bool = True)[source]

Bases: object

Grid of points in \(x\)-space.

This object represents a suitable grid of momentum fractions to be used to evaluate the PDF over.

property raw: ndarray

Untransformed grid.

property size: int

Number of pointrs.

tolist() list[source]

Raw grid as Python list.

dump() dict[source]

Representation as dictionary.

classmethod load(doc: dict)[source]

Create object from dictinary.

classmethod fromcard(value: list, log: bool)[source]

Create object from theory card config.

The config can invoke other grid generation methods.

class eko.interpolation.InterpolatorDispatcher(xgrid: XGrid | Sequence | ndarray[Any, dtype[_ScalarType_co]], polynomial_degree, mode_N=True)[source]

Bases: object

Setup 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 (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

get_interpolation(targetgrid: ndarray[Any, dtype[_ScalarType_co]] | Sequence)[source]

Compute interpolation matrix between targetgrid and xgrid.

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

targetgrid (array) – grid to interpolate to

Returns:

R – interpolation matrix $R_{ij}$, where $i$ is the index over targetgrid, and $j$ is the index on the internal basis (the one stored in the InterpolatorDispatcher instance)

Return type:

array

to_dict()[source]

Return the configuration for the underlying xgrid.

An instance can be constructed again with just the information returned by this method.

Note

The output dictionary contains a numpy array.

Returns:

ret – full grid configuration

Return type:

dict

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]

Create 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.lambertgrid(n_pts, x_min=1e-07, x_max=1.0)[source]

Create 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 implemented 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.matchings module

Holds the classes that define the FNS.

class eko.matchings.Segment(origin: float, target: float, nf: int)[source]

Bases: object

Oriented path in the threshold landscape.

origin: float

Starting point.

target: float

Final point.

nf: int

Number of active flavors.

property is_downward: bool

Return True if origin is bigger than target.

class eko.matchings.Matching(scale: float, hq: int, inverse: bool)[source]

Bases: object

Matching between two different segments.

The meaning of the flavor index hq is the PID of the corresponding heavy quark.

scale: float
hq: int
inverse: bool
class eko.matchings.Atlas(matching_scales: HeavyQuarks[float], origin: Tuple[float, int])[source]

Bases: object

Holds information about the matching scales.

These scales are the \(Q^2\) has to pass in order to get there from a given \(Q^2_{ref}\).

normalize(target: Tuple[float, int]) Tuple[float, int][source]

Fill number of flavors if needed.

classmethod ffns(nf: int, mu2: float)[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.

The origin is set with that number of flavors.

path(target: Tuple[float, int]) List[Segment][source]

Determine the path to the target evolution point.

Essentially, the path is always monotonic in the number of flavors, increasing or decreasing the active flavors by one unit every time a matching happens at the suitable scale.

Examples

Since this can result in a counter-intuitive behavior, let’s walk through some examples.

Starting with the intuitive one: >>> Atlas([10, 20, 30], (5, 3)).path((25, 5)) [Segment(5, 10, 3), Segment(10, 20, 4), Segment(20, 25, 5)]

If the number of flavor has been reached, it will continue walking without matchin again. >>> Atlas([10, 20, 30], (5, 3)).path((25, 4)) [Segment(5, 10, 3), Segment(10, 25, 4)]

It is irrelevant the scale you start from, to step from 3 to 4 you have to cross the charm matching scale, whether this means walking upward or downward. >>> Atlas([10, 20, 30], (15, 3)).path((25, 5)) [Segment(15, 10, 3), Segment(10, 20, 4), Segment(20, 25, 5)]

An actual backward evolution is defined by lowering the number of flavors going from origin to target. >>> Atlas([10, 20, 30], (25, 5)).path((5, 3)) [Segment(25, 20, 5), Segment(20, 10, 4), Segment(10, 5, 3)]

But the only difference is in the matching between two segments, since a single segment is always happening in a fixed number of flavors, and it is completely analogue for upward or downward evolution.

Note

Since the only task required to determine a path is interleaving the correct matching scales, this is done by slicing the walls.

matched_path(target: Tuple[float, int]) List[Segment | Matching][source]

Determine the path to the target, including matchings.

In practice, just a wrapper around path() adding the intermediate matchings.

eko.matchings.nf_default(mu2: float, atlas: Atlas) int[source]

Determine the number of active flavors in the default flow.

Default flow is defined by the natural sorting of the matching scales:

\[\mu_c < \mu_b < \mu_t\]

So, the flow is defined starting with 3 flavors below the charm matching, and increasing by one every time a matching scale is passed while increasing the scale.

eko.matchings.is_downward_path(path: List[Segment]) bool[source]

Determine if a path is downward.

Criterias are:

  • in the number of active flavors when the path list contains more than one Segment, note this can be different from each Segment.is_downward

  • in \(\mu^2\), when just one single Segment is given

eko.matchings.flavor_shift(is_downward: bool) int[source]

Determine the shift to number of light flavors.

eko.matchings.lepton_number(q2)[source]

Compute the number of leptons.

Note: muons and electrons are always massless as for up, down and strange.

Parameters:

q2 (float) – scale

Returns:

Number of leptons

Return type:

int

eko.mellin module

Tools for Mellin inversion.

We provide all necessary toold to deal with 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_path(t, r, o)[source]

Compute 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.Talbot_jac(t, r, _o)[source]

Compute 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.line_path(t, m, c)[source]

Compute 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.mellin.line_jac(_t, m, _c)[source]

Compute 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.edge_path(t, m, c, phi)[source]

Compute 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.edge_jac(t, m, _c, phi)[source]

Compute 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

class eko.mellin.Path(*args, **kwargs)[source]

Bases: Path

Mellin path dispatcher.

Parameters:
  • t (float) – way parameter

  • logx (float) – Mellin inversion point

  • axis_offset (bool) – add offset on the real axis

class_type = jitclass.Path#7f9dc23047f0<t:float64,r:float64,o:int8>

eko.member module

Atomic operator member.

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]

Copy implementation.

classmethod id_like(other)[source]

Create an identity operator.

Parameters:

other (OpMember) – reference member

Returns:

1 all spaces

Return type:

cls

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

Bases: object

Operator member name in operator evolution space.

Parameters:

name (str) – operator name

property target

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

property input

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

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

Bases: object

Abstract base class to hold a dictionary of interpolation matrices.

Parameters:
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

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

to_flavor_basis_tensor(qed: bool)[source]

Convert the computations into an rank 4 tensor.

A sparse tensor defined with dot-notation (e.g. S.g) is converted to a plain rank-4 array over flavor operator space and momentum fraction operator space.

If qed is passed, the unified intrinsic basis is used.

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

Bases: OperatorBase

Operator above space of real numbers.

eko.msbar_masses module

RGE for the \(\overline{MS}\) masses.

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

Provide exact \(\overline{MS}\) RGE kernel.

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

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

  • oreder (tuple(int,int)) – perturbative order

  • nf (int) – number of active flavours

Returns:

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

\[k_{exact} = e^{-\int_{a_s(\mu_{h,0}^2)}^{a_s(\mu^2)}\gamma_m(a_s)/ \beta(a_s)da_s}\]

Return type:

float

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

Provide expanded \(\overline{MS}\) RGE kernel.

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

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

  • order (tuple(int,int)) – perturbative order

  • nf (int) – number of active flavours

Returns:

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

\[\begin{split}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 ] \\ & + \frac{a_s^3}{6} [ -2 b_3 c_0 - b_1^3 c_0 (1 + c_0) (2 + c_0) - 2 b_2 c_1 \\ & - 3 b_2 c_0 c_1 + b_1^2 (2 + 3 c_0 (2 + c_0)) c_1 + c_1^3 + 3 c_1 c_2 \\ & + b_1 (b_2 c_0 (4 + 3 c_0) - 3 (1 + c_0) c_1^2 - (2 + 3 c_0) c_2) + 2 c_3 ]\end{split}\]

Return type:

float

eko.msbar_masses.ker_dispatcher(q2_to, q2m_ref, strong_coupling, xif2, nf)[source]

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

Parameters:
  • q2_to (float) – final scale

  • q2m_ref (float) – initial scale

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

  • xif2 (float) – factorization to renormalization scale ratio

  • nf (int) – number of active flavours

Returns:

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

Return type:

ker

eko.msbar_masses.compute_matching_coeffs_up(nf: int)[source]

Upward \(\overline{MS}\) matching coefficients.

Used at threshold when moving to a regime with more flavors [LS15].

Note [LS15] (eq 5.) reports the backward relation, multiplied by a factor of 4 (and 4^2 …)

Parameters:

nf (int) – number of active flavors in the lower patch

Returns:

downward matching coefficient matrix

Return type:

dict

eko.msbar_masses.compute_matching_coeffs_down(nf: int)[source]

Downward \(\overline{MS}\) matching coefficients.

Used at threshold when moving to a regime with less flavors [LS15] 5.

Parameters:

nf (int) – number of active flavors in the lower patch

Returns:

downward matching coefficient matrix

Return type:

dict

eko.msbar_masses.solve(m2_ref, q2m_ref, strong_coupling, nf_ref, xif2)[source]

Compute the \(\overline{MS}\) masses.

Solves the equation \(m_{\overline{MS}}(m) = m\) for a fixed number of nf.

Parameters:
  • m2_ref (float) – squared initial mass reference

  • q2m_ref (float) – squared initial scale

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

  • nf_ref (int) – number of active flavours at the scale q2m_ref, where the solution is searched

  • xif2 (float) – \(\mu_F^2/\mu_R^2\)

Returns:

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

Return type:

float

eko.msbar_masses.evolve(m2_ref, q2m_ref, strong_coupling, thresholds_ratios, xif2, q2_to, nf_ref=None, nf_to=None)[source]

Perform the \(\overline{MS}\) mass evolution up to given scale.

It allows for different number of active flavors.

Parameters:
  • m2_ref (float) – squared initial mass reference

  • q2m_ref (float) – squared initial scale

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

  • xif2 (float) – \(\mu_F^2/\mu_R^2\)

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

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

  • nf_to (int) – number of flavor active at the target scale

Returns:

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

Return type:

float

eko.msbar_masses.compute(masses_ref: HeavyQuarks[ReferenceRunning[float]], couplings: CouplingsInfo, order: Tuple[int, int], evmeth: CouplingEvolutionMethod, matching: List[float], xif2: float = 1.0)[source]

Compute the \(\overline{MS}\) masses.

Computation is performed solving the equation \(m_{\bar{MS}}(\mu) = \mu\) for each heavy quark and consistent boundary contitions.

Parameters:
  • masses_ref – reference scale squared (a.k.a. \(Q_{ref}\))

  • couplings – couplings configuration

  • order – perturbative order

  • evmeth – evolution method

  • matching – threshold matching scale ratios

  • xif2 – squared ratio of factorization to central scale

Returns:

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

Return type:

list

eko.version module

Define versions.

Both package and data versions are stored here.