Detailed package info

This is a comprehensive documentation of the code for developers and power users. It is automatically generated from documentation strings embedded in the code itself.

Note

For all code examples it is assumed that gepard package is first imported as:

import gepard as g

1. Theory framework

theory module

Definition of theory frameworks.

class gepard.theory.Theory(**kwargs)[source]

Bases: object

Class of theory frameworks for calculation of observables.

Parameters:
  • name – short unique model name

  • texname – TeX model name for (e.g. for plot annotation)

  • description – longer description of the model

This is a base class for a complete theory framework. It implements the generic routines for calculation of observables and uncertainty propagation.

__init__(**kwargs) None[source]
chisq_single(points: DataSet, asym: bool = False, **kwargs) float[source]

Return total chi-square.

Parameters:
  • points – measurements with uncertainties, observable is named in observable attribute of each DataPoint, value is val and uncertainty is err.

  • asym – if measurements provide asymmetric uncertainties errplus and errminus, this enables their usage

  • **kwargs – keyword arguments

Returns:

Total chi^square for points.

Notes

If the theory or model provide uncertainties, they are ignored - only experimental uncertainties are taken into account.

pull(pt: DataPoint)[source]

Return pull of a single Datapoint.

chisq(points: DataSet, asym: bool = False, **kwargs) float

Return total chi-square.

Parameters:
  • points – measurements with uncertainties, observable is named in observable attribute of each DataPoint, value is val and uncertainty is err.

  • asym – if measurements provide asymmetric uncertainties errplus and errminus, this enables their usage

  • **kwargs – keyword arguments

Returns:

Total chi^square for points.

Notes

If the theory or model provide uncertainties, they are ignored - only experimental uncertainties are taken into account.

predict(pt: DataPoint, uncertainty: bool = False, **kwargs) float[source]

Give prediction for DataPoint pt.

Parameters:
  • pt – instance of DataPoint

  • uncertainty – if available, produce tuple (mean, uncertainty)

  • **kwargs – keyword arguments

Keyword Arguments:
  • observable – string. Default is pt.observable. It is acceptable also to pass CFF or x-space GPD as observable, e.g., observable = ‘ImH’

  • parameters – dictionary which will temporarily update model’s parameters

  • orig_conventions – give prediction using original conventions of the given DataPoint (e.g. for plotting)

Returns:

Predicted value for observable. If uncertainty is requested, tuple (value, uncertainty) is returned.

_XGAMMA_int(t: ndarray, pt: DataPoint)[source]

Return gamma* DVCS or DVMP cross sections differential in t (array version).

Parameters:
  • t – array of Mandelstam t values

  • pt – datapoint with the rest of kinematics

Returns:

numpy array of cross-sections

XGAMMA(pt: DataPoint) float[source]

Return total gamma* DVCS or DVMP cross section.

Parameters:

pt – datapoint

Returns:

Total or differential cross-section. For DVMP, this calculates only longitudinal_gamma* part. If pt has no momentum transfer defined, it calculates total xs.

model module

Base classes for models of “soft” hadronic structure functions.

  • Elastic electromagnetic form factors

  • Generalized parton distribution functions (GPDs)

  • Compton form factors (for DVCS)

Actual implementations are in other files. Here only generic treatment of generalized parametrized model is coded.

Todo

  • Model defined by grid of numbers

  • Flavored models

class gepard.model.Model(**kwargs)[source]

Bases: Theory

Base class for all models.

Instance of this class specifies structure of relevant hadrons. Methods provided are typically GPDs, CFFs, elastic FFs, DVMP transition TFFs, DAs, etc. Main subclasses are:

  • ParameterModel which depends on real parameters (some of which can be provided by minimization routine in the fitting procedure).

  • NeuralModel where structure functions are represented as neural nets (not yet implemented)

  • GridModel where values of structure functions are represented as grids of numbers, which may be interpolated (not yet implemented)

These are then further subclassed to model actual structure functions.

__init__(**kwargs) None[source]
class gepard.model.ParameterModel(**kwargs)[source]

Bases: Model

Base class for all parametrized models.

parameters

dict {‘par0’: float, …}. Actual value of parameter.

parameters_fixed

dict {‘par0’: bool, …}. Is parameter value fixed? Considered False if non-existent. Once Fitter object is created, user should manipulate this dict only via Fitter methods.

parameters_limits

dict {‘par0: (float, float), …}. Allowed range for fitting. Considered (-inf, inf) if non-existent. Once Fitter object is created, user should manipulate this dict only via Fitter methods.

__init__(**kwargs) None[source]
add_parameters(newpars: dict)[source]

Append newpars to parameters.

add_parameters_limits(newlimits: dict)[source]

Append newlimits to parameters_limits.

_release_parameters(*pars: str)[source]

Release parameters for fitting.

Parameters:

*pars – Names of parameters to be released

Notes

User should relase and fix parameters using methods of Fitter instance. This then calls this private ParameterModel method.

_fix_parameters(*pars: str)[source]

Fix parameters so they are not fitting variables.

Parameters:

*pars – Names of parameters to be fixed. If first name is ‘ALL’ then fix all parameters.

Notes

User should relase and fix parameters using methods of Fitter instance. This then calls this private ParameterModel method.

free_parameters() List[str][source]

Return list of names of free fitting parameters.

print_parameters()[source]

Print fitting parameters and their errors.

dvcs module

Deeply Virtual Compton Scattering (DVCS) observables.

class gepard.dvcs.DVCS(**kwargs)[source]

Bases: Theory

DVCS observables base class.

Implements cross-section for electroproduction of polarized real photon, on polarized nucleon target, as well as derived cross-sections and asymmetries. Subclasses should, as their methods, implement expressions for squared BH, interference and squared DVCS amplitudes: - TBH2unp, TINTunp, TDVCS2unp (for unpolarized target) - TBH2LP, TINTLP, TDVCS2LP (for longitudinally polarized target) - TBH2TP, TINTTP, TDVCS2TP (for transversally polarized target)

Implemented subclasses are:

  • BMK - From Belitsky, Mueller and Kirchner arXiv:hep-ph/0112108

  • hotfixedBMK - simple valence xB region improvement by Dieter

  • BM10ex

  • BM10 - Belitsky, Mueller and Ji

  • BM10tw2 - BM10 with higher twists set to zero

PreFacSigma(pt: DataPoint) float[source]

Overall prefactor of 4-fold XS.

Parameters:

pt – instance of DataPoint

Returns:

Overal prefactor from Eq. (22) of BMK, times 2pi because of proton Phi integration and times y/Q2 because of dy -> dQ2. Converted to nanobarns.

XS(pt: DataPoint, **kwargs) float[source]

Differential 4-fold e p –> e p gamma cross section.

Parameters:
  • pt – DataPoint instance

  • **kwargs – keyword arguments

Keyword Arguments:
  • vars (dict) – dict which can override pt kinematics, e.g., vars={‘phi’: 0}

  • flip (string, list) – flip the sign of pt attribute, e.g., flip=’in1polarization’

Returns:

Cross section.

This is for calculation of fully beam and target polarized cross-section. Most of the usually measured cross-sections and asymmetries are derived from this one.

Todo

  • transversal target BMK/Trento switching

_XGAMMA_DVCS_t_Approx(pt)[source]

Partial DVCS (gamma* p -> gamma p) cross section differential in t.

Parameters:

pt – instance of DataPoint

Returns:

Cross section. Approx. formula used in NPB10 paper.

_XGAMMA_DVCS_t_Ex(pt)[source]

Partial DVCS (gamma* p -> gamma p) cross section differential in t.

Parameters:
  • pt – instance of DataPoint

  • **kwargs – keyword arguments

Returns:

Cross section

Todo

Htilde is missing here. Small at low-x, but bigger than some terms that are retained.

_XGAMMA_DVCS_t(pt)

Partial DVCS (gamma* p -> gamma p) cross section differential in t.

Parameters:
  • pt – instance of DataPoint

  • **kwargs – keyword arguments

Returns:

Cross section

Todo

Htilde is missing here. Small at low-x, but bigger than some terms that are retained.

XSintphi(pt: DataPoint, **kwargs) float[source]

XS integrated over azimuthal angle.

Parameters:
  • pt – instance of DataPoint

  • **kwargs – keyword arguments

Returns:

Cross section (unpolarized)

Todo

This is unnecessary function. Should be changed to XUU.

_phiharmonic(fun, pt, **kwargs)[source]

Return fun evaluated for phi=pt.phi, or proper harmonic.

_XUU(pt, **kwargs)[source]

4-fold beam helicity-independent cross section.

XUU(pt, **kwargs)[source]

4-fold beam helicity-independent cross section.

XUUw(pt, **kwargs)[source]

Weighted 4-fold beam helicity-independent cross section.

_XLU(pt, **kwargs)[source]

4-fold beam helicity-dependent cross section.

XLU(pt, **kwargs)[source]

4-fold beam helicity-dependent cross section.

XLUw(pt, **kwargs)[source]

Weighted 4-fold beam helicity-dependent cross section.

_XUD(pt, **kwargs)[source]

Cross-section difference on (long. or trans.) polarized target.

XUL(pt, **kwargs)[source]

Cross-section difference on longitudinally polarized target.

XUT(pt, **kwargs)[source]

Cross-section difference on transversally polarized target.

_XCLU(pt, **kwargs)[source]

4-fold beam charge-spin cross section difference.

Parameters:
  • pt – instance of DataPoint

  • **kwargs – keyword arguments

Returns:

D_{CS,U}/2, as measured by COMPASS.

XCLU(pt, **kwargs)[source]

4-fold beam charge-spin cross section difference.

Parameters:
  • pt – instance of DataPoint

  • **kwargs – keyword arguments

Returns:

D_{CS,U}/2, as measured by COMPASS. (NOT TESTED.)

_XCUU(pt, **kwargs)[source]

4-fold beam charge-spin cross section sum.

Parameters:
  • pt – instance of DataPoint

  • **kwargs – keyword arguments

Returns:

S_{CS,U}/2, as measured by COMPASS.

XCUU(pt, **kwargs)[source]

4-fold beam charge-spin cross section sum.

Parameters:
  • pt – instance of DataPoint

  • **kwargs – keyword arguments

Returns:

S_{CS,U}/2, as measured by COMPASS. (NOT TESTED.)

_AC(pt, **kwargs)[source]

Calculate beam charge asymmtery.

AC(pt, **kwargs)[source]

Calculate beam charge asymmetry.

_ALU(pt, **kwargs)[source]

Calculate beam spin asymmetry.

_ALUapprox(pt, **kwargs)[source]

Calculate beam spin asymmetry.

Parameters:
  • pt – instance of DataPoint

  • **kwargs – keyword arguments

Returns:

Beam spin asymmetry A_LU. For harmonics, fast approximate formula is used.

_ALUexact(pt, **kwargs)[source]

Calculate beam spin asymmetry.

ALU(pt, **kwargs)

Calculate beam spin asymmetry.

_TSA(pt, **kwargs)[source]

Target spin asymmetry (transversal or longitudinal).

TSA(pt, **kwargs)[source]

Calculate target spin asymmetry (AUL or AUT).

Parameters:
  • pt – instance of DataPoint

  • **kwargs – keyword arguments

Returns:

Target spin asymmetry. Whether AUL or AUT is calculated is determined by the pt.in2polarizationvector being ‘L’ or ‘T’.

AUL(pt, **kwargs)[source]

Calculate longitudinal target spin asymmetry.

AUT(pt, **kwargs)[source]

Calculate transversal target spin asymmetry.

_BTSA(pt, **kwargs)[source]

Calculate beam-target spin asymmetry (BTSA).

BTSA(pt, **kwargs)[source]

Calculate beam-target spin asymmetry (ALL or ALT).

Parameters:
  • pt – instance of DataPoint

  • **kwargs – keyword arguments

Returns:

Beam-target spin asymmetry. Whether ALL or ALT is calculated is determined by the pt.in2polarizationvector being ‘L’ or ‘T’.

ALL(pt, **kwargs)[source]

Calculate beam and longitudinal target spin asymmetry.

ALT(pt, **kwargs)[source]

Calculate beam and transversal target spin asymmetry.

_CBTSA(pt, chargepar=-1, **kwargs)[source]

Calculate double charge and spin asymmetries.

ALTI(pt, **kwargs)[source]

Calculate {A_LT,I} as defined by HERMES 1106.2990 Eq. (19).

ALTBHDVCS(pt, **kwargs)[source]

Calculate {A_LT,BHDVCS} as defined by HERMES 1106.2990 Eq. (18).

_ALUI(pt, **kwargs)[source]

Beam spin asymmetry, interference part.

ALUI(pt, **kwargs)[source]

Calculate beam spin asymmetry, interference part.

_ALUDVCS(pt, **kwargs)[source]

Calculate ALU as defined by HERMES 0909.3587 Eq. (2.3).

ALUDVCS(pt, **kwargs)[source]

Calculate beam spin asymmetry, DVCS part.

_AUTI(pt, **kwargs)[source]

Calculate trans. target asymmetry, defined by HERMES 0802.2499 Eq. (15).

AUTI(pt, **kwargs)[source]

Calculate trans. target asymmetry, defined by HERMES 0802.2499 Eq. (15).

_AUTDVCS(pt, **kwargs)[source]

Calculate trans. target asymmetry, defined by HERMES 0802.2499 Eq. (14).

AUTDVCS(pt, **kwargs)[source]

Calculate trans. target asymmetry, defined by HERMES 0802.2499 Eq. (14).

AC0minusr1(pt)[source]

Return some combination liked by Dieter.

XLUw2C(pt)[source]

Im(C^I) as defined by HALL A.

XUUw2C(pt)[source]

Re(C^I) or Re(C^I + Del C^I) as defined by HALL A.

Parameters:

pt – instance of DataPoint

Returns:

Harmonic coefficients.

Todo

Although it is attributed to FTn=0, Re(C^I + Del C^I) is only a part of zeroth harmonic.

XwA(pt)[source]

Ratio of first two cos harmonics of w-weighted cross section.

BCSA(pt, **kwargs)[source]

Beam charge-spin asymmetry as measured by COMPASS.

cff module

Definitions of Compton Form Factors.

gepard.cff._flatten(T)[source]

Flatten the tuple.

class gepard.cff.CFF(**kwargs)[source]

Bases: ParameterModel

Base class for CFFs.

CFFs are set to be zero here. Actual models are built by subclassing this.

allCFFs = ['ImH', 'ReH', 'ImE', 'ReE', 'ImHt', 'ReHt', 'ImEt', 'ReEt']
allCFFsb = ['ImH', 'ReH', 'ImE', 'ReE', 'ImHt', 'ReHt', 'ImEt', 'ReEb']
allCFFeffs = ['ImHeff', 'ReHeff', 'ImEeff', 'ReEeff', 'ImHteff', 'ReHteff', 'ImEteff', 'ReEteff']
__init__(**kwargs)[source]
print_CFFs(pt: DataPoint, format: str = None)[source]

Print values of CFFs at given kinematic point.

Parameters:
  • pt – DataPoint instance

  • format – ‘mma’ to get output pastable into Mathematica

ReEb(pt: DataPoint)[source]

Return Re(E-bar).

Parameters:

pt – DataPoint instance

Returns:

xi*ReEt

cff(pt: DataPoint) ndarray[source]

Return array(ReH, ImH, ReE, …) for kinematic point.

Parameters:

pt – DataPoint instance

Returns:

Eight CFFs (ReH, ImH, ReE, …, ImEt).

This is generic top-class method that only gathers results of dedicated ReH(), ImH(), etc. functions. The idea is that this method can be implemented more efficiently in subclases where components of CFFs can be calculated at the same time or in parallel.

is_within_model_kinematics(pt: DataPoint)[source]

Check if kinematics is sensible.

Parameters:

pt – DataPoint instance

Returns:

True if pt kinematics is inside the model validity range.

ImE(pt)
ImEeff(pt)
ImEt(pt)
ImEteff(pt)
ImH(pt)
ImHeff(pt)
ImHt(pt)
ImHteff(pt)
ReE(pt)
ReEeff(pt)
ReEt(pt)
ReEteff(pt)
ReH(pt)
ReHeff(pt)
ReHt(pt)
ReHteff(pt)
name = 'ReEteff'
class gepard.cff.MellinBarnesCFF(**kwargs)[source]

Bases: CFF

Class of models built by Mellin-Barnes integration.

__init__(**kwargs) None[source]
cff(pt: DataPoint) ndarray[source]

Return array(ReH, ImH, ReE, …) for kinematic point.

ReH(pt: DataPoint) float[source]

Return Re(CFF H) for kinematic point.

ImH(pt: DataPoint) float[source]

Return Im(CFF H) for kinematic point.

ReE(pt: DataPoint) float[source]

Return Re(CFF E) for kinematic point.

ImE(pt: DataPoint) float[source]

Return Im(CFF E) for kinematic point.

class gepard.cff.DispersionCFF(**kwargs)[source]

Bases: CFF

Use dispersion relations for real parts of CFFs.

methods: ReH, ReE, ReHt, ReEt, subtraction Subclass should implement ansaetze for ImH, ImE, ImHt, ImEt and subtraction. This class implements just dispersion integrals.

__init__(**kwargs) None[source]
dispargV(x, fun, pt)[source]

Integrand of the dispersion integral (vector case).

Parameters:
  • x – integration variable, longitudinal momentum fraction

  • fun – function providing imaginary part of CFF

  • pt – DataPoint instance

Returns:

Integrand of dispersion integral. With variable change x->x^(1/(1-ga))=u in order to tame the singularity at x=0.

dispargA(x, fun, pt)[source]

Integrand of the dispersion integral (axial-vector case).

Parameters:
  • x – integration variable, longitudinal momentum fraction

  • fun – function providing imaginary part of CFF

  • pt – DataPoint instance

Returns:

Integrand of dispersion integral. With variable change x->x^(1/(1-ga))=u in order to tame the singularity at x=0.

subtraction(pt)[source]

Subtraction constant.

_ReV(pt, imfun, subsign)[source]

Real part of vector CFFs H and E.

Parameters:
  • pt – DataPoint instance

  • imfun – function providing imaginary part of CFF

  • subsign – sign of subtraction constant (-1 for H, +1 for E)

Returns:

Dispersion integral over imfun minus/plus subtraction constant.

_ReA(pt, imfun)[source]

Real part of axial vector CFFs Ht and Et.

Parameters:
  • pt – DataPoint instance

  • imfun – function providing imaginary part of CFF

Returns:

Dispersion integral over imfun.

ReH(pt, imfun=None)[source]

Real part of CFF H.

ReE(pt, imfun=None)[source]

Real part of CFF E.

ReHt(pt, imfun=None)[source]

Real part of CFF Ht.

ReEt(pt, imfun=None)[source]

Real part of CFF Et.

class gepard.cff.PionPole[source]

Bases: object

Various options for pion-pole contribution.

DMfixpole(pt)[source]

Return fixed pion-pole as used by Dieter.

DMfreepole(pt)[source]

Return free pion-pole as proposed by Dieter.

class gepard.cff.DispersionFixedPoleCFF(**kwargs)[source]

Bases: DispersionCFF, PionPole

Model for CFFs as in arXiv:0904.0458.

__init__(**kwargs)[source]
subtraction(pt)[source]

Dispersion relations subtraction constant.

ImH(pt: DataPoint, xi: float | ndarray = 0)[source]

Imaginary part of CFF H.

ImHt(pt: DataPoint, xi: float | ndarray = 0)[source]

Imaginary part of CFF Ht i.e. tilde{H}.

ImE(pt: DataPoint, xi: float | ndarray = 0)[source]

Imaginary part of CFF E.

ReEt(pt: DataPoint)[source]

Instead of disp. rel. use pole formula.

class gepard.cff.DispersionFreePoleCFF(**kwargs)[source]

Bases: DispersionFixedPoleCFF

Model for CFFs as in arXiv:0904.0458. + free pion pole.

__init__(**kwargs)[source]
ReEt(pt)[source]

Instead of disp. rel. use pole formula.

class gepard.cff.HybridCFF(**kwargs)[source]

Bases: MellinBarnesCFF, DispersionFixedPoleCFF, CFF

Combines MB model for small xB and DR model for valence xB.

Todo

Check general consistency for xi != pt.xi for imag CFFs.

is_within_model_kinematics(pt)[source]

Test kinematics of datapoint.

ImH(pt, xi=0)[source]

Imaginary part of CFF H.

ReH(pt)[source]

Real part of CFF H.

ImE(pt, xi=0)[source]

Imaginary part of CFF E.

ReE(pt)[source]

Real part of CFF E.

ImHt(pt, xi=0)[source]

Imaginary part of CFF Ht.

ReHt(pt)[source]

Real part of CFF Ht.

ImEt(pt)[source]

Imaginary part of CFF Et.

ReEt(pt)[source]

Real part of CFF Et. Provided by subclasses.

class gepard.cff.HybridFixedPoleCFF(**kwargs)[source]

Bases: HybridCFF

HybridCFF model with fixed pion pole model for ReEt.

ReEt(pt)[source]

Instead of disp. rel. use fixed pion pole formula.

class gepard.cff.HybridFreePoleCFF(**kwargs)[source]

Bases: HybridCFF

HybridCFF model with free pion pole model for ReEt.

ReEt(pt)[source]

Instead of disp. rel. use free pion pole formula.

gpd module

GPD models in conformal moment j-space.

param j:

complex-space conformal moment

param t:

momentum transfer squared

par: dictionary {‘parmeter1_name’: value1, …}

returns:

conformal moments of GPD as npts x 4 array, where npts is number of points on Mellin-Barnes contour and 4 flavors are: (1) – singlet quark (2) – gluon (3) – u_valence (4) – d_valence

Valence here means “valence-like GPD”. See hep-ph/0703179 for description.

gepard.gpd.qj(j: ndarray, t: float, poch: int, norm: float, al0: float, alp: float, alpf: float = 0, val: int = 0) ndarray[source]

GPD building block Q_j with reggeized t-dependence.

Parameters:
  • j – complex conformal moment

  • t – momentum transfer

  • poch – pochhammer (beta+1)

  • norm – overall normalization

  • al0 – Regge intercept

  • alp – Regge slope (leading pole)

  • alpf – Regge slope (full)

  • val – 0 for sea, 1 for valence partons

Returns:

conformal moment of GPD with Reggeized t-dependence,

Examples

>>> qj(0.5+0.3j, 0.,  4, 3, 0.5, -0.8)  
(5.668107685...-4.002820...j)

Notes

There are two implementations of Regge t dependence, one uses alpf \(\alpha'\) parameter (alp should be set to zero)

\[Q_j = N \frac{B(1-\alpha(t)+j, \beta+1)}{B(2-\alpha_0, \beta+1)}\]

This is then equal to Eq. (41) of arXiv:hep-ph/0605237 without residual \(F(\Delta^2)\).

The default uses alp parameter, and expression is from Eq. (40) or (41) of arXiv:0904.0458, which takes into account only the leading pole:

\[Q_j = N \frac{B(1-\alpha_0+j, \beta+1)}{B(2-\alpha_0, \beta+1)} \frac{1+j-\alpha_0}{1+j-\alpha(t)}\]

and is again defined without residual beta(t) from Eq. (19).

gepard.gpd.betadip(j: ndarray, t: float, m02: float, delm2: float, pp: int) ndarray[source]

GPD residual dipole t-dependence function.

Parameters:
  • j – conformal moment

  • t – momentum transfer

  • m02 – mass param

  • delm2 – mass param - interplay with j

  • pp – exponent (=2 for dipole)

Returns:

Dipole residual t-dependence (beta from Eq. (19) of 0904.0458)

gepard.gpd.betaexp(t: float, m02: float) float[source]

GPD residual exponential t-dependence function.

Parameters:
  • t – momentum transfer

  • m02 – mass parameter

Returns:

Exponential residual t-dependence (beta from Eq. (19) of 0904.0458)

gepard.gpd.toy(j: complex, *args) Tuple[complex, complex, complex, complex][source]

Return ‘toy’ (no params) singlet GPD ansatz.

gepard.gpd.test(j: complex, t: float, par: dict) Tuple[complex, complex, complex, complex][source]

Return simple testing singlet GPD ansatz.

gepard.gpd.singlet_ng_constrained(j: ndarray, t: float, par: dict, residualt: str = 'dipole') ndarray[source]

Singlet-only GPD ansatz, with ng parameter constrained by sum-rule.

Parameters:
  • j – conformal moment

  • t – momentum transfer

  • par – parameters dict

  • residualt – residual t-dependence type (‘dipole’ or ‘exp’)

Returns:

GPD ansatz, as numpy array

Notes

This ansatz is used for all published KM fits, for the sea parton part.

gepard.gpd.singlet_ng_constrained_E(j: ndarray, t: float, par: dict, residualt: str = 'dipole') ndarray[source]

Singlet-only GPD ansatz, with ng parameter constrained by sum-rule.

Parameters:
  • j – conformal moment

  • t – momentum transfer

  • par – parameters dict

  • residualt – residual t-dependence type (‘dipole’ or ‘exp’)

Returns:

GPD ansatz, as numpy array

Notes

This ansatz is used for all published KM fits, for the sea parton part.

gepard.gpd.ansatz07(j: ndarray, t: float, par: dict) ndarray[source]

GPD ansatz from paper hep-ph/0703179.

gepard.gpd.ansatz07_fixed(j: ndarray, t: float, type: str) ndarray[source]

GPD ansatz from hep-ph/0703179 with fixed parameters.

Parameters:
  • j – conformal moment

  • t – momentum transfer

  • type – ‘soft’, ‘hard’, ‘softNS’, ‘hardNS’

Returns:

GPD ansatz, as numpy array

Notes

This is the same as ansatz07, only instead of passing parameter dict, user passes type string choosing particular fixed parameter choices from the paper above.

class gepard.gpd.GPD(**kwargs)[source]

Bases: ParameterModel

Base class of all GPD models.

Parameters:
  • p – pQCD order (0 = LO, 1 = NLO, 2 = NNLO)

  • scheme – pQCD scheme (‘msbar’ or ‘csbar’)

  • nf – number of active quark flavors

  • Q02 – Initial Q0^2 for GPD evolution.

  • r20 – Initial mu0^2 for alpha_strong definition.

  • asp – alpha_strong/(2*pi) at scale r20 for (LO, NLO, NNLO)

  • residualt – residual t dependence (‘dipole’ or ‘exp’)

__init__(**kwargs) None[source]
class gepard.gpd.ConformalSpaceGPD(**kwargs)[source]

Bases: GPD, MellinBarnes

Base class of GPD models built in conformal moment space.

Parameters:
  • c – intersection of Mellin-Barnes curve with real axis

  • phi – angle of Mellin-Barnes curve with real axis

Notes

This just takes care of initialization of Mellin-Barnes contour points, and evolved Wilson coeffs, and MB Gauss integration weights. Actual choice and code for GPDs is provided by subclasses.

__init__(**kwargs) None[source]
pw_strengths()[source]

Strengths of SO(3) partial waves.

pw_strengths_E()[source]

Strengths of SO(3) partial waves for gpd E.

H(eta: float, t: float) ndarray[source]

Return (npts, 4) array H_j^a for all j-points and 4 flavors.

E(eta: float, t: float) ndarray[source]

Return (npts, 4) array E_j^a for all j-points and 4 flavors.

Hx(pt: DataPoint) ndarray[source]

Return x-space GPD H.

Parameters:

pt – datapoint with kinematics info

Returns:

x-space GPD. 3-dim vector (singlet quark, gluon, non-singlet quark) is returned by transforming original conformal moment space (j-space) model.

Todo

Non-singlet component is set to zero. We need to check normalization/symmetrization first.

skewness_Hx(pt: DataPoint) ndarray[source]

Return skewness of GPD H.

Parameters:

pt – datapoint with kinematics info

Returns:

2-dim vector (singlet quark, gluon) of ratio GPD(x, x, t) / GPD(x, 0, t).

Todo

Non-singlet component is missing. See Todo for Hx.

Ex(pt: DataPoint) ndarray[source]

Return x-space GPD E.

Parameters:

pt – datapoint with kinematics info

Returns:

x-space GPD. 3-dim vector (singlet quark, gluon, non-singlet quark) is returned by transforming original conformal moment space (j-space) model.

Todo

Non-singlet component is set to zero. We need to check normalization/symmetrization first.

class gepard.gpd.TestGPD(**kwargs)[source]

Bases: ConformalSpaceGPD

Simple testing ansatz for GPDs.

__init__(**kwargs) None[source]
H(eta: float, t: float) ndarray[source]

Return (npts, 4) array H_j^a for all j-points and 4 flavors.

class gepard.gpd.PWNormGPD(**kwargs)[source]

Bases: ConformalSpaceGPD

Singlet-only model for GPDs with three SO(3) partial waves.

Parameters:
  • p – pQCD order (0 = LO, 1 = NLO, 2 = NNLO)

  • scheme – pQCD scheme (‘msbar’ or ‘csbar’)

  • nf – number of active quark flavors

  • Q02 – Initial Q0^2 for GPD evolution.

  • r20 – Initial mu0^2 for alpha_strong definition.

  • asp – alpha_strong/(2*pi) at scale r20 for (LO, NLO, NNLO)

  • residualt – residual t dependence (‘dipole’ or ‘exp’)

  • c – intersection of Mellin-Barnes curve with real axis

  • phi – angle of Mellin-Barnes curve with real axis

Notes

Subleading PWs are proportional to the leading one, and only their norms are fitting parameters. Norms of second PWs is given by parameers ‘secs’ (quarks) and ‘secg’ gluons, and norm of third PWs is given by ‘this’ and ‘thig’.

This is used for modelling sea partons in KM10-KM20 models.

__init__(**kwargs) None[source]
H(eta: float, t: float) ndarray[source]

Return (npts, 4) array H_j^a for all j-points and 4 flavors.

E(eta: float, t: float) ndarray[source]

Return (npts, 4) array E_j^a for all j-points and 4 flavors.

dvmp module

Deeply Virtual Meson Production (DVMP) observables.

class gepard.dvmp.DVMP(**kwargs)[source]

Bases: Theory

DVMP observables.

Implements cross-section for electroproduction of meson.

_XGAMMA_DVMP_t_Approx(pt: DataPoint) float[source]

Partial (longitudinal) gamma* p -> V p cross section w.r.t. Mandelstam t.

Parameters:

pt – instance of DataPoint

Returns:

Cross-section differential in t. Approximate formula valid for small xB.

_XGAMMA_DVMP_t(pt: DataPoint) float

Partial (longitudinal) gamma* p -> V p cross section w.r.t. Mandelstam t.

Parameters:

pt – instance of DataPoint

Returns:

Cross-section differential in t. Approximate formula valid for small xB.

class gepard.dvmp.MellinBarnesTFF(**kwargs)[source]

Bases: ParameterModel

DVMP Transition Form Factors modelled as Mellin-Barnes integral.

__init__(**kwargs)[source]
tff(xi: float, t: float, Q2: float, meson: str = 'rho0') ndarray[source]

Return array(ReH_V, ImH_V, ReE_V, …) of DVMP transition FFs.

Only vector mesons are implemented, thus _V.

ImH_V(pt: DataPoint) ndarray[source]

Return Im(TFF H) for kinematic point.

ReH_V(pt: DataPoint) ndarray[source]

Return Re(TFF H) for kinematic point.

eff module

Models for Elastic (Electromagnetic) Form Factors.

class gepard.eff.EFF[source]

Bases: object

Dirac and Pauli elastic (electromagnetic) form factors F_1 and F_2.

class gepard.eff.DipoleEFF[source]

Bases: EFF

Dipole approximation from Dieter’s notebook.

F1(pt)[source]

Dirac elastic proton form factor - dipole parametrization.

F2(pt)[source]

Pauli elastic proton form factor - dipole parametrization.

class gepard.eff.KellyEFF[source]

Bases: EFF

Kelly’s approximation from DM’s notebook.

F1(pt)[source]

Dirac elastic nucleon form factor - Kelly’s parametrization.

F2(pt)[source]

Dirac elastic nucleon form factor - Kelly’s parametrization.

_pF1(t)[source]

Dirac elastic proton form factor - Kelly’s parametrization.

_pF2(t)[source]

Pauli elastic proton form factor - Kelly’s parametrization.

_nF1(t)[source]

Dirac elastic neutron form factor - Kelly’s parametrization.

_nF2(t)[source]

Pauli elastic neutron form factor - Kelly’s parametrization.

class gepard.eff.ZeroEFF[source]

Bases: EFF

Set F1=F2=0 to get just DVCS^2.

F1(pt)[source]

Elastic em Dirac form factor F1 set to zero.

F2(pt)[source]

Elastic em Pauli form factor F2 set to zero.

2. Perturbative QCD

adim module

QCD anomalous dimensions.

Note

Functions in this module have as first argument Mellin moment n=j+1, where j is conformal moment used everywhere else. Thus, they should always be called as f(j+1, …).

gepard.adim.non_singlet_LO(n: complex, nf: int, prty: int = 1) complex[source]

Non-singlet LO anomalous dimension.

Parameters:
  • n (complex) – which moment (= Mellin moment for integer n)

  • nf (int) – number of active quark flavors

  • prty (int) – 1 for NS^{+}, -1 for NS^{-}, irrelevant at LO

Returns:

Non-singlet LO anomalous dimension.

gepard.adim.singlet_LO(n: complex, nf: int, prty: int = 1) ndarray[source]

Singlet LO anomalous dimensions.

Parameters:
  • n (complex) – which moment (= Mellin moment for integer n)

  • nf (int) – number of active quark flavors

  • prty (int) – C parity, irrelevant at LO

Returns:

2x2 complex matrix ((QQ, QG),

(GQ, GG))

gepard.adim.non_singlet_NLO(n: complex, nf: int, prty: int) complex[source]

Non-singlet anomalous dimension.

Parameters:
  • n (complex) – which moment (= Mellin moment for integer n)

  • nf (int) – number of active quark flavors

  • prty (int) – 1 for NS^{+}, -1 for NS^{-}

Returns:

Non-singlet NLO anomalous dimension.

gepard.adim.singlet_NLO(n: complex, nf: int, prty: int = 1) ndarray[source]

Singlet NLO anomalous dimensions matrix.

Parameters:
  • n (complex) – which moment (= Mellin moment for integer n)

  • nf (int) – number of active quark flavors

  • prty (int) – C parity

Returns:

Matrix (LO, NLO) where each is in turn 2x2 complex matrix ((QQ, QG), (GQ, GG))

gepard.adim.block(n: complex, nf: int) ndarray[source]

LO+NLO block anomalous dimensions matrix.

Parameters:
  • n (complex) – which moment (= Mellin moment for integer n)

  • nf (int) – number of active quark flavors

Returns:

block[k, p, i, j]
  • k is index if n is array, otherwise 0

  • p is pQCD order (0=LO, 1=NLO)

  • i, j in [Q, G, NS+, NS-], forming 4x4 complex matrix ((QQ, QG, 0, 0), (GQ, GG, 0, 0), ( 0, 0, NS+, 0), ( 0, 0, 0, NS-))

qcd module

Perturbative QCD functions.

Coefficients of beta function of QCD up to order \(\alpha_s^5\) (NNNLO), normalized by

\[\frac{d a_s}{d \ln \mu_r^2} = \beta_0 a_s^2 + \beta_1 a_s^3 + \cdots\]

with

\[a_s = \frac{\alpha_{s}}{4\pi}.\]

Beyond NLO the QCD colour factors are hard-wired in this routine, and the numerical coefficients are truncated to six digits.

QCD coupling alpha strong is obtained by integrating the evolution equation for a fixed number of massless flavours NF. Except at leading order (LO), result is obtained using a fourth-order Runge-Kutta integration.

gepard.qcd.beta(p: int, nf: int) float[source]

QCD beta function coefficient.

Parameters:
  • p (int) – pQCD order, 0=LO, 1=NLO, 2=NNLO

  • nf (int) – number of active quark flavors

Returns:

QCD beta function coefficient

Return type:

float

Examples

>>> beta(0, 3)
-9.0
gepard.qcd._fbeta1(a: float, nf: int) float[source]
gepard.qcd.as2pf(p: int, nf: int, r2: float, as0: float, r20: float) float[source]

QCD beta function coefficient.

Parameters:
  • p – pQCD order, 0=LO, 1=NLO, 2=NNLO

  • nf – number of active quark flavors

  • r2 – final momentum scale squared

  • as0 – initial value for a_strong/(2 Pi)

  • r20 – initial momentum scale squared

Returns:

final value for a_strong/(2 Pi)

Return type:

float

Examples

>>> as2pf(0, 3, 8., 0.3, 4.)
0.15497879500975464

evolution module

GPD evolution operator mathcal{E}.

returns:

Numpy array evol[p, k, j] where flavor index j takes values in the evolution basis: (1) – singlet quark (2) – gluon (3) – NS(+) (4) – NS(-) (not tested!)

Notes

GPD models may be defined in different basis and should provide appropriate transformation matrix

Todo

  • Array indices ordering is a bit of a mess

gepard.evolution.lambdaf(gam0) ndarray[source]

Eigenvalues of the LO singlet anomalous dimensions matrix.

Parameters:

gam0 – matrix of LO anomalous dimensions

Returns:

lam[a, k] a in [+, -] and k is MB contour point index

gepard.evolution.projectors(gam0) Tuple[ndarray, ndarray][source]

Projectors on evolution quark-gluon singlet eigenaxes.

Parameters:

gam0 – LO anomalous dimension

Returns:

eigenvalues of LO an. dimm matrix lam[a, k] # Eq. (123)
pr: Projector pr[k, a, i, j] # Eq. (122)

k is MB contour point index a in [+, -] i,j in {Q, G}

Return type:

lam

gepard.evolution.rnlof(m, j) Tuple[ndarray, ndarray, ndarray][source]

Return projected NLO mu-independent P(bet*gam0 - gam1)P.

Parameters:
  • m – instance of the model

  • j – MB contour points (overrides m.jpoints)

Returns:

lam: eigenvalues of LO an. dimm matrix lam[a, k] # Eq. (123) pr: Projector pr[k, a, i, j] # Eq. (122) r1proj: r1proj[a,b] = sum_ij pr[a,i,j] R1[i,j] pr[b,i,j] # Eq. (124) a,b in {+,-}; i,j in {Q, G}

Return type:

Tuple of three arrays, namely

gepard.evolution.rnlonsf(m, j, prty) ndarray[source]

Return NLO mu-independent part of evolution operator.

Parameters:
  • m – instance of the model

  • j – MB contour points (overrides m.jpoints)

  • prty – 1 for NS^{+}, -1 for NS^{-}

Returns:

beta/gama ratio from Eq. (117)

Return type:

r1

gepard.evolution.erfunc(m, lamj, lamk, R) ndarray[source]

Mu-dep. part of NLO evolution operator. Eq. (126).

gepard.evolution.erfunc_nd(m, lamj, lamk, R) ndarray[source]

Mu-dep. part of NLO evolution operator. Eq. (126).

gepard.evolution.cb1(m, Q2, zn, zk, NS: bool = False)[source]

Non-diagonal part of NLO evol op.

Parameters:
  • m – instance of the model

  • Q2 – evolution point

  • zn – non-diagonal evolution Mellin-Barnes integration point (array)

  • zk – COPE Mellin-Barnes integration point (not array! - FIXME)

  • NS – do we want non-singlet?

Returns:

non-diagonal part of evol. op. from Eq. (140)

Return type:

B_jk

Note

It’s multiplied by GAMMA(3/2) GAMMA(K+3) / (2^(K+1) GAMMA(K+5/2)) so it’s ready to be combined with diagonally evolved C_K = 1 + … where in the end everything will be multiplied by

(2^(K+1) GAMMA(K+5/2)) / ( GAMMA(3/2) GAMMA(K+3) )

gepard.evolution.evolop(m, j, Q2: float, process_class: str) ndarray[source]

GPD evolution operator.

Parameters:
  • m – instance of the model

  • j – MB contour points (overrides m.jpoints)

  • Q2 – final evolution momentum squared

  • process_class – DIS, DVCS or DVMP

Returns:

Array corresponding Eq. (121) of Towards DVCS paper. evolop[k, p, i, j] - k is index of point on MB contour, - p is pQCD order (0=LO, 1=NLO) - i, j in [Q, G]

Todo

Argument should not be a process class but GPD vs PDF, or we should avoid it altogether somehow. This serves here only to get the correct choice of evolution scheme (msbar vs csbar).

gepard.evolution.evolopns(m, j, Q2: float, process_class: str) ndarray[source]

GPD evolution operator (NSP case only atm).

Parameters:
  • m – instance of the model

  • j – MB contour points (overrides m.jpoints)

  • Q2 – final evolution momentum squared

  • process_class – DIS, DVCS or DVMP

Returns:

Array corresponding Eq. (116) of Towards DVCS paper. evolopns[k, p] - k is index of point on MB contour, - p is pQCD order (0=LO, 1=NLO)

Notes

Code duplication, should be merged with evolop function

wilson module

Wilson coefficients and evolved Wilson coefficients.

gepard.wilson._fshu(j: ndarray) ndarray[source]

Shuvaev factor.

gepard.wilson.calc_wc(m: Theory, j: ndarray, process_class: str)[source]

Calculate Wilson coeffs.

Parameters:
  • m – instance of the Theory class

  • j – MB contour point(s) (overrides m.jpoints)

  • process_class – ‘DIS’, ‘DVCS’ or ‘DVMP’

Returns:

k in range(npts), p in [LO, NLO], f in [Q,G,NSP]

Return type:

wc[k, p, f]

gepard.wilson.calc_wce(m: Theory, Q2: float, process_class: str)[source]

Calculate evolved Wilson coeffs for given Q2, for all PWs.

Parameters:
  • Q2 – final evolution scale

  • m – instance of the Theory

  • process_class – ‘DIS’, ‘DVCS’ or ‘DVMP’

Returns:

s in range(npwmax), k in range(npts), j in [Q,G,NSP]

Return type:

wce[s,k,j]

gepard.wilson.calc_j2x(m: Theory, x: float, eta: float, Q2: float)[source]

Calculate j2x coeffs, combined with evolution operator.

Parameters:
  • m – instance of the Theory

  • x – long. momentum fraction argument of GPD

  • eta – skewness

  • Q2 – final evolution scale

Returns:

s in range(npwmax); k in range(npts); i,j in [Q,G,NSP]

Return type:

wce[s,k,i,j]

Todo

Implement general (eta != x) j to x transform.

3. Fitting

fitter module

Implementation of fitting algorithms.

class gepard.fitter.Fitter(**kwargs)[source]

Bases: object

Superclass for fitting procedures/algorithms.

__init__(**kwargs) None[source]
class gepard.fitter.MinuitFitter(fitpoints: DataSet, theory: Theory, **kwargs)[source]

Bases: Fitter

Fits using iminuit Python frontend to MINUIT2 C++ library.

Parameters:
__init__(fitpoints: DataSet, theory: Theory, **kwargs) None[source]
covsync()[source]

Synchronize covariance and theory parameter errors with iminuit.

fit()[source]

Perform simple fit.

For better control, use iminuit’s functions.

fix_parameters(*args)[source]

fix_parameters(‘p1’, ‘p2’, …).

release_parameters(*args)[source]

release_parameters(‘p1’, ‘p2’, …).

limit_parameters(dct)[source]

limit_parameters({‘p1’: (lo, hi), …}.

free_parameters() List[str][source]

Return list of names of free fitting parameters.

print_parameters()[source]

Values and errors for free parameters.

4. Utilites

data module

Representation of experimental data.

DataPoint – class for kinematical points, possibly representing measurements DataSet – container for DataPoint instances

dset – dictionary with public datasets

gepard.data.loaddata(resource)[source]

Return dictionary {id : DataSet, …} out of files in resource package.

exception gepard.data.KinematicsError[source]

Bases: Exception

Exception thrown when kinematics is unphysical.

class gepard.data.DataPoint(kindict: dict = None, **kwargs)[source]

Bases: dict

Kinematical point. May correspond to actual measurement.

Parameters:
  • xB (float) – Bjorken x_B

  • t (float) – momentum transfer to target squared

  • Q2 (float) – Q^2

  • phi (float) – azimuthal angle

  • FTn (int) – harmonic of azimuthal angle, e.g. -1 for sin(phi)

  • observable (str) – name of the measured observable

  • val (float) – measurement value

  • err (float) – total uncertainty of val

  • units (dict) – pysical units of variables

  • kindict (dict) – for old, alternative passing of kinematics values, like g.DataPoint({'xB': 0.1, 't': -0.2, 'Q2': 4.0})

There are other args as well.

Examples

>>> pt = g.DataPoint(xB=0.1, t=-0.2, Q2=4.0)

Todo

kindict is temporarily kept for backward compatibility. You should not rely on it.

__init__(kindict: dict = None, **kwargs)[source]
prepare()[source]

Pre-calculate some functions of kinematics.

copy()[source]

Copy the DataPoint object.

update_from_grid(gridline, dataset)[source]

Take data gridline, and update point atributes.

Parameters:
  • gridline (list) – list constructed from one row of data grid in data file.

  • dataset (DataSet) – container that is to contain this DataPoint

Note

It is assumed that data gridline is of the form:

x1  x2 ....   y1  y1stat y1syst

where y1syst need not be present, or can have one or two values, syst+ and syst- (Further elements are ignored.)

Todo

This passing of higher-level DataSet as argument sounds wrong! (See comment about aquisition in class docstring.)

to_conventions()[source]

Transform datapoint into gepard’s conventions.

from_conventions()[source]

Transform stuff from Approach’s conventions into original data’s.

orig_conventions(val)[source]

Like from_conventions, but for the prediction val.

class gepard.data.DataSet(datapoints=None, datafile=None)[source]

Bases: list

A container for DataPoint instances.

Parameters:
  • datapoints (list) – list of DataPoint instances

  • datafile (str) – data file to be parsed, represented as string

Either take explicit list of DataPoints or get them by parsing datafile.

Information that is common to all data of a given dataset (i.e. which is contained in a preamble of datafile is accessible via attributes:

observable

name of observable measured

Type:

str

collaboration

name of experimenatal collaboration

Type:

str

units

pysical units of variables

Type:

dict

newunits

internal pysical units of variables

Type:

dict

There are other attributes as well.

__init__(datapoints=None, datafile=None)[source]
parse(dataFile)[source]

Parse dataresource string and return tuple (preamble, data).

Parameters:

dataFile (str) – datafile to be parsed, as string

Return type:

(preamble, data) (tuple)

preamble is dictionary obtained by converting datafile preamble items into dictionary items like this:

y1 = AC from datafile goes into   {'y1' : 'AC', ...}

data is actual numerical grid of experimental data converted into list of lists.

df()[source]

Return pandas DataFrame of a DataSet.

gepard.data._str2num(s: str) int | float[source]

Convert string to number, taking care if it should be int or float.

Parameters:

s – string to be converted

Returns:

number Union[int, float]

http://mail.python.org/pipermail/tutor/2003-November/026136.html

gepard.data._complete_xBWQ2(kin)[source]

Make trio {xB, W, Q2} complete if two of them are given in ‘kin’.

gepard.data._complete_tmt(kin)[source]

Make duo {t, tm} complete if one of them is given in ‘kin’.

gepard.data._fill_kinematics(kin, old={})[source]

Update kinematics in place.

Complete set of kinematical variables is {xB, t, Q2, W, s, xi, tm, phi}. Using standard identities, missing values are calculated, if possible, first solely from values given in ‘kin’, and then, second, using values in ‘old’, if provided.

Parameters:
gepard.data.select(dataset: DataSet, criteria: List[str] = [], logic='AND')[source]

Filter point of dataset satisfying criteria.

Parameters:
  • dataset – DataSet or any collection of DataPoints

  • criteria – list of strings describing selection criteria

  • logic – ‘AND’ or ‘OR’: how to logically combine criteria

Returns:

DataSet of DataPoints satisfying criterion.

Example

criteria=[‘xB > 0.1’, ‘observable == ALU’]

gepard.data.list_data(ids)[source]

List basic info about datasets specified by id numbers.

gepard.data.describe_data(pts)[source]

Print observables and where they come from.

gepard.data.FTF(data, testdata=None, cosmax=None, sinmax=None, inverse=False)[source]

Fourier series fit to data (takes and returns pandas dataframe).

cosmax, sinmax - index of highest harmonics inverse = True - data values ARE harmonics, calculate function(phi)

gepard.data.FTFMC(data, nsamples=100, cosmax=None, sinmax=None, inverse=False)[source]

Fourier series fit to data with MC error propagation (takes and returns pandas dataframe). :param cosmax: index of highest cosine harmonic :param sinmax: index of highest sine harmonics :param inverse: if True data values ARE harmonics, calculate function(phi) (default: False)

gepard.data.cvsets(df_in, nfolds=3, shuffle=True)[source]

Return n-fold cross-validation sets [(train1, valid1), (train2, valid2), …].

Copies data.

gepard.data.FTanalyse(bins, HMAX=2, nf=4, Nrep=1)[source]

Determine the number of harmonics present in bins according to n-fold cross-validation.

Parameters:
  • bins – dictionary with bins

  • HMAX – highest harmonic used in searches

  • NS – number of replicas for MC error propagation

  • Nrep – number of CV repetitions (probably wrong to make > 1)

  • nf – number of folds for cross-validation

special module

Special functions of complex arguments.

poch: Pochhammer symbol dpsi: polygamma function Sn: harmonic sums

gepard.special.parity(n: int) int[source]

Parity of n, i. e. (-1)^n.

gepard.special.pochhammer(z: complex | ndarray, m: int) complex | ndarray[source]

Pochhammer symbol.

Parameters:
  • z – complex argument

  • m – integer index

Returns:

pochhammer(z,m)

Return type:

complex

gepard.special.poch(z: complex | ndarray, m: int) complex | ndarray

Pochhammer symbol.

Parameters:
  • z – complex argument

  • m – integer index

Returns:

pochhammer(z,m)

Return type:

complex

gepard.special.dpsi_one(z: complex, m: int) complex[source]

Polygamma - m’th derivative of Euler gamma at z.

gepard.special.S1(z: complex | ndarray) complex | ndarray[source]

Harmonic sum S_1.

gepard.special.S2(z: complex | ndarray) complex | ndarray[source]

Harmonic sum S_2.

gepard.special.S3(z: complex | ndarray) complex | ndarray[source]

Harmonic sum S_3.

gepard.special.S4(z: complex | ndarray) complex | ndarray[source]

Harmonic sum S_4.

gepard.special.S2_prime(z: complex | ndarray, prty: int) complex | ndarray[source]

Curci et al Eq. (5.25).

gepard.special.S3_prime(z: complex | ndarray, prty: int) complex | ndarray[source]

Curci et al Eq. (5.25).

gepard.special.delS2(z: complex | ndarray) complex | ndarray[source]

Return Harmonic sum S_2 difference.

Parameters:

z – complex argument

Returns:

delS2((z+1)/2) From Eq. (4.13) of 1310.5394. Note halving of the argument.

gepard.special.deldelS2(j: complex | ndarray, k: int) complex | ndarray[source]

Return diference of harmonic sum S_2 differences.

Parameters:
  • j – complex argument

  • k – integer index

Returns:

Equal to delS2((j+1)/2, (k+1)/2) From Eq. (4.38) of 1310.5394. Note halving of the argument.

gepard.special.Sm1(z: complex | ndarray, k: int) complex | ndarray[source]

Aux fun. FIXME: not tested.

gepard.special.MellinF2(n: complex | ndarray) complex | ndarray[source]

Return Mellin transform i.e. x^(N-1) moment of Li2(x)/(1+x).

Parameters:

n – complex argument

Returns:

According to Eq. (33) in Bluemlein and Kurth, hep-ph/9708388

gepard.special.SB3(j: complex | ndarray) complex | ndarray[source]

Eq. (4.44e) of arXiv:1310.5394.

gepard.special.S2_tilde(n: complex | ndarray, prty: int) complex | ndarray[source]

Eq. (30) of Bluemlein and Kurth, hep-ph/9708388.

mellin module

Mellin-Barnes integrations.

Todo

Various MB routines should be integrated into one.

class gepard.mellin.MellinBarnes(**kwargs)[source]

Bases: object

Base class for models built by Mellin-Barnes integration.

__init__(**kwargs) None[source]
_mellin_barnes_integral(xi, wce, gpd)[source]

Return convolution of evolved Wilson coefs and GPDs.

_mellin_barnes_integral_HE(xi, wce, H, E)[source]

Return convolution of evolved Wilson coefs and GPDs H and E.

_dis_mellin_barnes_integral(xi, wce, pdf)[source]

Return convolution of evolved Wilson coefs and PDFs.

_j2x_mellin_barnes_integral(x, eta, wce, gpd)[source]

Return convolution of j->x coef, evolution operator and GPD.

_j2x_mellin_barnes_integral_E(x, eta, wce, gpd)[source]

Return convolution of j->x coef, evolution operator and GPD E.

Notes

Not tested at all!

quadrature module

Quadrature formulas, abscissas and weights.

Also numerical derivative routine.

gepard.quadrature.mellin_barnes(c, phi, accuracy: int = 3) Tuple[ndarray, ndarray][source]

Construct basic MB array.

Parameters:
  • c – point of crossing the real axis

  • phi – angle with the real axis

  • accuracy – 2*accuracy points per interval

Returns:

(complex coordinates of MB contour, integration weighths)

gepard.quadrature.nd_mellin_barnes(accuracy: int = 3) Tuple[ndarray, ndarray][source]

Construct MB array for non-diagonal evolution operator.

Parameters:

accuracy – 2*accuracy points per interval

Returns:

(complex coordinates of MB contour, integration weighths)

Note

Code should be merged with standard mellin_barnes.

gepard.quadrature.quadSciPy81(func, a, b, args=())[source]

Compute a definite integral using 81-order Gaussian quadrature.

gepard.quadrature.quadSciPy35(func, a, b, args=())[source]

Compute a definite integral using 35-order Gaussian quadrature.

gepard.quadrature.quadSciPy18transposed(func, a, b, args=())[source]

Compute a definite integral using 18-order Gaussian quadrature.

gepard.quadrature.quadSciPy10(func, a, b, args=())[source]

Compute a definite integral using tenth-order Gaussian quadrature.

gepard.quadrature.quadSciPy10transposed(func, a, b, args=())[source]

Compute a definite integral using fifth-order Gaussian quadrature.

gepard.quadrature.quad10(func, a, b, args=())[source]

Compute a definite integral using tenth-order Gaussian quadrature.

gepard.quadrature.quadSciPy5(func, a, b, args=())[source]

Compute a definite integral using fifth-order Gaussian quadrature.

gepard.quadrature.quadSciPy5transposed(func, a, b, args=())[source]

Compute a definite integral using fifth-order Gaussian quadrature.

gepard.quadrature.quadSciPy4(func, a, b, args=())[source]

Compute a definite integral using fifth-order Gaussian quadrature.

gepard.quadrature.PVquadrature(func, a, b, args=())

Compute a definite integral using 18-order Gaussian quadrature.

gepard.quadrature.Hquadrature(func, a, b, args=())

Compute a definite integral using fifth-order Gaussian quadrature.

gepard.quadrature.tquadrature(func, a, b, args=())

Compute a definite integral using fifth-order Gaussian quadrature.

gepard.quadrature.bquadrature(func, a, b, args=())

Compute a definite integral using tenth-order Gaussian quadrature.

gepard.quadrature.rthtquadrature(func, a, b, args=())

Compute a definite integral using tenth-order Gaussian quadrature.

gepard.quadrature.deriv(func, x, h, nevals)[source]

Return derivative func’(x) using Ridders-Neville algorithm.

constants module

Constants, with their numerical values, not to be changed by user.

Also, fixed dictionaries.