LPAIR’s \(\ggll\)#
The \(pp \rightarrow p^{(\ast)}(\ggll)p^{(\ast)}\) process as previously implemented in the early 1990s in LPAIR
[BDSV91] can be reached through the lpair
process in CepGen.
More generally, this implementation was designed to compute the cross-section and to generate events for the \(\ggll\) process for ee, ep, and pp collisions.
A phenomenological review of both this process and its first implementation in LPAIR
may be found in [Ver83].
Process-specific options#
mode
(integer): kinematic regime to generate and the size of the phase space to perform the integration. It can take the following values:ProcessMode.ElasticElastic := 1
: elastic emission of photons from the incoming protons (default value if unspecified),ProcessMode.ElasticInelastic := 2 / ProcessMode.InelasticElastic := 3
: elastic scattering of one photon and an inelastic/semi-exclusive emission of the other photon, resulting in the excitation/fragmentation of the outgoing proton state,ProcessMode.InelasticInelastic := 4
: both protons fragmented in the final state.
pair
(integer/PDG): PDG identifier of the lepton to be produced in the final state. It can hence take the following values:PDG.electron := 11
for the \(e^+e^-\) pair production,PDG.muon := 13
for the \(\mu^+\mu^-\) pair production, andPDG.tau := 15
for the \(\tau^+\tau^-\) pair production.
Full object reference#
-
class LPAIR : public Process#
Matrix element for the \(\gamma\gamma\to\ell^{+}\ell^{-}\) process as defined in [28].
Private Functions
-
bool orient()#
Calculate energies and momenta of full event content, in the CM system.
-
inline double periPP() const#
Compute the squared matrix element squared for the \(\gamma\gamma\rightarrow\ell^{+}\ell^{-}\) process.
Note
Its expression is of the form:
\[M = \frac{1}{4bt_1 t_2}\sum_{i=1}^2\sum_{j=1}^2 u_i v_j t_{ij} = \frac{1}{4}\frac{u_1 v_1 t_{11}+u_2 v_1 t_{21}+u_1 v_2 t_{12}+u_2 v_2 t_{22}}{t_1 t_2 b}\]where \(b\) = bb_ is defined in computeWeight as:\[b = t_1 t_2+\left(w_{\gamma\gamma}\sin^2{\theta^{\rm CM}_6}+4m_\ell\cos^2{\theta^{\rm CM}_6}\right) p_g^2\]- Returns:
Convolution of the form factor or structure functions with the squared central two-photons matrix element (for a pair of spin \(-\frac{1}{2}-\)point particles)
-
double pickin()#
Describe the kinematics of the process \(p_1+p_2\to p_3+p_4+p_5\) in terms of Lorentz-invariant variables.
Note
These variables (along with others) will then be fed into the periPP method (thus are essential for the evaluation of the full matrix element).
- Returns:
Value of the Jacobian after the operation
Private Members
-
const ParticleProperties pair_#
-
const bool symmetrise_#
-
const bool randomise_charge_#
-
double ml_ = {0.}#
mass of the outgoing leptons
-
double ml2_ = {0.}#
squared mass of the outgoing leptons
-
double charge_factor_ = {0.}#
-
mode::Kinematics beams_mode_#
-
double re_ = {0.}#
-
double ep1_ = {0.}#
energy of the first proton-like incoming particle
-
double ep2_ = {0.}#
energy of the second proton-like incoming particle
-
double w12_ = {0.}#
\(\delta_2=m_1^2-m_2^2\) as defined in [28]
-
double ss_ = {0.}#
-
double p12_ = {0.}#
\(p_{12} = \frac{1}{2}\left(s-m_{p_1}^2-m_{p_2}^2\right)\)
-
double sl1_ = {0.}#
-
double e1mp1_ = {0.}#
-
double p_cm_ = {0.}#
-
double mom_prefactor_ = {0.}#
-
double gamma_cm_ = {0.}#
-
double beta_gamma_cm_ = {0.}#
-
std::unique_ptr<formfac::Parameterisation> formfac1_#
-
std::unique_ptr<formfac::Parameterisation> formfac2_#
-
std::unique_ptr<strfun::Parameterisation> strfun_#
-
bool is_strfun_sy_ = {false}#
-
double m_u_t1_ = {0.}#
\(t_1\), first parton normalised virtuality
-
double m_u_t2_ = {0.}#
\(t_2\), second parton normalised virtuality
-
double m_u_s2_ = {0.}#
\(s_2\)
-
double m_w4_ = {0.}#
\(w_4\), squared invariant mass of the two-parton system
-
double m_theta4_ = {0.}#
polar angle of the two-photon system
-
double m_phi6_cm_ = {0.}#
\(\phi_6^{\rm CM}\), azimutal angle of the first outgoing lepton xx6 = \(\frac{1}{2}\left(1-\cos\theta^{\rm CM}_6\right)\) definition (3D rotation of the first outgoing lepton with respect to the two-photon centre-of-mass system).
Note
If the nm_ optimisation flag is set this angle coefficient value becomes
\[\frac{1}{2}\left(\frac{a_{\rm map}}{b_{\rm map}}\frac{\beta-1}{\beta+1}+1\right)\]with \(a_{\rm map}=\frac{1}{2}\left(w_4-t_1-t_2\right)\), \(b_{\rm map}=\frac{1}{2}\sqrt{\left(\left(w_4-t_1-t_2\right)^2-4t_1t_2\right)\left(1-4\frac{w_6}{w_4}\right)}\), and \(\beta=\left(\frac{a_{\rm map}+b_{\rm map}}{a_{\rm map}-b_{\rm map}}\right)^{2x_5-1}\) and the Jacobian element is scaled by a factor \(\frac{1}{2}\frac{\left(a_{\rm map}^2-b_{\rm map}^2\cos^2\theta^{\rm CM}_6\right)}{a_{\rm map}b_{\rm map}}\log\left(\frac{a_{\rm map}+b_{\rm map}}{a_{\rm map}-b_{\rm map}}\right)\)
-
double m_x6_ = {0.}#
-
double s1_ = {0.}#
-
double s2_ = {0.}#
-
double sa1_ = {0.}#
-
double sa2_ = {0.}#
-
double p1k2_ = {0.}#
-
double p2k1_ = {0.}#
-
double ec4_ = {0.}#
central system energy
-
double pc4_ = {0.}#
central system 3-momentum norm
-
double pt4_ = {0.}#
central system transverse momentum
-
double mc4_ = {0.}#
central system invariant mass
-
double cos_theta4_ = {0.}#
central system polar angle cosine
-
double sin_theta4_ = {0.}#
central system polar angle sine
-
double q2dq_ = {0.}#
-
double epsilon_ = {0.}#
-
double alpha4_ = {0.}#
-
double beta4_ = {0.}#
-
double gamma4_ = {0.}#
-
double alpha5_ = {0.}#
-
double gamma5_ = {0.}#
-
double alpha6_ = {0.}#
-
double gamma6_ = {0.}#
-
double bb_ = {0.}#
-
double gram_ = {0.}#
-
double dd5_ = {0.}#
-
std::array<double, 2> deltas1_#
-
std::array<double, 2> deltas2_#
-
double delta_ = {0.}#
Invariant used to tame divergences in the matrix element computation.
Note
Defined as
\[\Delta = \left(p_1\cdot p_2\right)\left(q_1\cdot q_2\right)-\left(p_1\cdot q_2\right)\left(p_2\cdot q_1\right)\]with \(p_i, q_i\) the 4-momenta associated to the incoming proton-like particle and to the photon emitted from it.
-
double eph1_ = {0.}#
-
double eph2_ = {0.}#
Private Static Attributes
-
static constexpr double constb_ = 0.5 * M_1_PI * M_1_PI * M_1_PI#
-
bool orient()#