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