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