PYTHIA  8.311
PythiaCascade.h
1 // PythiaCascade.h is a part of the PYTHIA event generator.
2 // Copyright (C) 2024 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 // Author: Torbjorn Sjostrand.
6 
7 #ifndef Pythia8_PythiaCascade_H
8 #define Pythia8_PythiaCascade_H
9 
10 #include "Pythia8/Pythia.h"
11 
12 namespace Pythia8 {
13 
14 //==========================================================================
15 
16 // Wrapper class for Pythia usage in cosmic ray or detector cascades.
17 // Here it is assumed that the user wants to control the full
18 // evolution. The complete event record is then kept somewhere
19 // separate from PYTHIA, and specific input on production and decay
20 // generates the next subevent.
21 
22 // Intended flow:
23 
24 // - init sets up all generation, given a maximum energy. This may be
25 // time-consuming, less so if MPI initialization data can be reused.
26 
27 // - sigmaSetuphN prepares a possible collision for a given hadron by
28 // calculating the hadron-nucleon cross section (if possible).
29 // Moderately time-consuming.
30 
31 // - sigmahA uses the hN cross section calculated by sigaSetuphN and
32 // returns hadron-ion cross section for a collision with a specified
33 // nucleus. It can be called several times to cover a mix of target
34 // nuclei, with minimal time usage.
35 
36 // - nextColl performs the hadron-nucleus collision, as a sequences of
37 // hadron-nucleon ones. Can be quite time-consuming.
38 
39 // - nextDecay can be used anytime to decay a particle. Each
40 // individual decay is rather fast, but there may be many of them.
41 
42 // - stat can be used at end of run to give a summary of error
43 // messages. Negligible time usage.
44 
45 // - references to particleData() and rndm() can be used in the main
46 // - program.
47 
48 // Note: when a hadron interacts with a medium particle, the latter is
49 // added to the event record. This program uses these additional
50 // status codes for target particles in the medium:
51 // 181: the first (or only) target nucleon in a collision.
52 // 182: nucleons from further subcollisions with the same target nucleus.
53 
54 // Warning: the large boosts involved at the higher cosmic ray
55 // energies are not perfectly managed, which leads to non-negligible
56 // energy-momentum non-conservation. The worst (sub)events are caught
57 // and regenerated.
58 
59 //--------------------------------------------------------------------------
60 
62 
63 public:
64 
65  // Default constructor, all setup is done in init().
66  PythiaCascade() = default;
67 
68  //--------------------------------------------------------------------------
69 
70  // Initialize PythiaCascade for a given maximal incoming energy.
71 
72  // Set eMax to the maximal incoming energy (in GeV) on fixed target.
73 
74  // Keep listFinal = false if you want to get back the full event
75  // record, else only the "final" particles of the collision or decay
76  // are returned.
77 
78  // Keep rapidDecays = false if you want to do every decay yourself,
79  // else all particle types with a tau0 below smallTau0 will be
80  // decayed immediately. The tau0 units are mm/c, so default gives c
81  // * tau0 = 1e-10 mm = 100 fm. Note that time dilation effects are
82  // not included, and they can be large.
83 
84  // Set reuseMPI = false in case you cannot reuse the
85  // pythiaCascade.mpi file, or any other relevant stored file. The
86  // saved initialization must have been done with the same or a
87  // larger eMax than the one now being used.
88 
89  void init(double eMaxIn = 1e9, bool listFinalIn = false,
90  bool rapidDecaysIn = false, double smallTau0In = 1e-10,
91  bool reuseMPI = true, string initFile = "pythiaCascade.mpi") {
92 
93  // Store input for future usage.
94  eMax = eMaxIn;
95  listFinal = listFinalIn;
96  rapidDecays = rapidDecaysIn;
97  smallTau0 = smallTau0In;
98 
99  // Proton mass.
100  mp = pythiaMain.particleData.m0(2212);
101 
102  // Main Pythia object for managing the cascade evolution in a
103  // nucleus. Can also do decays, but no hard processes.
104  pythiaMain.readString("ProcessLevel:all = off");
105  pythiaMain.readString("13:mayDecay = on");
106  pythiaMain.readString("211:mayDecay = on");
107  pythiaMain.readString("321:mayDecay = on");
108  pythiaMain.readString("130:mayDecay = on");
109  pythiaMain.settings.flag("ParticleDecays:limitTau0", rapidDecays);
110  pythiaMain.settings.parm("ParticleDecays:tau0Max", smallTau0);
111 
112  // Redure statistics printout to relevant ones.
113  pythiaMain.readString("Stat:showProcessLevel = off");
114  pythiaMain.readString("Stat:showPartonLevel = off");
115 
116  // Initialize.
117  pythiaMain.init();
118 
119  // Secondary Pythia object for performing individual collisions,
120  // or decays. Variable incoming beam type and energy.
121  pythiaColl.readString("Beams:allowVariableEnergy = on");
122  pythiaColl.readString("Beams:allowIDAswitch = on");
123 
124  // Set up for fixed-target collisions.
125  pythiaColl.readString("Beams:frameType = 3");
126  pythiaColl.settings.parm("Beams:pzA", eMax);
127  pythiaColl.readString("Beams:pzB = 0.");
128 
129  // Must use the soft and low-energy QCD processes.
130  pythiaColl.readString("SoftQCD:all = on");
131  pythiaColl.readString("LowEnergyQCD:all = on");
132 
133  // Primary (single) decay to be done by pythiaColl, to circumvent
134  // limitTau0.
135  pythiaColl.readString("13:mayDecay = on");
136  pythiaColl.readString("211:mayDecay = on");
137  pythiaColl.readString("321:mayDecay = on");
138  pythiaColl.readString("130:mayDecay = on");
139 
140  // Secondary decays to be done by pythiaMain, respecting limitTau0.
141  pythiaColl.readString("HadronLevel:Decay = off");
142 
143  // Reduce printout and relax energy-momentum conservation.
144  // (Unusually large errors unfortunate consequence of large boosts.)
145  pythiaColl.readString("Print:quiet = on");
146  pythiaColl.readString("Check:epTolErr = 0.01");
147  pythiaColl.readString("Check:epTolWarn = 0.0001");
148  pythiaColl.readString("Check:mTolErr = 0.01");
149 
150  // Redure statistics printout to relevant ones.
151  pythiaColl.readString("Stat:showProcessLevel = off");
152  pythiaColl.readString("Stat:showPartonLevel = off");
153 
154  // Reuse MPI initialization file if it exists; else create a new one.
155  if (reuseMPI)
156  pythiaColl.readString("MultipartonInteractions:reuseInit = 3");
157  else pythiaColl.readString("MultipartonInteractions:reuseInit = 1");
158  pythiaColl.settings.word("MultipartonInteractions:initFile", initFile);
159 
160  // Initialize.
161  pythiaColl.init();
162 
163  }
164 
165  //--------------------------------------------------------------------------
166 
167  // Calculate the average number of collisions. Average number of
168  // hadron-nucleon collisions in a hadron-nucleus one. If not in
169  // table then interpolate, knowing that <n> - 1 propto A^{2/3}.
170 
171  double nCollAvg(int A) {
172  for (int i = 0; i < nA; ++i) {
173  if (A == tabA[i]) {
174  return (sigmaNow < tabBorder) ? 1. + tabSlopeLo[i] * sigmaNow
175  : 1. + tabOffset[i] + tabSlope[i] * sigmaNow;
176  } else if (A < tabA[i]) {
177  double nColl1 = (sigmaNow < tabBorder) ? tabSlopeLo[i - 1] * sigmaNow
178  : tabOffset[i - 1] + tabSlope[i - 1] * sigmaNow;
179  double nColl2 = (sigmaNow < tabBorder) ? tabSlopeLo[i] * sigmaNow
180  : tabOffset[i] + tabSlope[i] * sigmaNow;
181  double wt1 = (tabA[i] - A) / (tabA[i] - tabA[i - 1]);
182  return 1. + wt1 * pow( A / tabA[i - 1], 2./3.) * nColl1
183  + (1. - wt1) * pow( A / tabA[i], 2./3.) * nColl2;
184  }
185  }
186 
187  return numeric_limits<double>::quiet_NaN();
188  }
189 
190  //--------------------------------------------------------------------------
191 
192  // Calculate the hadron-proton collision cross section.
193  // Return false if not possible to find.
194 
195  bool sigmaSetuphN(int idNowIn, Vec4 pNowIn, double mNowIn) {
196 
197  // Cannot handle low-energy hadrons.
198  if (pNowIn.e() - mNowIn < eKinMin) return false;
199 
200  // Cannot handle hadrons above maximum energy set at initialization.
201  if (pNowIn.e() > eMax) {
202  logger.ERROR_MSG("too high energy");
203  return false;
204  }
205 
206  // Save incoming quantities for reuse in later methods.
207  idNow = idNowIn;
208  pNow = pNowIn;
209  mNow = mNowIn;
210 
211  // Calculate hadron-nucleon cross section. Check if cross section
212  // vanishes.
213  eCMNow = (pNow + Vec4(0, 0, 0, mp)).mCalc();
214  sigmaNow = pythiaColl.getSigmaTotal(idNow, 2212, eCMNow, mNow, mp);
215  if (sigmaNow <= 0.) {
216  if (eCMNow - mNow - mp > eKinMin)
217  logger.ERROR_MSG("vanishing cross section");
218  return false;
219  }
220 
221  // Done.
222  return true;
223 
224  }
225 
226  //--------------------------------------------------------------------------
227 
228  // Calculate the hadron-nucleus cross section for a given nucleon
229  // number A, using hN cross section from sigmaSetuphN. Interpolate
230  // where not (yet) available.
231 
232  double sigmahA(int A) {
233 
234  // Restrict to allowed range 1 <= A <= 208.
235  if (A < 1 || A > 208) {
236  logger.ERROR_MSG("A is outside of valid range (1 <= A <= 208)");
237  return 0.;
238  }
239 
240  // Correction factor for number of h-nucleon collisions per
241  // h-nucleus one.
242  double sigmahA = A * sigmaNow / nCollAvg(A);
243 
244  // Done.
245  return sigmahA;
246 
247  }
248 
249  //--------------------------------------------------------------------------
250 
251  // Generate a collision, and return the event record.
252  // Input (Z, A) of nucleus, and optionally collision vertex.
253 
254  Event& nextColl(int Znow, int Anow, Vec4 vNow = Vec4() ) {
255 
256  // References to the two event records. Clear main event record.
257  Event& eventMain = pythiaMain.event;
258  Event& eventColl = pythiaColl.event;
259  eventMain.clear();
260 
261  // Restrict to allowed range 1 <= A <= 208.
262  if (Anow < 1 || Anow > 208) {
263  logger.ERROR_MSG("A is outside of valid range (1 <= A <= 208)");
264  return eventMain;
265  }
266 
267  // Insert incoming particle in cleared main event record.
268  eventMain.append(90, -11, 0, 0, 1, 1, 0, 0, pNow, mNow);
269  int iHad = eventMain.append(idNow, 12, 0, 0, 0, 0, 0, 0, pNow, mNow);
270  eventMain[iHad].vProd(vNow);
271 
272  // Set up for collisions on a nucleus.
273  int np = Znow;
274  int nn = Anow - Znow;
275  int sizeOld = 0;
276  int sizeNew = 0;
277  Vec4 dirNow = pNow / pNow.pAbs();
278  Rndm& rndm = pythiaMain.rndm;
279 
280  // Drop rate of geometric series. (Deuterium has corrected nCollAvg.)
281  double probMore = 1. - 1. / nCollAvg(Anow);
282 
283  // Loop over varying number of hit nucleons in target nucleus.
284  for (int iColl = 1; iColl <= Anow; ++iColl) {
285  if (iColl > 1 && rndm.flat() > probMore) break;
286 
287  // Pick incoming projectile: trivial for first subcollision, else ...
288  int iProj = iHad;
289  int procType = 0;
290 
291  // ... find highest-pLongitudinal particle from latest subcollision.
292  if (iColl > 1) {
293  iProj = 0;
294  double pMax = 0.;
295  for (int i = sizeOld; i < sizeNew; ++i)
296  if ( eventMain[i].isFinal() && eventMain[i].isHadron()) {
297  double pp = dot3(dirNow, eventMain[i].p());
298  if (pp > pMax) {
299  iProj = i;
300  pMax = pp;
301  }
302  }
303 
304  // No further subcollision if no particle with enough energy.
305  if ( iProj == 0 || eventMain[iProj].e() - eventMain[iProj].m()
306  < eKinMin) break;
307 
308  // Choose process; only SD or ND at perturbative energies.
309  double eCMSub = (eventMain[iProj].p() + Vec4(0, 0, 0, mp)).mCalc();
310  if (eCMSub > 10.) procType = (rndm.flat() < probSD) ? 4 : 1;
311  }
312 
313  // Pick one p or n from target.
314  int idProj = eventMain[iProj].id();
315  bool doProton = rndm.flat() < (np / double(np + nn));
316  if (doProton) np -= 1;
317  else nn -= 1;
318  int idNuc = (doProton) ? 2212 : 2112;
319 
320  // Do a projectile-nucleon subcollision. Return empty event if failure.
321  pythiaColl.setBeamIDs(idProj, idNuc);
322  pythiaColl.setKinematics(eventMain[iProj].p(), Vec4());
323  if (!pythiaColl.next(procType)) {
324  eventMain.clear();
325  return eventMain;
326  }
327 
328  // Insert target nucleon. Mothers are (0,iProj) to mark who it
329  // interacted with. Always use proton mass for simplicity.
330  int statusNuc = (iColl == 1) ? -181 : -182;
331  int iNuc = eventMain.append( idNuc, statusNuc, 0, iProj, 0, 0, 0, 0,
332  0., 0., 0., mp, mp);
333  eventMain[iNuc].vProdAdd(vNow);
334 
335  // Update full energy of the event with the proton mass.
336  eventMain[0].e( eventMain[0].e() + mp);
337  eventMain[0].m( eventMain[0].p().mCalc() );
338 
339  // Insert secondary produced particles (but skip intermediate partons)
340  // into main event record and shift to correct production vertex.
341  sizeOld = eventMain.size();
342  for (int iSub = 3; iSub < eventColl.size(); ++iSub) {
343  if (!eventColl[iSub].isFinal()) continue;
344  int iNew = eventMain.append(eventColl[iSub]);
345  eventMain[iNew].mothers(iNuc, iProj);
346  eventMain[iNew].vProdAdd(vNow);
347  }
348  sizeNew = eventMain.size();
349 
350  // Update daughters of colliding hadrons and other history.
351  eventMain[iProj].daughters(sizeOld, sizeNew - 1);
352  eventMain[iNuc].daughters(sizeOld, sizeNew - 1);
353  eventMain[iProj].statusNeg();
354  eventMain[iProj].tau(0.);
355 
356  // End of loop over interactions in a nucleus.
357  }
358 
359  // Optionally do decays of short-lived particles.
360  if (rapidDecays) pythiaMain.moreDecays();
361 
362  // Optionally compress event record.
363  if (listFinal) compress();
364 
365  // Return generated collision.
366  return eventMain;
367 
368  }
369 
370  //--------------------------------------------------------------------------
371 
372  // Generate a particle decay, and return the event record.
373  // You can allow sequential decays, if they occur rapidly enough.
374 
375  Event& nextDecay(int idNowIn, Vec4 pNowIn, double mNowIn,
376  Vec4 vNow = Vec4() ) {
377 
378  // Save incoming quantities. (Not needed, but by analogy with
379  // collisions.)
380  idNow = idNowIn;
381  pNow = pNowIn;
382  mNow = mNowIn;
383 
384  // References to the event records. Clear them.
385  Event& eventMain = pythiaMain.event;
386  Event& eventColl = pythiaColl.event;
387  eventMain.clear();
388  eventColl.clear();
389 
390  // Insert incoming particle in cleared collision event record.
391  eventColl.append(90, -11, 0, 0, 1, 1, 0, 0, pNow, mNow);
392  int iHad = eventColl.append(idNow, 12, 0, 0, 0, 0, 0, 0, pNow, mNow);
393  eventColl[iHad].vProd(vNow);
394 
395  // Decay incoming particle. Return empty event if fail. Copy event record.
396  if (!pythiaColl.moreDecays(iHad)) return eventMain;
397  eventMain = eventColl;
398 
399  // Optionally do secondary decays of short-lived particles.
400  if (rapidDecays) pythiaMain.moreDecays();
401 
402  // Optionally compress event record.
403  if (listFinal) compress();
404 
405  // Return generated collision.
406  return eventMain;
407 
408  }
409 
410  //--------------------------------------------------------------------------
411 
412  // Compress the event record by removing initial and intermediate
413  // particles. Keep line 0, since the += operator for event records
414  // only adds from 1 on.
415 
416  void compress() {
417 
418  // Reference to the main event record. Original and new size.
419  Event& eventMain = pythiaMain.event;
420  int sizeOld = eventMain.size();
421  int sizeNew = 1;
422 
423  // Loop through all particles and move up the final ones to the top.
424  // Remove history information. Update event record size.
425  for (int i = 1; i < sizeOld; ++i) if (eventMain[i].isFinal()) {
426  eventMain[sizeNew] = eventMain[i];
427  eventMain[sizeNew].mothers( 0, 0);
428  eventMain[sizeNew].daughters( 0, 0);
429  ++sizeNew;
430  }
431 
432  // Shrink event record to new size.
433  eventMain.popBack( sizeOld - sizeNew);
434 
435  }
436 
437  //--------------------------------------------------------------------------
438 
439  // Summary of aborts, errors and warnings.
440 
441  void stat() {
442  pythiaMain.stat();
443  pythiaColl.stat();
444  logger.errorStatistics();
445  }
446 
447  //--------------------------------------------------------------------------
448 
449  // Possibility to access particle data and random numbers from
450  // pythiaMain.
451 
452  ParticleData& particleData() {return pythiaMain.particleData;}
453 
454  Rndm& rndm() {return pythiaMain.rndm;}
455 
456 //--------------------------------------------------------------------------
457 
458 private:
459 
460  // The Pythia instances used for decays and for collisions. Could
461  // be made public, but cleaner to allow only limited access as
462  // above.
463  Pythia pythiaMain, pythiaColl;
464 
465  // Logger instance for errors in this class.
466  Logger logger;
467 
468  // Save quantities.
469  bool listFinal, rapidDecays;
470  int idNow;
471  double eMax, smallTau0, mp, mNow, eCMNow, sigmaNow;
472  Vec4 pNow;
473 
474  // Fraction of single diffractive events beyond first collision in nucleus.
475  static constexpr double probSD = 0.3;
476 
477  // Hadrons below the kinetic energy threshold cannot interact with medium.
478  // (At least not using Pythia.) Used both in lab and in CM frame.
479  static constexpr double eKinMin = 0.2;
480 
481  // Studied nuclei by A number, with offset and slope of <nColl>(sigma):
482  // 1H, 2H, 4He, 9Be, 12C, 14N, 16O, 27Al, 40Ar, 56Fe, 63Cu, 84Kr,
483  // 107Ag, 129Xe, 197Au, 208Pb.
484  static constexpr int nA = 16;
485  static constexpr double tabBorder = 20.;
486  static const double tabA[nA];
487  static const double tabOffset[nA];
488  static const double tabSlope[nA];
489  static const double tabSlopeLo[nA];
490 
491 };
492 
493 // Constant definitions.
494 const double PythiaCascade::tabA[] = {1, 2, 4, 9, 12, 14, 16, 27, 40, 56,
495  63, 84, 107, 129, 197, 208};
496 const double PythiaCascade::tabOffset[] = {0., 0.03, 0.08, 0.15, 0.20, 0.20,
497  0.20, 0.26, 0.30, 0.34, 0.40, 0.40, 0.40, 0.50, 0.50, 0.60};
498 const double PythiaCascade::tabSlope[] = {0., 0.0016, 0.0033, 0.0075, 0.0092,
499  0.0105, 0.012, 0.017, 0.022, 0.027, 0.028, 0.034, 0.040, 0.044, 0.055,
500  0.055};
501 const double PythiaCascade::tabSlopeLo[] = {0., 0.0031, 0.0073, 0.015, 0.0192,
502  0.0205, 0.022, 0.03, 0.037, 0.044, 0.048, 0.054, 0.06, 0.069, 0.08, 0.085};
503 
504 //==========================================================================
505 
506 } // end namespace Pythia8
507 
508 #endif // end Pythia8_PythiaCascade_H
double nCollAvg(int A)
Definition: PythiaCascade.h:171
bool flag(string keyIn)
Give back current value, with check that key exists.
Definition: Settings.cc:1473
PythiaCascade()=default
Default constructor, all setup is done in init().
Rndm rndm
Random number generator.
Definition: Pythia.h:344
ParticleData particleData
ParticleData: the particle data table/database.
Definition: Pythia.h:341
The Event class holds all info on the generated event.
Definition: Event.h:453
Settings settings
Settings: databases of flags/modes/parms/words to control run.
Definition: Pythia.h:338
void clear()
Clear event record.
Definition: Event.h:473
bool setKinematics(double eCMIn)
Switch beam kinematics.
Definition: Pythia.cc:1398
Event event
The event record for the complete event history.
Definition: Pythia.h:329
bool init()
Initialize.
Definition: Pythia.cc:491
double dot3(const Vec4 &v1, const Vec4 &v2)
Scalar and cross product of 3-vector parts.
Definition: Basics.cc:606
void stat()
Summary of aborts, errors and warnings.
Definition: PythiaCascade.h:441
Definition: Logger.h:23
double sigmahA(int A)
Definition: PythiaCascade.h:232
bool setBeamIDs(int idAin, int idBin=0)
Switch to new beam particle identities; for similar hadrons only.
Definition: Pythia.h:233
Definition: Basics.h:385
bool moreDecays()
Special routine to allow more decays if on/off switches changed.
Definition: Pythia.h:259
int append(Particle entryIn)
Put a new particle at the end of the event record; return index.
Definition: Event.h:507
Event & nextDecay(int idNowIn, Vec4 pNowIn, double mNowIn, Vec4 vNow=Vec4())
Definition: PythiaCascade.h:375
bool sigmaSetuphN(int idNowIn, Vec4 pNowIn, double mNowIn)
Definition: PythiaCascade.h:195
void init(double eMaxIn=1e9, bool listFinalIn=false, bool rapidDecaysIn=false, double smallTau0In=1e-10, bool reuseMPI=true, string initFile="pythiaCascade.mpi")
Initialize PythiaCascade for a given maximal incoming energy.
Definition: PythiaCascade.h:89
Intended flow:
Definition: PythiaCascade.h:61
void stat()
Main routine to provide final statistics on generation.
Definition: Pythia.cc:1726
void popBack(int nRemove=1)
Remove last n entries.
Definition: Event.h:561
int size() const
Event record size.
Definition: Event.h:504
double getSigmaTotal()
Get total cross section for two hadrons in the event record or standalone.
Definition: Pythia.h:272
Header for classes to set beam momentum and interaction vertex spread.
Definition: Analysis.h:20
The Pythia class contains the top-level routines to generate an event.
Definition: Pythia.h:71
double m(const Vec4 &v1)
Invariant mass and its square.
Definition: Basics.cc:576
void compress()
Definition: PythiaCascade.h:416
ParticleData & particleData()
Definition: PythiaCascade.h:452
bool next()
Generate the next event.
Definition: Pythia.h:229
bool readString(string, bool warn=true, int subrun=SUBRUNDEFAULT)
Read in one update for a setting or particle data from a single line.
Definition: Pythia.cc:365
double flat()
Generate next random number uniformly between 0 and 1.
Definition: Basics.cc:189
void errorStatistics() const
Print error statistics.
Definition: Logger.h:113
This class holds a map of all ParticleDataEntries.
Definition: ParticleData.h:422
Definition: Basics.h:32
Event & nextColl(int Znow, int Anow, Vec4 vNow=Vec4())
Definition: PythiaCascade.h:254