PYTHIA  8.312
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  bool 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. Return if failure.
117  if (!pythiaMain.init()) return false;
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  if (!pythiaColl.init()) return false;
162  return true;
163 
164  }
165 
166  //--------------------------------------------------------------------------
167 
168  // Calculate the average number of collisions. Average number of
169  // hadron-nucleon collisions in a hadron-nucleus one. If not in
170  // table then interpolate, knowing that <n> - 1 propto A^{2/3}.
171 
172  double nCollAvg(int A) {
173  for (int i = 0; i < nA; ++i) {
174  if (A == tabA[i]) {
175  return (sigmaNow < tabBorder) ? 1. + tabSlopeLo[i] * sigmaNow
176  : 1. + tabOffset[i] + tabSlope[i] * sigmaNow;
177  } else if (A < tabA[i]) {
178  double nColl1 = (sigmaNow < tabBorder) ? tabSlopeLo[i - 1] * sigmaNow
179  : tabOffset[i - 1] + tabSlope[i - 1] * sigmaNow;
180  double nColl2 = (sigmaNow < tabBorder) ? tabSlopeLo[i] * sigmaNow
181  : tabOffset[i] + tabSlope[i] * sigmaNow;
182  double wt1 = (tabA[i] - A) / (tabA[i] - tabA[i - 1]);
183  return 1. + wt1 * pow( A / tabA[i - 1], 2./3.) * nColl1
184  + (1. - wt1) * pow( A / tabA[i], 2./3.) * nColl2;
185  }
186  }
187 
188  return numeric_limits<double>::quiet_NaN();
189  }
190 
191  //--------------------------------------------------------------------------
192 
193  // Calculate the hadron-proton collision cross section.
194  // Return false if not possible to find.
195 
196  bool sigmaSetuphN(int idNowIn, Vec4 pNowIn, double mNowIn) {
197 
198  // Cannot handle low-energy hadrons.
199  if (pNowIn.e() - mNowIn < eKinMin) return false;
200 
201  // Cannot handle hadrons above maximum energy set at initialization.
202  if (pNowIn.e() > eMax) {
203  logger.ERROR_MSG("too high energy");
204  return false;
205  }
206 
207  // Save incoming quantities for reuse in later methods.
208  idNow = idNowIn;
209  pNow = pNowIn;
210  mNow = mNowIn;
211 
212  // Calculate hadron-nucleon cross section. Check if cross section
213  // vanishes.
214  eCMNow = (pNow + Vec4(0, 0, 0, mp)).mCalc();
215  sigmaNow = pythiaColl.getSigmaTotal(idNow, 2212, eCMNow, mNow, mp);
216  if (sigmaNow <= 0.) {
217  if (eCMNow - mNow - mp > eKinMin)
218  logger.ERROR_MSG("vanishing cross section");
219  return false;
220  }
221 
222  // Done.
223  return true;
224 
225  }
226 
227  //--------------------------------------------------------------------------
228 
229  // Calculate the hadron-nucleus cross section for a given nucleon
230  // number A, using hN cross section from sigmaSetuphN. Interpolate
231  // where not (yet) available.
232 
233  double sigmahA(int A) {
234 
235  // Restrict to allowed range 1 <= A <= 208.
236  if (A < 1 || A > 208) {
237  logger.ERROR_MSG("A is outside of valid range (1 <= A <= 208)");
238  return 0.;
239  }
240 
241  // Correction factor for number of h-nucleon collisions per
242  // h-nucleus one.
243  double sigmahA = A * sigmaNow / nCollAvg(A);
244 
245  // Done.
246  return sigmahA;
247 
248  }
249 
250  //--------------------------------------------------------------------------
251 
252  // Generate a collision, and return the event record.
253  // Input (Z, A) of nucleus, and optionally collision vertex.
254 
255  Event& nextColl(int Znow, int Anow, Vec4 vNow = Vec4() ) {
256 
257  // References to the two event records. Clear main event record.
258  Event& eventMain = pythiaMain.event;
259  Event& eventColl = pythiaColl.event;
260  eventMain.clear();
261 
262  // Restrict to allowed range 1 <= A <= 208.
263  if (Anow < 1 || Anow > 208) {
264  logger.ERROR_MSG("A is outside of valid range (1 <= A <= 208)");
265  return eventMain;
266  }
267 
268  // Insert incoming particle in cleared main event record.
269  eventMain.append(90, -11, 0, 0, 1, 1, 0, 0, pNow, mNow);
270  int iHad = eventMain.append(idNow, 12, 0, 0, 0, 0, 0, 0, pNow, mNow);
271  eventMain[iHad].vProd(vNow);
272 
273  // Set up for collisions on a nucleus.
274  int np = Znow;
275  int nn = Anow - Znow;
276  int sizeOld = 0;
277  int sizeNew = 0;
278  Vec4 dirNow = pNow / pNow.pAbs();
279  Rndm& rndm = pythiaMain.rndm;
280 
281  // Drop rate of geometric series. (Deuterium has corrected nCollAvg.)
282  double probMore = 1. - 1. / nCollAvg(Anow);
283 
284  // Loop over varying number of hit nucleons in target nucleus.
285  for (int iColl = 1; iColl <= Anow; ++iColl) {
286  if (iColl > 1 && rndm.flat() > probMore) break;
287 
288  // Pick incoming projectile: trivial for first subcollision, else ...
289  int iProj = iHad;
290  int procType = 0;
291 
292  // ... find highest-pLongitudinal particle from latest subcollision.
293  if (iColl > 1) {
294  iProj = 0;
295  double pMax = 0.;
296  for (int i = sizeOld; i < sizeNew; ++i)
297  if ( eventMain[i].isFinal() && eventMain[i].isHadron()) {
298  double pp = dot3(dirNow, eventMain[i].p());
299  if (pp > pMax) {
300  iProj = i;
301  pMax = pp;
302  }
303  }
304 
305  // No further subcollision if no particle with enough energy.
306  if ( iProj == 0 || eventMain[iProj].e() - eventMain[iProj].m()
307  < eKinMin) break;
308 
309  // Choose process; only SD or ND at perturbative energies.
310  double eCMSub = (eventMain[iProj].p() + Vec4(0, 0, 0, mp)).mCalc();
311  if (eCMSub > 10.) procType = (rndm.flat() < probSD) ? 4 : 1;
312  }
313 
314  // Pick one p or n from target.
315  int idProj = eventMain[iProj].id();
316  bool doProton = rndm.flat() < (np / double(np + nn));
317  if (doProton) np -= 1;
318  else nn -= 1;
319  int idNuc = (doProton) ? 2212 : 2112;
320 
321  // Do a projectile-nucleon subcollision. Return empty event if failure.
322  pythiaColl.setBeamIDs(idProj, idNuc);
323  pythiaColl.setKinematics(eventMain[iProj].p(), Vec4());
324  if (!pythiaColl.next(procType)) {
325  eventMain.clear();
326  return eventMain;
327  }
328 
329  // Insert target nucleon. Mothers are (0,iProj) to mark who it
330  // interacted with. Always use proton mass for simplicity.
331  int statusNuc = (iColl == 1) ? -181 : -182;
332  int iNuc = eventMain.append( idNuc, statusNuc, 0, iProj, 0, 0, 0, 0,
333  0., 0., 0., mp, mp);
334  eventMain[iNuc].vProdAdd(vNow);
335 
336  // Update full energy of the event with the proton mass.
337  eventMain[0].e( eventMain[0].e() + mp);
338  eventMain[0].m( eventMain[0].p().mCalc() );
339 
340  // Insert secondary produced particles (but skip intermediate partons)
341  // into main event record and shift to correct production vertex.
342  sizeOld = eventMain.size();
343  for (int iSub = 3; iSub < eventColl.size(); ++iSub) {
344  if (!eventColl[iSub].isFinal()) continue;
345  int iNew = eventMain.append(eventColl[iSub]);
346  eventMain[iNew].mothers(iNuc, iProj);
347  eventMain[iNew].vProdAdd(vNow);
348  }
349  sizeNew = eventMain.size();
350 
351  // Update daughters of colliding hadrons and other history.
352  eventMain[iProj].daughters(sizeOld, sizeNew - 1);
353  eventMain[iNuc].daughters(sizeOld, sizeNew - 1);
354  eventMain[iProj].statusNeg();
355  eventMain[iProj].tau(0.);
356 
357  // End of loop over interactions in a nucleus.
358  }
359 
360  // Optionally do decays of short-lived particles.
361  if (rapidDecays) pythiaMain.moreDecays();
362 
363  // Optionally compress event record.
364  if (listFinal) compress();
365 
366  // Return generated collision.
367  return eventMain;
368 
369  }
370 
371  //--------------------------------------------------------------------------
372 
373  // Generate a particle decay, and return the event record.
374  // You can allow sequential decays, if they occur rapidly enough.
375 
376  Event& nextDecay(int idNowIn, Vec4 pNowIn, double mNowIn,
377  Vec4 vNow = Vec4() ) {
378 
379  // Save incoming quantities. (Not needed, but by analogy with
380  // collisions.)
381  idNow = idNowIn;
382  pNow = pNowIn;
383  mNow = mNowIn;
384 
385  // References to the event records. Clear them.
386  Event& eventMain = pythiaMain.event;
387  Event& eventColl = pythiaColl.event;
388  eventMain.clear();
389  eventColl.clear();
390 
391  // Insert incoming particle in cleared collision event record.
392  eventColl.append(90, -11, 0, 0, 1, 1, 0, 0, pNow, mNow);
393  int iHad = eventColl.append(idNow, 12, 0, 0, 0, 0, 0, 0, pNow, mNow);
394  eventColl[iHad].vProd(vNow);
395 
396  // Decay incoming particle. Return empty event if fail. Copy event record.
397  if (!pythiaColl.moreDecays(iHad)) return eventMain;
398  eventMain = eventColl;
399 
400  // Optionally do secondary decays of short-lived particles.
401  if (rapidDecays) pythiaMain.moreDecays();
402 
403  // Optionally compress event record.
404  if (listFinal) compress();
405 
406  // Return generated collision.
407  return eventMain;
408 
409  }
410 
411  //--------------------------------------------------------------------------
412 
413  // Compress the event record by removing initial and intermediate
414  // particles. Keep line 0, since the += operator for event records
415  // only adds from 1 on.
416 
417  void compress() {
418 
419  // Reference to the main event record. Original and new size.
420  Event& eventMain = pythiaMain.event;
421  int sizeOld = eventMain.size();
422  int sizeNew = 1;
423 
424  // Loop through all particles and move up the final ones to the top.
425  // Remove history information. Update event record size.
426  for (int i = 1; i < sizeOld; ++i) if (eventMain[i].isFinal()) {
427  eventMain[sizeNew] = eventMain[i];
428  eventMain[sizeNew].mothers( 0, 0);
429  eventMain[sizeNew].daughters( 0, 0);
430  ++sizeNew;
431  }
432 
433  // Shrink event record to new size.
434  eventMain.popBack( sizeOld - sizeNew);
435 
436  }
437 
438  //--------------------------------------------------------------------------
439 
440  // Summary of aborts, errors and warnings.
441 
442  void stat() {
443  pythiaMain.stat();
444  pythiaColl.stat();
445  logger.errorStatistics();
446  }
447 
448  //--------------------------------------------------------------------------
449 
450  // Possibility to access particle data and random numbers from
451  // pythiaMain.
452 
453  ParticleData& particleData() {return pythiaMain.particleData;}
454 
455  Rndm& rndm() {return pythiaMain.rndm;}
456 
457 //--------------------------------------------------------------------------
458 
459 private:
460 
461  // The Pythia instances used for decays and for collisions. Could
462  // be made public, but cleaner to allow only limited access as
463  // above.
464  Pythia pythiaMain, pythiaColl;
465 
466  // Logger instance for errors in this class.
467  Logger logger;
468 
469  // Save quantities.
470  bool listFinal, rapidDecays;
471  int idNow;
472  double eMax, smallTau0, mp, mNow, eCMNow, sigmaNow;
473  Vec4 pNow;
474 
475  // Fraction of single diffractive events beyond first collision in nucleus.
476  static constexpr double probSD = 0.3;
477 
478  // Hadrons below the kinetic energy threshold cannot interact with medium.
479  // (At least not using Pythia.) Used both in lab and in CM frame.
480  static constexpr double eKinMin = 0.2;
481 
482  // Studied nuclei by A number, with offset and slope of <nColl>(sigma):
483  // 1H, 2H, 4He, 9Be, 12C, 14N, 16O, 27Al, 40Ar, 56Fe, 63Cu, 84Kr,
484  // 107Ag, 129Xe, 197Au, 208Pb.
485  static constexpr int nA = 16;
486  static constexpr double tabBorder = 20.;
487  static const double tabA[nA];
488  static const double tabOffset[nA];
489  static const double tabSlope[nA];
490  static const double tabSlopeLo[nA];
491 
492 };
493 
494 // Constant definitions.
495 const double PythiaCascade::tabA[] = {1, 2, 4, 9, 12, 14, 16, 27, 40, 56,
496  63, 84, 107, 129, 197, 208};
497 const double PythiaCascade::tabOffset[] = {0., 0.03, 0.08, 0.15, 0.20, 0.20,
498  0.20, 0.26, 0.30, 0.34, 0.40, 0.40, 0.40, 0.50, 0.50, 0.60};
499 const double PythiaCascade::tabSlope[] = {0., 0.0016, 0.0033, 0.0075, 0.0092,
500  0.0105, 0.012, 0.017, 0.022, 0.027, 0.028, 0.034, 0.040, 0.044, 0.055,
501  0.055};
502 const double PythiaCascade::tabSlopeLo[] = {0., 0.0031, 0.0073, 0.015, 0.0192,
503  0.0205, 0.022, 0.03, 0.037, 0.044, 0.048, 0.054, 0.06, 0.069, 0.08, 0.085};
504 
505 //==========================================================================
506 
507 } // end namespace Pythia8
508 
509 #endif // end Pythia8_PythiaCascade_H
double nCollAvg(int A)
Definition: PythiaCascade.h:172
bool flag(string keyIn)
Give back current value, with check that key exists.
Definition: Settings.cc:1478
PythiaCascade()=default
Default constructor, all setup is done in init().
Rndm rndm
Random number generator.
Definition: Pythia.h:338
ParticleData particleData
ParticleData: the particle data table/database.
Definition: Pythia.h:335
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:332
void clear()
Clear event record.
Definition: Event.h:473
bool setKinematics(double eCMIn)
Switch beam kinematics.
Definition: Pythia.cc:1420
Event event
The event record for the complete event history.
Definition: Pythia.h:323
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:625
void stat()
Summary of aborts, errors and warnings.
Definition: PythiaCascade.h:442
Definition: Logger.h:23
double sigmahA(int A)
Definition: PythiaCascade.h:233
bool setBeamIDs(int idAin, int idBin=0)
Switch to new beam particle identities; for similar hadrons only.
Definition: Pythia.cc:1398
Definition: Basics.h:386
bool moreDecays()
Special routine to allow more decays if on/off switches changed.
Definition: Pythia.h:253
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:376
bool sigmaSetuphN(int idNowIn, Vec4 pNowIn, double mNowIn)
Definition: PythiaCascade.h:196
Intended flow:
Definition: PythiaCascade.h:61
void stat()
Main routine to provide final statistics on generation.
Definition: Pythia.cc:1748
bool 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
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:266
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:595
void compress()
Definition: PythiaCascade.h:417
ParticleData & particleData()
Definition: PythiaCascade.h:453
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:255