libepa Programmer reference

Contents

Index

Introduction

A charged particle at rest is surrounded by electric field which can be represented as a bunch of virtual photons with zero time component of the photons momenta. If the particle is moving with the speed close to the speed of light, momenta of the photons of its field acquire longitudinal component, so that the photons virtuality Q 2 - q 2 ω 2 , where q is the photon 4-momentum, ω is the photon energy. Under the equivalent photon approximation (EPA), these photons are treated as real and are assumed to be distributed according to a known spectrum.

In many experiments in high energy physics, collisions of two ultrarelativistic charged particles are studied. Particles often do not collide head-on, but rather pass at some distance from each other and collide with their electromagnetic fields. If both particles survive in the collision, such collision is called ultraperipheral. Ultraperipheral collisions (UPC) can be considered as photon-photon collisions with the photons meeting the requirements of the equivalent photon approximation.

libepa is a library for calculation of cross sections of ultraperipheral collisions under the equivalent photon approximation. The cross sections are calculated as integrals over the convolution of the photon EPA spectra with the photon-photon cross section provided by the user. By default, the integration is performed with the help of GNU Scientific Library GSL. The user can supply their own integration routines.

This document is the programmer reference. See the paper 2311.01353 for the description of ultraperipheral collisions and examples.

Assumptions

Notation

Overview

libepa is build around anonymous (lambda) functions and their composition. Consider calculation of the cross section for the production of a pair of muons with invariant mass 100 GeV in ultraperipheral collisions of two protons with the energy 13 TeV:

    const double muon_mass = 105.6583745e-3; // GeV
    const double collision_energy = 13e3; // GeV
    const double invariant_mass = 100; // GeV
    auto luminosity = pp_luminosity(collision_energy);
    auto photons_to_muons = photons_to_fermions(muon_mass);
    auto cross_section = xsection(photons_to_muons, luminosity);
    double result = cross_section(invariant_mass); // barn

Here pp_luminosity returns an anonymous function with the signature double (double rs). The function calculates photon-photon luminosity in the proton-proton collision with the energy 13 TeV for the photons invariant mass rs ( s , where s is the Mandelstam variable). The collision energy is captured by the function. photons_to_fermions returns a function with the same signature, but it calculates the cross section for the production of a pair of muons in barns taking the square of the muons invariant mass as the argument. These functions are passed to xsection which returns a function with the same signature again, computing the wanted cross section (xsection's function simply computes the product of the luminosity and the photon-photon cross section; other cross section functions are more complicated).

Many functions calculate numerical integrals. For example, pp_luminosity calculates the integral of the product of photon spectra of the protons. Integration routines are decoupled from EPA calculations and are passed in forms of integrators — functions with the signature double (const std::function<double (double)>& f, double a, double b). Integrator computes the integral of f from a to b. To denote (semi-)infinite intervals, a can be set to -infinity, and b can be set to infinity. libepa uses GNU Scientific Library (GSL) for integration and provides convenient functions producing integrators based on GSL algorithms, see qag_integrator and cquad_integrator. GSL routines use estimations of absolute and relative errors to control the accuracy of calculation. Calculation continues until either the absolute or the relative change in the estimation of the integral value between succeeding iterations is less than the parameters epsabs or epsrel, see the GSL documentation on numerical integration . In libepa, these are referred to as relative_error and absolute_error. They are captured by integrators. Their default values used by the libepa wrappers are controlled by the variables default_absolute_error and default_relative_error. Note that to correctly utilize absolute_error, one has to know the magnitude of the integral in advance. This is rather inconvenient, so by default this value is not used (default_absolute_error is set to 0). The default value for default_relative_error is 10-3 meaning that the default accuracy of numerical integration should not be much worse than 0.1%.

Many functions calculate multiple integrals. For example, calculation of the differential fiducial cross section taking into account the probability of strong interaction between the colliding particles (the survival probability also known as the survival factor) provided by the function xsection_fid_b requires calculation of quintuple integrals of the phase space parameters, assuming that the photon spectra are known analytically. Add an extra level of integration if the integrated cross section is desired. In this case, to avoid accumulation of numerical errors, each inner integral has to be calculated with significantly greater accuracy than the outer integral. The recommended approach in this case is to reduce the relative_error parameter for each inner integral by an order of magnitude. In libepa functions, this reduction is controlled by the default_error_step parameter. Default values of these parameters sometimes result in computations taking unbearably long time or in numerical errors spoiling the calculation. Increasing default_relative_error and default_error_step might help in this case.

In libepa language, integration level of an integral is the nesting number of the integral in the calculation. For example, when computing dx dy f ( x , y ) , the integral with respect to dx is of integration level 0, and the integral with respect to dy is of integration level 1. Integrator generator is a function producing integrators for different integration levels. The difference between these integrators when the default GSL-based algorithms are used is in their captured relative_error parameter which is calculated as follows:
relative_error = default_relative_error × default_error_stepintegration level.
Available integrator generators based on GSL routines are qag_integrator and cquad_integrator

GSL integration routines require workspace which is allocated while constructing the integrator. The workspace is wrapped in std::shared_ptr, so copying the integrator is cheap, and the workspace is deallocated along with the integrator. Note, however, that two copies of such integrator cannot be used in parallel. To add a compile-time constraint on that, copying of the workspace will probably be forbidden in the future.

API Reference

EPA

    #include <epa/epa.hpp>

All libepa symbols are defined in namespace epa. All values are in units of powers of GeV except cross sections which are in barns. Calculations assume ħ = c = 1 where ħ is the Planck constant, and c is the speed of light.

Version

Current libepa can be determined through cpp macros EPA_VERSION_MAJOR, EPA_VERSION_MINOR and EPA_VERSION_PATCH which expand to integer literals.

Error handling

When encountering an error, libepa throws an exception. For integration errors, these are exceptions of class gsl::Error. Standard exceptions might be thrown by the functions from the standard library. To help with debugging, the variable bool print_backtrace can be set to true. In this case libepa will print to the standard error output the backtrace of the calculation.

Constants
The following constants are used in the functions (some of the constants are only used to define other constants):
Type Constant Value Units Comment
double infinity INF Double precision positive infinity
double pi M_PI The value of π
double alpha 7.2973525693e-3 Fine structure constant α *
double planck 0.1973269804 GeV fm Planck constant ħ *
double barn ħ 2 × 10 -2 GeV2 barn Coefficient converting GeV−2 to barns
double fm 1 / ( 10 barn ) (GeV fm)−1 Coefficient converting fm to GeV−1
double electron_charge 1.602176634e-19 Cl Electron charge in SI units (e) *
double light_speed 299792458 m/s Speed of light in SI units (c) *
double amu ( 1.66053906660 · 10 -27 · c 2 / e ) × 10 -9 GeV Atomic mass unit in GeV *
Constants denoted with * were taken from E. Tiesinga, P. J. Mohr, D. B. Newell, B. N. Taylor, CODATA recommended values of the fundamental physical constants: 2018, Reviews of Modern Physics 93, 025010 (2021).
Integrators
Integrator
    typedef std::function<
      double (
        const std::function<double (double)>&, // f(x)
        double, // a
        double  // b
      )
    > Integrator;
  
Integrator is a function taking as a parameter a function of one argument (f(x)) and computing an integral of this function on a given interval [a, b].
Parameters governing GSL integration routines:
Type Variable Default value
double default_absolute_error 0
double default_relative_error 1e-3
double default_error_step 1e-1
size_t default_integration_limit 1000
size_t default_cquad_integration_limit 100
gsl::integration::QAGMethod default_integration_method gsl::integration::GAUSS41
Integrator default_integrator qag_integrator

GSL integration routines perform integration by dividing the integration region into successively smaller intervals and estimating the integral value across these intervals using different approaches. At each iteration, the integral value is computed and compared to that of the previous iteration. Iterations stop when the absolute or the relative difference between these values is smaller than the provided parameters epsabs and epsrel respectively. In the following routines epsabs is referred to as absolute_error and epsrel is referred to as relative_error. Default values of these parameters are provided by default_absolute_error and default_relative_error. Correct application of absolute_error requires estimation of the order of the magnitude of the integral, so its usage is discouraged and its default value is set to zero.

Integration of inner integrals requires higher accuracy than that of outer integrals. For this reason integrator generators reduce relative_error by multiplying it by error_step at each integration level. default_error_step provides the default value of error_step.

GSL integration routines use preallocated workspace. default_integration_limit and default_cquad_integration_limit define the default size of the workspace for the QAG and CQUAD algorithms. When the integration runs out of workspace, an exception of type gsl::Error is thrown.

GSL QAG integration routine divides the integration region into subintervals by calculating the function value at n points, where n is defined by the integration_method parameter. Supported values are those defined in enum gsl::integration::QAGMethod: GAUSS15, GAUSS21, GAUSS31, GAUSS41, GAUSS51, GAUSS61. The default value is provided by default_integration_method.

default_integrator defines the default integrator generator.

qag_integrator
    Integrator qag_integrator(
      double absolute_error = default_absolute_error,
      double relative_error = default_relative_error,
      gsl::integration::QAGMethod integration_method = default_integration_method,
      std::shared_ptr<gsl::integration::QAGWorkspace> workspace = nullptr
    );
    
Returns an integrator based on quadratic adaptive integration routines for general (used-defined) integrand implemented in GSL (QAG). See above for the description of the parameters. If workspace is nullptr, a new workspace sufficient for default_integration_limit is allocated.
    struct qag_integrator_keys {
      double absolute_error = default_absolute_error;
      double relative_error = default_relative_error;
      gsl::integration::QAGMethod method = default_integration_method;
      std::shared_ptr<gsl::integration::QAGWorkspace> workspace;
    };

    Integrator qag_integrator(const qag_integrator_keys);
    
A keyword variant of qag_integrator. Use as follows:
    auto integrate = qag_integrator({ .method = gsl::integration::GAUSS21 });
    
    Integrator qag_integrator(unsigned level)
    
Returns a QAG integrator with relative_error = default_relative_error × default_error_steplevel. Other parameters are initialized with defaults.
qag_integrator_generator
    std::function<Integrator (unsigned /* integration_level */)> qag_integrator_generator(
      double relative_error = default_relative_error,
      double error_step     = default_error_step
    );
  
Returns a generator of QAG integrators with relative_error = relative_error × error_steplevel. Other parameters are initialized with defaults.
cquad_integrator
    Integrator cquad_integrator(
      double absolute_error = default_absolute_error,
      double relative_error = default_relative_error,
      std::shared_ptr<gsl::integration::CQuadWorkspace> workspace = nullptr
    );
    
Returns an integrator based on CQUAD doubly-adaptive algorithm implemented in GSL. See above for the description of the parameters. If workspace is nullptr, a new workspace sufficient for default_cquad_integration_limit is allocated.
    struct cquad_integrator_keys {
      double absolute_error = default_absolute_error;
      double relative_error = default_relative_error;
      std::shared_ptr<gsl::integration::CQuadWorkspace> workspace;
    };

    Integrator cquad_integrator(const cquad_integrator_keys&);
    
A keyword variant of cquad_integrator. Use as follows:
    auto integrate = cquad_integrator({ .relative_error = 1e-2 });
    
    Integrator cquad_integrator(unsigned level = 0);
    
Returns a CQUAD integrator with relative_error = relative_error × error_steplevel. Other parameters are initialized with defaults.
Form factors
FormFactor
    typedef std::function<double (double /* Q^2 */)> FormFactor;
  
Electromagnetic form factor is used to describe the structure of elementary particle as seen in electromagnetic interactions. For the theory, refer to S. Pacetti, R. B. Ferroli, E. Tomasi-Gustaffson, Proton electromagnetic form factors: Basic notions, present achievements and future perspectives, Phys. Rep. 550-1 (2015). Form factor of a pointlike particle equals 1. Form factor of a finite size particle is a function mapping the photon 3-momentum squared Q 2 to [ 0 ; 1 ] .
form_factor_monopole
    FormFactor form_factor_monopole(double lambda2);
  
Monopole form factor: F 1 ( Q 2 ) = 1 1 + Q 2 Λ 2 . lambda2 is Λ 2 in GeV2.
form_factor_dipole
    FormFactor form_factor_dipole(double lambda2);
  
Dipole form factor: F 2 ( Q 2 ) = 1 ( 1 + Q 2 Λ 2 ) 2 . lambda2 is Λ 2 in GeV2.
Equivalent photon spectra
Spectrum
    typedef std::function<double (double /* ω */)> Spectrum;
  
Equivalent photon spectrum describes the distribution of photons moving along with the particle with respect to the photons energy ω.
spectrum
    Spectrum spectrum(
      unsigned Z,
      double gamma,
      FormFactor, // F
      Integrator = default_integrator(0)
    );
  
Equivalent photon spectrum defined as an integral of the particle form factor: n ( ω ) = 2 Z 2 α π ω 0 [ F ( q 2 + ( ω / γ ) 2 ) q 2 + ( ω / γ ) 2 ] 2 q 3 d q .
Note that the equivalent photon spectrum of a pointlike particle is divergent.
spectrum_monopole
    Spectrum spectrum_monopole(unsigned Z, double gamma, double lambda2);
  
Equivalent photon spectrum for the monopole form factor: n 1 ( ω ) = Z 2 α π ω [ ( 2 a + 1 ) ln ( 1 + 1 a ) - 2 ] , where a = ( ω / Λ γ ) 2 .
spectrum_dipole
    Spectrum spectrum_dipole(unsigned Z, double gamma, double lambda2);
  
Equivalent photon spectrum for the dipole form factor: n 2 ( ω ) = Z 2 α π ω [ ( 4 a + 1 ) ln ( 1 + 1 a ) - 24 a 2 + 42 a + 17 6 ( a + 1 ) 2 ], where a = ( ω / Λ γ ) 2 .
Spectrum_b
    typedef std::function<double (double /* b */, double /* ω */)> Spectrum_b;
  
Distribution of photons moving along with the particle at distance b in the transverse plane with respect to photons energy ω. It is also referred to as the equivalent photon spectrum throughout this documentation.
spectrum_b
    Spectrum_b spectrum_b(
      unsigned Z,
      double gamma,
      FormFactor, // F
      Integrator = default_integrator(0)
    );
  
Equivalent photon spectrum defined as an integral of the particle form factor: n ( b , ω ) = Z 2 α π 2 ω [ 0 F ( q 2 + ( ω / γ ) 2 ) q 2 + ( ω / γ ) 2 J 1 ( b q ) q 2 d q ] 2 . Note that the form factor is weighted with an oscillating function (the Bessel function J 1) and integrated over semi-infinite interval. For small values of b this integration is very difficult to perform numerically. Take care with the default quadrature integrator.
spectrum_b_point
    Spectrum_b spectrum_b_point(unsigned Z, double gamma);
  
Equivalent photon spectrum for pointlike particle: n 0 ( b , ω ) = Z 2 α ω π 2 γ 2 K 1 2 ( b ω γ ) .
spectrum_b_monopole
    Spectrum_b spectrum_b_monopole(unsigned Z, double gamma, double lambda2);
  
Equivalent photon spectrum for the monopole form factor: n 1 ( b , ω ) = Z 2 α π 2 ω [ ω γ K 1 ( b ω γ ) - Λ 2 + ( ω γ ) 2 K 1 ( b Λ 2 + ( ω γ ) 2 ) ] 2 .
spectrum_b_dipole
    Spectrum_b spectrum_b_dipole(unsigned Z, double gamma, double lambda2);
  
Equivalent photon spectrum for the dipole form factor: n 2 ( b , ω ) = Z 2 α π 2 ω [ ω γ K 1 ( b ω γ ) - Λ 2 + ( ω γ ) 2 K 1 ( b Λ 2 + ( ω γ ) 2 ) - b Λ 2 2 K 0 ( b Λ 2 + ω γ 2 ) ] 2 .
Luminosities
Luminosity
    typedef std::function<double (double /*  s  */)> Luminosity;
  
Luminosity L is the number of photons coming from both of the colliding particles during the collision. Luminosity is the luminosity differentiated with respect to the invariant mass s of the photon-photon system ( d L / d s ) .
Luminosity_y
    typedef std::function<double (double /* s */, double /* y */)> Luminosity_y;
  
Luminosity_y is the luminosity differentiated with respect to the invariant mass s and rapidity y of the photon-photon system ( d L / d y d s ) .
Luminosity_fid
    typedef std::function<double (double /* s */, double /* ymin */, double /* ymax */)> Luminosity_fid;
  
Luminosity_fid is used to calculate fiducial cross section. It is the luminosity differentiated with respect to invariant mass s of the photon-photon system ( d L / d s ) and integrated with respect to rapidity of the photon-photon system in the range [ ymin, ymax ] instead of [−∞; ∞].
luminosity
    Luminosity luminosity(Spectrum nA, Spectrum nB, Integrator = default_integrator(0));
    
Photon-photon luminosity in ultraperipheral collision of particles A and B differentiated with respect to the invariant mass s of the photon-photon system: d L A B d s ( s ) = s 2 - n A ( s 2 e y ) n B ( s 2 e - y ) d y.
    Luminosity luminosity(Spectrum n, Integrator = default_integrator(0));
    
Convenience function for the collision of identical particles. It takes into account the symmetry of the integrand, improving on convergence and calculation time: d L d s ( s ) = s 0 n ( s 2 e y ) n ( s 2 e - y ) d y.
luminosity_y
    Luminosity_y luminosity_y(Spectrum nA, Spectrum nB);
    
Photon-photon luminosity in ultraperipheral collision of particles A and B differentiated with respect to both invariant mass of the photons s and rapidity of the system y: d L A B d s d y ( s , y ) = s 2 n A ( s 2 e y ) n B ( s 2 e - y )
    Luminosity_y luminosity_y(Spectrum);
    
Convenience function for the collision of identical particles.
luminosity_fid
    Luminosity_fid luminosity_fid(Spectrum nA, Spectrum nB, Integrator = default_integrator(0));
    
Photon-photon luminosity in ultraperipheral collisions of particles A and B differentiated with respect to the invariant mass of the photons s with parameterized bounds for the rapidity integral: d L A B d s ( s , y min , y max ) = s 2 y min y max n A ( s 2 e y ) n B ( s 2 e - y ) d y.
    Luminosity_fid luminosity_fid(Spectrum, Integrator = default_integrator(0));
    
Convenience function for the collision of identical particles. When y min = - y max , it takes into account the symmetry of the integrand, improving on convergence and calculation time.
Polarization
    struct Polarization {
      double parallel;
      double perpendicular;
    };
  
Describes the relative polarization of colliding photons.
Luminosity_b, Luminosity_y_b, Luminosity_fid_b
    typedef std::function<double (double /* s */, Polarization /* p */)> Luminosity_b

    typedef std::function<double (double /* s */, double /* y */, Polarization /* p */)> Luminosity_y_b;

    typedef std::function<double (double /* s */, Polarization /* p */, double /* ymin */, double /* ymax */)> Luminosity_fid_b;
  
When non-electromagnetic interactions between the colliding particles are taken into account, part of the collision plane is excluded, and in this case polarization of the colliding photons is not summed over. Luminosity functions taking into account non-electromagnetic interactions compute two values, one for the photons with parallel polarization, another for the photons with perpendicular polarization. The p parameter provides the weights to multiply these values. Pass p.parallel = 1, p.perpendicular = 0 to take only the parallel photon polarizations into account, p.parallel = 1, p.perpendicular = 1 to take both polarizations into account. Use the values calculated by XSection_b functions to calculate the cross section.
luminosity_b
      Luminosity_b luminosity_b(
        Spectrum_b nA,
        Spectrum_b nB,
        std::function<double (double /* b */)> upc_probability,
        const std::function<Integrator (unsigned)>& = default_integrator,
        unsigned integration_level = 0
      );
    
Photon-photon luminosity in ultraperipheral collisions of particles A and B differentiated with respect to the invariant mass s of the photon-photon system: d L A B d s ( s , { p , p } ) = p d L d s ( s ) + p d L d s ( s ) , where d L d s ( s ) = s 2 - d y d 2 b 1 d 2 b 2 n A ( b 1 , s 2 e y ) n B ( b 2 , s 2 e - y ) P ( b ) cos 2 φ, d L d s ( s ) = s 2 - d y d 2 b 1 d 2 b 2 n A ( b 1 , s 2 e y ) n B ( b 2 , s 2 e - y ) P ( b ) sin 2 φ, d L / d s is the luminosity of photons with parallel polarizations, d L / d s is the luminosity of photons with perpendicular polarizations, P ( b ) is upc_probability, the probability for the particles to survive in an ultraperipheral collision with the impact parameter b (the probability to avoid electromagnetic interactions), b = | b 1 - b 2 |. integration_level is the integration level of the outermost integral to be calculated by the returned function.
      Luminosity_b luminosity_b(
        Spectrum_b n,
        std::function<double (double /* b */)> upc_probability,
        const std::function<Integrator (unsigned)>& = default_integrator,
        unsigned integration_level = 0
      );
    
Convenience function for the collision of identical particles. It takes into account the symmetry of the integrands, improving on convergence and calculation time: d L d s ( s ) = s 0 d y d 2 b 1 d 2 b 2 n ( b 1 , s 2 e y ) n ( b 2 , s 2 e - y ) P ( b ) cos 2 φ, d L d s ( s ) = s 0 d y d 2 b 1 d 2 b 2 n ( b 1 , s 2 e y ) n ( b 2 , s 2 e - y ) P ( b ) sin 2 φ,
luminosity_y_b
      Luminosity_y_b luminosity_y_b(
        Spectrum_b nA,
        Spectrum_b nB,
        std::function<double (double /* b */)> upc_probability,
        const std::function<Integrator (unsigned)>& = default_integrator,
        unsigned integration_level = 0
      );
    
Photon-photon luminosity in ultraperipheral collisions of particles A and B differentiated with respect to the photons invariant mass s and rapidity y: d L A B d s d y ( s , y , { p , p } ) = p d L d s d y ( s , y ) + p d L d s d y ( s , y ) where d L d s d y ( s , y ) = s 2 d 2 b 1 d 2 b 2 n A ( b 1 , s 2 e y ) n B ( b 2 , s 2 e - y ) P ( b ) cos 2 φ, d L d s d y ( s , y ) = s 2 d 2 b 1 d 2 b 2 n A ( b 1 , s 2 e y ) n B ( b 2 , s 2 e - y ) P ( b ) sin 2 φ, d L / d s is the luminosity of photons with parallel polarizations, d L / d s is the luminosity of photons with perpendicular polarizations, P ( b ) is upc_probability, the probability for the particles to survive in an ultraperipheral collision with the impact parameter b (the probability to avoid electromagnetic interactions), b = | b 1 - b 2 |. integration_level is the integration level of the outermost integral to be calculated by the returned function.
      Luminosity_y_b luminosity_y_b(
        Spectrum_b,
        std::function<double (double /* b */)> upc_probability,
        const std::function<Integrator (unsigned)>& = default_integrator,
        unsigned integration_level = 0
      );
    
Convenience function for the collision of identical particles.
luminosity_fid_b
      Luminosity_fid_b luminosity_fid_b(
        Spectrum_b nA,
        Spectrum_b nB,
        std::function<double (double /* b */)> upc_probability,
        const std::function<Integrator (unsigned)>& = default_integrator,
        unsigned integration_level = 0
      );
    
Photon-photon luminosity in ultraperipheral collisions of particles A and B differentiated with respect to the photons invariant mass s with parameterized bounds on the rapidity integral: d L A B d s ( s , { p , p } , y min , y max ) = p d L d s ( s ) + p d L d s ( s ) , where d L d s ( s ) = s 2 y min y max d y d 2 b 1 d 2 b 2 n A ( b 1 , s 2 e y ) n B ( b 2 , s 2 e - y ) P ( b ) cos 2 φ, d L d s ( s ) = s 2 y min y max d y d 2 b 1 d 2 b 2 n A ( b 1 , s 2 e y ) n B ( b 2 , s 2 e - y ) P ( b ) sin 2 φ, d L / d s is the luminosity of photons with parallel polarizations, d L / d s is the luminosity of photons with perpendicular polarizations, P ( b ) is upc_probability, the probability for the particles to survive in an ultraperipheral collision with the impact parameter b (the probability to avoid electromagnetic interactions), b = | b 1 - b 2 |. integration_level is the integration level of the outermost integral to be calculated by the returned function.
      Luminosity_fid_b luminosity_fid_b(
        Spectrum_b,
        std::function<double (double /* b */)> upc_probability,
        const std::function<Integrator (unsigned)>& = default_integrator,
        unsigned integration_level = 0
      );
    
Convenience function for the collision of identical particles. When y min = - y max , it takes into account the symmetry of the integrand, improving on convergence and calculation time.
Cross sections of ultraperipheral collisions
XSection
    typedef std::function<double (double /* s */)> XSection;
  
Cross section for the invariant mass s (s is the Mandelstam parameter). Some functions use this type to denote cross section differentiated with respect to s: d σ / d s .
XSection_b
    typedef std::function<Polarization (double /* s */)> XSection_b;
  
Photon fusion cross section taking into account photons polarization. Returns two values: cross sections for the parallel and perpendicular photons polarizations respectively.
xsection
    XSection xsection(XSection /* photon_photon */, Luminosity  /* L */);
  
Differential cross section for the production of system X in ultraperipheral collisions with photon-photon luminosity L with non-electromagnetic interactions neglected: d σ A B A B X d s ( s ) = σ γ γ X ( s ) · d L d s ( s ) , where σ γ γ X is photon_photon, the cross section for the fusion of two (real) photons.
xsection_fid
    XSection xsection_fid(
      std::function<double (double /* s */, double /* pT */)> xsection_pT,
      Luminosity_fid,
      double mass,   // of the produced particle
      double pT_min,
      double eta_max,
      double ω1_min, // use 0 for no limit
      double ω1_max, // use infinity for no limit
      double ω2_min, // use 0 for no limit
      double ω2_max, // use infinity for no limit
      Integrator = default_integrator(0)
    )
    
Differential fiducial cross section for the production of a pair of particles χ + χ - in an ultraperipheral collision of particles A and B ( A B A B χ + χ - ) with non-electromagnetic interactions neglected and with the constraints on phase space p T ± > p ^ T , | η ± | < η ^ , ω 1, min < ω 1 < ω 1, max , ω 2, min < ω 2 < ω 2, max , where p ^ T is pT_min, η ^ is eta_max, p T ± and η ± are the transverse momentum and pseudorapidity of χ ± , ω 1 and ω 2 are the photons energies: d σ A B A B χ + χ - fid. d s = max ( p ^ T p ~ T ) s 2 1 - 4 m χ 2 s d p T d σ γ γ χ + χ - d p T ( s , p T ) d L A B d s ( s , max ( - y ^ , y ~ ) , min ( y ^ , Y ~ ) ) , p ~ T = { max ( P T ( 0 ) , P T ( max ( y ~ , - Y ~ ) ) ) if 4 m χ 2 s ( 1 + sinh 2 ( max ( y ~ , - Y ~ ) ) cosh 2 η ^ ) < 1, P T ( 0 ) otherwise, P T ( Y ) = s 4 [ 1 + 1 - 4 m χ 2 s ( 1 + sinh 2 Y cosh 2 η ^ ) cosh ( Y - η ^ ) - 1 - 1 - 4 m χ 2 s ( 1 + sinh 2 Y cosh 2 η ^ ) cosh ( Y + η ^ ) ], y ^ = ln ( 2 p T s · sinh η ^ + cosh 2 η ^ + m χ 2 p T 2 1 + 1 - p T 2 + m χ 2 s / 4 ), y ~ = ln max ( 2 ω 1 , min s , s 2 ω 2 , max ) , Y ~ = ln min ( 2 ω 1 , max s , s 2 ω 2 , min ) .
    XSection xsection_fid(
      std::function<double (double /* s */, double /* pT */)> xsection_pT,
      Luminosity_fid luminosity,
      double mass,
      double pT_min,
      double eta_max,
      double ω_min,
      double ω_max,
      Integrator integrator = default_integrator(0)
    )
    
Convenience function for the case of symmetric bounds on energy losses of the colliding particles. Equivalent to xsection_fid(xsection_pT, luminosity, mass, pT_min, eta_max, ω_min, ω_max, ω_min, ω_max, integrator).
    XSection xsection_fid(
      std::function<double (double /* s */, double /* pT */)> xsection_pT,
      Luminosity_fid luminosity,
      double mass,
      double pT_min,
      double eta_max,
      Integrator integrator = default_integrator(0)
    )
    
Convenience function for no bounds on the energy losses of the colliding particles. Equivalent to xsection_fid(xsection_pT, luminosity, mass, pT_min, eta_max, 0, infinity, 0, infinity, integrator).
xsection_fid_b
    XSection xsection_fid_b(
      std::function<Polarization (double /* s */, double /* pT */)> xsection_pT,
      Luminosity_fid_b,
      double mass,   // of the produced particle
      double pT_min,
      double eta_max,
      double ω1_min, // use 0 for no limit
      double ω1_max, // use infinity for no limit
      double ω2_min, // use 0 for no limit
      double ω2_max, // use infinity for no limit
      Integrator = cquad_integrator(0u)
    )
    
Differential fiducial cross section for the production of a pair of particles χ + χ - in an ultraperipheral collision of particles A and B ( A B A B χ + χ - ) with non-electromagnetic interactions taken into account and with the constraints on phase space p T ± > p ^ T , | η ± | < η ^ , ω 1, min < ω 1 < ω 1, max , ω 2, min < ω 2 < ω 2, max , where p ^ T is pT_min, η ^ is eta_max, p T ± and η ± are the transverse momentum and pseudorapidity of χ ± , ω 1 and ω 2 are the photons energies: d σ A B A B χ + χ - fid. d s ( s ) = max ( p ^ T , p ~ T ) s 2 1 - 4 m χ 2 s d p T d L A B d s ( s , { d σ γ γ χ + χ - d p T ( s ) , d σ γ γ χ + χ - d p T ( s ) } , max ( - y ^ , y ~ ) , min ( y ^ , Y ~ ) ) , where d σ γ γ χ + χ - d p T ( s ) and d σ γ γ χ + χ - d p T ( s ) are the cross sections for the fusion of photons with parallel and perpendicular polarizations respectively provided by the function xsection_pT. See xsection_fid for the definitions of p ~ T , y ^ , y ~ and Y ~ .
    XSection xsection_fid_b(
      std::function<Polarization (double /* s */, double /* pT */)> xsection_pT,
      Luminosity_fid luminosity,
      double mass,
      double pT_min,
      double eta_max,
      double ω_min,
      double ω_max,
      Integrator integrator = cquad_integrator(0u)
    )
    
Convenience function for symmetric bounds on energy losses of the colliding particles. Equivalent to xsection_fid_b(xsection_pT, luminosity, mass, pT_min, eta_max, ω_min, ω_max, ω_min, ω_max, integrator).
    XSection xsection_fid_b(
      std::function<Polarization (double /* s */, double /* pT */)> xsection_pT,
      Luminosity_fid luminosity,
      double mass,
      double pT_min,
      double eta_max,
      Integrator integrator = cquad_integrator(0u)
    )
    
Convenience function for no bounds on energy losses of the colliding particles. Equivalent to xsection_fid_b(xsection_pT, luminosity, mass, pT_min, eta_max, 0, infinity, 0, infinity, integrator).
Cross sections of photons fusion
photons_to_fermions
    std::function<double (double /* s */)> photons_to_fermions(double mass, double charge = 1);
  
Cross section for the reaction γ γ χ + χ - : σ γ γ χ + χ - ( s ) = 4 π α 2 Z 4 s [ ( 1 + 4 m χ 2 s - 8 m χ 4 s ) ln 1 + 1 - 4 m χ 2 / s 1 - 1 - 4 m χ 2 / s - ( 1 + 4 m χ 2 s ) 1 - 4 m χ 2 s ], where χ ± is a fermion with the mass m χ (mass) and charge Z (charge), s is the invariant mass of the photons.
photons_to_fermions_pT
    std::function<double (double /* s */, double /* pT */)> photons_to_fermions_pT(double mass, double charge = 1);
  
Differential cross section d σ γ γ χ + χ - d p T ( s , p T ) = 8 π α 2 Z 4 p T s ( p T 2 + m χ 2 ) · 1 - 2 ( p T 4 + m χ 4 ) s ( p T 2 + m χ 2 ) 1 - 4 ( p T 2 + m χ 2 ) s , where χ ± is a fermion with the mass m χ (mass) and charge Z (charge), s is the invariant mass of the photons, pT is the transverse momentum of each of the fermions in the center of mass system.
photons_to_fermions_b
    std::function<Polarization (double /* s */)> photons_to_fermions_b(double mass, double charge = 1);
  
Cross section for the reaction γ γ χ + χ - with polarized photons. Returns two values: σ γ γ χ + χ - ( s ) = 4 π α 2 Z 4 s [ ( 1 + 4 m χ 2 s - 12 m χ 4 s ) ln 1 + 1 - 4 m χ 2 / s 1 - 1 - 4 m χ 2 / s - ( 1 + 6 m χ 2 s ) 1 - 4 m χ 2 s ], σ γ γ χ + χ - ( s ) = 4 π α 2 Z 4 s [ ( 1 + 4 m χ 2 s - 4 m χ 4 s ) ln 1 + 1 - 4 m χ 2 / s 1 - 1 - 4 m χ 2 / s - ( 1 + 2 m χ 2 s ) 1 - 4 m χ 2 s ]. Here σ γ γ χ + χ - ( s ) is the cross section for the photons with parallel polarization, σ γ γ χ + χ - ( s ) is the cross section for the photons with perpendicular polarization, χ ± is a fermion with the mass m χ (mass) and charge Z (charge), s is the invariant mass of the photons.
photons_to_fermions_pT_b
    std::function<Polarization (double /* s */, double /* pT */)> photons_to_fermions_pT_b(double mass, double charge = 1);
  
Differential cross section returning two values: d σ γ γ χ + χ - d p T ( s , p T ) = 8 π α 2 Z 4 p T s ( p T 2 + m χ 2 ) · 1 - 2 ( p T 4 + 2 m χ 4 ) s ( p T 2 + m χ 2 ) 1 - 4 ( p T 2 + m χ 2 ) s , d σ γ γ χ + χ - d p T ( s , p T ) = 8 π α 2 Z 4 p T s ( p T 2 + m χ 2 ) · 1 - 2 p T 4 s ( p T 2 + m χ 2 ) 1 - 4 ( p T 2 + m χ 2 ) s , Here σ γ γ χ + χ - ( s ) is the cross section for the photons with parallel polarization, σ γ γ χ + χ - ( s ) is the cross section for the photons with perpendicular polarization, χ ± is a fermion with the mass m χ (mass) and charge Z (charge), s is the invariant mass of the photons, pT is the transverse momentum of each of the fermions in the center of mass system.

Proton

  #include <epa/proton.hpp>
This section is devoted to functions specific for proton and proton-proton collisions. They are defined in the epa namespace as well.
Constants
Type Constant Value Units Comment
double proton_mass 0.93827208816 GeV Proton mass
double proton_radius 0.8414 * fm GeV−1 Proton radius
double proton_magnetic_moment 2.79284734463 μB Proton magnetic moment in units of Bohr magneton
double proton_dipole_form_factor_lambda2 12.0 / proton_radius2 GeV2 The Λ 2 parameter of the proton form factor
Numerical values are taken from E. Tiesinga, P. J. Mohr, D. B. Newell, B. N. Taylor, CODATA recommended values of the fundamental physical constants: 2018, Reviews of Modern Physics 93, 025010 (2021).
Form factor
proton_dipole_form_factor
    FormFactor proton_dipole_form_factor(double lambda2 = proton_dipole_form_factor_lambda2);
  
Proton Dirac form factor in dipole approximation: F p ( Q 2 ) = 1 ( 1 + Q 2 Λ 2 ) 2 · 1 + μ p Q 2 4 m p 2 1 + Q 2 4 m p 2 , where μ p is the proton magnetic moment, m p is the proton mass, Λ 2 is lambda2.
Equivalent photon spectrum
proton_dipole_spectrum
    Spectrum proton_dipole_spectrum(double energy, double lambda2 = proton_dipole_form_factor_lambda2);
  
Equivalent photon spectrum for a proton of given energy with both Dirac and Pauli form factors taken in the dipole approximation: n p ( ω ) = α π ω { ( 1 + 4 u - ( μ p 2 - 1 ) u v ) ln ( 1 + 1 u ) - 24 u 2 + 42 u + 17 6 ( u + 1 ) 2 - μ p 2 - 1 ( v - 1 ) 3 [ 1 + u / v v - 1 ln u + v u + 1 - 6 u 2 ( v 2 - 3 v + 3 ) + 3 u ( 3 v 2 - 9 v + 10 ) + 2 v 2 - 7 v + 11 6 ( u + 1 ) 2 ] } , where u = ( ω Λ γ ) 2 , v = ( 2 m p Λ ) 2 , m p is the proton mass, μ p is the proton magnetic moment, Λ 2 is lambda2. This spectrum was obtained by integrating the function n p ( ω ) = 2 α π ω 0 D ( q 2 + ( ω / γ ) 2 ) ( q 2 + ( ω / γ ) 2 ) 2 q 3 d q , where D ( Q 2 ) = G E 2 ( Q 2 ) + Q 2 4 m p 2 G M 2 ( Q 2 ) 1 + Q 2 4 m p 2 , G E ( Q 2 ) and G M ( Q 2 ) are the Sachs electric and magnetic form factors: G E ( Q 2 ) = 1 ( 1 + Q 2 Λ 2 ) 2 , G M ( Q 2 ) = μ p ( 1 + Q 2 Λ 2 ) 2 , This is the correct spectrum, but its _b counterpart has not been derived yet.
proton_dipole_spectrum_Dirac
    Spectrum proton_dipole_spectrum_Dirac(double energy, double lambda2 = proton_dipole_form_factor_lambda2);
  
Equivalent photon spectrum for a proton of given energy with the Dirac form factor taken in dipole approximation and the Pauli form factor neglected: n p, Dirac ( ω ) = α π ω { ( 1 + 4 u - 2 ( μ p - 1 ) u v ) ln ( 1 + 1 u ) + μ p - 1 ( v - 1 ) 4 ( μ p - 1 v - 1 ( 1 + 4 u + 3 v ) - 2 ( 1 + u v ) ) ln u + v u + 1 - ( 4 + 1 u + 1 + 1 6 ( u + 1 ) 2 ) + 2 ( μ p - 1 ) [ ( u v + 1 + u / v ( v - 1 ) 3 ) 1 u + 1 + 1 2 ( u v - 1 + u / v ( v - 1 ) 2 ) 1 ( u + 1 ) 2 + 1 3 ( u v + 1 + u / v v - 1 ) 1 ( u + 1 ) 3 ] - ( μ p - 1 ) 2 [ 1 ( v - 1 ) 4 + 1 + 3 u + 2 v ( u + 1 ) ( v - 1 ) 4 - 1 + 2 u + v 2 ( u + 1 ) 2 ( v - 1 ) 3 + 1 2 ( u + 1 ) 2 ( v - 1 ) 2 ] } or, equivalently, n p, Dirac ( ω ) = α π ω { ( 1 + 4 u - 2 ( μ p - 1 ) u v ) ln ( 1 + 1 u ) + μ p - 1 ( v - 1 ) 4 [ μ p - 1 v - 1 ( 1 + 4 u + 3 v ) - 2 ( 1 + u v ) ] ln u + v u + 1 - 24 u 2 + 42 u + 17 6 ( u + 1 ) 2 + ( μ p - 1 ) 6 u 2 ( v 2 - 3 v + 3 ) + 3 u ( 3 v 2 - 9 v + 10 ) + 2 v 2 - 7 v + 11 3 ( u + 1 ) 2 ( v - 1 ) 3 - ( μ p - 1 ) 2 24 u 2 + 6 u ( v + 7 ) - v 2 + 8 v + 17 6 ( u + 1 ) 2 ( v - 1 ) 4 }. Here u = ( ω Λ γ ) 2 , v = ( 2 m p Λ ) 2 , m p is the proton mass, μ p is the proton magnetic moment, Λ 2 is lambda2.
proton_dipole_spectrum_b_Dirac
    Spectrum_b proton_dipole_spectrum_b_Dirac(
      double energy,
      double lambda2 = proton_dipole_form_factor_lambda2
    )
  
Equivalent photon spectrum for a proton of given energy with the Dirac form factor taken in dipole approximation and the Pauli form factor neglected: n p , Dirac ( b , ω ) = α π 2 ω [ ω γ K 1 ( b ω γ ) - ( 1 + ( μ p - 1 ) Λ 4 16 m p 4 ( 1 - Λ 2 4 m p 2 ) 2 ) Λ 2 + ω 2 γ 2 K 1 ( b Λ 2 + ω 2 γ 2 ) + ( μ p - 1 ) Λ 4 16 m p 4 ( 1 - Λ 2 4 m p 4 ) 2 4 m p 2 + ω 2 γ 2 K 1 ( b 4 m p 2 + ω 2 γ 2 ) - 1 - μ p Λ 2 4 m p 2 1 - Λ 2 4 m p 2 · b Λ 2 2 K 0 ( b Λ 2 + ω 2 γ 2 ) ] 2 . Here m p is the proton mass, μ p is the proton magnetic moment, Λ 2 is lambda2.
Luminosity
pp_elastic_slope
    double pp_elastic_slope(double collision_energy);
  
The slope of the cross section for elastic scattering of two protons. This is the parameter B used in arXiv:1112.3243: B = B 0 + 2 B 1 ln ( E / E 0 ) + 4 B 2 ln 2 ( E / E 0 ) , where B 0 = 12 GeV - 2 , B 1 = - 0.22 ± 0.17 GeV - 2 , B 2 = 0.037 ± 0.006 GeV - 2 , E is collision_energy, E 0 = 1 GeV.
pp_upc_probability
    std::function<double (double /* b */)> pp_upc_probability(double collision_energy);
  
The probability to avoid non-electromagnetic interactions in an ultraperipheral collision of two protons with the impact parameter b: P ( b ) = ( 1 - e - b 2 2 B ) 2 . Uses pp_elastic_slope to obtain the parameter B.
pp_luminosity
    Luminosity pp_luminosity(double collision_energy, Integrator = default_integrator(0));
  
Photon-photon luminosity in ultraperipheral proton-proton collisions with non-electromagnetic interactions neglected. This function calls proton_dipole_spectrum and passes the result to luminosity.
pp_luminosity_y
    Luminosity_y pp_luminosity_y(double collision_energy);
  
Photon-photon luminosity in ultraperipheral proton-proton collisions differentiated with respect to both invariant mass of the photons s and rapidity of the system y . This function calls proton_dipole_spectrum and passes the result to luminosity_y.
pp_luminosity_fid
    Luminosity_fid pp_luminosity_fid(double collision_energy, Integrator = default_integrator(0));
  
Photon-photon luminosity in ultraperipheral proton-proton collisions with parameterized bounds for the rapidity integral. This function calls proton_dipole_spectrum and passes the result to luminosity_fid.
ppx_luminosity_b
    Luminosity_b ppx_luminosity_b(
      Spectrum   spectrum, // optional --- when provided might speed up calculations
      Spectrum_b spectrum_b,
      double B, // = pp_elastic_slope(collision_energy)
      const std::function<Integrator (unsigned)>& integrator_generator = default_integrator,
      unsigned integration_level = 0
    )
  
Photon-photon luminosity in ultraperipheral collisions of two identical particles with the same probability to avoid strong interactions as protons but with arbitrary spectrum. Use this function if you want to work with a custom proton spectrum. d L A B d s ( s , { p , p } ) = p d L d s ( s ) + p d L d s ( s ) , d L d s ( s ) = π 2 s 0 b 1 d b 1 0 b 2 d b 2 - d y n ( b 1 , s 2 e y ) n ( b 2 , s 2 e - y ) { 1 - 2 e - b 1 2 + b 2 2 2 B [ I 0 ( b 1 b 2 B ) + I 2 ( b 1 b 2 B ) ] + e - b 1 2 + b 2 2 B [ I 0 ( 2 b 1 b 2 B ) + I 2 ( 2 b 1 b 2 B ) ] } , d L d s ( s ) = π 2 s 0 b 1 d b 1 0 b 2 d b 2 - d y n ( b 1 , s 2 e y ) n ( b 2 , s 2 e - y ) { 1 - 2 e - b 1 2 + b 2 2 2 B [ I 0 ( b 1 b 2 B ) - I 2 ( b 1 b 2 B ) ] + e - b 1 2 + b 2 2 B [ I 0 ( 2 b 1 b 2 B ) - I 2 ( 2 b 1 b 2 B ) ] } . The spectrum argument is optional. If provided, it it is used to calculate the 1 in the braces. This should improve the convergence of the integral and might also result in better accuracy. B is the elastic slope of the cross section, see pp_elastic_slope.
ppx_luminosity_y_b
    Luminosity_y_b ppx_luminosity_y_b(
      Spectrum_b spectrum_b,
      double B, // = pp_elastic_slope(collision_energy)
      const std::function<Integrator (unsigned)>& integrator_generator = default_integrator,
      unsigned integration_level = 0
    )
  
Photon-photon luminosity differentiated with respect to the rapidity of the system y in ultraperipheral collisions of two identical particles with the same probability to avoid strong interactions as protons but with arbitrary spectrum. Use this function if you want to work with a custom proton spectrum. d L A B d s d y ( s , y , { p , p } ) = p d L d s d y ( s , y ) + p d L d s d y ( s , y ) , where d L d s d y ( s , y ) = π 2 s 0 b 1 d b 1 0 b 2 d b 2 n ( b 1 , s 2 e y ) n ( b 2 , s 2 e - y ) { 1 - 2 e - b 1 2 + b 2 2 2 B [ I 0 ( b 1 b 2 B ) + I 2 ( b 1 b 2 B ) ] + e - b 1 2 + b 2 2 B [ I 0 ( 2 b 1 b 2 B ) + I 2 ( 2 b 1 b 2 B ) ] } , d L d s d y ( s , y ) = π 2 s 0 b 1 d b 1 0 b 2 d b 2 n ( b 1 , s 2 e y ) n ( b 2 , s 2 e - y ) { 1 - 2 e - b 1 2 + b 2 2 2 B [ I 0 ( b 1 b 2 B ) - I 2 ( b 1 b 2 B ) ] + e - b 1 2 + b 2 2 B [ I 0 ( 2 b 1 b 2 B ) - I 2 ( 2 b 1 b 2 B ) ] } . Here B is the elastic slope of the cross section. See pp_elastic_slope.
ppx_luminosity_fid_b
    Luminosity_y_b ppx_luminosity_fid_b(
      Spectrum   spectrum, // optional --- when provided might speed up calculations
      Spectrum_b spectrum_b,
      double B, // = pp_elastic_slope(collision_energy)
      const std::function<Integrator (unsigned)>& integrator_generator = default_integrator,
      unsigned integration_level = 0
    )
  
Photon-photon luminosity with parameterized bounds for the rapidity integral in ultraperipheral collisions of two identical particles with the same probability to avoid strong interactions as protons but with arbitrary spectrum. Use this function if you want to work with a custom proton spectrum. d L A B d s ( s , { p , p } , y min , y max ) = p d L d s ( s ) + p d L d s ( s , y min , y max ) , d L d s ( s , y min , y max ) = π 2 s 0 b 1 d b 1 0 b 2 d b 2 - y min y max d y n ( b 1 , s 2 e y ) n ( b 2 , s 2 e - y ) { 1 - 2 e - b 1 2 + b 2 2 2 B [ I 0 ( b 1 b 2 B ) + I 2 ( b 1 b 2 B ) ] + e - b 1 2 + b 2 2 B [ I 0 ( 2 b 1 b 2 B ) + I 2 ( 2 b 1 b 2 B ) ] } , d L d s ( s , y min , y max ) = π 2 s 0 b 1 d b 1 0 b 2 d b 2 - y min y max d y n ( b 1 , s 2 e y ) n ( b 2 , s 2 e - y ) { 1 - 2 e - b 1 2 + b 2 2 2 B [ I 0 ( b 1 b 2 B ) - I 2 ( b 1 b 2 B ) ] + e - b 1 2 + b 2 2 B [ I 0 ( 2 b 1 b 2 B ) - I 2 ( 2 b 1 b 2 B ) ] } . Here B is the elastic slope of the cross section. See pp_elastic_slope.
pp_luminosity_b
  Luminosity_b pp_luminosity_b(
    double collision_energy,
    const std::function<Integrator (unsigned)>& integrator_generator = default_integrator,
    unsigned integration_level = 0
  )
  
Photon-photon luminosity in ultraperipheral collisions of two protons. This function calls ppx_luminosity_b passing it the results of proton_dipole_spectrum_Dirac and proton_dipole_spectrum_b_Dirac.
pp_luminosity_y_b
  Luminosity_y_b pp_luminosity_y_b(
    double collision_energy,
    const std::function<Integrator (unsigned)>& integrator_generator = default_integrator,
    unsigned integration_level = 0
  )
  
Photon-photon luminosity differentiated with respect to the system rapidity y in ultraperipheral collisions of two protons. This function calls ppx_luminosity_y_b passing it the results of proton_dipole_spectrum_Dirac and proton_dipole_spectrum_b_Dirac.
pp_luminosity_fid_b
    Luminosity_fid_b pp_luminosity_fid_b(
      double collision_energy,
      const std::function<Integrator (unsigned)>& integrator_generator = default_integrator,
      unsigned integration_level = 0
    )
  
Photon-photon luminosity with parameterized bounds of the rapidity integral in ultraperipheral collisions of two protons. This function calls ppx_luminosity_fid_b passing it the results of proton_dipole_spectrum_Dirac and proton_dipole_spectrum_b_Dirac.
Cross sections
In the following, l denotes the produced fermion and it stands for "lepton": it is assumed that the produced particles do not rescatter on the protons and, in particular, do not participate in strong interactions, hence the name.
pp_to_ppll
    XSection pp_to_ppll(
      double collision_energy,
      double mass,   // of the produced particle
      double charge, // of the produced partice
      Integrator = default_integrator(0)
    )
    
Differential cross section d σ / d s for the production of a fermion pair in ultraperipheral collisions of protons. Non-electromagnetic interactions of protons are neglected.
    XSection pp_to_ppll(
      double collision_energy,
      double mass,   // of the produced particle
      Integrator = default_integrator(0)
    )
    
Equivalent to the above with charge = 1.
    XSection pp_to_ppll(
      double collision_energy,
      double mass,    // of the produced particle
      double charge   // of the produced particle
      double pT_min   // minimal transverse momentum of the produced particle
      double eta_max  // maximal pseudorapidity of the produced particle
      const std::function<Integrator (unsigned)>& integrator_generator = default_integrator,
      unsigned integration_level
    )
    
Differential fiducial cross section d σ / d s for the production of a fermion pair in ultraperipheral collisions of protons with the constraints p T > p T min , | η | < η max where p T and η are the transverse momentum and pseudorapidity of the produced particle.
    XSection pp_to_ppll(
      double collision_energy,
      double mass,    // of the produced particle
      double pT_min   // minimal transverse momentum of the produced particle
      double eta_max  // maximal pseudorapidity of the produced particle
      const std::function<Integrator (unsigned)>& integrator_generator = default_integrator,
      unsigned integration_level
    )
    
Equivalent to the above with charge = 1.
pp_to_ppll_b
    XSection pp_to_ppll_b(
      double collision_energy,
      double mass,   // of the produced particle
      double charge, // of the produced partice
      Integrator = default_integrator(0)
    )
    
Differential cross section d σ / d s for the production of a fermion pair in ultraperipheral collisions of protons. Strong and electromagnetic interactions of protons are taken into account.
    XSection pp_to_ppll_b(
      double collision_energy,
      double mass,   // of the produced particle
      Integrator = default_integrator(0)
    )
    
Equivalent to the above with charge = 1.
    XSection pp_to_ppll_b(
      double collision_energy,
      double mass,    // of the produced particle
      double charge,  // of the produced partice
      double pT_min,  // minimal transverse momentum of the produced particle
      double eta_max, // maximal pseudorapidity of the produced particle
      const std::function<Integrator (unsigned)>& integrator_generator = pp_to_ppll_b_integrator(0),
      unsigned integration_level = 0
    )
    
Fiducial cross section d σ / d s for the production of a fermion pair in ultraperipheral collisions of protons. Strong and electromagnetic interactions of protons are taken into account.
    XSection pp_to_ppll_b(
      double collision_energy,
      double mass,    // of the produced particle
      double pT_min,  // minimal transverse momentum of the produced particle
      double eta_max, // maximal pseudorapidity of the produced particle
      const std::function<Integrator (unsigned)>& integrator_generator = pp_to_ppll_b_integrator(0),
      unsigned integration_level = 0
    )
    
Equivalent to the above with charge = 1.
pp_to_ppll_b_integrator
    std::function<Integrator (unsigned)> pp_to_ppll_b_integrator(unsigned level = 0);
  
Default integrator for pp_to_ppll_b. Uses CQUAD for integration level level and QAG for other levels.