libepa Programmer referencealphaamubarncquad_integratordefault_absolute_error
default_cquad_integration_limit
default_error_stepdefault_integration_limit
default_integration_method
default_integratordefault_relative_error
electron_chargefmFormFactorform_factor_dipoleform_factor_monopoleinfinityIntegratorlight_speedLuminosityLuminosity_bLuminosity_fidLuminosity_fid_bLuminosity_yLuminosity_y_bluminosityluminosity_bluminosity_fidluminosity_fid_bluminosity_yluminosity_y_bpiphotons_to_fermionsphotons_to_fermions_b
photons_to_fermions_pT
photons_to_fermions_pT_b
planckPolarizationpp_elastic_slopepp_luminositypp_luminosity_bpp_luminosity_fidpp_luminosity_fid_bpp_luminosity_ypp_luminosity_y_bpp_upc_probabilitypp_to_ppllpp_to_ppll_bpp_to_ppll_b_integrator
ppx_luminosity_bppx_luminosity_y_bppx_luminosity_fid_bprint_backtraceproton_dipole_form_factor
proton_dipole_form_factor_lambda2
proton_dipole_spectrum
proton_dipole_spectrum_b_Dirac
proton_dipole_spectrum_Dirac
proton_magnetic_moment
proton_massproton_radiusqag_integratorqag_integrator_generator
SpectrumSpectrum_bspectrumspectrum_bspectrum_b_dipolespectrum_b_monopolespectrum_b_pointspectrum_dipolespectrum_monopoleXSectionXSection_bxsectionxsection_fidxsection_fid_bA 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 , where 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.
libepa — a C++/Python library for calculations of cross
sections of ultraperipheral collisions.
Computer Physics Communications 305, 109347 (2024);
arXiv:2311.01353
— a paper with a more detailed description of the physics of
ultraperipheral collisions, also discussing the examples and comparing
libepa to other programs.
*_b family of functions deal with the case when the source
particles interact non-electromagnetically (strongly or weakly). Other
functions assume only electromagnetic interactions.
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
(, where 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.
libepa language, integration level of an integral is
the nesting number of the integral in the calculation. For example, when
computing
,
the integral with respect to is of integration level
0, and the integral with respect to 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.
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.
#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
where is the Planck constant, and
is the speed of light.
Current libepa can be determined through cpp macros
EPA_VERSION_MAJOR, EPA_VERSION_MINOR and
EPA_VERSION_PATCH which expand to integer literals.
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.
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].
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 points, where
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.
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
to
.
form_factor_monopole
FormFactor form_factor_monopole(double lambda2);
Monopole form factor:
lambda2 is in
GeV2.
form_factor_dipole
FormFactor form_factor_dipole(double lambda2);
Dipole form factor:
lambda2 is in
GeV2.
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:
spectrum_monopole
Spectrum spectrum_monopole(unsigned Z, double gamma, double lambda2);
Equivalent photon spectrum for the
monopole form factor:
where
spectrum_dipole
Spectrum spectrum_dipole(unsigned Z, double gamma, double lambda2);
Equivalent photon spectrum for the
dipole form factor:
where
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:
Note that the form factor is weighted with an oscillating function (the
Bessel function ) and
integrated over semi-infinite interval. For small values of
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:
spectrum_b_monopole
Spectrum_b spectrum_b_monopole(unsigned Z, double gamma, double lambda2);
Equivalent photon spectrum for the
monopole form factor:
spectrum_b_dipole
Spectrum_b spectrum_b_dipole(unsigned Z, double gamma, double lambda2);
Equivalent photon spectrum for the
dipole form factor:
Luminosity
typedef std::function<double (double /* */)> Luminosity;
Luminosity 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
of the photon-photon system
.
Luminosity_y
typedef std::function<double (double /* */, double /* y */)> Luminosity_y;
Luminosity_y is the luminosity differentiated with respect to
the invariant mass and rapidity
y of the photon-photon system
.
Luminosity_fid
typedef std::function<double (double /* */, 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
of the photon-photon system
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
and differentiated with
respect to the invariant mass of the
photon-photon system:
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:
luminosity_y
Luminosity_y luminosity_y(Spectrum nA, Spectrum nB);
Photon-photon luminosity in ultraperipheral collision of particles
and differentiated with
respect to both invariant mass of the photons
and rapidity of the system
:
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
and differentiated with
respect to the invariant mass of the photons
with parameterized bounds for the
rapidity integral:
Luminosity_fid luminosity_fid(Spectrum, Integrator = default_integrator(0));
Convenience function for the collision of identical particles. When
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 /* */, Polarization /* p */)> Luminosity_b
typedef std::function<double (double /* */, double /* y */, Polarization /* p */)> Luminosity_y_b;
typedef std::function<double (double /* */, 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
and differentiated with
respect to the invariant mass of the
photon-photon system:
where
is the luminosity of photons with parallel polarizations,
is the luminosity of photons with perpendicular polarizations,
is upc_probability, the probability for the particles to survive
in an ultraperipheral collision with the impact parameter
(the probability to avoid electromagnetic
interactions),
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:
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
and differentiated with
respect to the photons invariant mass
and rapidity :
where
is the luminosity of photons with parallel polarizations,
is the luminosity of photons with perpendicular polarizations,
is upc_probability, the probability for the particles to survive
in an ultraperipheral collision with the impact parameter
(the probability to avoid electromagnetic
interactions),
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
and differentiated with
respect to the photons invariant mass
with parameterized bounds on the rapidity integral:
where
is the luminosity of photons with parallel polarizations,
is the luminosity of photons with perpendicular polarizations,
is upc_probability, the probability for the particles to survive
in an ultraperipheral collision with the impact parameter
(the probability to avoid electromagnetic
interactions),
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
it takes into account the symmetry of the integrand, improving on
convergence and calculation time.
XSection
typedef std::function<double (double /* */)> XSection;
Cross section for the invariant mass
( is the Mandelstam parameter).
Some functions use this type to denote cross section differentiated with respect
to :
XSection_b
typedef std::function<Polarization (double /* */)> 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
in ultraperipheral collisions with photon-photon
luminosity with non-electromagnetic interactions
neglected:
where
is photon_photon, the cross section for the fusion of two (real)
photons.
xsection_fid
XSection xsection_fid(
std::function<double (double /* */, 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 and
with non-electromagnetic interactions neglected and with the constraints on phase space
,
,
,
,
where
is pT_min,
is eta_max,
and
are the transverse momentum and pseudorapidity of
,
and
are the photons energies:
XSection xsection_fid(
std::function<double (double /* */, 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 /* */, 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 /* */, 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 and
with non-electromagnetic interactions taken into account and with the
constraints on phase space
,
,
,
,
where
is pT_min,
is eta_max,
and
are the transverse momentum and pseudorapidity of
,
and
are the photons energies:
where
and
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
,
,
and
.
XSection xsection_fid_b(
std::function<Polarization (double /* */, 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 /* */, 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).
photons_to_fermions
std::function<double (double /* */)> photons_to_fermions(double mass, double charge = 1);
Cross section for the reaction
where is a fermion with
the mass
(mass) and charge (charge),
is the invariant mass of the photons.
photons_to_fermions_pT
std::function<double (double /* */, double /* pT */)> photons_to_fermions_pT(double mass, double charge = 1);
Differential cross section
where is a fermion with
the mass
(mass) and charge (charge),
is the invariant mass of the photons,
is the transverse momentum of
each of the fermions in the center of mass system.
photons_to_fermions_b
std::function<Polarization (double /* */)> photons_to_fermions_b(double mass, double charge = 1);
Cross section for the reaction
with polarized photons. Returns two values:
Here
is the cross section for the photons with parallel polarization,
is the cross section for the photons with perpendicular polarization, is a fermion with the mass
(mass) and
charge (charge),
is the invariant mass of the photons.
photons_to_fermions_pT_b
std::function<Polarization (double /* */, double /* pT */)> photons_to_fermions_pT_b(double mass, double charge = 1);
Differential cross section returning two values:
Here
is the cross section for the photons with parallel polarization,
is the cross section for the photons with perpendicular polarization, is a fermion with the mass
(mass) and
charge (charge),
is the invariant mass of the photons,
is the transverse momentum of
each of the fermions in the center of mass system.
#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.
proton_dipole_form_factor
FormFactor proton_dipole_form_factor(double lambda2 = proton_dipole_form_factor_lambda2);
Proton Dirac form factor in dipole approximation:
where is the proton
magnetic moment, is the
proton mass, is
lambda2.
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:
where
is the proton mass,
is the proton magnetic
moment, is
lambda2. This spectrum was obtained by integrating the function
where
and
are the Sachs electric and magnetic form factors:
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:
or, equivalently,
Here
is the proton mass,
is the proton magnetic
moment, 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:
Here
is the proton mass,
is the proton magnetic moment,
is lambda2.
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 used in arXiv:1112.3243:
where
is collision_energy,
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:
Uses pp_elastic_slope to obtain
the parameter .
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
and rapidity of the system
. 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.
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 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.
where
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.
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 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.
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
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
for the production of a fermion pair in ultraperipheral collisions of
protons with the constraints
where 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
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
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.