StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
SigmaTotal.cc
1 // SigmaTotal.cc is a part of the PYTHIA event generator.
2 // Copyright (C) 2018 Torbjorn Sjostrand.
3 // PYTHIA is licenced under the GNU GPL v2 or later, see COPYING for details.
4 // Please respect the MCnet Guidelines, see GUIDELINES for details.
5 
6 // Function definitions (not found in the header) for the SigmaTotal class,
7 // and the auxiliary classes SigmaTotAux, SigmaTotOwn, SigmaSaSDL, SigmaMBR,
8 // SigmaABMST and SigmaRPP.
9 
10 #include "Pythia8/SigmaTotal.h"
11 
12 namespace Pythia8 {
13 
14 //==========================================================================
15 
16 // The SigmaTotAux class.
17 // Base class for the individual implementations.
18 
19 //--------------------------------------------------------------------------
20 
21 // alpha_em(0).
22 const double SigmaTotAux::ALPHAEM = 0.00729353;
23 
24 // Conversion coefficients = 1/(16pi) * (mb <-> GeV^2).
25 const double SigmaTotAux::HBARC2 = 0.38938;
26 const double SigmaTotAux::CONVERTEL = 0.0510925;
27 
28 // Proton and pion masses, and their squares. Euler's constant.
29 const double SigmaTotAux::MPROTON = 0.9382720;
30 const double SigmaTotAux::SPROTON = 0.8803544;
31 const double SigmaTotAux::MPION = 0.1349766;
32 const double SigmaTotAux::SPION = 0.0182187;
33 const double SigmaTotAux::GAMMAEUL = 0.577215665;
34 
35 // Reference scale for nominal definition of t slope.
36 const double SigmaTotAux::TABSREF = 2e-3;
37 
38 // Numerical integration parameters.
39 const int SigmaTotAux::NPOINTS = 1000;
40 const double SigmaTotAux::TABSMAX = 1.;
41 const double SigmaTotAux::MINSLOPEEL = 10.;
42 
43 //--------------------------------------------------------------------------
44 
45 // Initialize Coulomb correction parameters.
46 
47 bool SigmaTotAux::initCoulomb(Settings& settings,
48  ParticleData* particleDataPtrIn) {
49 
50  // Save pointer to particle database.
51  particleDataPtr = particleDataPtrIn,
52 
53  // User-set values for total and elastic cross section.
54  tryCoulomb = settings.flag("SigmaElastic:Coulomb");
55  rhoOwn = settings.parm("SigmaElastic:rho");
56  tAbsMin = settings.parm("SigmaElastic:tAbsMin");
57  lambda = settings.parm("SigmaElastic:lambda");
58  phaseCst = settings.parm("SigmaElastic:phaseConst");
59 
60  return true;
61 
62 }
63 
64 //--------------------------------------------------------------------------
65 
66 // Add Coulomb corrections to the elastic and total cross sections.
67 
68 bool SigmaTotAux::addCoulomb() {
69 
70  // Trivial case when there should be no Coulomb contribution.
71  hasCou = false;
72  sigTotCou = sigTot;
73  sigElCou = sigEl;
74 
75  // Relative sign (or zero) for Coulomb term in elastic scattering.
76  int iChA = particleDataPtr->chargeType(idA);
77  int iChB = particleDataPtr->chargeType(idB);
78  chgSgn = 0.;
79  if (iChA * iChB > 0) chgSgn = 1.;
80  if (iChA * iChB < 0) chgSgn = -1.;
81 
82  // Done if no Coulomb corrections.
83  if (!tryCoulomb || iChA * iChB == 0) return false;
84 
85  // Reduce hadronic part of elastic cross section by tMin cut.
86  sigElCou = sigEl * exp( - bEl * tAbsMin);
87  if (tAbsMin < 0.9 * TABSMAX) {
88 
89  // Loop through t range according to dt/t^2.
90  double sigCou = 0.;
91  double sigInt = 0.;
92  double xRel, tAbs, form2, phase;
93  for (int i = 0; i < NPOINTS; ++i) {
94  xRel = (i + 0.5) / NPOINTS;
95  tAbs = tAbsMin * TABSMAX / (tAbsMin + xRel * (TABSMAX - tAbsMin));
96 
97  // Evaluate Coulomb and interference terms.
98  form2 = pow4(lambda/(lambda + tAbs));
99  sigCou += pow2(form2);
100  phase = chgSgn * ALPHAEM * (-phaseCst - log(0.5 * bEl * tAbs));
101  sigInt += form2 * exp(-0.5 * bEl * tAbs) * tAbs
102  * (rhoOwn * cos(phase) + sin(phase));
103  }
104 
105  // Include common factors to give new elastic and total cross sections.
106  sigCou *= pow2(ALPHAEM) / (4. * CONVERTEL * tAbsMin) ;
107  sigInt *= - chgSgn * ALPHAEM * sigTot / tAbsMin;
108  sigElCou += (sigCou + sigInt) / NPOINTS;
109  hasCou = true;
110  }
111  sigTotCou = sigTot - sigEl + sigElCou;
112 
113  // Done.
114  return true;
115 
116 }
117 
118 //--------------------------------------------------------------------------
119 
120 // Coulomb contribution to the differential elastic cross sections.
121 
122 double SigmaTotAux::dsigmaElCoulomb( double t) {
123 
124  // Coulomb contribution and interference.
125  double form2 = pow4(lambda/(lambda - t));
126  double phase = chgSgn * ALPHAEM * (-phaseCst - log(-0.5 * bEl * t));
127  return pow2(chgSgn * ALPHAEM * form2) / (4. * CONVERTEL * t*t)
128  - chgSgn * ALPHAEM * form2 * sigTot * exp(0.5 * bEl * t)
129  * (rhoOwn * cos(phase) + sin(phase)) / (-t);
130 
131 }
132 
133 //==========================================================================
134 
135 // The SigmaTotal class.
136 // Parametrizes total, elastic and diffractive cross sections,
137 // making use of the other clases below.
138 
139 //--------------------------------------------------------------------------
140 
141 // Definitions of static variables.
142 
143 // Minimum threshold below which no cross sections will be defined.
144 const double SigmaTotal::MMIN = 2.;
145 
146 //--------------------------------------------------------------------------
147 
148 // Store pointer to Info and initialize data members.
149 
150 void SigmaTotal::init( Info* infoPtrIn, Settings& settings,
151  ParticleData* particleDataPtrIn, Rndm* rndmPtrIn) {
152 
153  // Store pointers.
154  infoPtr = infoPtrIn;
155  particleDataPtr = particleDataPtrIn;
156  settingsPtr = &settings;
157  rndmPtr = rndmPtrIn;
158 
159  // Choice of mode.
160  modeTotEl = settings.mode("SigmaTotal:mode");
161  modeDiff = settings.mode("SigmaDiffractive:mode");
162 
163 }
164 
165 //--------------------------------------------------------------------------
166 
167 // Calculate, or recalculate for new beams or new energy.
168 
169 bool SigmaTotal::calc(int idA, int idB, double eCM) {
170 
171  // Initial values.
172  isCalc = ispp = false;
173  s = eCM * eCM;
174 
175  // Find hadron masses and check that energy is enough.
176  // For mesons use the corresponding vector meson masses.
177  idAbsA = abs(idA);
178  idAbsB = abs(idB);
179  int idModA = (idAbsA < 100 || idAbsA > 1000) ? idAbsA : 10 * (idAbsA/10) + 3;
180  int idModB = (idAbsB < 100 || idAbsB > 1000) ? idAbsB : 10 * (idAbsB/10) + 3;
181  if (idAbsA == 22) idModA = 113;
182  if (idAbsB == 22) idModB = 113;
183  if (idAbsA == 990) idModA = idAbsA;
184  if (idAbsB == 990) idModB = idAbsB;
185  double mA = particleDataPtr->m0(idModA);
186  double mB = particleDataPtr->m0(idModB);
187  if (eCM < mA + mB + MMIN) {
188  infoPtr->errorMsg("Error in SigmaTotal::calc: too low energy");
189  return false;
190  }
191 
192  // Most options only work for pp/ppbar, so may need to modify choice.
193  // Treat a neutron like a proton (except no Coulomb term).
194  modeTotElNow = modeTotEl;
195  modeDiffNow = modeDiff;
196  if (idAbsA == 2112) idAbsA = 2212;
197  if (idAbsB == 2112) idAbsB = 2212;
198  if (idAbsA != 2212 || idAbsB != 2212) {
199  modeTotElNow = min(1, modeTotEl);
200  modeDiffNow = min(1, modeDiff);
201  }
202  ispp = (idAbsA == 2212 && idAbsB == 2212 && idA * idB > 0);
203 
204  // Set up pointer to class that handles total and elastic cross sections.
205  if (sigTotElPtr) delete sigTotElPtr;
206  if (modeTotElNow == 0) sigTotElPtr = new SigmaTotOwn();
207  else if (modeTotElNow == 1) sigTotElPtr = new SigmaSaSDL();
208  else if (modeTotElNow == 2) sigTotElPtr = new SigmaMBR();
209  else if (modeTotElNow == 3) sigTotElPtr = new SigmaABMST();
210  else sigTotElPtr = new SigmaRPP();
211 
212  // Initialize and calculate for selected total/elastic class.
213  sigTotElPtr->init( infoPtr, *settingsPtr, particleDataPtr, rndmPtr);
214  if ( !sigTotElPtr->calcTotEl( idA, idB, s, mA, mB) ) return false;
215 
216  // Set up pointer to class that handles diffractive cross sections.
217  if (sigDiffPtr) delete sigDiffPtr;
218  if (modeDiffNow == 0) sigDiffPtr = new SigmaTotOwn();
219  else if (modeDiffNow == 1) sigDiffPtr = new SigmaSaSDL();
220  else if (modeDiffNow == 2) sigDiffPtr = new SigmaMBR();
221  else sigDiffPtr = new SigmaABMST();
222 
223  // Initialize and calculate for selected diffractive class.
224  if (sigDiffPtr != sigTotElPtr)
225  sigDiffPtr->init( infoPtr, *settingsPtr, particleDataPtr, rndmPtr);
226  if ( !sigDiffPtr->calcDiff( idA, idB, s, mA, mB) ) return false;
227 
228  // Inelastic nondiffractive by unitarity.
229  double sigTotTmp = sigTotElPtr->sigTot;
230  sigND = sigTotTmp - sigTotElPtr->sigEl - sigDiffPtr->sigXB
231  - sigDiffPtr->sigAX - sigDiffPtr->sigXX - sigDiffPtr->sigAXB;
232  if (sigND < 0.) {
233  infoPtr->errorMsg("Error in SigmaTotal::init: sigND < 0");
234  return false;
235  } else if (sigND < 0.4 * sigTotTmp) infoPtr->errorMsg("Warning in "
236  "SigmaTotal::init: sigND suspiciously low");
237 
238  // Done.
239  isCalc = true;
240  return true;
241 
242 }
243 
244 //==========================================================================
245 
246 // The SigmaTotOwn class.
247 // Parametrizes total, elastic and diffractive cross sections
248 // by user settings.
249 
250 //--------------------------------------------------------------------------
251 
252 // Initialize data members.
253 
254 void SigmaTotOwn::init(Info* , Settings& settings,
255  ParticleData* particleDataPtrIn, Rndm*) {
256 
257  // Main user-set values for total and elastic cross sections.
258  sigTot = settings.parm("SigmaTotal:sigmaTot");
259  sigEl = settings.parm("SigmaTotal:sigmaEl");
260  bEl = settings.parm("SigmaElastic:bSlope");
261 
262  // Initialize parameters for Coulomb corrections to elastic scattering.
263  initCoulomb( settings, particleDataPtrIn);
264 
265  // User-set values for diffractive cross sections.
266  sigXB = settings.parm("SigmaTotal:sigmaXB");
267  sigAX = settings.parm("SigmaTotal:sigmaAX");
268  sigXX = settings.parm("SigmaTotal:sigmaXX");
269  sigAXB = settings.parm("SigmaTotal:sigmaAXB");
270 
271  // Set diffraction parameters.
272  pomFlux = settings.mode("SigmaDiffractive:PomFlux");
273 
274  // Set up Pomeron flux constants, see HardDiffraction::init.
275  a0 = 1. + settings.parm("SigmaDiffractive:PomFluxEpsilon");
276  ap = settings.parm("SigmaDiffractive:PomFluxAlphaPrime");
277 
278  // Schuler-Sjöstrand.
279  if (pomFlux == 1) {
280  b0 = 2.3;
281 
282  // Bruni-Ingelman.
283  } else if (pomFlux == 2) {
284  A1 = 6.38;
285  A2 = 0.424;
286  a1 = 8.;
287  a2 = 3.;
288 
289  // Streng-Berger.
290  } else if (pomFlux == 3) {
291  a1 = 4.7;
292 
293  // Donnachie-Landshoff.
294  } else if (pomFlux == 4) {
295  A1 = 0.27;
296  a1 = 8.38;
297  A2 = 0.56;
298  a2 = 3.78;
299  A3 = 0.18;
300  a3 = 1.36;
301 
302  // MBR.
303  } else if (pomFlux == 5) {
304  A1 = 0.9;
305  a1 = 4.6;
306  A2 = 0.1;
307  a2 = 0.6;
308  a0 = 1. + settings.parm("SigmaDiffractive:MBRepsilon");
309  ap = settings.parm("SigmaDiffractive:MBRalpha");
310 
311  // H1 Fit A, B.
312  } else if (pomFlux == 6 || pomFlux == 7) {
313  ap = 0.06;
314  b0 = 5.5;
315  a0 = (pomFlux == 6) ? 1.1182 : 1.1110;
316  }
317 
318  // b_min for double diffraction, suppression of small gaps, minimal CD mass.
319  bMinDD = settings.parm("SigmaDiffractive:OwnbMinDD");
320  dampenGap = settings.flag("SigmaDiffractive:OwndampenGap");
321  ygap = settings.parm("SigmaDiffractive:Ownygap");
322  ypow = settings.parm("SigmaDiffractive:Ownypow");
323  expPygap = exp(ypow * ygap);
324  mMinCDnow = settings.parm("SigmaDiffractive:OwnmMinCD");
325 
326 }
327 
328 //--------------------------------------------------------------------------
329 
330 // With total and elastic cross section already set, only add Coulomb.
331 
332 bool SigmaTotOwn::calcTotEl( int idAin, int idBin, double , double , double) {
333 
334  // Save some input.
335  idA = idAin;
336  idB = idBin;
337  isExpEl = true;
338 
339  // Possibly allow Coulomb correction + interference.
340  addCoulomb();
341 
342  // Done.
343  return true;
344 
345 }
346 
347 //--------------------------------------------------------------------------
348 
349 // Differential elastic cross section.
350 
351 double SigmaTotOwn::dsigmaEl( double t, bool useCoulomb, bool ) {
352 
353  // Hadronic contribution: simple exponential.
354  double dsig = CONVERTEL * pow2(sigTot) * (1. + pow2(rhoOwn)) * exp(bEl * t);
355 
356  // Possibly add Coulomb contribution and interference.
357  if (useCoulomb && hasCou) dsig += dsigmaElCoulomb(t);
358 
359  // Done.
360  return dsig;
361 
362 }
363 
364 //--------------------------------------------------------------------------
365 
366 // Differential single diffractive cross sections as provided by the
367 // x-weighted Pomeron flux, see HardDiffraction::xfPom for details.
368 
369 double SigmaTotOwn::dsigmaSD(double xi, double t, bool , int ) {
370 
371  // Common setup.
372  wtNow = 1.;
373  yNow = -log(xi);
374 
375  // Schuler-Sjöstrand.
376  if (pomFlux == 1) {
377  b = 2. * b0 + 2. * ap * yNow;
378  wtNow = exp(b * t);
379 
380  // Bruni-Ingelman.
381  } else if (pomFlux == 2) {
382  wtNow = A1 * exp(a1 * t) + A2 * exp(a2 * t);
383 
384  // Streng-Berger.
385  } else if (pomFlux == 3) {
386  b = a1 + 2. * ap * yNow;
387  wtNow = pow( xi, 2. - 2. * a0) * exp(b * t);
388 
389  // Donnachie-Landshoff.
390  } else if (pomFlux == 4) {
391  Q = 2. * ap * yNow;
392  wtNow = pow( xi, 2. - 2. * a0) * ( A1 * exp((Q + a1) * t)
393  + A2 * exp((Q + a2) * t) + A3 * exp((Q + a3) * t) );
394 
395  // MBR.
396  } else if (pomFlux == 5) {
397  Q = 2. * ap * yNow;
398  wtNow = pow( xi, 2. - 2. * a0)
399  * (A1 * exp((Q + a1) * t) + A2 * exp((Q + a2) * t) );
400 
401  // H1 Fit A, B.
402  } else if (pomFlux == 6 || pomFlux == 7) {
403  b = b0 + 2. * ap * yNow;
404  wtNow = pow( xi, 2. - 2. * a0) * exp(b * t);
405  }
406  // Optionally dampen with 1 / (1 + exp( -p * (y - y_gap))).
407  if (dampenGap) wtNow /= 1. + expPygap * pow( xi, ypow);
408 
409  // Done.
410  return wtNow;
411 
412 }
413 
414 //--------------------------------------------------------------------------
415 
416 // Differential double diffractive cross sections.
417 // Naive non-authorized generalization of the simple Pomeron fluxes above.
418 
419 double SigmaTotOwn::dsigmaDD(double xi1, double xi2, double t, int ) {
420 
421  // Common setup.
422  wtNow = 1.;
423  yNow = - log(xi1 * xi2 * s / SPROTON);
424 
425  // Schuler-Sjöstrand.
426  if (pomFlux == 1) {
427  b = max( bMinDD, 2. * ap * yNow);
428  wtNow = exp(b * t);
429 
430  // Bruni-Ingelman.
431  } else if (pomFlux == 2) {
432  wtNow = A1 * exp(a1 * t) + A2 * exp(a2 * t);
433 
434  // Streng-Berger.
435  } else if (pomFlux == 3) {
436  b = max( bMinDD, 2. * ap * yNow);
437  wtNow = pow( xi1 * xi2, 2. - 2. * a0) * exp(b * t);
438 
439  // Donnachie-Landshoff.
440  } else if (pomFlux == 4) {
441  Q = max( bMinDD, 2. * ap * yNow);
442  wtNow = pow( xi1 * xi2, 2. - 2. * a0) * exp(Q * t);
443 
444  // MBR.
445  } else if (pomFlux == 5) {
446  Q = max( bMinDD, 2. * ap * yNow);
447  wtNow = pow( xi1 * xi2, 2. - 2. * a0) * exp(Q * t);
448 
449  // H1 Fit A, B.
450  } else if (pomFlux == 6 || pomFlux == 7) {
451  b = max( bMinDD, 2. * ap * yNow);
452  wtNow = pow( xi1 * xi2, 2. - 2. * a0) * exp(b * t);
453  }
454 
455  // Optionally dampen with 1 / (1 + exp( -p * (y - y_gap))).
456  if (dampenGap) wtNow /= 1. + expPygap * pow( xi1 * xi2 * s / SPROTON, ypow);
457 
458  // Done.
459  return wtNow;
460 
461 }
462 
463 //--------------------------------------------------------------------------
464 
465 // Differential central diffractive cross section.
466 // Naive non-authorized generalization of the simple Pomeron fluxes above.
467 
468 double SigmaTotOwn::dsigmaCD( double xi1, double xi2, double t1, double t2,
469  int ) {
470 
471  // Common setup.
472  wtNow = 1.;
473  yNow1 = -log(xi1);
474  yNow2 = -log(xi2);
475 
476  // Schuler-Sjöstrand.
477  if (pomFlux == 1) {
478  b1 = 2. * b0 + 2. * ap * yNow1;
479  b2 = 2. * b0 + 2. * ap * yNow2;
480  wtNow = exp(b1 * t1 + b2 * t2);
481 
482  // Bruni-Ingelman.
483  } else if (pomFlux == 2) {
484  wtNow = (A1 * exp(a1 * t1) + A2 * exp(a2 * t1))
485  * (A1 * exp(a1 * t2) + A2 * exp(a2 * t2));
486 
487  // Streng-Berger.
488  } else if (pomFlux == 3) {
489  b1 = a1 + 2. * ap * yNow1;
490  b2 = a1 + 2. * ap * yNow2;
491  wtNow = pow( xi1 * xi2, 2. - 2. * a0) * exp(b1 * t1 + b2 * t2);
492 
493  // Donnachie-Landshoff.
494  } else if (pomFlux == 4) {
495  Q1 = 2. * ap * yNow1;
496  Q2 = 2. * ap * yNow2;
497  wtNow = pow( xi1 * xi2, 2. - 2. * a0) * ( A1 * exp((Q1 + a1) * t1)
498  + A2 * exp((Q1 + a2) * t1) + A3 * exp((Q1 + a3) * t1) )
499  * ( A1 * exp((Q2 + a1) * t2) + A2 * exp((Q2 + a2) * t2)
500  + A3 * exp((Q2 + a3) * t2) );
501 
502  // MBR.
503  } else if (pomFlux == 5) {
504  Q1 = 2. * ap * yNow1;
505  Q2 = 2. * ap * yNow2;
506  wtNow = pow( xi1 * xi2, 2. - 2. * a0)
507  * (A1 * exp((Q1 + a1) * t1) + A2 * exp((Q1 + a2) * t1) )
508  * (A1 * exp((Q2 + a1) * t2) + A2 * exp((Q2 + a2) * t2) );
509 
510  // H1 Fit A, B.
511  } else if (pomFlux == 6 || pomFlux == 7) {
512  b1 = b0 + 2. * ap * yNow1;
513  b2 = b0 + 2. * ap * yNow2;
514  wtNow = pow( xi1 * xi2, 2. - 2. * a0) * exp(b1 * t1 + b2 * t2);
515  }
516 
517  // Optionally dampen with 1 / (1 + exp( -p * (y - y_gap))) for both gaps.
518  if (dampenGap) wtNow /= (1. + expPygap * pow( xi1, ypow))
519  * (1. + expPygap * pow( xi2, ypow));
520 
521  // Done.
522  return wtNow;
523 
524 }
525 
526 //==========================================================================
527 
528 // The SigmaSaSDL class.
529 // Formulae are taken from:
530 // G.A. Schuler and T. Sjostrand, Phys. Rev. D49 (1994) 2257,
531 // Z. Phys. C73 (1997) 677
532 // which borrows some total cross sections from
533 // A. Donnachie and P.V. Landshoff, Phys. Lett. B296 (1992) 227.
534 
535 // Implemented processes with their process number iProc:
536 // = 0 : p + p; = 1 : pbar + p;
537 // = 2 : pi+ + p; = 3 : pi- + p; = 4 : pi0/rho0 + p;
538 // = 5 : phi + p; = 6 : J/psi + p;
539 // = 7 : rho + rho; = 8 : rho + phi; = 9 : rho + J/psi;
540 // = 10 : phi + phi; = 11 : phi + J/psi; = 12 : J/psi + J/psi.
541 // = 13 : gamma + p (preliminary); 14 : gamma + gamma (preliminary);
542 // = 15 : Pom + p (preliminary).
543 // For now a neutron is treated like a proton.
544 
545 //--------------------------------------------------------------------------
546 
547 // Definitions of static variables.
548 // Note that a lot of parameters are hardcoded as const here, rather
549 // than being interfaced for public change, since any changes would
550 // have to be done in a globally consistent manner. Which basically
551 // means a rewrite/replacement of the whole class.
552 
553 // General constants in total cross section parametrization:
554 // sigmaTot = X * s^epsilon + Y * s^eta (pomeron + reggeon).
555 // Christine added Pom + p, gamma + gamma and gamma + p.
556 const double SigmaSaSDL::EPSILON = 0.0808;
557 const double SigmaSaSDL::ETA = -0.4525;
558 const double SigmaSaSDL::X[] = { 21.70, 21.70, 13.63, 13.63, 13.63,
559  10.01, 0.970, 8.56, 6.29, 0.609, 4.62, 0.447, 0.0434, 67.7e-3, 211e-6};
560 const double SigmaSaSDL::Y[] = { 56.08, 98.39, 27.56, 36.02, 31.79,
561  1.51, -0.146, 13.08, -0.62, -0.060, 0.030, -0.0028, 0.00028, 129e-3,
562  215e-6};
563 
564 // Type of the two incoming hadrons as function of the process number:
565 // = 0 : p/n ; = 1 : pi/rho/omega; = 2 : phi; = 3 : J/psi.
566 const int SigmaSaSDL::IHADATABLE[] = { 0, 0, 1, 1, 1, 2, 3, 1, 1,
567  1, 2, 2, 3};
568 const int SigmaSaSDL::IHADBTABLE[] = { 0, 0, 0, 0, 0, 0, 0, 1, 2,
569  3, 2, 3, 3};
570 
571 // Hadron-Pomeron coupling beta(t) = beta(0) * exp(b*t).
572 // 0 = Pom + p, 1 = Pom + pi/rho/omega, 2 = Pom + phi, 3 = Pom + J/psi.
573 const double SigmaSaSDL::BETA0[] = { 4.658, 2.926, 2.149, 0.208};
574 const double SigmaSaSDL::BHAD[] = { 2.3, 1.4, 1.4, 0.23};
575 
576 // Pomeron trajectory alpha(t) = 1 + epsilon + alpha' * t
577 const double SigmaSaSDL::ALPHAPRIME = 0.25;
578 
579 // Conversion coefficients = 1/(16pi) * (mb <-> GeV^2) * (G_3P)^n,
580 // with n = 0 elastic, n = 1 single and n = 2 double diffractive.
581 const double SigmaSaSDL::CONVERTSD = 0.0336;
582 const double SigmaSaSDL::CONVERTDD = 0.0084;
583 
584 const double SigmaSaSDL::GAMMAFAC[3] = {2.2, 23.6, 18.4};
585 const double SigmaSaSDL::VMDMASS[3] = {0.77549, 0.78265, 1.0146};
586 
587 // Parameters and coefficients for single diffractive scattering.
588 const int SigmaSaSDL::ISDTABLE[] = { 0, 0, 1, 1, 1, 2, 3, 4, 5,
589  6, 7, 8, 9};
590 const double SigmaSaSDL::CSD[10][8] = {
591  { 0.213, 0.0, -0.47, 150., 0.213, 0.0, -0.47, 150., } ,
592  { 0.213, 0.0, -0.47, 150., 0.267, 0.0, -0.47, 100., } ,
593  { 0.213, 0.0, -0.47, 150., 0.232, 0.0, -0.47, 110., } ,
594  { 0.213, 7.0, -0.55, 800., 0.115, 0.0, -0.47, 110., } ,
595  { 0.267, 0.0, -0.46, 75., 0.267, 0.0, -0.46, 75., } ,
596  { 0.232, 0.0, -0.46, 85., 0.267, 0.0, -0.48, 100., } ,
597  { 0.115, 0.0, -0.50, 90., 0.267, 6.0, -0.56, 420., } ,
598  { 0.232, 0.0, -0.48, 110., 0.232, 0.0, -0.48, 110., } ,
599  { 0.115, 0.0, -0.52, 120., 0.232, 6.0, -0.56, 470., } ,
600  { 0.115, 5.5, -0.58, 570., 0.115, 5.5, -0.58, 570. } };
601 
602 // Parameters and coefficients for double diffractive scattering.
603 const int SigmaSaSDL::IDDTABLE[] = { 0, 0, 1, 1, 1, 2, 3, 4, 5,
604  6, 7, 8, 9};
605 const double SigmaSaSDL::CDD[10][9] = {
606  { 3.11, -7.34, 9.71, 0.068, -0.42, 1.31, -1.37, 35.0, 118., } ,
607  { 3.11, -7.10, 10.6, 0.073, -0.41, 1.17, -1.41, 31.6, 95., } ,
608  { 3.12, -7.43, 9.21, 0.067, -0.44, 1.41, -1.35, 36.5, 132., } ,
609  { 3.13, -8.18, -4.20, 0.056, -0.71, 3.12, -1.12, 55.2, 1298., } ,
610  { 3.11, -6.90, 11.4, 0.078, -0.40, 1.05, -1.40, 28.4, 78., } ,
611  { 3.11, -7.13, 10.0, 0.071, -0.41, 1.23, -1.34, 33.1, 105., } ,
612  { 3.12, -7.90, -1.49, 0.054, -0.64, 2.72, -1.13, 53.1, 995., } ,
613  { 3.11, -7.39, 8.22, 0.065, -0.44, 1.45, -1.36, 38.1, 148., } ,
614  { 3.18, -8.95, -3.37, 0.057, -0.76, 3.32, -1.12, 55.6, 1472., } ,
615  { 4.18, -29.2, 56.2, 0.074, -1.36, 6.67, -1.14, 116.2, 6532. } };
616 
617 //--------------------------------------------------------------------------
618 
619 // Store pointer to Info and initialize data members.
620 
621 void SigmaSaSDL::init( Info* infoPtrIn, Settings& settings,
622  ParticleData* particleDataPtrIn, Rndm*) {
623 
624  // Store pointer.
625  infoPtr = infoPtrIn;
626 
627  // Initialize parameters for Coulomb corrections to elastic scattering.
628  initCoulomb( settings, particleDataPtrIn);
629 
630  // User-set values to dampen diffractive cross sections.
631  doDampen = settings.flag("SigmaDiffractive:dampen");
632  maxXBOwn = settings.parm("SigmaDiffractive:maxXB");
633  maxAXOwn = settings.parm("SigmaDiffractive:maxAX");
634  maxXXOwn = settings.parm("SigmaDiffractive:maxXX");
635  maxAXBOwn = settings.parm("SigmaDiffractive:maxAXB");
636 
637  // Parameters for diffractive systems.
638  epsSaS = settings.parm("SigmaDiffractive:SaSepsilon");
639  sigmaPomP = settings.parm("Diffraction:sigmaRefPomP");
640  mPomP = settings.parm("Diffraction:mRefPomP");
641  pPomP = settings.parm("Diffraction:mPowPomP");
642  zeroAXB = settings.flag("SigmaTotal:zeroAXB");
643  sigAXB2TeV = settings.parm("SigmaTotal:sigmaAXB2TeV");
644 
645  // Diffractive mass spectrum starts at m + mMin0 and has a low-mass
646  // enhancement, factor cRes, up to around m + mRes0. Special for CD.
647  mMin0 = settings.parm("SigmaDiffractive:mMin");
648  cRes = settings.parm("SigmaDiffractive:lowMEnhance");
649  mRes0 = settings.parm("SigmaDiffractive:mResMax");
650  mMinCDnow = settings.parm("SigmaDiffractive:mMinCD");
651 
652  // Derived quantities.
653  alP2 = 2. * ALPHAPRIME;
654  s0 = 1. / ALPHAPRIME;
655 
656 }
657 
658 //--------------------------------------------------------------------------
659 
660 // Calculate total and (integrated) elastic cross sections.
661 
662 bool SigmaSaSDL::calcTotEl( int idAin, int idBin, double sIn, double mAin,
663  double mBin) {
664 
665  // Find appropriate combination of incoming beams.
666  idA = idAin;
667  idB = idBin;
668  s = sIn;
669  isExpEl = true;
670  if (!findBeamComb( idA, idB, mAin, mBin)) return false;
671  double sEps = pow( s, EPSILON);
672  double sEta = pow( s, ETA);
673 
674  // Ordinary hadron-hadron collisions. Total cross section.
675  if (iProc < 13) {
676  sigTot = X[iProc] * sEps + Y[iProc] * sEta;
677 
678  // Elastic slope parameter and cross section.
679  bEl = 2. * bA + 2. * bB + 4. * sEps - 4.2;
680  sigEl = CONVERTEL * pow2(sigTot) * (1. + pow2(rhoOwn)) / bEl;
681 
682  // gamma + p and p + gamma.
683  } else if (iProc == 13){
684  sigTot = X[iProc] * sEps + Y[iProc] * sEta;
685  double sigElNow = 0.;
686 
687  // Loop over VMD states on side A for elastic cross section.
688  for (int iA = 0; iA < 3; ++iA){
689  double bANow = BHAD[iHadAtmp[iA]];
690  double bBNow = BHAD[iHadBtmp[iA]];
691  double bElNow = 2. * bANow + 2. * bBNow + 4. * sEps - 4.2;
692  double tmpTot = X[iProcVP[iA]] * sEps
693  + Y[iProcVP[iA]] * sEta;
694  sigElNow += multVP[iA] * CONVERTEL * pow2(tmpTot)
695  * (1. + pow2(rhoOwn)) / bElNow;
696  }
697  sigEl = sigElNow;
698 
699  // gamma + gamma
700  } else if (iProc == 14){
701  sigTot = X[iProc] * sEps + Y[iProc] * sEta;
702  double sigElNow = 0.;
703 
704  // Loop over VMD states on side A and B for elastic cross section.
705  for (int iA = 0; iA < 3; ++iA)
706  for (int iB = 0; iB < 3; ++iB) {
707  double bANow = BHAD[iHadAtmp[iA]];
708  double bBNow = BHAD[iHadBtmp[iB]];
709  double bElNow = 2. * bANow + 2. * bBNow + 4. * sEps - 4.2;
710  double tmpTot = X[iProcVV[iA][iB]] * sEps
711  + Y[iProcVV[iA][iB]] * sEta;
712  sigElNow += multVV[iA][iB] * CONVERTEL * pow2(tmpTot)
713  * (1. + pow2(rhoOwn)) / bElNow;
714  }
715  sigEl = sigElNow;
716 
717  // Primitive implementation of Pomeron + p. No elastic scattering.
718  } else if (iProc == 15) {
719  sigTot = sigmaPomP * pow( sqrt(s) / mPomP, pPomP);
720  sigEl = 0.;
721  }
722 
723  // Possibly add Coulomb correction and interference.
724  addCoulomb();
725 
726  // Done.
727  return true;
728 
729 }
730 
731 //--------------------------------------------------------------------------
732 
733 // Differential elastic cross section.
734 
735 double SigmaSaSDL::dsigmaEl( double t, bool useCoulomb, bool ) {
736 
737  // Ordinary hadron-hadron collisions.
738  double dsig = 0.;
739  if (iProc < 13) {
740  dsig = CONVERTEL * pow2(sigTot) * (1. + pow2(rhoOwn)) * exp(bEl * t);
741 
742  // gamma + p and p + gamma: loop over VMD states on side A.
743  } else if (iProc == 13) {
744  double sEps = pow( s, EPSILON);
745  double sEta = pow( s, ETA);
746  double dsigNow = 0.;
747  for (int iA = 0; iA < 3; ++iA){
748  // Elastic slope parameter and cross section.
749  double bANow = BHAD[iHadAtmp[iA]];
750  double bBNow = BHAD[iHadBtmp[iA]];
751  double bElNow = 2. * bANow + 2. * bBNow + 4. * sEps - 4.2;
752  double tmpTot = X[iProcVP[iA]] * sEps
753  + Y[iProcVP[iA]] * sEta;
754  // Hadronic contribution: simple exponential.
755  dsigNow += multVP[iA] * CONVERTEL * pow2(tmpTot)
756  * (1. + pow2(rhoOwn)) * exp(bElNow * t);
757  }
758  dsig = dsigNow;
759 
760  // gamma + gamma: loop over VMD states on side A and B.
761  } else if (iProc == 14) {
762  double sEps = pow( s, EPSILON);
763  double sEta = pow( s, ETA);
764  double dsigNow = 0.;
765  for (int iA = 0; iA < 3; ++iA)
766  for (int iB = 0; iB < 3; ++iB){
767  // Elastic slope parameter and cross section.
768  double bANow = BHAD[iHadAtmp[iA]];
769  double bBNow = BHAD[iHadBtmp[iB]];
770  double bElNow = 2. * bANow + 2. * bBNow + 4. * sEps - 4.2;
771  double tmpTot = X[iProcVV[iA][iB]] * sEps
772  + Y[iProcVV[iA][iB]] * sEta;
773  // Hadronic contribution: simple exponential.
774  dsigNow += multVV[iA][iB] * CONVERTEL * pow2(tmpTot)
775  * (1. + pow2(rhoOwn)) * exp(bElNow * t);
776  }
777  dsig = dsigNow;
778 
779  // No elastic scattering for Pomeron + p.
780  } else dsig = 0.;
781 
782  // Possibly add Coulomb contribution and interference.
783  if (useCoulomb && hasCou) dsig += dsigmaElCoulomb(t);
784 
785  // Done.
786  return dsig;
787 
788 }
789 
790 //--------------------------------------------------------------------------
791 
792 // Diffractive cross sections.
793 
794 bool SigmaSaSDL::calcDiff( int idAin, int idBin, double sIn, double mAin,
795  double mBin) {
796 
797  // Find appropriate combination of incoming beams.
798  s = sIn;
799  idA = idAin;
800  idB = idBin;
801  mA = mAin;
802  mB = mBin;
803  if (!findBeamComb( idAin, idBin, mAin, mBin)) return false;
804  sigXB = sigAX = sigXX = sigAXB = 0.;
805 
806  // Ordinary hadron-hadron collisions.
807  if (iProc < 13) {
808  // Lookup coefficients for single and double diffraction.
809  int iSD = ISDTABLE[iProc];
810  int iDD = IDDTABLE[iProc];
811  double eCM = sqrt(s);
812  double sum1, sum2, sum3, sum4;
813 
814  // Single diffractive scattering A + B -> X + B cross section.
815  mMinXB = mA + mMin0;
816  double sMinXB = pow2(mMinXB);
817  mResXB = mA + mRes0;
818  sResXB = pow2(mResXB);
819  double sRMavgXB = mResXB * mMinXB;
820  double sRMlogXB = log(1. + sResXB/sMinXB);
821  double sMaxXB = CSD[iSD][0] * s + CSD[iSD][1];
822  double BcorrXB = CSD[iSD][2] + CSD[iSD][3] / s;
823  sum1 = log( (2.*bB + alP2 * log(s/sMinXB))
824  / (2.*bB + alP2 * log(s/sMaxXB)) ) / alP2;
825  sum2 = cRes * sRMlogXB / (2.*bB + alP2 * log(s/sRMavgXB) + BcorrXB) ;
826  sigXB = CONVERTSD * X[iProc] * BETA0[iHadB] * max( 0., sum1 + sum2);
827 
828  // Single diffractive scattering A + B -> A + X cross section.
829  mMinAX = mB + mMin0;
830  double sMinAX = pow2(mMinAX);
831  mResAX = mB + mRes0;
832  sResAX = pow2(mResAX);
833  double sRMavgAX = mResAX * mMinAX;
834  double sRMlogAX = log(1. + sResAX/sMinAX);
835  double sMaxAX = CSD[iSD][4] * s + CSD[iSD][5];
836  double BcorrAX = CSD[iSD][6] + CSD[iSD][7] / s;
837  sum1 = log( (2.*bA + alP2 * log(s/sMinAX))
838  / (2.*bA + alP2 * log(s/sMaxAX)) ) / alP2;
839  sum2 = cRes * sRMlogAX / (2.*bA + alP2 * log(s/sRMavgAX) + BcorrAX) ;
840  sigAX = CONVERTSD * X[iProc] * BETA0[iHadA] * max( 0., sum1 + sum2);
841 
842  // Order single diffractive correctly.
843  if (swapped) {
844  swap( bB, bA);
845  swap( sigXB, sigAX);
846  swap( mMinXB, mMinAX);
847  swap( mResXB, mResAX);
848  swap( iHadA, iHadB);
849  }
850 
851  // Double diffractive scattering A + B -> X1 + X2 cross section.
852  double y0min = log( s * SPROTON / (sMinXB * sMinAX) ) ;
853  double sLog = log(s);
854  double Delta0 = CDD[iDD][0] + CDD[iDD][1] / sLog
855  + CDD[iDD][2] / pow2(sLog);
856  sum1 = (y0min * (log( max( 1e-10, y0min/Delta0) ) - 1.) + Delta0)/ alP2;
857  if (y0min < 0.) sum1 = 0.;
858  double sMaxXX = s * ( CDD[iDD][3] + CDD[iDD][4] / sLog
859  + CDD[iDD][5] / pow2(sLog) );
860  double sLogUp = log( max( 1.1, s * s0 / (sMinXB * sRMavgAX) ));
861  double sLogDn = log( max( 1.1, s * s0 / (sMaxXX * sRMavgAX) ));
862  sum2 = cRes * log( sLogUp / sLogDn ) * sRMlogAX / alP2;
863  sLogUp = log( max( 1.1, s * s0 / (sMinAX * sRMavgXB) ));
864  sLogDn = log( max( 1.1, s * s0 / (sMaxXX * sRMavgXB) ));
865  sum3 = cRes * log(sLogUp / sLogDn) * sRMlogXB / alP2;
866  double BcorrXX = CDD[iDD][6] + CDD[iDD][7] / eCM + CDD[iDD][8] / s;
867  sum4 = pow2(cRes) * sRMlogAX * sRMlogXB
868  / max( 0.1, alP2 * log( s * s0 / (sRMavgAX * sRMavgXB) ) + BcorrXX);
869  sigXX = CONVERTDD * X[iProc] * max( 0., sum1 + sum2 + sum3 + sum4);
870 
871  // Central diffractive scattering A + B -> A + X + B, only p and pbar.
872  mMinAXB = 1.;
873  if ( (idAbsA == 2212 || idAbsA == 2112)
874  && (idAbsB == 2212 || idAbsB == 2112) && !zeroAXB) {
875  double sMinAXB = pow2(mMinAXB);
876  double sRefAXB = pow2(2000.);
877  sigAXB = sigAXB2TeV * pow( log(0.06 * s / sMinAXB), 1.5 )
878  / pow( log(0.06 * sRefAXB / sMinAXB), 1.5 );
879  }
880 
881  // Option with user-requested damping of diffractive cross sections.
882  if (doDampen) {
883  sigXB = sigXB * maxXBOwn / (sigXB + maxXBOwn);
884  sigAX = sigAX * maxAXOwn / (sigAX + maxAXOwn);
885  sigXX = sigXX * maxXXOwn / (sigXX + maxXXOwn);
886  sigAXB = (maxAXBOwn > 0.) ? sigAXB * maxAXBOwn / (sigAXB + maxAXBOwn)
887  : 0.;
888  }
889 
890  // gamma + p: loop over VMD states on side A.
891  } else if (iProc == 13) {
892  double sigAXNow = 0.;
893  double sigXBNow = 0.;
894  double sigXXNow = 0.;
895 
896  for (int iA = 0; iA < 3; ++iA){
897 
898  // Lookup coefficients for single and double diffraction.
899  int iSD = ISDTABLE[iProcVP[iA]];
900  int iDD = IDDTABLE[iProcVP[iA]];
901  double eCM = sqrt(s);
902  double sum1, sum2, sum3, sum4;
903 
904  // Single diffractive scattering A + B -> X + B cross section.
905  mMinXB = mAtmp[iA] + mMin0;
906  double sMinXB = pow2(mMinXB);
907  mResXB = mAtmp[iA] + mRes0;
908  sResXB = pow2(mResXB);
909  double sRMavgXB = mResXB * mMinXB;
910  double sRMlogXB = log(1. + sResXB/sMinXB);
911  double sMaxXB = CSD[iSD][0] * s + CSD[iSD][1];
912  double BcorrXB = CSD[iSD][2] + CSD[iSD][3] / s;
913  double bBNow = BHAD[iHadBtmp[iA]];
914  sum1 = log( (2.*bBNow + alP2 * log(s/sMinXB))
915  / (2.*bBNow + alP2 * log(s/sMaxXB)) ) / alP2;
916  sum2 = cRes * sRMlogXB
917  / (2.*bBNow + alP2 * log(s/sRMavgXB) + BcorrXB) ;
918  sigXBNow += multVP[iA] * CONVERTSD * X[iProcVP[iA]]
919  * BETA0[iHadBtmp[iA]] * max( 0., sum1 + sum2);
920 
921  // Single diffractive scattering A + B -> A + X cross section.
922  mMinAX = mBtmp[iA] + mMin0;
923  double sMinAX = pow2(mMinAX);
924  mResAX = mBtmp[iA] + mRes0;
925  sResAX = pow2(mResAX);
926  double sRMavgAX = mResAX * mMinAX;
927  double sRMlogAX = log(1. + sResAX/sMinAX);
928  double sMaxAX = CSD[iSD][4] * s + CSD[iSD][5];
929  double BcorrAX = CSD[iSD][6] + CSD[iSD][7] / s;
930  double bANow = BHAD[iHadAtmp[iA]];
931  sum1 = log( (2.*bANow + alP2 * log(s/sMinAX))
932  / (2.*bANow + alP2 * log(s/sMaxAX)) ) / alP2;
933  sum2 = cRes * sRMlogAX
934  / (2.*bANow + alP2 * log(s/sRMavgAX) + BcorrAX) ;
935  sigAXNow += multVP[iA] * CONVERTSD * X[iProcVP[iA]]
936  * BETA0[iHadAtmp[iA]] * max( 0., sum1 + sum2);
937 
938  // Double diffractive scattering A + B -> X1 + X2 cross section.
939  double y0min = log( s * SPROTON / (sMinXB * sMinAX) ) ;
940  double sLog = log(s);
941  double Delta0 = CDD[iDD][0] + CDD[iDD][1] / sLog
942  + CDD[iDD][2] / pow2(sLog);
943  sum1 = (y0min * (log( max( 1e-10, y0min/Delta0) ) - 1.) + Delta0)/ alP2;
944  if (y0min < 0.) sum1 = 0.;
945  double sMaxXX = s * ( CDD[iDD][3] + CDD[iDD][4] / sLog
946  + CDD[iDD][5] / pow2(sLog) );
947  double sLogUp = log( max( 1.1, s * s0 / (sMinXB * sRMavgAX) ));
948  double sLogDn = log( max( 1.1, s * s0 / (sMaxXX * sRMavgAX) ));
949  sum2 = cRes * log( sLogUp / sLogDn ) * sRMlogAX / alP2;
950  sLogUp = log( max( 1.1, s * s0 / (sMinAX * sRMavgXB) ));
951  sLogDn = log( max( 1.1, s * s0 / (sMaxXX * sRMavgXB) ));
952  sum3 = cRes * log(sLogUp / sLogDn) * sRMlogXB / alP2;
953  double BcorrXX = CDD[iDD][6] + CDD[iDD][7] / eCM + CDD[iDD][8] / s;
954  sum4 = pow2(cRes) * sRMlogAX * sRMlogXB
955  /max( 0.1, alP2 * log( s * s0 / (sRMavgAX * sRMavgXB) ) + BcorrXX);
956  sigXXNow += multVP[iA] * CONVERTDD * X[iProcVP[iA]]
957  * max( 0., sum1 + sum2 + sum3 + sum4);
958 
959  // End of VMD loop.
960  }
961 
962  // Order single diffractive correctly.
963  if (swapped) {
964  swap( bB, bA);
965  swap( sigXB, sigAX);
966  swap( mMinXB, mMinAX);
967  swap( mResXB, mResAX);
968  swap( iHadA, iHadB);
969  swap( sigXBNow, sigAXNow);
970  for (int i = 0; i < 3; ++i){
971  swap( iHadAtmp[i], iHadBtmp[i]);
972  swap( mAtmp[i], mBtmp[i]);
973  }
974  }
975 
976  // Option with user-requested damping of diffractive cross sections.
977  if (doDampen) {
978  sigXBNow = sigXBNow * maxXBOwn / (sigXBNow + maxXBOwn);
979  sigAXNow = sigAXNow * maxAXOwn / (sigAXNow + maxAXOwn);
980  sigXXNow = sigXXNow * maxXXOwn / (sigXXNow + maxXXOwn);
981  }
982  sigAX = sigAXNow;
983  sigXB = sigXBNow;
984  sigXX = sigXXNow;
985  sigAXB = 0.;
986 
987  // gamma + gamma: loop over VMD states on sides A and B.
988  } else if (iProc == 14) {
989  double sigAXNow = 0.;
990  double sigXBNow = 0.;
991  double sigXXNow = 0.;
992  for (int iA = 0; iA < 3; ++iA)
993  for (int iB = 0; iB < 3; ++iB){
994 
995  // Lookup coefficients for single and double diffraction.
996  int iSD = ISDTABLE[iProcVV[iA][iB]];
997  int iDD = IDDTABLE[iProcVV[iA][iB]];
998  double eCM = sqrt(s);
999  double sum1, sum2, sum3, sum4;
1000 
1001  // Single diffractive scattering A + B -> X + B cross section.
1002  mMinXB = mAtmp[iA] + mMin0;
1003  double sMinXB = pow2(mMinXB);
1004  mResXB = mAtmp[iA] + mRes0;
1005  sResXB = pow2(mResXB);
1006  double sRMavgXB = mResXB * mMinXB;
1007  double sRMlogXB = log(1. + sResXB/sMinXB);
1008  double sMaxXB = CSD[iSD][0] * s + CSD[iSD][1];
1009  double BcorrXB = CSD[iSD][2] + CSD[iSD][3] / s;
1010  double bBNow = BHAD[iHadBtmp[iB]];
1011  sum1 = log( (2.*bBNow + alP2 * log(s/sMinXB))
1012  / (2.*bBNow + alP2 * log(s/sMaxXB)) ) / alP2;
1013  sum2 = cRes * sRMlogXB
1014  / (2.*bBNow + alP2 * log(s/sRMavgXB) + BcorrXB) ;
1015  sigXBNow += multVV[iA][iB] * CONVERTSD * X[iProcVV[iA][iB]]
1016  * BETA0[iHadBtmp[iB]] * max( 0., sum1 + sum2);
1017 
1018  // Single diffractive scattering A + B -> A + X cross section.
1019  mMinAX = mBtmp[iB] + mMin0;
1020  double sMinAX = pow2(mMinAX);
1021  mResAX = mBtmp[iB] + mRes0;
1022  sResAX = pow2(mResAX);
1023  double sRMavgAX = mResAX * mMinAX;
1024  double sRMlogAX = log(1. + sResAX/sMinAX);
1025  double sMaxAX = CSD[iSD][4] * s + CSD[iSD][5];
1026  double BcorrAX = CSD[iSD][6] + CSD[iSD][7] / s;
1027  double bANow = BHAD[iHadAtmp[iA]];
1028  sum1 = log( (2.*bANow + alP2 * log(s/sMinAX))
1029  / (2.*bANow + alP2 * log(s/sMaxAX)) ) / alP2;
1030  sum2 = cRes * sRMlogAX
1031  / (2.*bANow + alP2 * log(s/sRMavgAX) + BcorrAX) ;
1032  sigAXNow += multVV[iA][iB] * CONVERTSD * X[iProcVV[iA][iB]]
1033  * BETA0[iHadAtmp[iA]] * max( 0., sum1 + sum2);
1034 
1035  // Double diffractive scattering A + B -> X1 + X2 cross section.
1036  double y0min = log( s * SPROTON / (sMinXB * sMinAX) ) ;
1037  double sLog = log(s);
1038  double Delta0 = CDD[iDD][0] + CDD[iDD][1] / sLog
1039  + CDD[iDD][2] / pow2(sLog);
1040  sum1 = (y0min * (log( max( 1e-10, y0min/Delta0) ) - 1.) + Delta0)/ alP2;
1041  if (y0min < 0.) sum1 = 0.;
1042  double sMaxXX = s * ( CDD[iDD][3] + CDD[iDD][4] / sLog
1043  + CDD[iDD][5] / pow2(sLog) );
1044  double sLogUp = log( max( 1.1, s * s0 / (sMinXB * sRMavgAX) ));
1045  double sLogDn = log( max( 1.1, s * s0 / (sMaxXX * sRMavgAX) ));
1046  sum2 = cRes * log( sLogUp / sLogDn ) * sRMlogAX / alP2;
1047  sLogUp = log( max( 1.1, s * s0 / (sMinAX * sRMavgXB) ));
1048  sLogDn = log( max( 1.1, s * s0 / (sMaxXX * sRMavgXB) ));
1049  sum3 = cRes * log(sLogUp / sLogDn) * sRMlogXB / alP2;
1050  double BcorrXX = CDD[iDD][6] + CDD[iDD][7] / eCM + CDD[iDD][8] / s;
1051  sum4 = pow2(cRes) * sRMlogAX * sRMlogXB
1052  /max( 0.1, alP2 * log( s * s0 / (sRMavgAX * sRMavgXB) ) + BcorrXX);
1053  sigXXNow += multVV[iA][iB] * CONVERTDD * X[iProcVV[iA][iB]]
1054  * max( 0., sum1 + sum2 + sum3 + sum4);
1055 
1056  // End of VMD loop.
1057  }
1058 
1059  // Option with user-requested damping of diffractive cross sections.
1060  if (doDampen) {
1061  sigXBNow = sigXBNow * maxXBOwn / (sigXBNow + maxXBOwn);
1062  sigAXNow = sigAXNow * maxAXOwn / (sigAXNow + maxAXOwn);
1063  sigXXNow = sigXXNow * maxXXOwn / (sigXXNow + maxXXOwn);
1064  }
1065  sigAX = sigAXNow;
1066  sigXB = sigXBNow;
1067  sigXX = sigXXNow;
1068  sigAXB = 0.;
1069 
1070  // No diffractive scattering for Pomeron + p.
1071  } else return false;
1072 
1073  // Done.
1074  return true;
1075 
1076 }
1077 
1078 //--------------------------------------------------------------------------
1079 
1080 // Differential single diffractive cross sections.
1081 
1082 double SigmaSaSDL::dsigmaSD(double xi, double t, bool isXB, int ) {
1083 
1084  // Calculate mass; later checked to be above threshold. Optional weight.
1085  double m2X = xi * s;
1086  double mX = sqrt(m2X);
1087  double epsWt = pow( m2X, -epsSaS);
1088 
1089  // Ordinary hadron-hadron collisions.
1090  if (iProc < 13) {
1091 
1092  // Separate XB and AX cases since A and B masses/couplings may differ.
1093  if (isXB) {
1094  if (mX < mMinXB || pow2(mX + mMinAX) > s) return 0.;
1095 
1096  // Return differential AB -> XB cross section value.
1097  double bXB = 2. * bB + alP2 * log(1. / xi) ;
1098  return CONVERTSD * X[iProc] * BETA0[iHadB] * exp(bXB * t) * (1. - xi)
1099  * ( 1. + cRes * sResXB / (sResXB + m2X) ) * epsWt;
1100 
1101  } else {
1102  if (mX < mMinAX || pow2(mX + mMinXB) > s) return 0.;
1103 
1104  // Return differential AB -> AX cross section value.
1105  double bAX = 2. * bA + alP2 * log(1. / xi) ;
1106  return CONVERTSD * X[iProc] * BETA0[iHadA] * exp(bAX * t) * (1. - xi)
1107  * ( 1. + cRes * sResAX / (sResAX + m2X) ) * epsWt;
1108  }
1109 
1110  // gamma + p: loop over VMD states on side A.
1111  } else if (iProc == 13) {
1112  double sigNow = 0.;
1113  for (int iA = 0; iA < 3; ++iA){
1114 
1115  // Mass thresholds depend on VMD state.
1116  mMinXB = mAtmp[iA] + mMin0;
1117  mResXB = mAtmp[iA] + mRes0;
1118  sResXB = pow2(mResXB);
1119  mMinAX = mBtmp[iA] + mMin0;
1120  mResAX = mBtmp[iA] + mRes0;
1121  sResAX = pow2(mResAX);
1122 
1123  // Calculate and sum differential AB -> XB cross section value.
1124  if (isXB && mX > mMinXB && pow2(mX + mMinAX) < s) {
1125  double bBNow = BHAD[iHadBtmp[iA]];
1126  double bXBNow = 2. * bBNow + alP2 * log(1. / xi) ;
1127  sigNow += multVP[iA] * CONVERTSD * X[iProcVP[iA]]
1128  * BETA0[iHadBtmp[iA]] * exp(bXBNow * t) * (1. - xi)
1129  * ( 1. + cRes * sResXB / (sResXB + m2X) );
1130 
1131  // Calculate and sum differential AB -> AX cross section value.
1132  } else if (!isXB && mX > mMinAX && pow2(mX + mMinXB) < s) {
1133  double bANow = BHAD[iHadAtmp[iA]];
1134  double bAXNow = 2. * bANow + alP2 * log(1. / xi) ;
1135  sigNow += multVP[iA] * CONVERTSD * X[iProcVP[iA]]
1136  * BETA0[iHadAtmp[iA]] * exp(bAXNow * t) * (1. - xi)
1137  * ( 1. + cRes * sResAX / (sResAX + m2X) );
1138  }
1139  }
1140  return sigNow * epsWt;
1141 
1142  // gamma + gamma: loop over VMD states on side A and B.
1143  } else if (iProc == 14) {
1144  double sigNow = 0.;
1145  for (int iA = 0; iA < 3; ++iA)
1146  for (int iB = 0; iB < 3; ++iB){
1147 
1148  // Mass thresholds depend on VMD state.
1149  mMinXB = mAtmp[iA] + mMin0;
1150  mResXB = mAtmp[iA] + mRes0;
1151  sResXB = pow2(mResXB);
1152  mMinAX = mBtmp[iB] + mMin0;
1153  mResAX = mBtmp[iB] + mRes0;
1154  sResAX = pow2(mResAX);
1155 
1156  // Calculate and sum differential AB -> XB cross section value.
1157  if (isXB && mX > mMinXB && pow2(mX + mMinAX) < s) {
1158  double bBNow = BHAD[iHadBtmp[iB]];
1159  double bXBNow = 2. * bBNow + alP2 * log(1. / xi) ;
1160  sigNow += multVV[iA][iB] * CONVERTSD * X[iProcVV[iA][iB]]
1161  * BETA0[iHadBtmp[iB]] * exp(bXBNow * t) * (1. - xi)
1162  * ( 1. + cRes * sResXB / (sResXB + m2X) );
1163 
1164  // Calculate and sum differential AB -> AX cross section value.
1165  } else if (!isXB && mX > mMinAX && pow2(mX + mMinXB) < s) {
1166  double bANow = BHAD[iHadAtmp[iA]];
1167  double bAXNow = 2. * bANow + alP2 * log(1. / xi) ;
1168  sigNow += multVV[iA][iB] * CONVERTSD * X[iProcVV[iA][iB]]
1169  * BETA0[iHadAtmp[iA]] * exp(bAXNow * t) * (1. - xi)
1170  * ( 1. + cRes * sResAX / (sResAX + m2X) );
1171  }
1172  }
1173  return sigNow * epsWt;
1174 
1175  // No diffractive scattering for Pomeron + p.
1176  } else return 0.;
1177 
1178 }
1179 
1180 //--------------------------------------------------------------------------
1181 
1182 // Differential double diffractive cross section.
1183 
1184 double SigmaSaSDL::dsigmaDD(double xi1, double xi2, double t, int ) {
1185 
1186  // Calculate masses and check that they are above threshold. Optional weight.
1187  double m2X1 = xi1 * s;
1188  double mX1 = sqrt(m2X1);
1189  double m2X2 = xi2 * s;
1190  double mX2 = sqrt(m2X2);
1191  double epsWt = pow( m2X1 * m2X2, -epsSaS);
1192 
1193  // Ordinary hadron-hadron collisions.
1194  if (iProc < 13) {
1195  if (mX1 < mMinXB || mX2 < mMinAX) return 0.;
1196 
1197  // Return differential AB -> X1X2 cross section value.
1198  double bXX = alP2 * log( exp(4.) + s * s0 / (m2X1 * m2X2) );
1199  return CONVERTDD * BETA0[iHadA] * BETA0[iHadB] * exp(bXX * t)
1200  * (1. - pow2(mX1 + mX2) / s)
1201  * (s * SPROTON / (s * SPROTON + m2X1 * m2X2))
1202  * ( 1. + cRes * sResXB / (sResXB + m2X1) )
1203  * ( 1. + cRes * sResAX / (sResAX + m2X2) ) * epsWt;
1204 
1205  // gamma + p: loop over VMD states on side A.
1206  } else if (iProc == 13){
1207  double sigNow = 0.;
1208  for (int iA = 0; iA < 3; ++iA){
1209 
1210  // Mass thresholds depend on VMD state.
1211  mMinXB = mAtmp[iA] + mMin0;
1212  mResXB = mAtmp[iA] + mRes0;
1213  sResXB = pow2(mResXB);
1214  mMinAX = mBtmp[iA] + mMin0;
1215  mResAX = mBtmp[iA] + mRes0;
1216  sResAX = pow2(mResAX);
1217 
1218  // Calculate and sum differential AB -> X1X2 cross section value.
1219  if (mX1 > mMinXB && mX2 > mMinAX) {
1220  double bXX = alP2 * log( exp(4.) + s * s0 / (m2X1 * m2X2) );
1221  sigNow += multVP[iA] * CONVERTDD * BETA0[iHadAtmp[iA]]
1222  * BETA0[iHadBtmp[iA]] * exp(bXX * t)
1223  * (1. - pow2(mX1 + mX2) / s)
1224  * (s * SPROTON / (s * SPROTON + m2X1 * m2X2))
1225  * ( 1. + cRes * sResXB / (sResXB + m2X1) )
1226  * ( 1. + cRes * sResAX / (sResAX + m2X2) );
1227  }
1228  }
1229  return sigNow * epsWt;
1230 
1231  // gamma + gamma: loop over VMD states on side A and B.
1232  } else if (iProc == 14){
1233  double sigNow = 0.;
1234  for (int iA = 0; iA < 3; ++iA)
1235  for (int iB = 0; iB < 3; ++iB){
1236 
1237  // Mass thresholds depend on VMD state.
1238  mMinXB = mAtmp[iA] + mMin0;
1239  mResXB = mAtmp[iA] + mRes0;
1240  sResXB = pow2(mResXB);
1241  mMinAX = mBtmp[iB] + mMin0;
1242  mResAX = mBtmp[iB] + mRes0;
1243  sResAX = pow2(mResAX);
1244 
1245  // Calculate and sum differential AB -> X1X2 cross section value.
1246  if (mX1 > mMinXB && mX2 > mMinAX) {
1247  double bXX = alP2 * log( exp(4.) + s * s0 / (m2X1 * m2X2) );
1248  sigNow += multVV[iA][iB] * CONVERTDD * BETA0[iHadAtmp[iA]]
1249  * BETA0[iHadBtmp[iB]] * exp(bXX * t)
1250  * (1. - pow2(mX1 + mX2) / s)
1251  * (s * SPROTON / (s * SPROTON + m2X1 * m2X2))
1252  * ( 1. + cRes * sResXB / (sResXB + m2X1) )
1253  * ( 1. + cRes * sResAX / (sResAX + m2X2) );
1254  }
1255  }
1256  return sigNow * epsWt;
1257 
1258  // No diffractive scattering for Pomeron + p.
1259  } else return 0.;
1260 
1261 }
1262 
1263 //--------------------------------------------------------------------------
1264 
1265 // Differential central diffractive cross sections.
1266 // Simple extension of single diffraction, but without normalization.
1267 
1268 double SigmaSaSDL::dsigmaCD( double xi1, double xi2, double t1, double t2,
1269  int ) {
1270 
1271  // Ordinary hadron-hadron collisions.
1272  if (iProc < 13) {
1273 
1274  // Calculate mass; checked to be above threshold.
1275  double m2X = xi1 * xi2 * s;
1276  double mX = sqrt(m2X);
1277  if (mX < mMinCDnow || pow2(mX + mA + mB) > s) return 0.;
1278  wtNow = 1.;
1279 
1280  // Weight associated with differential AB -> AX cross section.
1281  double bAX = 2. * bA + alP2 * log(1. / xi1) ;
1282  wtNow *= CONVERTSD * X[iProc] * BETA0[iHadA] * exp(bAX * t1)
1283  * (1. - xi1);
1284 
1285  // Weight associated with differential AB -> XB cross section.
1286  double bXB = 2. * bB + alP2 * log(1. / xi2) ;
1287  wtNow *= CONVERTSD * X[iProc] * BETA0[iHadB] * exp(bXB * t2)
1288  * (1. - xi2);
1289 
1290  // Optional weight to tilt the mass spectrum.
1291  wtNow *= pow( m2X, -epsSaS);
1292 
1293  // Done.
1294  return wtNow;
1295 
1296  // No central diffraction for gamma + p, gamma + gamma or Pomeron + p.
1297  } else return 0.;
1298 
1299 }
1300 
1301 //--------------------------------------------------------------------------
1302 
1303 // Find beam combination.
1304 
1305 bool SigmaSaSDL::findBeamComb( int idAin, int idBin, double mAin,
1306  double mBin) {
1307 
1308  // Order flavour of incoming hadrons: idAbsA < idAbsB (restore later).
1309  idAbsA = abs(idAin);
1310  idAbsB = abs(idBin);
1311  mA = mAin;
1312  mB = mBin;
1313  swapped = false;
1314  if (idAbsA > idAbsB) {
1315  swap( idAbsA, idAbsB);
1316  swap( mA, mB);
1317  swapped = true;
1318  }
1319  sameSign = (idAin * idBin > 0);
1320 
1321  // Find process number.
1322  iProc = -1;
1323  if (idAbsA > 1000) {
1324  iProc = (sameSign) ? 0 : 1;
1325  } else if (idAbsA > 100 && idAbsB > 1000) {
1326  iProc = (sameSign) ? 2 : 3;
1327  if (idAbsA/10 == 11 || idAbsA/10 == 22) iProc = 4;
1328  if (idAbsA > 300) iProc = 5;
1329  if (idAbsA > 400) iProc = 6;
1330  if (idAbsA > 900) iProc = 15;
1331  } else if (idAbsA > 100) {
1332  iProc = 7;
1333  if (idAbsB > 300) iProc = 8;
1334  if (idAbsB > 400) iProc = 9;
1335  if (idAbsA > 300) iProc = 10;
1336  if (idAbsA > 300 && idAbsB > 400) iProc = 11;
1337  if (idAbsA > 400) iProc = 12;
1338  } else if (idAbsA == 22 || idAbsB == 22) {
1339  if (idAbsA == idAbsB) iProc = 14;
1340  if (idAbsA > 1000 || idAbsB > 1000) iProc = 13;
1341  }
1342  if (iProc == -1) return false;
1343 
1344  // Set up global variables.
1345  iHadA = IHADATABLE[iProc];
1346  iHadB = IHADBTABLE[iProc];
1347  bA = BHAD[iHadA];
1348  bB = BHAD[iHadB];
1349 
1350  // Set up VMD global variables for gamma + p.
1351  if (iProc == 13){
1352  for (int i = 0; i < 3; ++i){
1353  // VMD always on side a
1354  mAtmp[i] = VMDMASS[i];
1355  mBtmp[i] = mB;
1356  iHadAtmp[i] = (i == 2) ? 2 : 1;
1357  iHadBtmp[i] = 0;
1358  multVP[i] = ALPHAEM / GAMMAFAC[i];
1359  iProcVP[i] = (i == 2) ? 5 : 4;
1360  }
1361  }
1362 
1363  // Set up VMD global variables for gamma + gamma.
1364  if (iProc == 14){
1365  for (int i = 0; i < 3; ++i){
1366  mAtmp[i] = VMDMASS[i];
1367  mBtmp[i] = VMDMASS[i];
1368  iHadAtmp[i] = (i == 2) ? 2 : 1;
1369  iHadBtmp[i] = (i == 2) ? 2 : 1;
1370  for (int j = 0; j < 3; ++j){
1371  multVV[i][j] = pow2(ALPHAEM) / (GAMMAFAC[i] * GAMMAFAC[j]);
1372  iProcVV[i][j] = (i == 2 || j == 2) ? 8 : 7;
1373  if (i == 2 && j == 2) iProcVV[i][j] = 10;
1374  }
1375  }
1376  }
1377 
1378  // Done.
1379  return true ;
1380 
1381 }
1382 
1383 //==========================================================================
1384 
1385 // The SigmaMBR class.
1386 // It parametrizes pp/ppbar total, elastic and diffractive cross sections
1387 // according to the Minimum Bias Rockefeller (MBR) model.
1388 
1389 // Integration of MBR cross sections.
1390 const int SigmaMBR::NINTEG = 1000;
1391 const int SigmaMBR::NINTEG2 = 40;
1392 
1393 // Proton form factor appoximation with two exponents, [FFB1,FFB2] = GeV^-2.
1394 // Used for quick t-integration, while full expression used for t selection.
1395 const double SigmaMBR::FFA1 = 0.9;
1396 const double SigmaMBR::FFA2 = 0.1;
1397 const double SigmaMBR::FFB1 = 4.6;
1398 const double SigmaMBR::FFB2 = 0.6;
1399 
1400 //--------------------------------------------------------------------------
1401 
1402 // Initialize data members.
1403 
1404 void SigmaMBR::init( Info*, Settings& settings,
1405  ParticleData* particleDataPtrIn, Rndm*) {
1406 
1407  // Parameters for MBR model.
1408  eps = settings.parm("SigmaDiffractive:MBRepsilon");
1409  alph = settings.parm("SigmaDiffractive:MBRalpha");
1410  beta0gev = settings.parm("SigmaDiffractive:MBRbeta0");
1411  beta0mb = beta0gev * sqrt(HBARC2);
1412  sigma0mb = settings.parm("SigmaDiffractive:MBRsigma0");
1413  sigma0gev = sigma0mb/HBARC2;
1414  m2min = settings.parm("SigmaDiffractive:MBRm2Min");
1415  dyminSDflux = settings.parm("SigmaDiffractive:MBRdyminSDflux");
1416  dyminDDflux = settings.parm("SigmaDiffractive:MBRdyminDDflux");
1417  dyminCDflux = settings.parm("SigmaDiffractive:MBRdyminCDflux");
1418  dyminSD = settings.parm("SigmaDiffractive:MBRdyminSD");
1419  dyminDD = settings.parm("SigmaDiffractive:MBRdyminDD");
1420  dyminCD = settings.parm("SigmaDiffractive:MBRdyminCD") / 2.;
1421  dyminSigSD = settings.parm("SigmaDiffractive:MBRdyminSigSD");
1422  dyminSigDD = settings.parm("SigmaDiffractive:MBRdyminSigDD");
1423  dyminSigCD = settings.parm("SigmaDiffractive:MBRdyminSigCD") / sqrt(2.);
1424  a1 = FFA1;
1425  a2 = FFA2;
1426  b1 = FFB1;
1427  b2 = FFB2;
1428 
1429  // Initialize parameters for Coulomb corrections to elastic scattering.
1430  initCoulomb( settings, particleDataPtrIn);
1431 
1432  // No rho parameter implemented.
1433  rhoOwn = 0.;
1434 
1435 }
1436 
1437 //--------------------------------------------------------------------------
1438 
1439 // Total and elastic cross section.
1440 
1441 bool SigmaMBR::calcTotEl( int idAin, int idBin, double sIn, double , double) {
1442 
1443  // Save some input.
1444  idA = idAin;
1445  idB = idBin;
1446  s = sIn;
1447  isExpEl = true;
1448 
1449  // Total cross section and elastic/total parametrizations.
1450  double sCDF = pow2(1800.);
1451  double ratio = 0.;
1452  if (s <= sCDF) {
1453  double sign = (idA * idB > 0) ? 1. : -1.;
1454  sigTot = 16.79 * pow(s, 0.104) + 60.81 * pow(s, -0.32)
1455  - sign * 31.68 * pow(s, -0.54);
1456  ratio = 0.100 * pow(s, 0.06) + 0.421 * pow(s, -0.52)
1457  + sign * 0.160 * pow(s, -0.6);
1458  } else {
1459  double sigCDF = 80.03;
1460  double sF = pow2(22.);
1461  sigTot = sigCDF + ( pow2( log(s / sF)) - pow2( log(sCDF / sF)) )
1462  * M_PI / (3.7 / HBARC2);
1463  ratio = 0.066 + 0.0119 * log(s);
1464  }
1465  sigEl = sigTot * ratio;
1466  bEl = CONVERTEL * pow2(sigTot) / sigEl;
1467 
1468  // Possibly add Coulomb correction and interference.
1469  addCoulomb();
1470 
1471  // Done.
1472  return true ;
1473 
1474 }
1475 
1476 //--------------------------------------------------------------------------
1477 
1478 // Differential elastic cross section.
1479 
1480 double SigmaMBR::dsigmaEl( double t, bool useCoulomb, bool ) {
1481 
1482  // Hadronic contribution: simple exponential.
1483  double dsig = sigEl * bEl * exp(bEl * t);
1484 
1485  // Possibly add Coulomb contribution and interference.
1486  if (useCoulomb && hasCou) dsig += dsigmaElCoulomb(t);
1487 
1488  // Done.
1489  return dsig;
1490 
1491 }
1492 
1493 //--------------------------------------------------------------------------
1494 
1495 // Total diffractive cross sections, obtained from the ratio of
1496 // two integrals: the Regge cross section and the renormalized flux.
1497 
1498 bool SigmaMBR::calcDiff( int , int , double sIn, double , double ) {
1499 
1500  // Common variables,
1501  s = sIn;
1502  double cflux, csig, c1, step, f;
1503  double dymin0 = 0.;
1504 
1505  // Calculate SD cross section.
1506  double dymaxSD = log(s / m2min);
1507  cflux = pow2(beta0gev) / (16. * M_PI);
1508  csig = cflux * sigma0mb;
1509 
1510  // SD flux.
1511  c1 = cflux;
1512  double fluxsd = 0.;
1513  step = (dymaxSD - dyminSDflux) / NINTEG;
1514  for (int i = 0; i < NINTEG; ++i) {
1515  double dy = dyminSDflux + (i + 0.5) * step;
1516  f = exp(2. * eps * dy) * ( (a1 / (b1 + 2. * alph * dy))
1517  + (a2 / (b2 + 2. * alph * dy)) );
1518  f *= 0.5 * (1. + erf( (dy - dyminSD) / dyminSigSD));
1519  fluxsd = fluxsd + step * c1 * f;
1520  }
1521  if (fluxsd < 1.) fluxsd = 1.;
1522 
1523  // Regge SD cross section.
1524  c1 = csig * pow(s, eps);
1525  sigSD = 0.;
1526  sdpmax = 0.;
1527  step = (dymaxSD - dymin0) / NINTEG;
1528  for (int i = 0; i < NINTEG; ++i) {
1529  double dy = dymin0 + (i + 0.5) * step;
1530  f = exp(eps * dy) * ( (a1 / (b1 + 2. * alph * dy))
1531  + (a2 / (b2 + 2. * alph * dy)) );
1532  f *= 0.5 * (1. + erf( (dy - dyminSD) / dyminSigSD));
1533  if (f > sdpmax) sdpmax = f;
1534  sigSD = sigSD + step * c1 * f;
1535  }
1536  sdpmax *= 1.01;
1537  sigSD /= fluxsd;
1538 
1539  // Calculate DD cross section.
1540  // Note: dymaxDD = ln(s * s0 /mMin^4) with s0 = 1 GeV^2.
1541  double dymaxDD = log(s / pow2(m2min));
1542  cflux = sigma0gev / (16. * M_PI);
1543  csig = cflux * sigma0mb;
1544 
1545  // DD flux.
1546  c1 = cflux / (2. * alph);
1547  double fluxdd = 0.;
1548  step = (dymaxDD - dyminDDflux) / NINTEG;
1549  for (int i = 0; i < NINTEG; ++i) {
1550  double dy = dyminDDflux + (i + 0.5) * step;
1551  f = (dymaxDD - dy) * exp(2. * eps * dy)
1552  * ( exp(-2. * alph * dy * exp(-dy))
1553  - exp(-2. * alph * dy * exp(dy)) ) / dy;
1554  f *= 0.5 * (1. + erf( (dy - dyminDD) / dyminSigDD));
1555  fluxdd = fluxdd + step * c1 * f;
1556  }
1557  if (fluxdd < 1.) fluxdd = 1.;
1558 
1559  // Regge DD cross section.
1560  c1 = csig * pow(s, eps) / (2. * alph);
1561  ddpmax = 0.;
1562  sigDD = 0.;
1563  step = (dymaxDD - dymin0) / NINTEG;
1564  for (int i = 0; i < NINTEG; ++i) {
1565  double dy = dymin0 + (i + 0.5) * step;
1566  f = (dymaxDD - dy) * exp(eps * dy)
1567  * ( exp(-2. * alph * dy * exp(-dy))
1568  - exp(-2. * alph * dy * exp(dy)) ) / dy;
1569  f *= 0.5 * (1. + erf( (dy - dyminDD) / dyminSigDD));
1570  if (f > ddpmax) ddpmax = f;
1571  sigDD = sigDD + step * c1 * f;
1572  }
1573  ddpmax *= 1.01;
1574  sigDD /= fluxdd;
1575 
1576  // Calculate DPE (CD) cross section.
1577  double dymaxCD = log(s / m2min);
1578  cflux = pow4(beta0gev) / pow2(16. * M_PI);
1579  csig = cflux * pow2(sigma0mb / beta0mb);
1580  double dy1, dy2, f1, f2, step2;
1581 
1582  // DPE flux.
1583  c1 = cflux;
1584  double fluxdpe = 0.;
1585  step = (dymaxCD - dyminCDflux) / NINTEG;
1586  for (int i = 0; i < NINTEG; ++i) {
1587  double dy = dyminCDflux + (i + 0.5) * step;
1588  f = 0.;
1589  step2 = (dy - dyminCDflux) / NINTEG2;
1590  for (int j = 0; j < NINTEG2; ++j) {
1591  double yc = -0.5 * (dy - dyminCDflux) + (j + 0.5) * step2;
1592  dy1 = 0.5 * dy - yc;
1593  dy2 = 0.5 * dy + yc;
1594  f1 = exp(2. * eps * dy1) * ( (a1 / (b1 + 2. * alph * dy1))
1595  + (a2 / (b2 + 2. * alph * dy1)) );
1596  f2 = exp(2. * eps * dy2) * ( (a1 / (b1 + 2. * alph * dy2))
1597  + (a2 / (b2 + 2. * alph * dy2)) );
1598  f1 *= 0.5 * (1. + erf( (dy1 - dyminCD) / dyminSigCD) );
1599  f2 *= 0.5 * (1. + erf( (dy2 - dyminCD) / dyminSigCD) );
1600  f += f1 * f2 * step2;
1601  }
1602  fluxdpe += step * c1 * f;
1603  }
1604  if (fluxdpe < 1.) fluxdpe = 1.;
1605 
1606  // Regge DPE (CD) cross section.
1607  c1 = csig * pow(s, eps);
1608  sigCD = 0.;
1609  dpepmax = 0;
1610  step = (dymaxCD - dymin0) / NINTEG;
1611  for (int i = 0; i < NINTEG; ++i) {
1612  double dy = dymin0 + (i + 0.5) * step;
1613  f = 0.;
1614  step2 = (dy - dymin0) / NINTEG2;
1615  for (int j = 0; j < NINTEG2; ++j) {
1616  double yc = -0.5 * (dy - dymin0) + (j + 0.5) * step2;
1617  dy1 = 0.5 * dy - yc;
1618  dy2 = 0.5 * dy + yc;
1619  f1 = exp(eps * dy1) * ( (a1 / (b1 + 2. * alph * dy1))
1620  + (a2 / (b2 + 2. * alph * dy1)) );
1621  f2 = exp(eps * dy2) * ( (a1 / (b1 + 2. * alph * dy2))
1622  + (a2 / (b2 + 2. * alph * dy2)) );
1623  f1 *= 0.5 * (1. + erf( (dy1 - dyminCD) / dyminSigCD) );
1624  f2 *= 0.5 * (1. + erf( (dy2 - dyminCD) / dyminSigCD) );
1625  f += f1 * f2 * step2;
1626  }
1627  sigCD += step * c1 * f;
1628  if ( f > dpepmax) dpepmax = f;
1629  }
1630  dpepmax *= 1.01;
1631  sigCD /= fluxdpe;
1632 
1633  // Output to standard names. Done.
1634  sigXB = sigSD;
1635  sigAX = sigSD;
1636  sigXX = sigDD;
1637  sigAXB = sigCD;
1638  return true;
1639 
1640 }
1641 
1642 //--------------------------------------------------------------------------
1643 
1644 // Differential single diffractive cross sections, separated in xi and t steps.
1645 
1646 double SigmaMBR::dsigmaSD(double xi, double t, bool , int step) {
1647 
1648  // Rapidity gap size.
1649  double dy = -log(xi);
1650 
1651  // Step 1: evaluate cross section in xi, integrated over t.
1652  if (step == 1) {
1653  if (xi * s < m2min) return 0.;
1654  return exp(eps * dy)
1655  * ( (a1 / (b1 + 2. * alph * dy)) + (a2 / (b2 + 2. * alph * dy)) )
1656  * 0.5 * (1. + erf( (dy - dyminSD) / dyminSigSD));
1657 
1658  // Step 2: evaluate cross section in t for fixed xi.
1659  } else if (step == 2) {
1660  return pow2(pFormFac(t)) * exp(2. * alph * dy * t);
1661  }
1662 
1663  // Done.
1664  return 0.;
1665 
1666 }
1667 
1668 //--------------------------------------------------------------------------
1669 
1670 // Differential double diffractive cross sections, separated in xi and t steps.
1671 
1672 double SigmaMBR::dsigmaDD(double xi1, double xi2, double t, int step) {
1673 
1674  // Rapidity gap size (with implicit scale s_0 = 1 GeV).
1675  double dy = - log(xi1 * xi2 * s);
1676 
1677  // Step 1: evaluate cross section in xi1 and xi2, integrated over t.
1678  if (step == 1) {
1679  if (xi1 * s < m2min || xi2 * s < m2min || dy < 0.) return 0.;
1680  return exp(eps * dy) * ( exp(-2. * alph * dy * exp(-dy))
1681  - exp(-2. * alph * dy * exp(dy)) ) / dy
1682  * 0.5 * (1. + erf( (dy - dyminDD) / dyminSigDD));
1683 
1684  // Step 2: evaluate cross section in t for fixed xi1 and xi2.
1685  } else if (step == 2) {
1686  if ( t < -exp(dy) || t > -exp(-dy) ) return 0.;
1687  return exp(2. * alph * dy * t);
1688  }
1689 
1690  // Done.
1691  return 0.;
1692 
1693 }
1694 
1695 //--------------------------------------------------------------------------
1696 
1697 // Differential central diffractive cross section, separated in xi and t steps.
1698 
1699 double SigmaMBR::dsigmaCD(double xi1, double xi2, double t1, double t2,
1700  int step) {
1701 
1702  // Rapidity gap sizes.
1703  double dy1 = -log(xi1);
1704  double dy2 = -log(xi2);
1705 
1706  // Step 1: evaluate cross section in xi1 and xi2, integrated over t1 and t2.
1707  if (step == 1) {
1708  if (xi1 * xi2 * s < m2min) return 0.;
1709  double wt1 = exp(eps * dy1)
1710  * ( (a1 / (b1 + 2. * alph * dy1)) + (a2 / (b2 + 2. * alph * dy1)) )
1711  * 0.5 * (1. + erf( (dy1 - dyminCD) / dyminSigCD));
1712  double wt2 = exp(eps * dy2)
1713  * ( (a1 / (b1 + 2. * alph * dy2)) + (a2 / (b2 + 2. * alph * dy2)) )
1714  * 0.5 * (1. + erf( (dy2 - dyminCD) / dyminSigCD));
1715  return wt1 * wt2;
1716 
1717  // Step 2: evaluate cross section in t1 and t2 for fixed xi1 and xi2.
1718  } else if (step == 2) {
1719  return pow2(pFormFac(t1) * pFormFac(t2))
1720  * exp(2. * alph * (dy1 * t1 + dy2 * t2));
1721  }
1722 
1723  // Done.
1724  return 0.;
1725 
1726 }
1727 
1728 //==========================================================================
1729 
1730 // The SigmaABMST class.
1731 // It parametrizes pp/ppbar total, elastic and diffractive cross sections
1732 // according to Appleby, Barlow, Molson, Serluca and Toader (ABMST).
1733 
1734 //--------------------------------------------------------------------------
1735 
1736 // Definitions of static variables.
1737 
1738 // Parameters of parametrization: total and elastic cross sections.
1739 const double SigmaABMST::EPSI[] = {0.106231, 0.0972043, -0.510662, -0.302082};
1740 const double SigmaABMST::ALPP[] = { 0.0449029, 0.278037, 0.821595, 0.904556};
1741 const double SigmaABMST::NORM[] = { 228.359, 193.811, 518.686, 10.7843};
1742 const double SigmaABMST::SLOPE[] = { 8.38, 3.78, 1.36};
1743 const double SigmaABMST::FRACS[] = { 0.26, 0.56, 0.18};
1744 const double SigmaABMST::TRIG[] = { 0.3, 5.03};
1745 const double SigmaABMST::LAM2P = 0.521223;
1746 const double SigmaABMST::BAPPR[] = { 8.5, 1.086};
1747 const double SigmaABMST::LAM2FF = 0.71;
1748 
1749 // Parameters of parametrization: diffractive cross section.
1750 const double SigmaABMST::MRES[4] = {1.44, 1.52, 1.68, 2.19};
1751 const double SigmaABMST::WRES[4] = {0.325, 0.13, 0.14, 0.45};
1752 const double SigmaABMST::CRES[4] = {3.07, 0.4149, 1.108, 0.9515};
1753 const double SigmaABMST::AFAC[4] = {0.624529, 3.09088, 4.0, 177.217};
1754 const double SigmaABMST::BFAC[4] = {2.5835, 4.51487, 3.03392, 5.86474};
1755 const double SigmaABMST::CFAC[4] = {0., 0.186211, 10.0, 21.0029};
1756 const double SigmaABMST::PPP[4] = {0.4, 0.5, 0.4597, 5.7575};
1757 const double SigmaABMST::EPS[2] = {0.08, -0.4525};
1758 const double SigmaABMST::APR[2] = {0.25, 0.93};
1759 const double SigmaABMST::CPI[6] = {13.63, 0.0808, 31.79, -0.4525, 0.93, 14.4};
1760 const double SigmaABMST::CNST[5] = {1., -0.05, -0.25, -1.15, 13.5};
1761 
1762 // Parameters for integration over t and xi for SD, DD and CD.
1763 const int SigmaABMST::NPOINTSTSD = 200;
1764 const double SigmaABMST::XIDIVSD = 0.1;
1765 const double SigmaABMST::DXIRAWSD = 0.01;
1766 const double SigmaABMST::DLNXIRAWSD = 0.1;
1767 // For DD either Monte Carlo integration or (slower) grid one.
1768 const bool SigmaABMST::MCINTDD = true;
1769 const int SigmaABMST::NPOINTSTDD = 20;
1770 const double SigmaABMST::XIDIVDD = 0.1;
1771 const double SigmaABMST::DXIRAWDD = 0.02;
1772 const double SigmaABMST::DLNXIRAWDD = 0.1;
1773 const double SigmaABMST::BMCINTDD = 2.;
1774 const int SigmaABMST::NPOINTMCDD = 200000;
1775 // For CD only Monte Carlo integration.
1776 const double SigmaABMST::BMCINTCD = 2.;
1777 const int SigmaABMST::NPOINTMCCD = 200000;
1778 
1779 //--------------------------------------------------------------------------
1780 
1781 // Initialize data members.
1782 
1783 void SigmaABMST::init( Info* , Settings& settings, ParticleData* ,
1784  Rndm* rndmPtrIn) {
1785 
1786  // Save pointer.
1787  rndmPtr = rndmPtrIn;
1788 
1789  // Common setup.
1790  m2minp = pow2(MPROTON + MPION);
1791  m2minm = pow2(MPROTON - MPION);
1792 
1793  // Allow Couplomb corrections for elastic scattering.
1794  tryCoulomb = settings.flag("SigmaElastic:Coulomb");
1795  tAbsMin = settings.parm("SigmaElastic:tAbsMin");
1796 
1797  // Setup for single diffraction.
1798  modeSD = settings.mode("SigmaDiffractive:ABMSTmodeSD");
1799  multSD = settings.parm("SigmaDiffractive:ABMSTmultSD");
1800  powSD = settings.parm("SigmaDiffractive:ABMSTpowSD");
1801  s0 = (modeSD % 2 == 0) ? 4000. : 100.;
1802  c0 = (modeSD % 2 == 0) ? 0.6 : 0.012;
1803 
1804  // Setup for double diffraction.
1805  modeDD = settings.mode("SigmaDiffractive:ABMSTmodeDD");
1806  multDD = settings.parm("SigmaDiffractive:ABMSTmultDD");
1807  powDD = settings.parm("SigmaDiffractive:ABMSTpowDD");
1808 
1809  // Setup for central diffraction.
1810  modeCD = settings.mode("SigmaDiffractive:ABMSTmodeCD");
1811  multCD = settings.parm("SigmaDiffractive:ABMSTmultCD");
1812  powCD = settings.parm("SigmaDiffractive:ABMSTpowCD");
1813  mMinCDnow = settings.parm("SigmaDiffractive:ABMSTmMinCD");
1814 
1815  // Setup to dampen diffractive events with small rapidity gaps.
1816  dampenGap = settings.flag("SigmaDiffractive:ABMSTdampenGap");
1817  ygap = settings.parm("SigmaDiffractive:ABMSTygap");
1818  ypow = settings.parm("SigmaDiffractive:ABMSTypow");
1819  expPygap = exp(ypow * ygap);
1820 
1821  // Setup to force minimal t fall-off like exp(b_min * t).
1822  useBMin = settings.flag("SigmaDiffractive:ABMSTuseBMin");
1823  bMinSD = settings.parm("SigmaDiffractive:ABMSTbMinSD");
1824  bMinDD = settings.parm("SigmaDiffractive:ABMSTbMinDD");
1825  bMinCD = settings.parm("SigmaDiffractive:ABMSTbMinCD");
1826 
1827 }
1828 
1829 //--------------------------------------------------------------------------
1830 
1831 // Calculate total and (integrated) elastic cross sections.
1832 
1833 bool SigmaABMST::calcTotEl( int idAin, int idBin, double sIn, double ,
1834  double ) {
1835 
1836  // Find appropriate combination of incoming beams.
1837  idA = idAin;
1838  idB = idBin;
1839  ispp = (idA * idB > 0);
1840  s = sIn;
1841  facEl = HBARC2 / (16. * M_PI);
1842  isExpEl = false;
1843 
1844  // Total cross section and the rho parameter using the full expression.
1845  complex amp = amplitude( 0., false, false);
1846  sigTot = HBARC2 * imag(amp);
1847  rhoOwn = real(amp) / imag(amp);
1848 
1849  // Total elastic cross section, by integration in exp( MINSLOPEEL * t).
1850  sigEl = 0.;
1851  for (int i = 0; i < NPOINTS; ++i) {
1852  double y = (i + 0.5) / NPOINTS;
1853  double t = log(y) / MINSLOPEEL;
1854  sigEl += dsigmaEl( t, false) / y;
1855  }
1856  sigEl /= NPOINTS * MINSLOPEEL;
1857 
1858  // Approximate exponential slope.
1859  bEl = log( dsigmaEl( -TABSREF, false) / dsigmaEl( 0., false) ) / (-TABSREF);
1860 
1861  // Done if no Coulomb corrections.
1862  hasCou = tryCoulomb;
1863  if (abs(idAin) == 2112 || abs(idBin) == 2112) hasCou = false;
1864  sigTotCou = sigTot;
1865  sigElCou = sigEl;
1866  if (!hasCou) return true;
1867 
1868  // Reduce hadronic part of elastic cross section by tMin cut.
1869  sigElCou = sigEl * exp( - bEl * tAbsMin);
1870  if (tAbsMin < 0.9 * TABSMAX) {
1871 
1872  // Loop through t range according to dt/t^2.
1873  double sumCou = 0.;
1874  double xRel, tAbs;
1875  for (int i = 0; i < NPOINTS; ++i) {
1876  xRel = (i + 0.5) / NPOINTS;
1877  tAbs = tAbsMin * TABSMAX / (tAbsMin + xRel * (TABSMAX - tAbsMin));
1878 
1879  // Evaluate cross section difference between with and without Coulomb.
1880  sumCou += pow2(tAbs) * (dsigmaEl( -tAbs, true)
1881  - dsigmaEl( -tAbs, false));
1882  }
1883 
1884  // Include common factors to give new elastic and total cross sections.
1885  sigElCou += sumCou * (TABSMAX - tAbsMin)/ (tAbsMin * TABSMAX * NPOINTS);
1886  }
1887  sigTotCou = sigTot - sigEl + sigElCou;
1888 
1889  // Done.
1890  return true;
1891 
1892 }
1893 
1894 //--------------------------------------------------------------------------
1895 
1896 // The scattering amplitude, from which cross sections are derived.
1897 
1898 complex SigmaABMST::amplitude( double t, bool useCoulomb,
1899  bool onlyPomerons) {
1900 
1901  // Common values.
1902  double snu = s - 2. * SPROTON + 0.5 * t;
1903  double ampt = FRACS[0] * exp(SLOPE[0] * t) + FRACS[1] * exp(SLOPE[1] * t)
1904  + FRACS[2] * exp(SLOPE[2] * t);
1905  complex amp[6], l2p[4], ll2p[4], d2p[4][3];
1906 
1907  // Two Pomeron and even and odd Reggeon exchange.
1908  for (int i = 0; i < 4; ++i)
1909  amp[i] = ((i < 3) ? complex(-NORM[i], 0.) : complex( 0., NORM[i]))
1910  * ampt * sModAlp( ALPP[i] * snu, 1. + EPSI[i] + ALPP[i] * t);
1911 
1912  // Two-pomeron exchange.
1913  amp[4] = complex(0., 0.);
1914  for (int i = 0; i < 4; ++i) {
1915  l2p[i] = ALPP[i] * complex( log(ALPP[i] * snu), -0.5 * M_PI);
1916  ll2p[i] = (1. + EPSI[i]) * l2p[i] / ALPP[i];
1917  for (int k = 0; k < 3; ++k) d2p[i][k] = SLOPE[k] + l2p[i];
1918  }
1919  for (int i = 0; i < 4; ++i)
1920  for (int j = 0; j < 4; ++j)
1921  for (int k = 0; k < 3; ++k)
1922  for (int l = 0; l < 3; ++l) {
1923  complex part = NORM[i] * NORM[j] * exp( ll2p[i] + ll2p[j] )
1924  * exp( t * d2p[i][k] * d2p[j][l] / (d2p[i][k] + d2p[j][l]) )
1925  * FRACS[k] * FRACS[l] / (d2p[i][k] + d2p[j][l]);
1926  if (i == 3) part *= complex( 0., 1.);
1927  if (j == 3) part *= complex( 0., 1.);
1928  amp[4] += part;
1929  }
1930  amp[4] *= LAM2P * complex( 0., 1.) / (16. * M_PI * snu);
1931 
1932  // Triple-gluon exchange.
1933  amp[5] = sqrt(16. * M_PI / HBARC2) * TRIG[0] * ((t < -TRIG[1])
1934  ? 1. / pow4(t) : exp(4. + 4. * t / TRIG[1]) / pow4(TRIG[1]));
1935 
1936  // Add up contributions.
1937  complex ampSum = 0.;
1938  if (onlyPomerons) ampSum = (amp[0] + amp[1]) / snu;
1939  else ampSum = (amp[0] + amp[1] + amp[2] + ((ispp) ? -amp[3] : amp[3])
1940  + amp[4]) / snu + ((ispp) ? amp[5] : -amp[5]);
1941 
1942  // Optional Coulomb term. Must not be used for t = 0.
1943  if (useCoulomb && t < 0.) {
1944  double bApp = BAPPR[0] + BAPPR[1] * 0.5 * log(s);
1945  double phase = (GAMMAEUL + log( -0.5 * t * (bApp + 8. / LAM2FF)) +
1946  - 4. * t / LAM2FF * log(- 4. * t / LAM2FF)
1947  - 2. * t / LAM2FF) * (ispp ? 1. : -1.);
1948  complex ampCou = exp( complex( 0., ALPHAEM * phase) ) * 8. * M_PI
1949  * ALPHAEM * ampt / t;
1950  ampSum += (ispp) ? -ampCou : +ampCou;
1951  }
1952 
1953  // Done.
1954  return ampSum;
1955 
1956 }
1957 
1958 //--------------------------------------------------------------------------
1959 
1960 // Diffractive cross sections.
1961 
1962 bool SigmaABMST::calcDiff( int idAin , int idBin, double sIn, double ,
1963  double ) {
1964 
1965  // Find appropriate combination of incoming beams.
1966  idA = idAin;
1967  idB = idBin;
1968  ispp = (idA * idB > 0);
1969  s = sIn;
1970  facEl = HBARC2 / (16. * M_PI);
1971 
1972  // Total cross section needed for central diffraction.
1973  complex amp = amplitude( 0., false, true);
1974  sigTot = HBARC2 * imag(amp);
1975 
1976  // Single diffractive cross sections by grid integration.
1977  sigXB = dsigmaSDintXiT( 0., 1., -100., 0.);
1978  sigAX = sigXB;
1979 
1980  // Double diffractive cross section by Monte Carlo or grid integration.
1981  if (MCINTDD) sigXX = dsigmaDDintMC();
1982  else sigXX = dsigmaDDintXi1Xi2T( 0., 1., 0., 1., -100., 0.);
1983 
1984  // Central diffractive cross section by Monte Carlo integration.
1985  sigAXB = dsigmaCDintMC();
1986 
1987  // Done.
1988  return true;
1989 
1990 }
1991 
1992 //--------------------------------------------------------------------------
1993 
1994 // Variations on differential SD cross sections xi * dsigma / dxi dt.
1995 
1996 double SigmaABMST::dsigmaSD(double xi, double t, bool, int) {
1997 
1998  // Calculate core SD cross section.
1999  double dSigSD = dsigmaSDcore( xi, t);
2000 
2001  // Optionally require falloff at least like exp(bMin * t);
2002  if (useBMin && bMinSD > 0.) {
2003  double dSigSDmx = dsigmaSDcore( xi, -SPION) * exp(bMinSD * t);
2004  if (dSigSD > dSigSDmx) dSigSD = dSigSDmx;
2005  }
2006 
2007  // Optionally dampen with 1 / (1 + exp( -p * (y - y_gap))).
2008  if (dampenGap) dSigSD /= 1. + expPygap * pow( xi, ypow);
2009 
2010  // Optionally multiply by s-dependent factor.
2011  if (modeSD > 1) dSigSD *= multSD * pow( s / SPROTON, powSD);
2012 
2013  // Done.
2014  return dSigSD;
2015 
2016 }
2017 
2018 //--------------------------------------------------------------------------
2019 
2020 // Core differential single diffractive cross sections xi * dsigma / dxi dt.
2021 
2022 double SigmaABMST::dsigmaSDcore(double xi, double t) {
2023 
2024  // Calculate mass and check it is above threshold. ABMST valid for |t| < 4.
2025  double m2X = xi * s;
2026  if (m2X < m2minp) return 0.;
2027  double tAbs = abs(t);
2028  if (modeSD % 2 == 0 && tAbs > 4.) return 0.;
2029 
2030  // Some basic parameters.
2031  double tmp = 3. + c0 * pow2( log( s / s0) );
2032  double scaleFac = (s < s0) ? 1. : 3. / tmp;
2033  double mCut = (s < s0) ? 3 : tmp;
2034  if (modeSD % 2 == 0) {
2035  scaleFac = 1.;
2036  mCut = (s < s0) ? 3. : 3. + c0 * log(s/s0);
2037  }
2038  double m2Cut = pow2(mCut);
2039  double xiCut = m2Cut / s;
2040  double xiThr = m2minp / s;
2041  bool isHighM = (m2X > m2Cut);
2042  double xiNow = (isHighM) ? xi : xiCut;
2043  double m2XNow = xiNow * s;
2044 
2045  // Trajectories.
2046  for (int i = 0; i < 2; ++i) {
2047  alp0[i] = 1. + EPS[i];
2048  alpt[i] = alp0[i] + APR[i] * t;
2049  }
2050  alpt[2] = CPI[4] * (t - SPION);
2051 
2052  // Individual amplitudes: PPP remains unmodified.
2053  double ampPPP = pow(xiNow, alp0[0] - 2. * alpt[0]) * pow(s/CNST[0], EPS[0]);
2054  if (t > CNST[2]) ampPPP *= (PPP[0] + PPP[1] * t);
2055  else ampPPP *= (AFAC[0] * exp(BFAC[0] * t) + CFAC[0]) * t / (t + CNST[1]);
2056  if (t < CNST[3]) ampPPP *= (1. + PPP[2] * (tAbs + CNST[3])
2057  + PPP[3] * pow2(tAbs + CNST[3]));
2058  double ampPPR = pow(xiNow, alp0[1] - 2. * alpt[0]) * pow(s/CNST[0], EPS[1]);
2059  double ampRRP = pow(xiNow, alp0[0] - 2. * alpt[1]) * pow(s/CNST[0], EPS[0]);
2060  double ampRRR = pow(xiNow, alp0[1] - 2. * alpt[1]) * pow(s/CNST[0], EPS[1]);
2061 
2062  // ABMST uses exponential slope plus constant term in t.
2063  if (modeSD % 2 == 0) {
2064  ampPPR *= (AFAC[1] * exp(BFAC[1] * t) + CFAC[1]);
2065  ampRRP *= (AFAC[2] * exp(BFAC[2] * t) + CFAC[2]);
2066  ampRRR *= (AFAC[3] * exp(BFAC[3] * t) + CFAC[3]);
2067 
2068  // Modified t form eliminates constant term and has broader exponential.
2069  } else {
2070  double modX[2], modX2[2], expX[2], sumX[2];
2071  double modY[3], modY2[3], expY[3], sumY[3];
2072  double den[3], modB[3], apC[3];
2073  for (int ia = 0; ia < 2; ++ia) {
2074  modX[ia] = -2. * APR[ia] * log(xiNow);
2075  modX2[ia] = pow2(modX[ia]);
2076  expX[ia] = exp(-4. * modX[ia]);
2077  sumX[ia] = 4. * modX[ia] + 1.;
2078  }
2079  for (int ib = 0; ib < 3; ++ib) {
2080  int ix = (ib == 0) ? 0 : 1;
2081  modY[ib] = BFAC[ib+1] + modX[ix];
2082  modY2[ib] = pow2(modY[ib]);
2083  expY[ib] = exp(-4. * modY[ib]);
2084  sumY[ib] = 4. * modY[ib] + 1.;
2085  den[ib] = AFAC[ib+1] * modX2[ix] * (1. - expY[ib] * sumY[ib])
2086  + CFAC[ib+1] * modY2[ib] * (1. - expX[ix] * sumX[ix]);
2087  modB[ib] = (AFAC[ib+1] * modX2[ix] * modY[ib] * (1. - expY[ib])
2088  + CFAC[ib+1] * modY2[ib] * modX[ix] * (1. - expX[ix]))
2089  / den[ib] - modX[ix];
2090  apC[ib] = pow2( AFAC[ib+1] * modX[ix] * (1. - expY[ib])
2091  + CFAC[ib+1] * modY[ib] * (1. - expX[ix]) ) / den[ib];
2092  }
2093  ampPPR *= apC[0] * exp(modB[0] * t);
2094  ampRRP *= apC[1] * exp(modB[1] * t);
2095  ampRRR *= apC[2] * exp(modB[2] * t);
2096  }
2097 
2098  // Form factor and pion amplitude, remains unmodified.
2099  double fFac = pFormFac(t);
2100  double cnstPi = CPI[5]/(4. * M_PI) * tAbs/pow2(t - SPION) * pow2(fFac);
2101  double sigPi = CPI[0] * pow(m2XNow, CPI[1])
2102  + CPI[2] * pow(m2XNow, CPI[3]);
2103  double ampPi = cnstPi * sigPi * pow(xiNow, 1. - 2. * alpt[2]);
2104 
2105  // Total high-mass contribution. Done if at high masses.
2106  double ampHM = scaleFac * (ampPPP + ampPPR + ampRRP + ampRRR + ampPi);
2107  if (isHighM) return xi * ampHM;
2108 
2109  // Begin handling of low-mass region.
2110  double ampBkg = 0.;
2111  double ampRes = 0.;
2112  double ampMatch = 0.;
2113  double ampLM = 0.;
2114 
2115  // Resonance contribution.
2116  double qRef = sqrt( (m2X - m2minp) * (m2X - m2minm) / (4. * m2X) );
2117  for (int i = 0; i < 4; ++i) {
2118  double m2Res = pow2(MRES[i]);
2119  double qRes = sqrt( (m2Res - m2minp) * (m2Res - m2minm) / (4. * m2Res) );
2120  double mGam = MRES[i] * WRES[i] * pow( qRef / qRes, 2. * i + 3.)
2121  * pow( (1. + 5. * qRes) / (1. + 5. * qRef), i + 1.);
2122  ampRes += CRES[i] * mGam / (pow2(m2X - m2Res) + pow2(mGam));
2123  ampMatch += CRES[i] * mGam / (pow2(m2Cut - m2Res) + pow2(mGam));
2124  }
2125  ampRes *= exp( CNST[4] * (t - CNST[1]) ) / xi;
2126  ampMatch *= exp( CNST[4] * (t - CNST[1]) ) / xiNow
2127  * (xi - xiThr) / (xiNow - xiThr);
2128 
2129  // Background contribution.
2130  double dAmpPPP = ampPPP * (alp0[0] - 2. * alpt[0]) / xiNow;
2131  double dAmpPPR = ampPPR * (alp0[1] - 2. * alpt[0]) / xiNow;
2132  double dAmpRRP = ampRRP * (alp0[0] - 2. * alpt[1]) / xiNow;
2133  double dAmpRRR = ampRRR * (alp0[1] - 2. * alpt[1]) / xiNow;
2134  double dSigPi = CPI[0] * CPI[1] * pow(m2XNow, CPI[1] - 1.)
2135  + CPI[2] * CPI[3] * pow(m2XNow, CPI[3] - 1.);
2136  double dAmpPi = cnstPi * (sigPi * (1. - 2. * alpt[2])
2137  * pow(xiNow, -2. * alpt[2])
2138  + dSigPi * pow(xiNow, 1. - 2. * alpt[2]) );
2139  double dAmpHM = scaleFac * (dAmpPPP + dAmpPPR + dAmpRRP + dAmpRRR
2140  + dAmpPi);
2141 
2142  // Original background which vanishes quadratically at threshold.
2143  if (modeSD % 2 == 0) {
2144  double coeff1 = (dAmpHM * (xiCut - xiThr) - ampHM) / pow2(xiCut - xiThr);
2145  double coeff2 = 2. * ampHM / (xiCut - xiThr) - dAmpHM;
2146  ampBkg = coeff1 * pow2(xi - xiThr) + coeff2 * (xi - xiThr);
2147 
2148  // Modified form with combination of linear and quadratic description.
2149  } else {
2150  double xiAtM3 = 9. / s;
2151  double coeff1 = dAmpHM;
2152  double coeff2 = ampHM - coeff1 * (xiCut - xiThr);
2153  double coeff3 = - coeff2 / pow2(xiAtM3 - xiThr);
2154  double coeff4 = (2. * coeff1 * (xiAtM3 - xiThr) + 2. * coeff2)
2155  / (xiAtM3 - xiThr) - coeff1;
2156  ampBkg = (xi < xiAtM3) ? coeff3 * pow2(xi - xiThr)
2157  + coeff4 * (xi - xiThr) : coeff1 * (xi - xiThr) + coeff2;
2158  }
2159 
2160  // Return total low-mass contribution.
2161  ampLM = ampRes - ampMatch + ampBkg;
2162  return xi * ampLM;
2163 
2164 }
2165 
2166 //--------------------------------------------------------------------------
2167 
2168 // Single diffractive cross sections integrated over [t_min, t_max].
2169 
2170 double SigmaABMST::dsigmaSDintT(double xi, double tMinIn, double tMaxIn) {
2171 
2172  // Calculate t range. Check if range closed.
2173  double mu1 = SPROTON / s;
2174  double mu3 = xi;
2175  double rootv = (1. - 4. * mu1) * (pow2(1. - mu1 - mu3) - 4. * mu1 * mu3);
2176  if (rootv <= 0.) return 0.;
2177  double tMin = -0.5 * s * (1. - 3. * mu1 - mu3 + sqrt(rootv));
2178  double tMax = s * s * mu1 * pow2(mu3 - mu1) / tMin;
2179  tMin = max( tMin, tMinIn);
2180  tMax = min( tMax, tMaxIn);
2181  if (tMin >= tMax) return 0.;
2182 
2183  // Prepare integration.
2184  double slope = -0.5 * log(xi);
2185  double etMin = exp(slope * tMin);
2186  double etMax = exp(slope * tMax);
2187 
2188  // Do integration by uniform steps in exp(slope * t).
2189  double dsig = 0.;
2190  double etNow, tNow;
2191  for (int i = 0; i < NPOINTSTSD; ++i) {
2192  etNow = etMin + (i + 0.5) * (etMax - etMin) / NPOINTSTSD;
2193  tNow = log(etNow) / slope;
2194  dsig += dsigmaSD( xi, tNow, true, 0) / etNow;
2195  }
2196 
2197  // Normalize and done.
2198  dsig *= (etMax - etMin) / (NPOINTSTSD * slope);
2199  return dsig;
2200 
2201 }
2202 
2203 //--------------------------------------------------------------------------
2204 
2205 // Single diffractive cross section integrated over
2206 // [xi_min, xi_max] * [t_min, t_max].
2207 
2208 double SigmaABMST::dsigmaSDintXiT( double xiMinIn, double xiMaxIn,
2209  double tMinIn, double tMaxIn) {
2210 
2211  // Restrictions on xi range. Check it is not closed.
2212  double sig = 0.;
2213  double xiMin = max(xiMinIn, m2minp / s);
2214  double xiMax = min(xiMaxIn, 1.);
2215  if (xiMin >= xiMax) return 0.;
2216  double xiNow;
2217 
2218  // Integration in xi: check size of affected range and adjust dxi.
2219  if (xiMax > XIDIVSD) {
2220  double xiMinRng = max( XIDIVSD, xiMin);
2221  int nxi = 2 + (xiMax - xiMinRng) / DXIRAWSD;
2222  double dxi = (xiMax - xiMinRng) / nxi;
2223  for (int ixi = 0; ixi < nxi; ++ixi) {
2224  xiNow = xiMinRng + (ixi + 0.5) * dxi;
2225  sig += dxi * dsigmaSDintT( xiNow, tMinIn, tMaxIn) / xiNow;
2226  }
2227  }
2228 
2229  // Integration in ln(xi): check size of affected range and adjust dlnxi.
2230  if (xiMin < XIDIVSD) {
2231  double xiMaxRng = min( XIDIVSD, xiMax);
2232  int nlnxi = 2 + log( xiMaxRng / xiMin) / DLNXIRAWSD;
2233  double dlnxi = log( xiMaxRng / xiMin) / nlnxi;
2234  for (int ilnxi = 0; ilnxi < nlnxi; ++ilnxi) {
2235  xiNow = xiMin * exp( dlnxi * (ilnxi + 0.5));
2236  sig += dlnxi * dsigmaSDintT( xiNow, tMinIn, tMaxIn);
2237  }
2238  }
2239 
2240  // Done.
2241  return sig;
2242 
2243 }
2244 
2245 //--------------------------------------------------------------------------
2246 
2247 // Differential double diffractive cross section.
2248 
2249 double SigmaABMST::dsigmaDD(double xi1, double xi2, double t, int ) {
2250 
2251  // Calculate mass and check it is above threshold. ABMST valid for |t| < 4.
2252  double m2X1 = xi1 * s;
2253  double m2X2 = xi2 * s;
2254  if (m2X1 < m2minp || m2X2 < m2minp) return 0.;
2255  if (modeSD % 2 == 0 && abs(t) > 4.) return 0.;
2256 
2257  // dSigma_DD(x1, x2, t) = dSigma_SD(x1, t) * dSigma_SD(x2, t) / dSigma_el(t).
2258  // Elastic using only Pomerons.
2259  double dSigDD = dsigmaSDcore( xi1, t) * dsigmaSDcore( xi2, t)
2260  / dsigmaEl( t, false, true);
2261 
2262  // Minimal fall-off relative to t = 0 value.
2263  if (useBMin && bMinDD > 0.) {
2264  double dSigDDmx = dsigmaSDcore( xi1, -SPION) * dsigmaSDcore( xi2, -SPION)
2265  * exp(bMinDD * t) / dsigmaEl( 0., false, true);
2266  if (dSigDD > dSigDDmx) dSigDD = dSigDDmx;
2267  }
2268 
2269  // Optionally dampen with 1 / (1 + exp( -p * (y - y_gap))).
2270  if (dampenGap) dSigDD /= 1. + expPygap * pow( xi1 * xi2 * s / SPROTON, ypow);
2271 
2272  // Optionally multiply by s-dependent factor.
2273  if (modeDD == 1) dSigDD *= multDD * pow( s / SPROTON, powDD);
2274 
2275  // Done.
2276  return dSigDD;
2277 
2278 }
2279 
2280 //--------------------------------------------------------------------------
2281 
2282 // Double diffractive cross section integrated by Monte Carlo.
2283 
2284 double SigmaABMST::dsigmaDDintMC() {
2285 
2286  // Set up parameters of integration.
2287  double sigSum = 0.;
2288  double xiMin = m2minp / s;
2289  double mu1 = SPROTON / s;
2290  double xi1, xi2, t;
2291 
2292  // Integrate flat in dln(xi1) * dln(xi2) * exp(b_min t) dt.
2293  for (int iPoint = 0; iPoint < NPOINTMCDD; ++iPoint) {
2294  xi1 = pow( xiMin, rndmPtr->flat() );
2295  xi2 = pow( xiMin, rndmPtr->flat() );
2296  t = log( rndmPtr->flat() ) / BMCINTDD;
2297 
2298  // Check that point is inside phase space.
2299  if (sqrt(xi1) + sqrt(xi2) > 1.) continue;
2300  if (!tInRange( t/s, 1., mu1, mu1, xi1, xi2)) continue;
2301 
2302  // Calculate and add cross section.
2303  sigSum += dsigmaDD( xi1, xi2, t) * exp(-BMCINTDD * t);
2304  }
2305 
2306  // Normalize and done.
2307  sigSum *= pow2(log(xiMin)) / (BMCINTDD * NPOINTMCDD);
2308  return sigSum;
2309 
2310 }
2311 
2312 //--------------------------------------------------------------------------
2313 
2314 // Double diffractive cross section integrated over [t_min, t_max].
2315 
2316 double SigmaABMST::dsigmaDDintT(double xi1, double xi2, double tMinIn,
2317  double tMaxIn) {
2318 
2319  // Calculate t range. Check if range closed.
2320  double mu1 = SPROTON / s;
2321  pair<double,double> tRng = tRange(1., mu1, mu1, xi1, xi2);
2322  double tMin = max( s * tRng.first, tMinIn);
2323  double tMax = min( s * tRng.second, tMaxIn);
2324  if (tMin >= tMax) return 0.;
2325 
2326  // Prepare integration.
2327  double slope = BMCINTDD;
2328  double etMin = exp(slope * tMin);
2329  double etMax = exp(slope * tMax);
2330 
2331  // Do integration by uniform steps in exp(slope * t).
2332  double dsig = 0.;
2333  double etNow, tNow;
2334  for (int i = 0; i < NPOINTSTDD; ++i) {
2335  etNow = etMin + (i + 0.5) * (etMax - etMin) / NPOINTSTDD;
2336  tNow = log(etNow) / slope;
2337  dsig += dsigmaDD( xi1, xi2, tNow) / etNow;
2338  }
2339 
2340  // Normalize and done.
2341  dsig *= (etMax - etMin) / (NPOINTSTDD * slope);
2342  return dsig;
2343 
2344 }
2345 
2346 //--------------------------------------------------------------------------
2347 
2348 // Double diffractive cross section integrated over
2349 // [xi2_min, xi2_max] * [t_min, t_max].
2350 
2351 double SigmaABMST::dsigmaDDintXi2T( double xi1, double xi2MinIn,
2352  double xi2MaxIn, double tMinIn, double tMaxIn) {
2353 
2354  // Restrictions on xi2 range. Check it is not closed.
2355  double dsig = 0.;
2356  double xi2Min = max( xi2MinIn, m2minp / s);
2357  double xi2Max = min( xi2MaxIn, 1. + xi1 - 2. * sqrt(xi1));
2358  if (xi2Min >= xi2Max) return 0.;
2359  double xi2Now;
2360 
2361  // Integration in xi2: check size of affected range and adjust dxi2.
2362  if (xi2Max > XIDIVDD) {
2363  double xi2MinRng = max( XIDIVDD, xi2Min);
2364  int nxi2 = 2 + (xi2Max - xi2MinRng) / DXIRAWDD;
2365  double dxi2 = (xi2Max - xi2MinRng) / nxi2;
2366  for (int ixi2 = 0; ixi2 < nxi2; ++ixi2) {
2367  xi2Now = xi2MinRng + (ixi2 + 0.5) * dxi2;
2368  dsig += dxi2 * dsigmaDDintT( xi1, xi2Now, tMinIn, tMaxIn)
2369  / xi2Now;
2370  }
2371  }
2372 
2373  // Integration in ln(xi2): check size of affected range and adjust dlnxi2.
2374  if (xi2Min < XIDIVDD) {
2375  double xi2MaxRng = min( XIDIVDD, xi2Max);
2376  int nlnxi2 = 2 + log( xi2MaxRng / xi2Min) / DLNXIRAWDD;
2377  double dlnxi2 = log( xi2MaxRng / xi2Min) / nlnxi2;
2378  for (int ilnxi2 = 0; ilnxi2 < nlnxi2; ++ilnxi2) {
2379  xi2Now = xi2Min * exp( dlnxi2 * (ilnxi2 + 0.5));
2380  dsig += dlnxi2 * dsigmaDDintT( xi1, xi2Now, tMinIn, tMaxIn);
2381  }
2382  }
2383 
2384  // Done.
2385  return dsig;
2386 
2387 }
2388 
2389 //--------------------------------------------------------------------------
2390 
2391 // Double diffractive cross section integrated over
2392 // [xi1_min, xi1_max] * [xi2_min, xi2_max] * [t_min, t_max].
2393 
2394 double SigmaABMST::dsigmaDDintXi1Xi2T( double xi1MinIn, double xi1MaxIn,
2395  double xi2MinIn, double xi2MaxIn, double tMinIn, double tMaxIn) {
2396 
2397  // Restrictions on xi1 range. Check it is not closed.
2398  double dsig = 0.;
2399  double xi1Min = max( xi1MinIn, m2minp / s);
2400  double xi1Max = min( xi1MaxIn, 1.);
2401  if (xi1Min >= xi1Max) return 0.;
2402  double xi1Now;
2403 
2404  // Integration in xi1: check size of affected range and adjust dxi1.
2405  if (xi1Max > XIDIVDD) {
2406  double xi1MinRng = max( XIDIVDD, xi1Min);
2407  int nxi1 = 2 + (xi1Max - xi1MinRng) / DXIRAWDD;
2408  double dxi1 = (xi1Max - xi1MinRng) / nxi1;
2409  for (int ixi1 = 0; ixi1 < nxi1; ++ixi1) {
2410  xi1Now = xi1MinRng + (ixi1 + 0.5) * dxi1;
2411  dsig += dxi1 * dsigmaDDintXi2T( xi1Now, xi2MinIn, xi2MaxIn,
2412  tMinIn, tMaxIn) / xi1Now;
2413  }
2414  }
2415 
2416  // Integration in ln(xi1): check size of affected range and adjust dlnxi1.
2417  if (xi1Min < XIDIVDD) {
2418  double xi1MaxRng = min( XIDIVDD, xi1Max);
2419  int nlnxi1 = 2 + log( xi1MaxRng / xi1Min) / DLNXIRAWDD;
2420  double dlnxi1 = log( xi1MaxRng / xi1Min) / nlnxi1;
2421  for (int ilnxi1 = 0; ilnxi1 < nlnxi1; ++ilnxi1) {
2422  xi1Now = xi1Min * exp( dlnxi1 * (ilnxi1 + 0.5));
2423  dsig += dlnxi1 * dsigmaDDintXi2T( xi1Now, xi2MinIn, xi2MaxIn,
2424  tMinIn, tMaxIn);
2425  }
2426  }
2427 
2428  // Done.
2429  return dsig;
2430 
2431 }
2432 
2433 //--------------------------------------------------------------------------
2434 
2435 // Differential central diffractive cross section.
2436 
2437 double SigmaABMST::dsigmaCD(double xi1, double xi2, double t1, double t2,
2438  int ) {
2439 
2440  // ABMST valid for |t| < 4.
2441  if (modeSD % 2 == 0 && max( abs(t1), abs(t2)) > 4.) return 0.;
2442 
2443  // dSigma_CD(x1, x2, t1, t2)
2444  // = dSigma_SD(x1, t1) * dSigma_SD(x2, t2) / sigma_tot.
2445  double dSigCD = dsigmaSDcore( xi1, t1) * dsigmaSDcore( xi2, t2) / sigTot;
2446 
2447  // Minimal fall-off relative to t = 0 value.
2448  if (useBMin && bMinCD > 0.) {
2449  double dSigCDmx = dsigmaSDcore( xi1, -SPION) * dsigmaSDcore( xi2, -SPION)
2450  * exp(bMinCD * (t1 + t2)) / sigTot;
2451  if (dSigCD > dSigCDmx) dSigCD = dSigCDmx;
2452  }
2453 
2454  // Optionally dampen with 1 / (1 + exp( -p * (y - y_gap))) for both gaps.
2455  if (dampenGap) dSigCD /= (1. + expPygap * pow( xi1, ypow))
2456  * (1. + expPygap * pow( xi2, ypow));
2457 
2458  // Optionally multiply by s-dependent factor.
2459  if (modeCD == 1) dSigCD *= multCD * pow( s / SPROTON, powCD);
2460 
2461  // Done.
2462  return dSigCD;
2463 
2464 }
2465 
2466 //--------------------------------------------------------------------------
2467 
2468 // Central diffractive cross section integrated by Monte Carlo.
2469 
2470 double SigmaABMST::dsigmaCDintMC() {
2471 
2472  // Set up parameters of integration.
2473  double sigSum = 0.;
2474  double xiMin = m2minp / s;
2475  double xi1, xi2, t1, t2;
2476 
2477  // Integrate flat in dln(xi1) * exp(b_min t1) dt1 * (same with xi2, t2).
2478  for (int iPoint = 0; iPoint < NPOINTMCCD; ++iPoint) {
2479  xi1 = pow( xiMin, rndmPtr->flat() );
2480  xi2 = pow( xiMin, rndmPtr->flat() );
2481  t1 = log( rndmPtr->flat() ) / BMCINTCD;
2482  t2 = log( rndmPtr->flat() ) / BMCINTCD;
2483 
2484  // Check that point is inside phase space.
2485  if (xi1 * xi2 < xiMin) continue;
2486  if (xi1 * xi2 + 2. * xiMin > 1.) continue;
2487  if (!tInRange( t1, s, SPROTON, SPROTON, SPROTON, SPROTON + xi1 * s))
2488  continue;
2489  if (!tInRange( t1, s, SPROTON, SPROTON, SPROTON, SPROTON + xi2 * s))
2490  continue;
2491 
2492  // Calculate and add cross section.
2493  sigSum += dsigmaCD( xi1, xi2, t1, t2) * exp(-BMCINTCD * (t1 + t2));
2494  }
2495 
2496  // Normalize and done.
2497  sigSum *= pow2(log(xiMin) / BMCINTCD) / NPOINTMCCD;
2498  return sigSum;
2499 
2500 }
2501 
2502 //==========================================================================
2503 
2504 // The SigmaRPP class.
2505 // It parametrizes pp/ppbar total and elastic cross sections according to
2506 // the fit in Review of Particle Physics 2016.
2507 
2508 //--------------------------------------------------------------------------
2509 
2510 // Definitions of static variables.
2511 
2512 // Parameters of parametrization.
2513 const double SigmaRPP::EPS1[] = { 1., 0.614, 0.444, 1., 1., 1.};
2514 const double SigmaRPP::ALPP[] = { 0.151, 0.8, 0.8, 0.947};
2515 const double SigmaRPP::NORM[] = { 0.2478, 0.0078, 11.22, -0.150, 148.4, -26.6,
2516  -1.5, -0.0441, 0., 0.686, -3.82, -8.60, 64.1, 99.1, -58.0, 9.5};
2517 const double SigmaRPP::BRPP[] = { 3.592, 0.622, 5.44, 0.205, 5.643, 1.92,
2518  0.41, 0., 0., 3.013, 2.572, 12.25, 2.611, 11.28, 1.27};
2519 const double SigmaRPP::KRPP[] = { 0.3076, 0.0998, 1.678, 0.190, -26.1};
2520 const double SigmaRPP::LAM2FF = 0.71;
2521 
2522 //--------------------------------------------------------------------------
2523 
2524 // Calculate total and (integrated) elastic cross sections.
2525 
2526 bool SigmaRPP::calcTotEl( int idAin, int idBin, double sIn, double ,
2527  double ) {
2528 
2529  // Find appropriate combination of incoming beams.
2530  idA = idAin;
2531  idB = idBin;
2532  ispp = (idA * idB > 0);
2533  s = sIn;
2534  facEl = CONVERTEL / (s * (s - 4. * SPROTON));
2535  isExpEl = false;
2536 
2537  // Total cross section and the rho parameter.
2538  complex amp = amplitude( 0., false);
2539  sigTot = imag(amp) / sqrt(s * ( s - 4. * SPROTON));
2540  rhoOwn = real(amp) / imag(amp);
2541 
2542  // Total elastic cross section, by integration in exp( MINSLOPEEL * t).
2543  sigEl = 0.;
2544  for (int i = 0; i < NPOINTS; ++i) {
2545  double y = (i + 0.5) / NPOINTS;
2546  double t = log(y) / MINSLOPEEL;
2547  sigEl += dsigmaEl( t, false) / y;
2548  }
2549  sigEl /= NPOINTS * MINSLOPEEL;
2550 
2551  // Approximate exponential slope.
2552  bEl = log( dsigmaEl( -TABSREF, false) / dsigmaEl( 0., false) ) / (-TABSREF);
2553 
2554  // Done if no Coulomb corrections.
2555  hasCou = tryCoulomb;
2556  if (abs(idAin) == 2112 || abs(idBin) == 2112) hasCou = false;
2557  sigTotCou = sigTot;
2558  sigElCou = sigEl;
2559  if (!hasCou) return true;
2560 
2561  // Reduce hadronic part of elastic cross section by tMin cut.
2562  sigElCou = sigEl * exp( - bEl * tAbsMin);
2563  if (tAbsMin < 0.9 * TABSMAX) {
2564 
2565  // Loop through t range according to dt/t^2.
2566  double sumCou = 0.;
2567  double xRel, tAbs;
2568  for (int i = 0; i < NPOINTS; ++i) {
2569  xRel = (i + 0.5) / NPOINTS;
2570  tAbs = tAbsMin * TABSMAX / (tAbsMin + xRel * (TABSMAX - tAbsMin));
2571 
2572  // Evaluate cross section difference between with and without Coulomb.
2573  sumCou += pow2(tAbs) * (dsigmaEl( -tAbs, true)
2574  - dsigmaEl( -tAbs, false));
2575  }
2576 
2577  // Include common factors to give new elastic and total cross sections.
2578  sigElCou += sumCou * (TABSMAX - tAbsMin)/ (tAbsMin * TABSMAX * NPOINTS);
2579  }
2580  sigTotCou = sigTot - sigEl + sigElCou;
2581 
2582  // Done.
2583  return true;
2584 
2585 }
2586 
2587 //--------------------------------------------------------------------------
2588 
2589 // Amplitude.
2590 
2591 complex SigmaRPP::amplitude( double t, bool useCoulomb) {
2592 
2593  // Modified s-related values.
2594  double shat = s - 2. * SPROTON + 0.5 * t;
2595  complex stlog = complex( log(shat), -0.5 * M_PI);
2596  complex taut = sqrt(abs(t)) * stlog;
2597 
2598  // Trajectories.
2599  double aP = EPS1[0] + ALPP[0] * t;
2600  double aRpos = EPS1[1] + ALPP[1] * t;
2601  double aRneg = EPS1[2] + ALPP[2] * t;
2602  double aO = EPS1[3] + ALPP[3] * t;
2603  double aOP = EPS1[4] + ALPP[0] * ALPP[3] * t / (ALPP[0] + ALPP[3]);
2604  double aPP = EPS1[5] + 0.5 * ALPP[0] * t;
2605  double aRPpos = EPS1[1] + ALPP[0] * ALPP[1] * t / (ALPP[0] + ALPP[1]);
2606  double aRPneg = EPS1[2] + ALPP[0] * ALPP[2] * t / (ALPP[0] + ALPP[2]);
2607 
2608  // Even terms.
2609  complex besArg = KRPP[0] * taut;
2610  complex besJ0n = besJ0(besArg);
2611  complex besJ1n = besJ1(besArg);
2612  complex besRat = (abs(besArg) < 0.01) ? 1. : 2. * besJ1n / besArg;
2613  complex fPosH = complex( 0., shat) * (NORM[0] * besRat
2614  * exp(BRPP[0] * t) * stlog * stlog
2615  + NORM[1] * besJ0n * exp(BRPP[1] * t) * stlog
2616  + NORM[2] * (besJ0n - besArg * besJ1n) * exp(BRPP[2] * t));
2617  complex fPosP = -NORM[3] * exp(BRPP[3] * t) * sModAlp( shat, aP);
2618  complex fPosPP = (-NORM[4] / stlog) * exp(BRPP[4] * t) * sModAlp( shat, aPP);
2619  complex fPosR = -NORM[5] * exp(BRPP[5] * t) * sModAlp( shat, aRpos);
2620  complex fPosRP = t * (NORM[6] / stlog) * exp(BRPP[6] * t)
2621  * sModAlp( shat, aRPpos);
2622  complex nPos = complex( 0., -shat) * NORM[7] * stlog * t
2623  * pow( 1. - t / KRPP[2], -5.);
2624  complex fPos = fPosH + fPosP + fPosPP + fPosR + fPosRP + nPos;
2625 
2626  // Odd terms.
2627  complex fNegMO = shat * (NORM[9] * cos(KRPP[1] * taut) * exp(BRPP[9] * t)
2628  * stlog + NORM[10] * exp(BRPP[10] * t));
2629  complex fNegO = complex( 0., NORM[11]) * exp(BRPP[11] * t)
2630  * sModAlp( shat, aO) * (1. + KRPP[4] * t);
2631  complex fNegOP = (complex( 0., NORM[12]) / stlog) * exp(BRPP[12] * t)
2632  * sModAlp( shat, aOP);
2633  complex fNegR = complex( 0., -NORM[13]) * exp(BRPP[13] * t)
2634  * sModAlp( shat, aRneg);
2635  complex fNegRP = t * (complex( 0., -NORM[14]) / stlog) * exp(BRPP[14] * t)
2636  * sModAlp( shat, aRPneg);
2637  complex nNeg = -shat * NORM[15] * stlog * t * pow( 1. - t / KRPP[3], -5.);
2638  complex fNeg = fNegMO + fNegO + fNegOP + fNegR + fNegRP + nNeg;
2639 
2640  // Combine nuclear part.
2641  complex ampSum = ispp ? fPos + fNeg : fPos - fNeg;
2642 
2643  // Optional Coulomb term. Must not be used for t = 0.
2644  complex ampCou = 0.;
2645  if (useCoulomb && t < 0.) {
2646  double bAppr = imag(ampSum) / ( sqrt(s * ( s - 4. * SPROTON))
2647  * 4. * M_PI * HBARC2 );
2648  double phase = (log( -0.5 * t * (bAppr + 8. / LAM2FF)) + GAMMAEUL
2649  - 4. * t / LAM2FF * log(- 4. * t / LAM2FF)
2650  - 2. * t / LAM2FF) * (ispp ? 1. : -1.);
2651  ampCou = exp( complex( 0., ALPHAEM * phase) ) * 8. * M_PI * HBARC2
2652  * ALPHAEM * s / t * pow(1 - t / LAM2FF, -4.);
2653  }
2654 
2655  // Combine and return.
2656  return ispp ? ampSum - ampCou : ampSum + ampCou;
2657 
2658 }
2659 
2660 //--------------------------------------------------------------------------
2661 
2662 // Complex Bessel functions J0 and J1.
2663 
2664 complex SigmaRPP::besJ0( complex x) {
2665  int mMax = 5. + 5. * abs(x);
2666  complex z = 0.25 * x * x;
2667  complex term = 1.;
2668  complex sum = term;
2669  for (int m = 1; m < mMax; ++m) {
2670  term *= - z / double(m * m);
2671  sum += term;
2672  }
2673  return sum;
2674 }
2675 
2676 complex SigmaRPP::besJ1( complex x) {
2677  int mMax = 5. + 5. * abs(x);
2678  complex z = 0.25 * x * x;
2679  complex term = 0.5 * x;
2680  complex sum = term;
2681  for (int m = 1; m < mMax; ++m) {
2682  term *= - z / double(m * (m+1));
2683  sum += term;
2684  }
2685  return sum;
2686 }
2687 
2688 //==========================================================================
2689 
2690 } // end namespace Pythia8