StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
SigmaTotal.h
1 // SigmaTotal.h is a part of the PYTHIA event generator.
2 // Copyright (C) 2020 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 // This file contains classes for cross section parametrizations.
7 // SigmaTotAux : base class for the different parametrizations.
8 // SigmaTotal : top-level class, making use of the classes below.
9 // SigmaTotOwn : user-set values.
10 // SigmaSaSDL : Schuler and Sjostrand, based on Donnachie and Landshoff.
11 // SigmaMBR : Min Bias Rockefeller.
12 // SigmaABMST : Appleby, Barlow, Molson, Serluca and Toader.
13 // SigmaRPP : Review of Particle Physics 2014.
14 
15 #ifndef Pythia8_SigmaTotal_H
16 #define Pythia8_SigmaTotal_H
17 
18 #include "Pythia8/Info.h"
19 #include "Pythia8/ParticleData.h"
20 #include "Pythia8/PhysicsBase.h"
21 #include "Pythia8/PythiaStdlib.h"
22 #include "Pythia8/PythiaComplex.h"
23 #include "Pythia8/Settings.h"
24 
25 namespace Pythia8 {
26 
27 //==========================================================================
28 
29 // The SigmaTotAux is base class for the different parametrizations
30 // made accessible via SigmaTotal. Many variables are public, since
31 // they can only be accessed via the SigmaTotal methods anyway.
32 
33 class SigmaTotAux {
34 
35 public:
36 
37  // Constructor.
38  SigmaTotAux() : isExpEl(), hasCou(), sigTot(), rhoOwn(), sigEl(), bEl(),
39  sigTotCou(), sigElCou(), sigXB(), sigAX(), sigXX(), sigAXB(), idA(),
40  idB(), tryCoulomb(), chgSgn(), tAbsMin(), lambda(), phaseCst(),
41  particleDataPtr(), rndmPtr() {}
42 
43  // Destructor.
44  virtual ~SigmaTotAux() {};
45 
46  // Store pointers and initialize data members.
47  virtual void init(Info*) = 0;
48 
49  // Calculate integrated total/elastic cross sections.
50  // Usage: calcTotEl( idAin, idBin, sIn, mAin, mBin).
51  virtual bool calcTotEl( int , int , double , double , double ) {return true;}
52 
53  // Store total and elastic cross section properties.
54  bool isExpEl, hasCou;
55  double sigTot, rhoOwn, sigEl, bEl, sigTotCou, sigElCou;
56 
57  // Differential elastic cross section, d(sigma_el) / dt.
58  virtual double dsigmaEl( double , bool = false, bool = true) {return 0.;}
59 
60  // Calculate integrated diffractive cross sections.
61  // Usage: calcDiff( idAin, idBin, sIn, mAin, mBin).
62  virtual bool calcDiff( int , int , double , double , double ) {return true;}
63 
64  // Store diffractive cross sections.
65  double sigXB, sigAX, sigXX, sigAXB;
66 
67  // Differential single diffractive cross section,
68  // xi * d(sigma_SD) / (dxi dt).
69  virtual double dsigmaSD( double , double, bool = true, int = 0) {return 0.;}
70 
71  // Possibility to separate xi and t choices for diffraction.
72  virtual bool splitDiff() {return false;}
73 
74  // Differential double diffractive cross section,
75  // xi1 * xi2 * d(sigma_DD) / (dxi1 dxi2 dt).
76  virtual double dsigmaDD( double , double , double , int = 0) {return 0.;}
77 
78  // Differential central diffractive cross section,
79  // xi1 * xi2 * d(sigma_CD) / (dxi1 dxi2 dt1 dt2).
80  virtual double dsigmaCD( double , double , double , double, int ) {
81  return 0.;}
82 
83  // Minimal central diffractive mass.
84  virtual double mMinCD() {return 1.;}
85 
86  // Standard methods to find t range of a 2 -> 2 process
87  // and to check whether a given t value is in that range.
88  pair<double,double> tRange( double sIn, double s1In, double s2In,
89  double s3In, double s4In) {
90  double lambda12 = pow2( sIn - s1In - s2In) - 4. * s1In * s2In;
91  double lambda34 = pow2( sIn - s3In - s4In) - 4. * s3In * s4In;
92  if (lambda12 < 0. || lambda34 < 0.) return make_pair( 0., 0.);
93  double tLow = -0.5 * (sIn - (s1In + s2In + s3In + s4In) + (s1In - s2In)
94  * (s3In - s4In) / sIn + sqrtpos(lambda12 * lambda34) / sIn);
95  double tUpp = ( (s3In - s1In) * (s4In - s2In) + (s1In + s4In - s2In - s3In)
96  * (s1In * s4In - s2In * s3In) / sIn ) / tLow;
97  return make_pair( tLow, tUpp); }
98  bool tInRange( double tIn, double sIn, double s1In, double s2In,
99  double s3In, double s4In) {
100  pair<double, double> tRng = tRange( sIn, s1In, s2In, s3In, s4In);
101  return (tIn > tRng.first && tIn < tRng.second); }
102 
103  // Commonly used proton form factor.
104  double pFormFac(double tIn) {return (4. * SPROTON - 2.79 * tIn)
105  / ((4. * SPROTON - tIn) * pow2(1. - tIn / 0.71)); }
106 
107 protected:
108 
109  // Constants: could only be changed in the code itself.
110  static const int NPOINTS;
111  static const double ALPHAEM, CONVERTEL, MPROTON, SPROTON, MPION, SPION,
112  GAMMAEUL, TABSREF, TABSMAX, MINSLOPEEL;
113 
114  // Initialization data, normally only set once.
115  int idA, idB;
116 
117  // Add Coulomb corrections to the elastic cross section.
118  bool tryCoulomb;
119  double chgSgn, tAbsMin, lambda, phaseCst;
120  virtual bool initCoulomb(Settings& settings,
121  ParticleData* particleDataPtrIn);
122  virtual bool addCoulomb();
123  virtual double dsigmaElCoulomb(double t);
124 
125  // Pointer to the particle data table.
126  ParticleData* particleDataPtr;
127 
128  // Pointer to the random number generator.
129  Rndm* rndmPtr;
130 
131 };
132 
133 //==========================================================================
134 
135 // The SigmaTotal class contains parametrizations of total, elastic and
136 // diffractive cross sections, and of the respective slope parameter.
137 
138 class SigmaTotal : public PhysicsBase {
139 
140 public:
141 
142  // Constructor.
143  SigmaTotal() : isCalc(false), ispp(), modeTotEl(), modeTotElNow(),
144  modeDiff(), modeDiffNow(), idAbsA(), idAbsB(), s(), sigND(),
145  sigTotElPtr(NULL), sigDiffPtr(NULL) {};
146 
147  // Destructor.
148  virtual ~SigmaTotal() { if (sigTotElPtr) delete sigTotElPtr;
149  if (sigDiffPtr) delete sigDiffPtr; }
150 
151  // Store pointers and initialize data members.
152  void init();
153 
154  // Calculate, or recalculate for new beams or new energy.
155  bool calc( int idA, int idB, double eCM);
156 
157  // Confirm that initialization worked.
158  bool hasSigmaTot() {return isCalc;}
159 
160  // Total and integrated elastic cross sections.
161  double sigmaTot() {return sigTotElPtr->sigTotCou;}
162  double rho() {return sigTotElPtr->rhoOwn;}
163  double sigmaEl() {return sigTotElPtr->sigElCou;}
164  bool bElIsExp() {return sigTotElPtr->isExpEl;}
165  double bSlopeEl() {return sigTotElPtr->bEl;}
166  bool hasCoulomb() {return sigTotElPtr->hasCou;}
167 
168  // Total elastic cross section.
169  bool calcTotEl( int idAin, int idBin, double sIn, double mAin, double mBin) {
170  return sigTotElPtr->calcTotEl( idAin, idBin, sIn, mAin, mBin); }
171 
172  // Differential elastic cross section.
173  double dsigmaEl( double t, bool useCoulomb = false,
174  bool onlyPomerons = false) {
175  return sigTotElPtr->dsigmaEl( t, useCoulomb, onlyPomerons); }
176 
177  // Integrated diffractive cross sections.
178  double sigmaXB() const {return sigDiffPtr->sigXB;}
179  double sigmaAX() const {return sigDiffPtr->sigAX;}
180  double sigmaXX() const {return sigDiffPtr->sigXX;}
181  double sigmaAXB() const {return sigDiffPtr->sigAXB;}
182  double sigmaND() const {return sigND;}
183 
184  // Differential single diffractive cross section.
185  double dsigmaSD( double xi, double t, bool isXB = true, int step = 0) {
186  return sigDiffPtr->dsigmaSD( xi, t, isXB, step); }
187 
188  // Possibility to separate xi and t choices for diffraction.
189  virtual bool splitDiff() {return sigDiffPtr->splitDiff();}
190 
191  // Differential double diffractive cross section.
192  double dsigmaDD( double xi1, double xi2, double t, int step = 0) {
193  return sigDiffPtr->dsigmaDD( xi1, xi2, t, step); }
194 
195  // Differential central diffractive cross section.
196  double dsigmaCD( double xi1, double xi2, double t1, double t2, int step = 0)
197  { return sigDiffPtr->dsigmaCD( xi1, xi2, t1, t2, step); }
198 
199  // Minimal central diffractive mass.
200  double mMinCD() {return sigDiffPtr->mMinCD();}
201 
202  // Sample the VMD states for resolved photons.
203  void chooseVMDstates(int idA, int idB, double eCM, int processCode);
204 
205  // Standard methods to find t range of a 2 -> 2 process
206  // and to check whether a given t value is in that range.
207  pair<double,double> tRange( double sIn, double s1In, double s2In,
208  double s3In, double s4In) {
209  return sigDiffPtr->tRange( sIn, s1In, s2In, s3In, s4In); }
210  bool tInRange( double tIn, double sIn, double s1In, double s2In,
211  double s3In, double s4In) {
212  return sigDiffPtr->tInRange( tIn, sIn, s1In, s2In, s3In, s4In); }
213 
214 private:
215 
216  // Constants: could only be changed in the code itself.
217  static const double MMIN;
218 
219  // Initialization data, normally only set once.
220  bool isCalc, ispp;
221 
222  int modeTotEl, modeTotElNow, modeDiff, modeDiffNow, idAbsA, idAbsB;
223  double s, sigND;
224 
225  // Pointer to class that handles total and elastic cross sections.
226  SigmaTotAux* sigTotElPtr;
227 
228  // Pointer to class that handles diffractive cross sections.
229  SigmaTotAux* sigDiffPtr;
230 
231 };
232 
233 //==========================================================================
234 
235 // The SigmaTotOwn class parametrizes total, elastic and diffractive
236 // cross sections by user settings.
237 
238 class SigmaTotOwn : public SigmaTotAux {
239 
240 public:
241 
242  // Constructor.
243  SigmaTotOwn() : dampenGap(), pomFlux(), s(), a0(), ap(), b0(), A1(), A2(),
244  A3(), a1(), a2(), a3(), bMinDD(), ygap(), ypow(), expPygap(), mMinCDnow(),
245  wtNow(), yNow(), yNow1(), yNow2(), b(), b1(), b2(), Q(), Q1(),
246  Q2() {};
247 
248  // Store pointers and initialize data members.
249  virtual void init(Info* infoPtrIn);
250 
251  // Calculate integrated total/elastic cross sections.
252  virtual bool calcTotEl( int idAin, int idBin, double , double , double);
253 
254  // Differential elastic cross section.
255  virtual double dsigmaEl( double t, bool useCoulomb = false, bool = true);
256 
257  // Calculate integrated diffractive cross sections.
258  virtual bool calcDiff( int , int , double sIn, double , double ) {
259  s = sIn; return true;}
260 
261  // Differential single diffractive cross section.
262  virtual double dsigmaSD( double xi, double t, bool = true, int = 0);
263 
264  // Differential double diffractive cross section.
265  virtual double dsigmaDD( double xi1, double xi2, double t, int = 0);
266 
267  // Differential central diffractive cross section.
268  virtual double dsigmaCD( double xi1, double xi2, double t1, double t2,
269  int = 0);
270 
271  // Minimal central diffractive mass.
272  virtual double mMinCD() {return mMinCDnow;}
273 
274 private:
275 
276  // Initialization data, normally only set once.
277  bool dampenGap;
278  int pomFlux;
279  double s, a0, ap, b0, A1, A2, A3, a1, a2, a3, bMinDD, ygap, ypow, expPygap,
280  mMinCDnow, wtNow, yNow, yNow1, yNow2, b, b1, b2, Q, Q1, Q2;
281 
282 };
283 
284 //==========================================================================
285 
286 // The SigmaSaSDL class parametrizes total, elastic and diffractive
287 // cross sections according to Schuler and Sjostrand, starting from
288 // Donnachie and Landshoff.
289 
290 class SigmaSaSDL : public SigmaTotAux {
291 
292 public:
293 
294  // Constructor.
295  SigmaSaSDL() : doDampen(), zeroAXB(), swapped(), sameSign(), idAbsA(),
296  idAbsB(), iProc(), iHadA(), iHadB(), iHadAtmp(), iHadBtmp(), iProcVP(),
297  iProcVV(), s(), mA(), mB(), bA(), bB(), maxXBOwn(), maxAXOwn(), maxXXOwn(),
298  maxAXBOwn(), epsSaS(), sigmaPomP(), mPomP(), pPomP(), sigAXB2TeV(),
299  mMin0(), cRes(), mRes0(), mMinCDnow(), alP2(), s0(), mMinXB(), mMinAX(),
300  mMinAXB(), mResXB(), mResAX(), sResXB(), sResAX(), wtNow(), mAtmp(),
301  mBtmp(), multVP(), multVV(), infoPtr() {};
302 
303  // Store pointers and initialize data members.
304  virtual void init(Info* infoPtrIn);
305 
306  // Calculate integrated total/elastic cross sections.
307  virtual bool calcTotEl( int idAin, int idBin, double sIn, double mAin,
308  double mBin);
309 
310  // Differential elastic cross section.
311  virtual double dsigmaEl( double t, bool useCoulomb = false, bool = true);
312 
313  // Calculate integrated diffractive cross sections.
314  virtual bool calcDiff( int idAin, int idBin, double sIn, double mAin,
315  double mBin);
316 
317  // Differential single diffractive cross section.
318  virtual double dsigmaSD( double xi, double t, bool isXB, int = 0);
319 
320  // Differential double diffractive cross section.
321  virtual double dsigmaDD( double xi1, double xi2, double t, int = 0);
322 
323  // Differential central diffractive cross section.
324  virtual double dsigmaCD( double xi1, double xi2, double t1, double t2,
325  int = 0);
326 
327  // Minimal central diffractive mass.
328  virtual double mMinCD() {return mMinCDnow;}
329 
330 private:
331 
332  // Constants: could only be changed in the code itself.
333  static const int IHADATABLE[], IHADBTABLE[], ISDTABLE[], IDDTABLE[], NVMD;
334  static const double EPSILON, ETA, X[], Y[], BETA0[], BHAD[], ALPHAPRIME,
335  CONVERTSD, CONVERTDD, VMDMASS[4], GAMMAFAC[4],
336  CSD[10][8], CDD[10][9];
337 
338  // Initialization data, normally only set once, and result of calculation.
339  bool doDampen, zeroAXB, swapped, sameSign;
340  int idAbsA, idAbsB, iProc, iHadA, iHadB, iHadAtmp[4],
341  iHadBtmp[4], iProcVP[4], iProcVV[4][4];
342  double s, mA, mB, bA, bB, maxXBOwn, maxAXOwn, maxXXOwn, maxAXBOwn, epsSaS,
343  sigmaPomP, mPomP, pPomP, sigAXB2TeV, mMin0, cRes, mRes0, mMinCDnow,
344  alP2, s0, mMinXB, mMinAX, mMinAXB, mResXB, mResAX, sResXB,
345  sResAX, wtNow, mAtmp[4], mBtmp[4], multVP[4], multVV[4][4];
346 
347  // Find combination of incoming beams.
348  bool findBeamComb( int idAin, int idBin, double mAin, double mBin);
349 
350  // Pointer to various information on the generation.
351  Info* infoPtr;
352 
353 };
354 
355 //==========================================================================
356 
357 // The SigmaMBR class parametrizes total and elastic cross sections
358 // according to the Minimum Bias Rockefeller (MBR) model.
359 
360 class SigmaMBR : public SigmaTotAux {
361 
362 public:
363 
364  // Constructor.
365  SigmaMBR() : s(), sigSD(), sigDD(), sigCD(), eps(), alph(), beta0gev(),
366  beta0mb(), sigma0mb(), sigma0gev(), m2min(), dyminSDflux(),
367  dyminDDflux(), dyminCDflux(), dyminSD(), dyminDD(), dyminCD(),
368  dyminSigSD(), dyminSigDD(), dyminSigCD(), a1(), a2(), b1(), b2(),
369  sdpmax(), ddpmax(), dpepmax() {};
370 
371  // Initialize data members.
372  virtual void init(Info* infoPtrIn);
373 
374  // Calculate integrated total/elastic cross sections.
375  virtual bool calcTotEl( int idAin, int idBin, double sIn, double , double );
376 
377  // Differential elastic cross section.
378  virtual double dsigmaEl(double t, bool useCoulomb = false, bool = true);
379 
380  // Calculate integrated diffractive cross sections.
381  virtual bool calcDiff( int , int , double sIn, double , double );
382 
383  // In MBR choice of xi and t values are separated.
384  virtual bool splitDiff() {return true;}
385 
386  // Differential single diffractive cross section.
387  virtual double dsigmaSD( double xi, double t, bool = true, int step = 0);
388 
389  // Differential double diffractive cross section.
390  virtual double dsigmaDD( double xi1, double xi2, double t, int step = 0);
391 
392  // Differential central diffractive cross section.
393  virtual double dsigmaCD( double xi1, double xi2, double t1, double t2,
394  int step = 0);
395 
396  // Minimal central diffractive mass.
397  virtual double mMinCD() {return sqrt(m2min);}
398 
399 private:
400 
401  // Integration of MBR cross sections and form factor approximation.
402  static const int NINTEG, NINTEG2;
403  static const double FFA1, FFA2, FFB1, FFB2;
404 
405  // Parameters of MBR model.
406  double s, sigSD, sigDD, sigCD, eps, alph, beta0gev,
407  beta0mb, sigma0mb, sigma0gev, m2min, dyminSDflux, dyminDDflux,
408  dyminCDflux, dyminSD, dyminDD, dyminCD, dyminSigSD, dyminSigDD,
409  dyminSigCD, a1, a2, b1, b2, sdpmax, ddpmax, dpepmax;
410 
411  // The error function erf(x) should normally be in your math library,
412  // but if not uncomment this simple parametrization by Sergei Winitzki.
413  //double erf(double x) { double x2 = x * x; double kx2 = 0.147 * x2;
414  // double tmp = sqrt(1. - exp(-x2 * (4./M_PI + kx2) / (1. + kx2)));
415  // return ((x >= 0.) ? tmp : -tmp); }
416 
417 };
418 
419 //==========================================================================
420 
421 // The SigmaABMST class parametrizes total and elastic cross sections
422 // according to Appleby, Barlow, Molson, Serluca and Toader (ABMST).
423 
424 class SigmaABMST : public SigmaTotAux {
425 
426 public:
427 
428  // Constructor.
429  SigmaABMST() : ispp(), dampenGap(), useBMin(), modeSD(), modeDD(), modeCD(),
430  s(), facEl(), m2minp(), m2minm(), alp0(), alpt(), s0(), c0(), ygap(),
431  ypow(), expPygap(), multSD(), powSD(), multDD(), powDD(), multCD(),
432  powCD(), mMinCDnow(), bMinSD(), bMinDD(), bMinCD() {};
433 
434  // Initialize data members.
435  virtual void init(Info* infoPtrIn);
436 
437  // Calculate integrated total/elastic cross sections.
438  virtual bool calcTotEl( int idAin, int idBin, double sIn, double , double );
439 
440  // Differential elastic cross section.
441  virtual double dsigmaEl( double t, bool useCoulomb = false,
442  bool onlyPomerons = false) {
443  return facEl * pow2(abs(amplitude( t, useCoulomb, onlyPomerons)));}
444 
445  // Calculate integrated diffractive cross sections.
446  virtual bool calcDiff( int idAin, int idBin, double sIn, double , double );
447 
448  // Differential single diffractive cross section.
449  virtual double dsigmaSD( double xi, double t, bool = true , int = 0);
450 
451  // Differential double diffractive cross section.
452  virtual double dsigmaDD( double xi1, double xi2, double t, int = 0);
453 
454  // Differential central diffractive cross section.
455  virtual double dsigmaCD( double xi1, double xi2, double t1, double t2,
456  int = 0);
457 
458  // Minimal central diffractive mass.
459  virtual double mMinCD() {return mMinCDnow;}
460 
461 private:
462 
463  // Constants: could only be changed in the code itself.
464  static const bool MCINTDD;
465  static const int NPOINTSTSD, NPOINTSTDD, NPOINTMCDD, NPOINTMCCD;
466  static const double EPSI[], ALPP[], NORM[], SLOPE[], FRACS[], TRIG[],
467  LAM2P, BAPPR[], LAM2FF, MRES[4], WRES[4], CRES[4],
468  AFAC[4], BFAC[4], CFAC[4], PPP[4], EPS[2], APR[2],
469  CPI[6], CNST[5], XIDIVSD, DXIRAWSD, DLNXIRAWSD,
470  XIDIVDD, DXIRAWDD, DLNXIRAWDD, BMCINTDD, BMCINTCD;
471 
472  // Initialization data, normally only set once.
473  bool ispp, dampenGap, useBMin;
474  int modeSD, modeDD, modeCD;
475  double s, facEl, m2minp, m2minm, alp0[2], alpt[3], s0, c0,
476  ygap, ypow, expPygap, multSD, powSD, multDD, powDD,
477  multCD, powCD, mMinCDnow, bMinSD, bMinDD, bMinCD;
478 
479  // The scattering amplitude, from which cross sections are derived.
480  complex amplitude( double t, bool useCoulomb = false,
481  bool onlyPomerons = false) ;
482 
483  // Auxiliary method for repetitive part of amplitude.
484  complex sModAlp( double sMod, double alpha) {
485  return exp(complex( 0., -0.5 * M_PI * alpha)) * pow( sMod, alpha); }
486 
487  // Core differential single diffractive cross section.
488  virtual double dsigmaSDcore(double xi, double t);
489 
490  // xi * d(sigma_SD) / (dxi dt) integrated over [t_min, t_max].
491  double dsigmaSDintT( double xi, double tMinIn, double tMaxIn);
492 
493  // d(sigma_SD) / (dxi dt) integrated over [xi_min, xi_max] * [t_min, t_max].
494  double dsigmaSDintXiT( double xiMinIn, double xiMaxIn, double tMinIn,
495  double tMaxIn );
496 
497  // d(sigma_DD) / (dxi1 dxi2 dt) integrated with Monte Carlo sampling.
498  double dsigmaDDintMC();
499 
500  // xi1 * xi2 * d(sigma_DD) / (dxi1 dxi2 dt) integrated over t.
501  double dsigmaDDintT( double xi1, double xi2, double tMinIn, double tMaxIn);
502 
503  // xi1 * d(sigma_DD) / (dxi1 dxi2 dt) integrated over xi2 and t.
504  double dsigmaDDintXi2T( double xi1, double xi2MinIn, double xi2MaxIn,
505  double tMinIn, double tMaxIn);
506 
507  // d(sigma_DD) / (dxi1 dxi2 dt) integrated over xi1, xi2 and t.
508  double dsigmaDDintXi1Xi2T( double xi1MinIn, double xi1MaxIn,
509  double xi2MinIn, double xi2MaxIn, double tMinIn, double tMaxIn);
510 
511  // d(sigma_CD) / (dxi1 dxi2 dt1 dt2) integrated with Monte Carlo sampling.
512  double dsigmaCDintMC();
513 
514 };
515 
516 //==========================================================================
517 
518 // The SigmaRPP class parametrizes total and elastic cross sections
519 // according to the fit in Review of Particle Physics 2014.
520 
521 class SigmaRPP : public SigmaTotAux {
522 
523 public:
524 
525  // Constructor.
526  SigmaRPP() : ispp(), s(), facEl() {};
527 
528  // Initialize data members.
529  virtual void init(Info* infoPtrIn) {
530  tryCoulomb = infoPtrIn->settingsPtr->flag("SigmaElastic:Coulomb");
531  tAbsMin = infoPtrIn->settingsPtr->parm("SigmaElastic:tAbsMin"); }
532 
533  // Calculate integrated total/elastic cross sections.
534  virtual bool calcTotEl( int idAin, int idBin, double sIn, double , double );
535 
536  // Differential elastic cross section.
537  virtual double dsigmaEl( double t, bool useCoulomb = false, bool = true) {
538  return facEl * pow2(abs(amplitude( t, useCoulomb))); }
539 
540 private:
541 
542  // Constants: could only be changed in the code itself.
543  static const double EPS1[], ALPP[], NORM[], BRPP[], KRPP[], LAM2FF;
544 
545  // Initialization data, normally only set once, and result of calculation.
546  bool ispp;
547  double s, facEl;
548 
549  // The scattering amplitude, from which cross sections are derived.
550  complex amplitude( double t, bool useCoulomb = false) ;
551 
552  // Auxiliary method for repetitive part of amplitude.
553  complex sModAlp( double sMod, double alpha) {
554  return exp(complex( 0., -0.5 * M_PI * alpha)) * pow( sMod, alpha); }
555 
556  // Auxiliary methods for Bessel J0 and J1 functions.
557  complex besJ0( complex x);
558  complex besJ1( complex x);
559 
560 };
561 
562 //==========================================================================
563 
564 } // end namespace Pythia8
565 
566 #endif // Pythia8_SigmaTotal_H