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