PYTHIA  8.313
HeavyIons.h
1 // HeavyIons.h is a part of the PYTHIA event generator.
2 // Copyright (C) 2025 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 the definition of the HeavyIons class which
7 // Provides Pythia with infrastructure to combine several nucleon
8 // collisions into a single heavy ion collision event. This file also
9 // includes the definition of the Angantyr class which implements the
10 // default model for heavy ion collisions in Pythia.
11 
12 #ifndef Pythia8_HeavyIons_H
13 #define Pythia8_HeavyIons_H
14 
15 #include "Pythia8/HIInfo.h"
16 #include "Pythia8/PhysicsBase.h"
17 
18 namespace Pythia8 {
19 
20 // Forward declaration of the Pythia class.
21 class Pythia;
22 
23 //==========================================================================
24 
25 // HeavyIons contains several standard Pythia objects to allow for
26 // the combination of different kinds of nucleon-nucleon events into
27 // a heavy ion collision event. The actual model for doing this must
28 // be implemented in a subclass overriding the the virtual'init()'
29 // and 'next()' functions.
30 
31 class HeavyIons : public PhysicsBase {
32 
33 public:
34 
35  // The constructor needs a reference to the main Pythia object to
36  // which it will belong. A HeavyIons object cannot belong to more
37  // than one main Pythia object.
38  HeavyIons(Pythia& mainPythiaIn)
39  : mainPythiaPtr(&mainPythiaIn), HIHooksPtr(0),
40  pythia(vector<Pythia*>(1, &mainPythiaIn)) {
41  }
42 
43  // Destructor.
44  virtual ~HeavyIons() {}
45 
46  // Virtual function to be implemented in a subclass. This will be
47  // called in the beginning of the Pythia::init function if the mode
48  // HeavyIon:mode is set non zero. The return value should be true
49  // if this object is able to handle the requested collision
50  // type. If false Pythia::init will set HeavyIon:mode to zero but
51  // will try to cope anyway.
52  virtual bool init() = 0;
53 
54  // Virtual function to be implemented in a subclass. This will be
55  // called in the beginning of the Pythia::next function if
56  // HeavyIon:mode is set non zero. After the call, Pythia::next will
57  // return immediately with the return value of this function.
58  virtual bool next() = 0;
59 
60  // Static function to allow Pythia to duplicate some setting names
61  // to be used for secondary Pythia objects.
62  static void addSpecialSettings(Settings& settings);
63 
64  // Return true if the beams in the Primary Pythia object contains
65  // heavy ions.
66  static bool isHeavyIon(Settings& settings);
67 
68  // Possibility to pass in pointer for special heavy ion user hooks.
69  bool setHIUserHooksPtr(HIUserHooksPtr userHooksPtrIn) {
70  HIHooksPtr = userHooksPtrIn; return true;
71  }
72 
73  // Set beam kinematics.
74  virtual bool setKinematics(double /*eCMIn*/) {
75  loggerPtr->ERROR_MSG("method not implemented for this heavy ion model");
76  return false; }
77  virtual bool setKinematics(double /*eAIn*/, double /*eBIn*/) {
78  loggerPtr->ERROR_MSG("method not implemented for this heavy ion model");
79  return false; }
80  virtual bool setKinematics(double, double, double, double, double, double) {
81  loggerPtr->ERROR_MSG("method not implemented for this heavy ion model");
82  return false; }
83  virtual bool setKinematics(Vec4, Vec4) {
84  loggerPtr->ERROR_MSG("method not implemented for this heavy ion model");
85  return false; }
86 
87  // Set beam particles.
88  virtual bool setBeamIDs(int /*idAIn*/, int /*idBIn*/ = 0) {
89  loggerPtr->ERROR_MSG("method not implemented for this heavy ion model");
90  return false; }
91 
92  // The HIInfo object contains information about the last generated heavy
93  // ion event as well as overall statistics of the generated events.
95 
96  // Print out statistics.
97  virtual void stat();
98 
99 protected:
100 
101  // Update the cross section in the main Pythia Info object using
102  // information in the hiInfo object.
103  void updateInfo();
104 
105  // If subclasses has additional Pythia objects for generating
106  // minimum bias nucleon collisions and the main Pythia object is
107  // set up to generated a hard signal process, this function can be
108  // used to clear all selected processes in a clone of the main
109  // Pythia object.
110  void clearProcessLevel(Pythia& pyt);
111 
112  // Duplicate setting on the form match: to settings on the form HImatch:
113  static void setupSpecials(Settings& settings, string match);
114 
115  // Copy settings on the form HImatch: to the corresponding match:
116  // in the given Pythia object.
117  static void setupSpecials(Pythia& p, string match);
118 
119  // Save current beam configuration.
120  int idProj, idTarg;
121 
122  // This is the pointer to the main Pythia object to which this
123  // object is assigned.
125 
126  // Object containing information on inclusive pp cross sections to
127  // be used in Glauber calculations in subclasses.
129 
130  // Optional HIUserHooks object able to modify the behavior of the
131  // HeavyIon model.
132  HIUserHooksPtr HIHooksPtr;
133 
134  // The internal Pythia objects. Index zero will always correspond
135  // to the mainPythiaPtr.
136  vector<Pythia *> pythia;
137 
138  // The names associated with the secondary pythia objects.
139  vector<string> pythiaNames;
140 
141  // The Info objects associated to the secondary the secondary
142  // pythia objects.
143  vector<Info*> info;
144 
145  // Helper class to gain access to the Info object in a pythia
146  // instance.
147  struct InfoGrabber : public UserHooks {
148 
149  // Only one function: return the info object.
150  Info * getInfo() {
151  return infoPtr;
152  }
153 
154  };
155 
156 };
157 
158 //==========================================================================
159 
160 // The default HeavyIon model in Pythia.
161 
162 class Angantyr : public HeavyIons {
163 
164 public:
165 
166  // Enumerate the different internal Pythia objects.
168  HADRON = 0, // For hadronization only.
169  MBIAS = 1, // Minimum Bias processed.
170  SASD = 2, // Single diffractive as one side of non-diffractive.
171  SIGPP = 3, // Optional object for signal processes (pp).
172  SIGPN = 4, // Optional object for signal processes (pn).
173  SIGNP = 5, // Optional object for signal processes (np).
174  SIGNN = 6, // Optional object for signal processes (nn).
175  ALL = 7 // Indicates all objects.
176  };
177 
178 public:
179 
180  // The constructor needs a reference to the main Pythia object to
181  // which it will belong. A Angantyr object cannot belong to more
182  // than one main Pythia object.
183  Angantyr(Pythia& mainPythiaIn);
184 
185  virtual ~Angantyr();
186 
187  // Initialize Angantyr.
188  virtual bool init() override;
189 
190  // Produce a collision involving heavy ions.
191  virtual bool next() override;
192 
193  // Set UserHooks for specific (or ALL) internal Pythia objects.
194  bool setUserHooksPtr(PythiaObject sel, UserHooksPtr userHooksPtrIn);
195 
196  // Set beam kinematics.
197  bool setKinematicsCM();
198  bool setKinematics(double eCMIn) override;
199  bool setKinematics(double eAIn, double eBIn) override;
200  bool setKinematics(double, double, double, double, double, double) override;
201  bool setKinematics(Vec4, Vec4) override;
202  bool setKinematics();
203 
204  // Set beam IDs.
205  bool setBeamIDs(int idAIn, int idBIn = 0) override;
206 
207  // Make sure the correct information is available irrespective of frame type.
208  void unifyFrames();
209 
210  // Print the Angantyr banner.
211  void banner() const;
212 
213  // Subcollisions for the current event.
215  return subColls; }
216 
217  // Get the underlying subcollision model.
219  return *collPtr.get(); }
220 
221  SubCollisionModel* subCollPtr() {
222  return collPtr.get();
223  }
224 
225  // Get the underlying impact parameter generator.
227  return *bGenPtr.get(); }
228 
229  // Projectile nucleus configuration for the current event.
230  const Nucleus& projectile() const {
231  return proj; }
232 
233  // Target nucleus configuration for the current event.
234  const Nucleus& target() const {
235  return targ; }
236 
237  // The underlying projectile nucleus model.
238  const NucleusModel& projectileModel() const {
239  return *projPtr.get(); }
240 
241  // The underlying target nucleus model.
242  const NucleusModel& targetModel() const {
243  return *targPtr.get(); }
244 
245  // Hadronic cross sections used by the subcollision model.
246  const SigmaTotal sigmaNN() const {
247  return sigTotNN; }
248 
249 protected:
250 
251  virtual void onInitInfoPtr() override {
253 
254  // Figure out what beams the user want.
255  void setBeamKinematics(int idA, int idB);
256 
257  // Initiaize a specific Pythia object and optionally run a number
258  // of events to get a handle of the cross section.
259  bool init(PythiaObject sel, string name, int n = 0);
260 
261  // Setup an EventInfo object from a Pythia instance.
262  EventInfo mkEventInfo(Pythia &, Info &, const SubCollision * coll = 0);
263 
264  // Generate events from the internal Pythia oblects;
265  EventInfo getSignal(const SubCollision& coll);
266  EventInfo getND() { return getMBIAS(0, 101); }
267  EventInfo getND(const SubCollision& coll) { return getMBIAS(&coll, 101); }
268  EventInfo getEl(const SubCollision& coll) { return getMBIAS(&coll, 102); }
269  EventInfo getSDP(const SubCollision& coll) { return getMBIAS(&coll, 103); }
270  EventInfo getSDT(const SubCollision& coll) { return getMBIAS(&coll, 104); }
271  EventInfo getDD(const SubCollision& coll) { return getMBIAS(&coll, 105); }
272  EventInfo getCD(const SubCollision& coll) { return getMBIAS(&coll, 106); }
273  EventInfo getSDabsP(const SubCollision& coll)
274  { return getSASD(&coll, 103); }
275  EventInfo getSDabsT(const SubCollision& coll)
276  { return getSASD(&coll, 104); }
277  EventInfo getMBIAS(const SubCollision * coll, int procid);
278  EventInfo getSASD(const SubCollision * coll, int procid);
279 
280  bool genAbs(SubCollisionSet& subCollsIn, list<EventInfo>& subEventsIn);
281  void addSASD(const SubCollisionSet& subCollsIn);
282  bool addDD(const SubCollisionSet& subCollsIn, list<EventInfo>& subEventsIn);
283  bool addSD(const SubCollisionSet& subCollsIn, list<EventInfo>& subEventsIn);
284  void addSDsecond(const SubCollisionSet& subCollsIn);
285  bool addCD(const SubCollisionSet& subCollsIn, list<EventInfo>& subEventsIn);
286  void addCDsecond(const SubCollisionSet& subCollsIn);
287  bool addEL(const SubCollisionSet& subCollsIn, list<EventInfo>& subEventsIn);
288  void addELsecond(const SubCollisionSet& subCollsIn);
289 
290  void resetEvent();
291  bool buildEvent(list<EventInfo>& subEventsIn);
292 
293  bool setupFullCollision(EventInfo& ei, const SubCollision& coll,
294  Nucleon::Status projStatus, Nucleon::Status targStatus);
295  bool isRemnant(const EventInfo& ei, int i, int past = 1 ) const {
296  int statNow = ei.event[i].status()*past;
297  if ( statNow == 63 ) return true;
298  if ( statNow > 70 && statNow < 80 )
299  return isRemnant(ei, ei.event[i].mother1(), -1);
300  return false;
301  }
302  bool fixIsoSpin(EventInfo& ei);
303  EventInfo& shiftEvent(EventInfo& ei);
304  static int getBeam(Event& ev, int i);
305 
306  // Generate a single diffractive
307  bool nextSASD(int proc);
308 
309  // Add a diffractive event to an exsisting one. Optionally connect
310  // the colours of the added event to the original.
311  bool addNucleonExcitation(EventInfo& orig, EventInfo& add,
312  bool colConnect = false);
313 
314  // Find the recoilers in the current event to conserve energy and
315  // momentum in addNucleonExcitation.
316  vector<int> findRecoilers(const Event& e, bool tside, int beam, int end,
317  const Vec4& pdiff, const Vec4& pbeam);
318 
319  // Add a sub-event to the final event record.
320  void addSubEvent(Event& evnt, Event& sub);
321  static void addJunctions(Event& evnt, Event& sub, int coloff);
322 
323  // Add a nucleus remnant to the given event. Possibly introducing
324  // a new particle type.
325  bool addNucleusRemnants();
326 
327 public:
328 
329  // Helper function to construct two transformations that would give
330  // the vectors p1 and p2 the total four momentum of p1p + p2p.
331  static bool
332  getTransforms(Vec4 p1, Vec4 p2, const Vec4& p1p,
333  pair<RotBstMatrix,RotBstMatrix>& R12);
334  static double mT2(const Vec4& p) { return p.pPos()*p.pNeg(); }
335  static double mT(const Vec4& p) { return sqrt(max(mT2(p), 0.0)); }
336 
337 private:
338 
339  // Private UserHooks class to select a specific process.
340  struct ProcessSelectorHook: public UserHooks {
341 
342  ProcessSelectorHook(): proc(0), b(-1.0) {}
343 
344  // Yes we can veto event after process-level selection.
345  virtual bool canVetoProcessLevel() {
346  return true;
347  }
348 
349  // Veto any unwanted process.
350  virtual bool doVetoProcessLevel(Event&) {
351  return proc > 0 && infoPtr->code() != proc;
352  }
353 
354  // Can set the overall impact parameter for the MPI treatment.
355  virtual bool canSetImpactParameter() const {
356  return b >= 0.0;
357  }
358 
359  // Set the overall impact parameter for the MPI treatment.
360  virtual double doSetImpactParameter() {
361  return b;
362  }
363 
364  // The wanted process;
365  int proc;
366 
367  // The selected b-value.
368  double b;
369 
370  };
371 
372  // Holder class to temporarily select a specific process
373  struct HoldProcess {
374 
375  // Set the given process for the given hook object.
376  HoldProcess(shared_ptr<ProcessSelectorHook> hook, int proc,
377  double b = -1.0) : saveHook(hook), saveProc(hook->proc), saveB(hook->b) {
378  hook->proc = proc;
379  hook->b = b;
380  }
381 
382  // Reset the process of the hook object given in the constructor.
383  ~HoldProcess() {
384  if ( saveHook ) {
385  saveHook->proc = saveProc;
386  saveHook->b = saveB;
387  }
388  }
389 
390  // The hook object.
391  shared_ptr<ProcessSelectorHook> saveHook;
392 
393  // The previous process of the hook object.
394  int saveProc;
395 
396  // The previous b-value of the hook object.
397  double saveB;
398 
399  };
400 
401  // The process selector for standard minimum bias processes.
402  shared_ptr<ProcessSelectorHook> selectMB;
403 
404  // The process selector for the SASD object.
405  shared_ptr<ProcessSelectorHook> selectSASD;
406 
407 private:
408 
409  static const int MAXTRY = 999;
410  static const int MAXEVSAVE = 999;
411 
412  // Flag set if there is a specific signal process specified beyond
413  // minimum bias.
414  bool hasSignal;
415 
416  // Whether to do hadronization and hadron level effects.
417  bool doHadronLevel;
418 
419  // Flag to determine whether to do single diffractive test.
420  bool doSDTest;
421 
422  // Flag to determine whether to do only Glauber modelling.
423  bool glauberOnly;
424 
425  // All subcollisions in current collision.
426  SubCollisionSet subColls;
427 
428  // The underlying SubCollisionModel for generating nucleon-nucleon
429  // subcollisions.
430  shared_ptr<SubCollisionModel> collPtr;
431 
432  // The impact parameter generator.
433  shared_ptr<ImpactParameterGenerator> bGenPtr;
434 
435  // The projectile and target nuclei in the current collision.
436  Nucleus proj;
437  Nucleus targ;
438 
439  // The underlying nucleus model for generating nuclons inside the
440  // projectile and target nucleus.
441  shared_ptr<NucleusModel> projPtr;
442  shared_ptr<NucleusModel> targPtr;
443 
444  // Flag to indicate whether variable energy is enabled.
445  bool doVarECM;
446 
447  // Different choices in choosing recoilers when adding
448  // diffractively excited nucleon.
449  int recoilerMode;
450 
451  // Different choices for handling impact parameters.
452  int bMode;
453 
454  // Critical internal error, abort the event.
455  bool doAbort;
456 
457 };
458 
459 //==========================================================================
460 
461 } // end namespace Pythia8
462 
463 #endif // Pythia8_HeavyIons_H
virtual bool init()=0
const SubCollisionModel & subCollisionModel() const
Get the underlying subcollision model.
Definition: HeavyIons.h:218
const NucleusModel & targetModel() const
The underlying target nucleus model.
Definition: HeavyIons.h:242
const SubCollisionSet & subCollisions() const
Subcollisions for the current event.
Definition: HeavyIons.h:214
Definition: PhysicsBase.h:27
void registerSubObject(PhysicsBase &pb)
Register a sub object that should have its information in sync with this.
Definition: PhysicsBase.cc:56
HeavyIons(Pythia &mainPythiaIn)
Definition: HeavyIons.h:38
Definition: Info.h:45
The Event class holds all info on the generated event.
Definition: Event.h:408
Definition: HINucleusModel.h:152
vector< Info * > info
Definition: HeavyIons.h:143
static bool isHeavyIon(Settings &settings)
Check the settings and return false of there are no heavy ion beams.
Definition: HeavyIons.cc:299
Definition: SigmaTotal.h:141
const NucleusModel & projectileModel() const
The underlying projectile nucleus model.
Definition: HeavyIons.h:238
static void setupSpecials(Settings &settings, string match)
Duplicate setting on the form match: to settings on the form HImatch:
Definition: HeavyIons.cc:33
Definition: HINucleusModel.h:187
const ImpactParameterGenerator impactParameterGenerator() const
Get the underlying impact parameter generator.
Definition: HeavyIons.h:226
Definition: HISubCollisionModel.h:502
const Nucleus & target() const
Target nucleus configuration for the current event.
Definition: HeavyIons.h:234
void clearProcessLevel(Pythia &pyt)
Definition: HeavyIons.cc:117
The SubCollisionSet gives a set of subcollisions between two nuclei.
Definition: HISubCollisionModel.h:88
int idProj
Save current beam configuration.
Definition: HeavyIons.h:120
HIUserHooksPtr HIHooksPtr
Definition: HeavyIons.h:132
virtual ~HeavyIons()
Destructor.
Definition: HeavyIons.h:44
Event event
The Event object.
Definition: HIBasics.h:33
virtual bool next()=0
const Nucleus & projectile() const
Projectile nucleus configuration for the current event.
Definition: HeavyIons.h:230
Definition: HISubCollisionModel.h:130
Definition: HeavyIons.h:147
virtual void onInitInfoPtr() override
Definition: HeavyIons.h:251
PythiaObject
Enumerate the different internal Pythia objects.
Definition: HeavyIons.h:167
vector< string > pythiaNames
The names associated with the secondary pythia objects.
Definition: HeavyIons.h:139
static void addSpecialSettings(Settings &settings)
The abstract HeavyIons class.
Definition: HeavyIons.cc:25
Pythia * mainPythiaPtr
Definition: HeavyIons.h:124
Definition: HISubCollisionModel.h:30
virtual bool canVetoProcessLevel()
Possibility to veto event after process-level selection.
Definition: UserHooks.h:60
virtual bool canSetImpactParameter() const
Can set the overall impact parameter for the MPI treatment.
Definition: UserHooks.h:235
bool setHIUserHooksPtr(HIUserHooksPtr userHooksPtrIn)
Possibility to pass in pointer for special heavy ion user hooks.
Definition: HeavyIons.h:69
Info * getInfo()
Only one function: return the info object.
Definition: HeavyIons.h:150
Status
Enum for specifying the status of a nucleon.
Definition: HINucleusModel.h:33
Class for storing Events and Info objects.
Definition: HIBasics.h:25
const SigmaTotal sigmaNN() const
Hadronic cross sections used by the subcollision model.
Definition: HeavyIons.h:246
Logger * loggerPtr
Pointer to logger.
Definition: PhysicsBase.h:87
SigmaTotal sigTotNN
Definition: HeavyIons.h:128
virtual bool setKinematics(double)
Set beam kinematics.
Definition: HeavyIons.h:74
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
HIInfo hiInfo
Definition: HeavyIons.h:94
Info * infoPtr
Definition: PhysicsBase.h:78
virtual bool setBeamIDs(int, int=0)
Set beam particles.
Definition: HeavyIons.h:88
UserHooks is base class for user access to program execution.
Definition: UserHooks.h:32
vector< Pythia * > pythia
Definition: HeavyIons.h:136
virtual bool doVetoProcessLevel(Event &)
Definition: UserHooks.h:64
Definition: HeavyIons.h:31
virtual void stat()
Print out statistics.
Definition: HeavyIons.cc:176
void updateInfo()
Update the Info object in the main Pythia object.
Definition: HeavyIons.cc:146
Definition: Basics.h:32
The default HeavyIon model in Pythia.
Definition: HeavyIons.h:162
Definition: HIInfo.h:27
Definition: Settings.h:196
virtual double doSetImpactParameter()
Set the overall impact parameter for the MPI treatment.
Definition: UserHooks.h:238