PYTHIA  8.312
HeavyIons.h
1 // HeavyIons.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 
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 setKinematics(double eCMIn) override;
198  bool setKinematics(double eAIn, double eBIn) override;
199  bool setKinematics(double, double, double, double, double, double) override;
200  bool setKinematics(Vec4, Vec4) override;
201  bool setKinematics();
202 
203  // Set beam IDs.
204  bool setBeamIDs(int idAIn, int idBIn = 0) override;
205 
206  // Make sure the correct information is available irrespective of frame type.
207  void unifyFrames();
208 
209  // Print the Angantyr banner.
210  void banner() const;
211 
212  // Subcollisions for the current event.
214  return subColls; }
215 
216  // Get the underlying subcollision model.
218  return *collPtr.get(); }
219 
220  // Get the underlying impact parameter generator.
222  return *bGenPtr.get(); }
223 
224  // Projectile nucleus configuration for the current event.
225  const Nucleus& projectile() const {
226  return proj; }
227 
228  // Target nucleus configuration for the current event.
229  const Nucleus& target() const {
230  return targ; }
231 
232  // The underlying projectile nucleus model.
233  const NucleusModel& projectileModel() const {
234  return *projPtr.get(); }
235 
236  // The underlying target nucleus model.
237  const NucleusModel& targetModel() const {
238  return *targPtr.get(); }
239 
240  // Hadronic cross sections used by the subcollision model.
241  const SigmaTotal sigmaNN() const {
242  return sigTotNN; }
243 
244 protected:
245 
246  virtual void onInitInfoPtr() override {
248 
249  // Figure out what beams the user want.
250  void setBeamKinematics(int idA, int idB);
251 
252  // Initiaize a specific Pythia object and optionally run a number
253  // of events to get a handle of the cross section.
254  bool init(PythiaObject sel, string name, int n = 0);
255 
256  // Setup an EventInfo object from a Pythia instance.
257  EventInfo mkEventInfo(Pythia &, Info &, const SubCollision * coll = 0);
258 
259  // Generate events from the internal Pythia oblects;
260  EventInfo getSignal(const SubCollision& coll);
261  EventInfo getND() { return getMBIAS(0, 101); }
262  EventInfo getND(const SubCollision& coll) { return getMBIAS(&coll, 101); }
263  EventInfo getEl(const SubCollision& coll) { return getMBIAS(&coll, 102); }
264  EventInfo getSDP(const SubCollision& coll) { return getMBIAS(&coll, 103); }
265  EventInfo getSDT(const SubCollision& coll) { return getMBIAS(&coll, 104); }
266  EventInfo getDD(const SubCollision& coll) { return getMBIAS(&coll, 105); }
267  EventInfo getCD(const SubCollision& coll) { return getMBIAS(&coll, 106); }
268  EventInfo getSDabsP(const SubCollision& coll)
269  { return getSASD(&coll, 103); }
270  EventInfo getSDabsT(const SubCollision& coll)
271  { return getSASD(&coll, 104); }
272  EventInfo getMBIAS(const SubCollision * coll, int procid);
273  EventInfo getSASD(const SubCollision * coll, int procid);
274 
275  bool genAbs(SubCollisionSet& subCollsIn, list<EventInfo>& subEventsIn);
276  void addSASD(const SubCollisionSet& subCollsIn);
277  bool addDD(const SubCollisionSet& subCollsIn, list<EventInfo>& subEventsIn);
278  bool addSD(const SubCollisionSet& subCollsIn, list<EventInfo>& subEventsIn);
279  void addSDsecond(const SubCollisionSet& subCollsIn);
280  bool addCD(const SubCollisionSet& subCollsIn, list<EventInfo>& subEventsIn);
281  void addCDsecond(const SubCollisionSet& subCollsIn);
282  bool addEL(const SubCollisionSet& subCollsIn, list<EventInfo>& subEventsIn);
283  void addELsecond(const SubCollisionSet& subCollsIn);
284 
285  bool buildEvent(list<EventInfo>& subEventsIn);
286 
287  bool setupFullCollision(EventInfo& ei, const SubCollision& coll,
288  Nucleon::Status projStatus, Nucleon::Status targStatus);
289  bool isRemnant(const EventInfo& ei, int i, int past = 1 ) const {
290  int statNow = ei.event[i].status()*past;
291  if ( statNow == 63 ) return true;
292  if ( statNow > 70 && statNow < 80 )
293  return isRemnant(ei, ei.event[i].mother1(), -1);
294  return false;
295  }
296  bool fixIsoSpin(EventInfo& ei);
297  EventInfo& shiftEvent(EventInfo& ei);
298  static int getBeam(Event& ev, int i);
299 
300  // Generate a single diffractive
301  bool nextSASD(int proc);
302 
303  // Add a diffractive event to an exsisting one. Optionally connect
304  // the colours of the added event to the original.
305  bool addNucleonExcitation(EventInfo& orig, EventInfo& add,
306  bool colConnect = false);
307 
308  // Find the recoilers in the current event to conserve energy and
309  // momentum in addNucleonExcitation.
310  vector<int> findRecoilers(const Event& e, bool tside, int beam, int end,
311  const Vec4& pdiff, const Vec4& pbeam);
312 
313  // Add a sub-event to the final event record.
314  void addSubEvent(Event& evnt, Event& sub);
315  static void addJunctions(Event& evnt, Event& sub, int coloff);
316 
317  // Add a nucleus remnant to the given event. Possibly introducing
318  // a new particle type.
319  bool addNucleusRemnants();
320 
321 public:
322 
323  // Helper function to construct two transformations that would give
324  // the vectors p1 and p2 the total four momentum of p1p + p2p.
325  static bool
326  getTransforms(Vec4 p1, Vec4 p2, const Vec4& p1p,
327  pair<RotBstMatrix,RotBstMatrix>& R12);
328  static double mT2(const Vec4& p) { return p.pPos()*p.pNeg(); }
329  static double mT(const Vec4& p) { return sqrt(max(mT2(p), 0.0)); }
330 
331 private:
332 
333  // Private UserHooks class to select a specific process.
334  struct ProcessSelectorHook: public UserHooks {
335 
336  ProcessSelectorHook(): proc(0), b(-1.0) {}
337 
338  // Yes we can veto event after process-level selection.
339  virtual bool canVetoProcessLevel() {
340  return true;
341  }
342 
343  // Veto any unwanted process.
344  virtual bool doVetoProcessLevel(Event&) {
345  return proc > 0 && infoPtr->code() != proc;
346  }
347 
348  // Can set the overall impact parameter for the MPI treatment.
349  virtual bool canSetImpactParameter() const {
350  return b >= 0.0;
351  }
352 
353  // Set the overall impact parameter for the MPI treatment.
354  virtual double doSetImpactParameter() {
355  return b;
356  }
357 
358  // The wanted process;
359  int proc;
360 
361  // The selected b-value.
362  double b;
363 
364  };
365 
366  // Holder class to temporarily select a specific process
367  struct HoldProcess {
368 
369  // Set the given process for the given hook object.
370  HoldProcess(shared_ptr<ProcessSelectorHook> hook, int proc,
371  double b = -1.0) : saveHook(hook), saveProc(hook->proc), saveB(hook->b) {
372  hook->proc = proc;
373  hook->b = b;
374  }
375 
376  // Reset the process of the hook object given in the constructor.
377  ~HoldProcess() {
378  if ( saveHook ) {
379  saveHook->proc = saveProc;
380  saveHook->b = saveB;
381  }
382  }
383 
384  // The hook object.
385  shared_ptr<ProcessSelectorHook> saveHook;
386 
387  // The previous process of the hook object.
388  int saveProc;
389 
390  // The previous b-value of the hook object.
391  double saveB;
392 
393  };
394 
395  // The process selector for standard minimum bias processes.
396  shared_ptr<ProcessSelectorHook> selectMB;
397 
398  // The process selector for the SASD object.
399  shared_ptr<ProcessSelectorHook> selectSASD;
400 
401 private:
402 
403  static const int MAXTRY = 999;
404  static const int MAXEVSAVE = 999;
405 
406  // Flag set if there is a specific signal process specified beyond
407  // minimum bias.
408  bool hasSignal;
409 
410  // Whether to do hadronization and hadron level effects.
411  bool doHadronLevel;
412 
413  // Flag to determine whether to do single diffractive test.
414  bool doSDTest;
415 
416  // Flag to determine whether to do only Glauber modelling.
417  bool glauberOnly;
418 
419  // All subcollisions in current collision.
420  SubCollisionSet subColls;
421 
422  // The underlying SubCollisionModel for generating nucleon-nucleon
423  // subcollisions.
424  shared_ptr<SubCollisionModel> collPtr;
425 
426  // The impact parameter generator.
427  shared_ptr<ImpactParameterGenerator> bGenPtr;
428 
429  // The projectile and target nuclei in the current collision.
430  Nucleus proj;
431  Nucleus targ;
432 
433  // The underlying nucleus model for generating nuclons inside the
434  // projectile and target nucleus.
435  shared_ptr<NucleusModel> projPtr;
436  shared_ptr<NucleusModel> targPtr;
437 
438  // Flag to indicate whether variable energy is enabled.
439  bool doVarECM;
440 
441  // Different choices in choosing recoilers when adding
442  // diffractively excited nucleon.
443  int recoilerMode;
444 
445  // Different choices for handling impact parameters.
446  int bMode;
447 
448  // Critical internal error, abort the event.
449  bool doAbort;
450 
451 };
452 
453 //==========================================================================
454 
455 } // end namespace Pythia8
456 
457 #endif // Pythia8_HeavyIons_H
virtual bool init()=0
const SubCollisionModel & subCollisionModel() const
Get the underlying subcollision model.
Definition: HeavyIons.h:217
const NucleusModel & targetModel() const
The underlying target nucleus model.
Definition: HeavyIons.h:237
const SubCollisionSet & subCollisions() const
Subcollisions for the current event.
Definition: HeavyIons.h:213
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:453
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:250
Definition: SigmaTotal.h:141
const NucleusModel & projectileModel() const
The underlying projectile nucleus model.
Definition: HeavyIons.h:233
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:221
Definition: HISubCollisionModel.h:491
const Nucleus & target() const
Target nucleus configuration for the current event.
Definition: HeavyIons.h:229
void clearProcessLevel(Pythia &pyt)
Definition: HeavyIons.cc:117
The SubCollisionSet gives a set of subcollisions between two nuclei.
Definition: HISubCollisionModel.h:84
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:55
virtual bool next()=0
const Nucleus & projectile() const
Projectile nucleus configuration for the current event.
Definition: HeavyIons.h:225
Definition: HISubCollisionModel.h:121
Definition: HeavyIons.h:147
virtual void onInitInfoPtr() override
Definition: HeavyIons.h:246
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:221
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:47
const SigmaTotal sigmaNN() const
Hadronic cross sections used by the subcollision model.
Definition: HeavyIons.h:241
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:195
virtual double doSetImpactParameter()
Set the overall impact parameter for the MPI treatment.
Definition: UserHooks.h:224