PYTHIA  8.317
HIInfo.h
1 // HIInfo.h is a part of the PYTHIA event generator.
2 // Copyright (C) 2026 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 HIInfo, EventInfo and
7 // HIUserHooks classes.
8 //
9 // EventInfo: stores full nucleon-nucleon events with corresponding Info.
10 // HIInfo: info about a Heavy Ion run and its produced events.
11 // HIUserHooks: User hooks for HeavyIons models.
12 
13 #ifndef Pythia8_HIInfo_H
14 #define Pythia8_HIInfo_H
15 
16 #include "Pythia8/Pythia.h"
17 #include "Pythia8/HIBasics.h"
18 #include "Pythia8/HISubCollisionModel.h"
19 
20 namespace Pythia8 {
21 
22 //==========================================================================
23 
24 // Class for collecting info about a Heavy Ion run and its produced
25 // events.
26 
27 class HIInfo {
28 
29 public:
30 
31  friend class HeavyIons;
32  friend class Angantyr;
33 
34  // Constructor.
35  HIInfo() = default;
36 
37  // The impact-parameter distance in the current event.
38  double b() const { return bSave; }
39 
40  // The impact parameter angle.
41  double phi() const { return phiSave; }
42 
43  // The summed elastic amplitude.
44  double T() { return TSave; }
45 
46  // The current sampling area in impact parameter (in millibarn).
47  double xSecArea() const {
48  return xSecScaleSave*FMSQ2MB;
49  }
50 
51  // The average NN non-diffractive impact parameter to be used to
52  // communicate to Pythia's MPI machinery.
53  double avNDb() const { return avNDbSave; }
54 
55  // The number of attempted impact parameter points.
56  long nAttempts() const { return NSave; }
57 
58  // The number of produced events.
59  long nAccepted() const { return NAccSave; }
60 
61  // The total number of sub-collisions in the event.
62  int nCollTot() const { return nCollSave[0]; }
63 
64  // The number of non-diffractive sub-collisions in the event
65  int nCollND() const { return nCollSave[1]; }
66 
67  // The number of single diffractive projectile excitation
68  // sub-collisions in the event
69  int nCollSDP() const { return nCollSave[2]; }
70 
71  // The number of single diffractive target excitation
72  // sub-collisions in the event
73  int nCollSDT() const { return nCollSave[3]; }
74 
75  // The number of double diffractive sub-collisions in the event
76  int nCollDD() const { return nCollSave[4]; }
77 
78  // The number of central diffractive sub-collisions in the event
79  int nCollCD() const { return nCollSave[5]; }
80 
81  // The number of elastic sub-collisions in the event
82  int nCollEl() const { return nCollSave[6]; }
83 
84  // The total number of interacting projectile nucleons in the event.
85  int nPartProj() const { return nProjSave[0]; }
86 
87  // The number of absorptively wounded projectile nucleons in the event.
88  int nAbsProj() const { return nProjSave[1]; }
89 
90  // The number of diffractively wounded projectile nucleons in the event.
91  int nDiffProj() const { return nProjSave[2]; }
92 
93  // The number of elastically scattered target nucleons in the event.
94  int nElProj() const { return nProjSave[3]; }
95 
96  // The total number of interacting target nucleons in the event.
97  int nPartTarg() const { return nTargSave[0]; }
98 
99  // The number of absorptively wounded target nucleons in the event.
100  int nAbsTarg() const { return nTargSave[1]; }
101 
102  // The number of diffractively wounded target nucleons in the event.
103  int nDiffTarg() const { return nTargSave[2]; }
104 
105  // The number of elastically scattered target nucleons in the event.
106  int nElTarg() const { return nTargSave[3]; }
107 
108  // The weight for this collision.
109  double weight() const { return weightSave; }
110 
111  // The sum of weights of the produced events.
112  double weightSum() const { return weightSumSave; }
113 
114  // The number of failed nucleon excitations in the current event.
115  int nFail() const {
116  return nFailSave;
117  }
118 
119  // Register a failed nucleon excitation.
120  void failedExcitation(const SubCollision& subColl) {
121  subColl.failed = true;
122  ++nFailSave;
123  }
124 
125  // Reset Glauber statistics.
126  void glauberReset();
127 
128  // The total cross section from the Glauber calculation.
129  double glauberTot() const {
130  return sigmaTotSave*FMSQ2MB;
131  }
132 
133  // The error on the total cross section from the Glauber
134  // calculation.
135  double glauberTotErr() const {
136  return sqrt(sigErr2TotSave/max(1.0, double(NSave)))*FMSQ2MB;
137  }
138 
139  // The non-diffractive cross section from the Glauber calculation.
140  double glauberND() const {
141  return sigmaNDSave*FMSQ2MB;
142  }
143 
144  // The error on the non-diffractive cross section from the Glauber
145  // calculation.
146  double glauberNDErr() const {
147  return sqrt(sigErr2NDSave/max(1.0, double(NSave)))*FMSQ2MB;
148  }
149 
150  // The total inelastic cross section from the Glauber calculation.
151  double glauberINEL() const {
152  return sigmaINELSave*FMSQ2MB;
153  }
154 
155  // The error on the total inelastic cross section from the Glauber
156  // calculation.
157  double glauberINELErr() const {
158  return sqrt(sigErr2INELSave/max(1.0, double(NSave)))*FMSQ2MB;
159  }
160 
161  // The elastic cross section from the Glauber calculation.
162  double glauberEL() const {
163  return sigmaELSave*FMSQ2MB;
164  }
165 
166  // The error on the elastic cross section from the Glauber
167  // calculation.
168  double glauberELErr() const {
169  return sqrt(sigErr2ELSave/max(1.0, double(NSave)))*FMSQ2MB;
170  }
171 
172  // The diffractive taget excitation cross section from the Glauber
173  // calculation.
174  double glauberDiffT() const {
175  return sigmaDiffTSave*FMSQ2MB;
176  }
177 
178  // The error on the diffractive taget excitation cross section from the
179  // Glauber calculation.
180  double glauberDiffTErr() const {
181  return sqrt(sigErr2DiffTSave/max(1.0, double(NSave)))*FMSQ2MB;
182  }
183 
184  // The diffractive projectile excitation cross section from the Glauber
185  // calculation.
186  double glauberDiffP() const {
187  return sigmaDiffPSave*FMSQ2MB;
188  }
189 
190  // The error on the diffractive projectile excitation cross section
191  // from the Glauber calculation.
192  double glauberDiffPErr() const {
193  return sqrt(sigErr2DiffPSave/max(1.0, double(NSave)))*FMSQ2MB;
194  }
195 
196  // The doubly diffractive excitation cross section from the Glauber
197  // calculation.
198  double glauberDDiff() const {
199  return sigmaDDiffSave*FMSQ2MB;
200  }
201 
202  // The error on the doubly diffractive excitation cross section from
203  // the Glauber calculation.
204  double glauberDDiffErr() const {
205  return sqrt(sigErr2DDiffSave/max(1.0, double(NSave)))*FMSQ2MB;
206  }
207 
208  // The elastic slope parameter.
209  double glauberBSlope() const {
210  return slopeSave/(sigmaTotSave*pow2(HBARC));
211  }
212 
213  // The error on the elastic slope parameter.
214  double glauberBSlopeErr() const {
215  return sqrt((slopeErrSave/pow2(slopeSave) +
216  sigErr2TotSave/pow2(sigmaTotSave))/max(1.0, double(NSave)))*
217  glauberBSlope();
218  }
219 
220  // Overestimates of cross sections of with meduim ions in a hadronic
221  // cascade.
222  const vector<double> & mediumXSecs() const {
223  return mediumXSecsSave;
224  }
225 
226 private:
227 
228  // Register a tried impact parameter point giving the total elastic
229  // amplitude, the impact parameter and impact parameter generation
230  // weight, together with the cross section scale.
231  void addAttempt(double T, double bin, double phiin,
232  double bweight, double xSecScale);
233 
234  // Reweight event for whatever reason.
235  void reweight(double w) {
236  weightSave *= w;
237  }
238 
239  // Select the primary process.
240  void select(Info & in) {
241  primInfo = in;
242  primInfo.hiInfo = this;
243  }
244 
245  // Accept an event and update statistics in info.
246  void accept();
247 
248  // Reject an attmpted event.
249  void reject() {}
250 
251  // Register Glauber statistics of the event.
252  void glauberStatistics();
253 
254  // Collect running averages for Glauber cross section:
255  static void runAvg(double & sig, double & sigErr2, double N,
256  double w) {
257  double delta = w - sig;
258  sig += delta/N;
259  sigErr2 += (delta*(w - sig) - sigErr2)/N;
260  }
261 
262  // Id of the colliding nuclei.
263  int idProjSave = 0, idTargSave = 0;
264 
265  // Impact parameter.
266  double bSave = 0.0;
267  double phiSave = 0.0;
268 
269  // Amplitude, attempts and weight.
270  long NSave = 0, NAccSave = 0;
271  double TSave = 0.0;
272  // OBS: Cross sections should be accessed through the main info class.
273  double sigmaTotSave = {}, sigmaNDSave = {},
274  sigmaELSave = {}, sigmaINELSave = {},
275  sigmaDiffPSave = {}, sigmaDiffTSave = {}, sigmaDDiffSave = {},
276  slopeSave ={};
277  double sigErr2TotSave = {}, sigErr2NDSave = {},
278  sigErr2ELSave = {}, sigErr2INELSave = {},
279  sigErr2DiffPSave = {}, sigErr2DiffTSave = {}, sigErr2DDiffSave = {},
280  slopeErrSave ={};
281  double avNDbSave = {};
282  double weightSave = {}, weightSumSave = {}, xSecScaleSave = {};
283 
284  // Number of collisions and paricipants. See accessor functions for
285  // indices.
286  vector<int> nCollSave{}, nProjSave{}, nTargSave{};
287 
288  // Map of primary processes and the number of events and the sum of
289  // weights.
290  map<int,double> sumPrimW{}, sumPrimW2{};
291  map<int,int> NPrim{};
292  map<int,string> NamePrim{};
293 
294  // The info object of the primary process.
295  Info primInfo{};
296 
297  // Number of failed nucleon excitations.
298  int nFailSave = 0;
299 
300  // Overestimates of cross sections of with meduim ions in a hadronic
301  // cascade.
302  vector<double> mediumXSecsSave = {};
303 
304 public:
305 
306  // Access to subcollision to be extracted by the user.
307  const SubCollisionSet* subCollisionsPtr() { return subCollisionsPtrSave; }
308 
309 private:
310 
311  // Set the subcollision pointer.
312  void setSubCollisions(const SubCollisionSet* subCollisionsPtrIn) {
313  subCollisionsPtrSave = subCollisionsPtrIn; }
314 
315  // Full information about the Glauber calculation, consisting of
316  // all subcollisions.
317  const SubCollisionSet* subCollisionsPtrSave{};
318 
319 };
320 
321 //==========================================================================
322 
323 // This is the heavy ion user hooks class which in the future may be
324 // used inside a Pythia object to generate heavy ion collisons. For
325 // now it is used outside Pythia and requires access to a number of
326 // Pythia objects.
327 
328 class HIUserHooks {
329 
330 public:
331 
332  // The default constructor is empty.
333  HIUserHooks(): idProjSave(0), idTargSave(0) {}
334 
335  // Virtual destructor.
336  virtual ~HIUserHooks() {}
337 
338  // Initialize this user hook.
339  virtual void init(int idProjIn, int idTargIn) {
340  idProjSave = idProjIn;
341  idTargSave = idTargIn;
342  }
343 
344  // A user-supplied impact parameter generator.
345  virtual bool hasImpactParameterGenerator() const { return false; }
346  virtual shared_ptr<ImpactParameterGenerator> impactParameterGenerator()
347  const { return nullptr; }
348 
349  // A user-supplied NucleusModel for the projectile and target.
350  virtual bool hasProjectileModel() const { return false; }
351  virtual shared_ptr<NucleusModel> projectileModel() const { return nullptr; }
352  virtual bool hasTargetModel() const { return false; }
353  virtual shared_ptr<NucleusModel> targetModel() const { return nullptr; }
354 
355  // A user-supplied SubCollisionModel for generating nucleon-nucleon
356  // subcollisions.
357  virtual bool hasSubCollisionModel() { return false; }
358  virtual shared_ptr<SubCollisionModel> subCollisionModel() { return nullptr; }
359 
360  // A user-supplied ordering of events in (inverse) hardness.
361  virtual bool hasEventOrdering() const { return false; }
362  virtual double eventOrdering(const Event &, const Info &) { return -1; }
363 
364  // A user-supplied method for fixing up proton-neutron mismatch in
365  // generated beams.
366  virtual bool canFixIsoSpin() const { return false; }
367  virtual bool fixIsoSpin(EventInfo &) { return false; }
368 
369  // A user-supplied method for shifting the event in impact parameter space.
370  virtual bool canShiftEvent() const { return false; }
371  virtual EventInfo & shiftEvent(EventInfo & ei) const { return ei; }
372 
373  // A user-supplied method of adding a diffractive excitation event
374  // to another event, optionally connecting their colours.
375  bool canAddNucleonExcitation() const { return false; }
376  bool addNucleonExcitation(EventInfo &, EventInfo &, bool) const {
377  return false; }
378 
379  // A user supplied wrapper around the Pythia::forceHadronLevel()
380  virtual bool canForceHadronLevel() const { return false; }
381  virtual bool forceHadronLevel(Pythia &) { return false; }
382 
383  // A user-supplied way of finding the remnants of an
384  // non-diffrcative pp collision (on the target side if tside is
385  // true) to be used to give momentum when adding.
386  virtual bool canFindRecoilers() const { return false; }
387  virtual vector<int>
388  findRecoilers(const Event &, bool /* tside */, int /* beam */, int /* end */,
389  const Vec4 & /* pdiff */, const Vec4 & /* pbeam */) const {
390  return vector<int>();
391  }
392 
393 protected:
394 
395  // Information set in the init() function.
396  // The PDG id of the projectile and target nuclei.
397  int idProjSave, idTargSave;
398 
399 };
400 
401 //==========================================================================
402 
403 } // end namespace Pythia8
404 
405 #endif // Pythia8_HIInfo_H
double glauberTot() const
The total cross section from the Glauber calculation.
Definition: HIInfo.h:129
int nCollSDP() const
Definition: HIInfo.h:69
constexpr double pow2(const double &x)
Powers of small integers - for balance speed/code clarity.
Definition: PythiaStdlib.h:175
virtual void init(int idProjIn, int idTargIn)
Initialize this user hook.
Definition: HIInfo.h:339
bool canAddNucleonExcitation() const
Definition: HIInfo.h:375
virtual bool canFindRecoilers() const
Definition: HIInfo.h:386
long nAttempts() const
The number of attempted impact parameter points.
Definition: HIInfo.h:56
double weight() const
The weight for this collision.
Definition: HIInfo.h:109
int nElProj() const
The number of elastically scattered target nucleons in the event.
Definition: HIInfo.h:94
bool failed
Whether the subcollision failed, i.e. has a failed excitation.
Definition: HISubCollisionModel.h:88
int nFail() const
The number of failed nucleon excitations in the current event.
Definition: HIInfo.h:115
virtual bool canShiftEvent() const
A user-supplied method for shifting the event in impact parameter space.
Definition: HIInfo.h:370
int nCollND() const
The number of non-diffractive sub-collisions in the event.
Definition: HIInfo.h:65
double glauberDiffTErr() const
Definition: HIInfo.h:180
Definition: Info.h:45
The Event class holds all info on the generated event.
Definition: Event.h:408
void glauberReset()
Reset Glauber statistics.
Definition: HIInfo.cc:132
double glauberBSlopeErr() const
The error on the elastic slope parameter.
Definition: HIInfo.h:214
double xSecArea() const
The current sampling area in impact parameter (in millibarn).
Definition: HIInfo.h:47
double glauberND() const
The non-diffractive cross section from the Glauber calculation.
Definition: HIInfo.h:140
double b() const
The impact-parameter distance in the current event.
Definition: HIInfo.h:38
int nCollCD() const
The number of central diffractive sub-collisions in the event.
Definition: HIInfo.h:79
int nElTarg() const
The number of elastically scattered target nucleons in the event.
Definition: HIInfo.h:106
double glauberINELErr() const
Definition: HIInfo.h:157
int nPartTarg() const
The total number of interacting target nucleons in the event.
Definition: HIInfo.h:97
double T()
The summed elastic amplitude.
Definition: HIInfo.h:44
virtual ~HIUserHooks()
Virtual destructor.
Definition: HIInfo.h:336
double glauberDiffP() const
Definition: HIInfo.h:186
The SubCollisionSet gives a set of subcollisions between two nuclei.
Definition: HISubCollisionModel.h:96
Definition: HIInfo.h:328
double weightSum() const
The sum of weights of the produced events.
Definition: HIInfo.h:112
virtual bool hasSubCollisionModel()
Definition: HIInfo.h:357
virtual bool hasImpactParameterGenerator() const
A user-supplied impact parameter generator.
Definition: HIInfo.h:345
constexpr double HBARC
Define conversion hbar * c = 0.2 GeV * fm = 1 and related.
Definition: PythiaStdlib.h:147
const SubCollisionSet * subCollisionsPtr()
Access to subcollision to be extracted by the user.
Definition: HIInfo.h:307
double glauberDiffT() const
Definition: HIInfo.h:174
int nPartProj() const
The total number of interacting projectile nucleons in the event.
Definition: HIInfo.h:85
HIInfo * hiInfo
Definition: Info.h:116
double glauberEL() const
The elastic cross section from the Glauber calculation.
Definition: HIInfo.h:162
int nCollEl() const
The number of elastic sub-collisions in the event.
Definition: HIInfo.h:82
double glauberDDiffErr() const
Definition: HIInfo.h:204
double glauberBSlope() const
The elastic slope parameter.
Definition: HIInfo.h:209
double glauberDiffPErr() const
Definition: HIInfo.h:192
int nCollDD() const
The number of double diffractive sub-collisions in the event.
Definition: HIInfo.h:76
Definition: HISubCollisionModel.h:30
double glauberINEL() const
The total inelastic cross section from the Glauber calculation.
Definition: HIInfo.h:151
virtual bool hasProjectileModel() const
A user-supplied NucleusModel for the projectile and target.
Definition: HIInfo.h:350
void failedExcitation(const SubCollision &subColl)
Register a failed nucleon excitation.
Definition: HIInfo.h:120
HIUserHooks()
The default constructor is empty.
Definition: HIInfo.h:333
double glauberNDErr() const
Definition: HIInfo.h:146
Class for storing Events and Info objects.
Definition: HIBasics.h:25
int idProjSave
Definition: HIInfo.h:397
int nAbsTarg() const
The number of absorptively wounded target nucleons in the event.
Definition: HIInfo.h:100
const vector< double > & mediumXSecs() const
Definition: HIInfo.h:222
int nAbsProj() const
The number of absorptively wounded projectile nucleons in the event.
Definition: HIInfo.h:88
HIInfo()=default
Constructor.
constexpr double FMSQ2MB
Define conversion between fm^2 and mb, in both directions.
Definition: PythiaStdlib.h:171
virtual bool canFixIsoSpin() const
Definition: HIInfo.h:366
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:72
int nCollSDT() const
Definition: HIInfo.h:73
int nDiffTarg() const
The number of diffractively wounded target nucleons in the event.
Definition: HIInfo.h:103
double glauberDDiff() const
Definition: HIInfo.h:198
double glauberELErr() const
Definition: HIInfo.h:168
double avNDb() const
Definition: HIInfo.h:53
virtual bool canForceHadronLevel() const
A user supplied wrapper around the Pythia::forceHadronLevel()
Definition: HIInfo.h:380
virtual bool hasEventOrdering() const
A user-supplied ordering of events in (inverse) hardness.
Definition: HIInfo.h:361
int nDiffProj() const
The number of diffractively wounded projectile nucleons in the event.
Definition: HIInfo.h:91
double glauberTotErr() const
Definition: HIInfo.h:135
double phi() const
The impact parameter angle.
Definition: HIInfo.h:41
Definition: HeavyIons.h:31
Definition: Basics.h:32
The default HeavyIon model in Pythia.
Definition: HeavyIons.h:169
Definition: HIInfo.h:27
long nAccepted() const
The number of produced events.
Definition: HIInfo.h:59
int nCollTot() const
The total number of sub-collisions in the event.
Definition: HIInfo.h:62