StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
Event.h
1 // Event.h is a part of the PYTHIA event generator.
2 // Copyright (C) 2018 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 // Header file for the Particle and Event classes.
7 // Particle: information on an instance of a particle.
8 // Junction: information on a junction between three colours.
9 // Event: list of particles in the current event.
10 
11 #ifndef Pythia8_Event_H
12 #define Pythia8_Event_H
13 
14 #include "Pythia8/Basics.h"
15 #include "Pythia8/ParticleData.h"
16 #include "Pythia8/PythiaStdlib.h"
17 
18 namespace Pythia8 {
19 
20 //==========================================================================
21 
22 // Forward references to ParticleDataEntry and ResonanceWidths classes.
23 class ParticleDataEntry;
24 class ResonanceWidths;
25 class Event;
26 
27 //==========================================================================
28 
29 // Particle class.
30 // This class holds info on a particle in general.
31 
32 class Particle {
33 
34 public:
35 
36  // Constructors.
37  Particle() : idSave(0), statusSave(0), mother1Save(0), mother2Save(0),
38  daughter1Save(0), daughter2Save(0), colSave(0), acolSave(0),
39  pSave(Vec4(0.,0.,0.,0.)), mSave(0.), scaleSave(0.), polSave(9.),
40  hasVertexSave(false), vProdSave(Vec4(0.,0.,0.,0.)), tauSave(0.),
41  pdePtr(0), evtPtr(0) { }
42  Particle(int idIn, int statusIn = 0, int mother1In = 0,
43  int mother2In = 0, int daughter1In = 0, int daughter2In = 0,
44  int colIn = 0, int acolIn = 0, double pxIn = 0., double pyIn = 0.,
45  double pzIn = 0., double eIn = 0., double mIn = 0.,
46  double scaleIn = 0., double polIn = 9.)
47  : idSave(idIn), statusSave(statusIn), mother1Save(mother1In),
48  mother2Save(mother2In), daughter1Save(daughter1In),
49  daughter2Save(daughter2In), colSave(colIn), acolSave(acolIn),
50  pSave(Vec4(pxIn, pyIn, pzIn, eIn)), mSave(mIn), scaleSave(scaleIn),
51  polSave(polIn), hasVertexSave(false), vProdSave(Vec4(0.,0.,0.,0.)),
52  tauSave(0.), pdePtr(0), evtPtr(0) { }
53  Particle(int idIn, int statusIn, int mother1In, int mother2In,
54  int daughter1In, int daughter2In, int colIn, int acolIn,
55  Vec4 pIn, double mIn = 0., double scaleIn = 0., double polIn = 9.)
56  : idSave(idIn), statusSave(statusIn), mother1Save(mother1In),
57  mother2Save(mother2In), daughter1Save(daughter1In),
58  daughter2Save(daughter2In), colSave(colIn), acolSave(acolIn),
59  pSave(pIn), mSave(mIn), scaleSave(scaleIn), polSave(polIn),
60  hasVertexSave(false), vProdSave(Vec4(0.,0.,0.,0.)), tauSave(0.),
61  pdePtr(0), evtPtr(0) { }
62  Particle(const Particle& pt) : idSave(pt.idSave),
63  statusSave(pt.statusSave), mother1Save(pt.mother1Save),
64  mother2Save(pt.mother2Save), daughter1Save(pt.daughter1Save),
65  daughter2Save(pt.daughter2Save), colSave(pt.colSave),
66  acolSave(pt.acolSave), pSave(pt.pSave), mSave(pt.mSave),
67  scaleSave(pt.scaleSave), polSave(pt.polSave),
68  hasVertexSave(pt.hasVertexSave), vProdSave(pt.vProdSave),
69  tauSave(pt.tauSave), pdePtr(pt.pdePtr), evtPtr(pt.evtPtr) { }
70  Particle& operator=(const Particle& pt) {if (this != &pt) {
71  idSave = pt.idSave; statusSave = pt.statusSave;
72  mother1Save = pt.mother1Save; mother2Save = pt.mother2Save;
73  daughter1Save = pt.daughter1Save; daughter2Save = pt.daughter2Save;
74  colSave = pt.colSave; acolSave = pt.acolSave; pSave = pt.pSave;
75  mSave = pt.mSave; scaleSave = pt.scaleSave; polSave = pt.polSave;
76  hasVertexSave = pt.hasVertexSave; vProdSave = pt.vProdSave;
77  tauSave = pt.tauSave; pdePtr = pt.pdePtr; evtPtr = pt.evtPtr; }
78  return *this; }
79 
80  // Member functions to set the Event and ParticleDataEntry pointers.
81  void setEvtPtr(Event* evtPtrIn) { evtPtr = evtPtrIn; setPDEPtr();}
82  void setPDEPtr(ParticleDataEntry* pdePtrIn = 0);
83 
84  // Member functions for input.
85  void id(int idIn) {idSave = idIn; setPDEPtr();}
86  void status(int statusIn) {statusSave = statusIn;}
87  void statusPos() {statusSave = abs(statusSave);}
88  void statusNeg() {statusSave = -abs(statusSave);}
89  void statusCode(int statusIn) {statusSave =
90  (statusSave > 0) ? abs(statusIn) : -abs(statusIn);}
91  void mother1(int mother1In) {mother1Save = mother1In;}
92  void mother2(int mother2In) {mother2Save = mother2In;}
93  void mothers(int mother1In = 0, int mother2In = 0)
94  {mother1Save = mother1In; mother2Save = mother2In;}
95  void daughter1(int daughter1In) {daughter1Save = daughter1In;}
96  void daughter2(int daughter2In) {daughter2Save = daughter2In;}
97  void daughters(int daughter1In = 0, int daughter2In = 0)
98  {daughter1Save = daughter1In; daughter2Save = daughter2In;}
99  void col(int colIn) {colSave = colIn;}
100  void acol(int acolIn) {acolSave = acolIn;}
101  void cols(int colIn = 0,int acolIn = 0) {colSave = colIn;
102  acolSave = acolIn;}
103  void p(Vec4 pIn) {pSave = pIn;}
104  void p(double pxIn, double pyIn, double pzIn, double eIn)
105  {pSave.p(pxIn, pyIn, pzIn, eIn);}
106  void px(double pxIn) {pSave.px(pxIn);}
107  void py(double pyIn) {pSave.py(pyIn);}
108  void pz(double pzIn) {pSave.pz(pzIn);}
109  void e(double eIn) {pSave.e(eIn);}
110  void m(double mIn) {mSave = mIn;}
111  void scale(double scaleIn) {scaleSave = scaleIn;}
112  void pol(double polIn) {polSave = polIn;}
113  void vProd(Vec4 vProdIn) {vProdSave = vProdIn; hasVertexSave = true;}
114  void vProd(double xProdIn, double yProdIn, double zProdIn, double tProdIn)
115  {vProdSave.p(xProdIn, yProdIn, zProdIn, tProdIn); hasVertexSave = true;}
116  void xProd(double xProdIn) {vProdSave.px(xProdIn); hasVertexSave = true;}
117  void yProd(double yProdIn) {vProdSave.py(yProdIn); hasVertexSave = true;}
118  void zProd(double zProdIn) {vProdSave.pz(zProdIn); hasVertexSave = true;}
119  void tProd(double tProdIn) {vProdSave.e(tProdIn); hasVertexSave = true;}
120  void tau(double tauIn) {tauSave = tauIn;}
121 
122  // Member functions for output.
123  int id() const {return idSave;}
124  int status() const {return statusSave;}
125  int mother1() const {return mother1Save;}
126  int mother2() const {return mother2Save;}
127  int daughter1() const {return daughter1Save;}
128  int daughter2() const {return daughter2Save;}
129  int col() const {return colSave;}
130  int acol() const {return acolSave;}
131  Vec4 p() const {return pSave;}
132  double px() const {return pSave.px();}
133  double py() const {return pSave.py();}
134  double pz() const {return pSave.pz();}
135  double e() const {return pSave.e();}
136  double m() const {return mSave;}
137  double scale() const {return scaleSave;}
138  double pol() const {return polSave;}
139  bool hasVertex() const {return hasVertexSave;}
140  Vec4 vProd() const {return vProdSave;}
141  double xProd() const {return vProdSave.px();}
142  double yProd() const {return vProdSave.py();}
143  double zProd() const {return vProdSave.pz();}
144  double tProd() const {return vProdSave.e();}
145  double tau() const {return tauSave;}
146 
147  // Member functions for output; derived int and bool quantities.
148  int idAbs() const {return abs(idSave);}
149  int statusAbs() const {return abs(statusSave);}
150  bool isFinal() const {return (statusSave > 0);}
151  int intPol() const;
152  bool isRescatteredIncoming() const {return
153  (statusSave == -34 || statusSave == -45 ||
154  statusSave == -46 || statusSave == -54);}
155 
156  // Member functions for output; derived double quantities.
157  double m2() const {return (mSave >= 0.) ? mSave*mSave
158  : -mSave*mSave;}
159  double mCalc() const {return pSave.mCalc();}
160  double m2Calc() const {return pSave.m2Calc();}
161  double eCalc() const {return sqrt(abs(m2() + pSave.pAbs2()));}
162  double pT() const {return pSave.pT();}
163  double pT2() const {return pSave.pT2();}
164  double mT() const {double temp = m2() + pSave.pT2();
165  return (temp >= 0.) ? sqrt(temp) : -sqrt(-temp);}
166  double mT2() const {return m2() + pSave.pT2();}
167  double pAbs() const {return pSave.pAbs();}
168  double pAbs2() const {return pSave.pAbs2();}
169  double eT() const {return pSave.eT();}
170  double eT2() const {return pSave.eT2();}
171  double theta() const {return pSave.theta();}
172  double phi() const {return pSave.phi();}
173  double thetaXZ() const {return pSave.thetaXZ();}
174  double pPos() const {return pSave.pPos();}
175  double pNeg() const {return pSave.pNeg();}
176  double y() const;
177  double eta() const;
178  double y(double mCut) const;
179  double y(double mCut, RotBstMatrix& M) const;
180  Vec4 vDec() const {return (tauSave > 0. && mSave > 0.)
181  ? vProdSave + tauSave * pSave / mSave : vProdSave;}
182  double xDec() const {return (tauSave > 0. && mSave > 0.)
183  ? vProdSave.px() + tauSave * pSave.px() / mSave : vProdSave.px();}
184  double yDec() const {return (tauSave > 0. && mSave > 0.)
185  ? vProdSave.py() + tauSave * pSave.py() / mSave : vProdSave.py();}
186  double zDec() const {return (tauSave > 0. && mSave > 0.)
187  ? vProdSave.pz() + tauSave * pSave.pz() / mSave : vProdSave.pz();}
188  double tDec() const {return (tauSave > 0. && mSave > 0.)
189  ? vProdSave.e() + tauSave * pSave.e() / mSave : vProdSave.e();}
190 
191  // Methods that can refer back to the event the particle belongs to.
192  virtual int index() const;
193  int iTopCopy() const;
194  int iBotCopy() const;
195  int iTopCopyId(bool simplify = false) const;
196  int iBotCopyId(bool simplify = false) const;
197  vector<int> motherList() const;
198  vector<int> daughterList() const;
199  vector<int> daughterListRecursive() const;
200  vector<int> sisterList(bool traceTopBot = false) const;
201  bool isAncestor(int iAncestor) const;
202  int statusHepMC() const;
203  bool isFinalPartonLevel() const;
204  bool undoDecay();
205 
206  // Further output, based on a pointer to a ParticleDataEntry object.
207  string name() const {
208  return (pdePtr != 0) ? pdePtr->name(idSave) : " ";}
209  string nameWithStatus(int maxLen = 20) const;
210  int spinType() const {
211  return (pdePtr != 0) ? pdePtr->spinType() : 0;}
212  int chargeType() const {
213  return (pdePtr != 0) ? pdePtr->chargeType(idSave) : 0;}
214  double charge() const {
215  return (pdePtr != 0) ? pdePtr->charge(idSave) : 0;}
216  bool isCharged() const {
217  return (pdePtr != 0) ? (pdePtr->chargeType(idSave) != 0) : false;}
218  bool isNeutral() const {
219  return (pdePtr != 0) ? (pdePtr->chargeType(idSave) == 0) : false;}
220  int colType() const {
221  return (pdePtr != 0) ? pdePtr->colType(idSave) : 0;}
222  double m0() const {
223  return (pdePtr != 0) ? pdePtr->m0() : 0.;}
224  double mWidth() const {
225  return (pdePtr != 0) ? pdePtr->mWidth() : 0.;}
226  double mMin() const {
227  return (pdePtr != 0) ? pdePtr->mMin() : 0.;}
228  double mMax() const {
229  return (pdePtr != 0) ? pdePtr->mMax() : 0.;}
230  double mSel() const {
231  return (pdePtr != 0) ? pdePtr->mSel() : 0.;}
232  double constituentMass() const {
233  return (pdePtr != 0) ? pdePtr->constituentMass() : 0.;}
234  double tau0() const {
235  return (pdePtr != 0) ? pdePtr->tau0() : 0.;}
236  bool mayDecay() const {
237  return (pdePtr != 0) ? pdePtr->mayDecay() : false;}
238  bool canDecay() const {
239  return (pdePtr != 0) ? pdePtr->canDecay() : false;}
240  bool doExternalDecay() const {
241  return (pdePtr != 0) ? pdePtr->doExternalDecay() : false;}
242  bool isResonance() const {
243  return (pdePtr != 0) ? pdePtr->isResonance() : false;}
244  bool isVisible() const {
245  return (pdePtr != 0) ? pdePtr->isVisible() : false;}
246  bool isLepton() const {
247  return (pdePtr != 0) ? pdePtr->isLepton() : false;}
248  bool isQuark() const {
249  return (pdePtr != 0) ? pdePtr->isQuark() : false;}
250  bool isGluon() const {
251  return (pdePtr != 0) ? pdePtr->isGluon() : false;}
252  bool isDiquark() const {
253  return (pdePtr != 0) ? pdePtr->isDiquark() : false;}
254  bool isParton() const {
255  return (pdePtr != 0) ? pdePtr->isParton() : false;}
256  bool isHadron() const {
257  return (pdePtr != 0) ? pdePtr->isHadron() : false;}
258  ParticleDataEntry& particleDataEntry() const {return *pdePtr;}
259 
260  // Member functions that perform operations.
261  void rescale3(double fac) {pSave.rescale3(fac);}
262  void rescale4(double fac) {pSave.rescale4(fac);}
263  void rescale5(double fac) {pSave.rescale4(fac); mSave *= fac;}
264  void rot(double thetaIn, double phiIn) {pSave.rot(thetaIn, phiIn);
265  if (hasVertexSave) vProdSave.rot(thetaIn, phiIn);}
266  void bst(double betaX, double betaY, double betaZ) {
267  pSave.bst(betaX, betaY, betaZ);
268  if (hasVertexSave) vProdSave.bst(betaX, betaY, betaZ);}
269  void bst(double betaX, double betaY, double betaZ, double gamma) {
270  pSave.bst(betaX, betaY, betaZ, gamma);
271  if (hasVertexSave) vProdSave.bst(betaX, betaY, betaZ, gamma);}
272  void bst(const Vec4& pBst) {pSave.bst(pBst);
273  if (hasVertexSave) vProdSave.bst(pBst);}
274  void bst(const Vec4& pBst, double mBst) {pSave.bst(pBst, mBst);
275  if (hasVertexSave) vProdSave.bst(pBst, mBst);}
276  void bstback(const Vec4& pBst) {pSave.bstback(pBst);
277  if (hasVertexSave) vProdSave.bstback(pBst);}
278  void bstback(const Vec4& pBst, double mBst) {pSave.bstback(pBst, mBst);
279  if (hasVertexSave) vProdSave.bstback(pBst, mBst);}
280  void rotbst(const RotBstMatrix& M) {pSave.rotbst(M);
281  if (hasVertexSave) vProdSave.rotbst(M);}
282  void offsetHistory( int minMother, int addMother, int minDaughter,
283  int addDaughter);
284  void offsetCol( int addCol);
285 
286 protected:
287 
288  // Constants: could only be changed in the code itself.
289  static const double TINY;
290 
291  // Properties of the current particle.
292  int idSave, statusSave, mother1Save, mother2Save, daughter1Save,
293  daughter2Save, colSave, acolSave;
294  Vec4 pSave;
295  double mSave, scaleSave, polSave;
296  bool hasVertexSave;
297  Vec4 vProdSave;
298  double tauSave;
299 
300  // Pointer to properties of the particle species.
301  // Should no be saved in a persistent copy of the event record.
302  // The //! below is ROOT notation that this member should not be saved.
303  // Event::restorePtrs() can be called to restore the missing information.
304  ParticleDataEntry* pdePtr;
305 
306  // Pointer to the whole event record to which the particle belongs (if any).
307  // As above it should not be saved.
308  Event* evtPtr;
309 
310 };
311 
312 // Invariant mass of a pair and its square.
313 // (Not part of class proper, but tightly linked.)
314 double m(const Particle&, const Particle&);
315 double m2(const Particle&, const Particle&);
316 
317 //==========================================================================
318 
319 // The junction class stores what kind of junction it is, the colour indices
320 // of the legs at the junction and as far out as legs have been traced,
321 // and the status codes assigned for fragmentation of each leg.
322 
323 class Junction {
324 
325 public:
326 
327  // Constructors.
328  Junction() : remainsSave(true), kindSave(0) {
329  for (int j = 0; j < 3; ++j) {
330  colSave[j] = 0; endColSave[j] = 0; statusSave[j] = 0; } }
331  Junction( int kindIn, int col0In, int col1In, int col2In)
332  : remainsSave(true), kindSave(kindIn) {colSave[0] = col0In;
333  colSave[1] = col1In; colSave[2] = col2In;
334  for (int j = 0; j < 3; ++j) {
335  endColSave[j] = colSave[j]; statusSave[j] = 0; } }
336  Junction(const Junction& ju) : remainsSave(ju.remainsSave),
337  kindSave(ju.kindSave) { for (int j = 0; j < 3; ++j) {
338  colSave[j] = ju.colSave[j]; endColSave[j] = ju.endColSave[j];
339  statusSave[j] = ju.statusSave[j]; } }
340  Junction& operator=(const Junction& ju) {if (this != &ju) {
341  remainsSave = ju.remainsSave; kindSave = ju.kindSave;
342  for (int j = 0; j < 3; ++j) { colSave[j] = ju.colSave[j];
343  endColSave[j] = ju.endColSave[j]; statusSave[j] = ju.statusSave[j]; } }
344  return *this; }
345 
346  // Set values.
347  void remains(bool remainsIn) {remainsSave = remainsIn;}
348  void col(int j, int colIn) {colSave[j] = colIn; endColSave[j] = colIn;}
349  void cols(int j, int colIn, int endColIn) {colSave[j] = colIn;
350  endColSave[j] = endColIn;}
351  void endCol(int j, int endColIn) {endColSave[j] = endColIn;}
352  void status(int j, int statusIn) {statusSave[j] = statusIn;}
353 
354  // Read out value.
355  bool remains() const {return remainsSave;}
356  int kind() const {return kindSave;}
357  int col(int j) const {return colSave[j];}
358  int endCol(int j) const {return endColSave[j];}
359  int status(int j) const {return statusSave[j];}
360 
361 private:
362 
363  // Kind, positions of the three ends and their status codes.
364  bool remainsSave;
365  int kindSave, colSave[3], endColSave[3], statusSave[3];
366 
367 };
368 
369 //==========================================================================
370 
371 // The Event class holds all info on the generated event.
372 
373 class Event {
374 
375 public:
376 
377  // Constructors.
378  Event(int capacity = 100) : startColTag(100), maxColTag(100),
379  savedSize(0), savedJunctionSize(0), savedPartonLevelSize(0),
380  scaleSave(0.), scaleSecondSave(0.),
381  headerList("----------------------------------------"),
382  particleDataPtr(0) { entry.reserve(capacity); }
383  Event& operator=(const Event& oldEvent);
384  Event(const Event& oldEvent) {*this = oldEvent;}
385 
386  // Initialize header for event listing, particle data table, and colour.
387  void init( string headerIn = "", ParticleData* particleDataPtrIn = 0,
388  int startColTagIn = 100) {
389  headerList.replace(0, headerIn.length() + 2, headerIn + " ");
390  particleDataPtr = particleDataPtrIn; startColTag = startColTagIn;}
391 
392  // Clear event record.
393  void clear() {entry.resize(0); maxColTag = startColTag;
394  savedPartonLevelSize = 0; scaleSave = 0.; scaleSecondSave = 0.;
395  clearJunctions();}
396 
397  // Clear event record, and set first particle empty.
398  void reset() {clear(); append(90, -11, 0, 0, 0., 0., 0., 0., 0.);}
399 
400  // Overload index operator to access element of event record.
401  Particle& operator[](int i) {return entry[i];}
402  const Particle& operator[](int i) const {return entry[i];}
403 
404  // Implement standard references to elements in the particle array.
405  Particle& front() {return entry.front();}
406  Particle& at(int i) {return entry.at(i);}
407  Particle& back() {return entry.back();}
408 
409  // Event record size.
410  int size() const {return entry.size();}
411 
412  // Put a new particle at the end of the event record; return index.
413  int append(Particle entryIn) {
414  entry.push_back(entryIn); setEvtPtr();
415  if (entryIn.col() > maxColTag) maxColTag = entryIn.col();
416  if (entryIn.acol() > maxColTag) maxColTag = entryIn.acol();
417  return entry.size() - 1;
418  }
419  int append(int id, int status, int mother1, int mother2, int daughter1,
420  int daughter2, int col, int acol, double px, double py, double pz,
421  double e, double m = 0., double scaleIn = 0., double polIn = 9.) {
422  entry.push_back( Particle(id, status, mother1, mother2, daughter1,
423  daughter2, col, acol, px, py, pz, e, m, scaleIn, polIn) ); setEvtPtr();
424  if (col > maxColTag) maxColTag = col;
425  if (acol > maxColTag) maxColTag = acol;
426  return entry.size() - 1;
427  }
428  int append(int id, int status, int mother1, int mother2, int daughter1,
429  int daughter2, int col, int acol, Vec4 p, double m = 0.,
430  double scaleIn = 0., double polIn = 9.) {
431  entry.push_back( Particle(id, status, mother1, mother2, daughter1,
432  daughter2, col, acol, p, m, scaleIn, polIn) ); setEvtPtr();
433  if (col > maxColTag) maxColTag = col;
434  if (acol > maxColTag) maxColTag = acol;
435  return entry.size() - 1;
436  }
437 
438  // Brief versions of append: no mothers and no daughters.
439  int append(int id, int status, int col, int acol, double px, double py,
440  double pz, double e, double m = 0., double scaleIn = 0.,
441  double polIn = 9.) { entry.push_back( Particle(id, status, 0, 0, 0, 0,
442  col, acol, px, py, pz, e, m, scaleIn, polIn) ); setEvtPtr();
443  if (col > maxColTag) maxColTag = col;
444  if (acol > maxColTag) maxColTag = acol;
445  return entry.size() - 1;
446  }
447  int append(int id, int status, int col, int acol, Vec4 p, double m = 0.,
448  double scaleIn = 0., double polIn = 9.) {entry.push_back( Particle(id,
449  status, 0, 0, 0, 0, col, acol, p, m, scaleIn, polIn) ); setEvtPtr();
450  if (col > maxColTag) maxColTag = col;
451  if (acol > maxColTag) maxColTag = acol;
452  return entry.size() - 1;
453  }
454 
455  // Set pointer to the event for a particle, by default latest one.
456  void setEvtPtr(int iSet = -1) {if (iSet < 0) iSet = entry.size() - 1;
457  entry[iSet].setEvtPtr( this);}
458 
459  // Add a copy of an existing particle at the end of the event record.
460  int copy(int iCopy, int newStatus = 0);
461 
462  // List the particles in an event.
463  void list(bool showScaleAndVertex = false,
464  bool showMothersAndDaughters = false, int precision = 3) const;
465 
466  // Remove last n entries.
467  void popBack(int nRemove = 1) { if (nRemove ==1) entry.pop_back();
468  else {int newSize = max( 0, size() - nRemove);
469  entry.resize(newSize);} }
470 
471  // Remove entries from iFirst to iLast, including endpoints, anf fix history.
472  // (To the extent possible; history pointers in removed range are zeroed.)
473  void remove(int iFirst, int iLast, bool shiftHistory = true);
474 
475  // Restore all ParticleDataEntry* pointers in the Particle vector.
476  // Useful when a persistent copy of the event record is read back in.
477  void restorePtrs() { for (int i = 0; i < size(); ++i) setEvtPtr(i); }
478 
479  // Save or restore the size of the event record (throwing at the end).
480  void saveSize() {savedSize = entry.size();}
481  void restoreSize() {entry.resize(savedSize);}
482  int savedSizeValue() {return savedSize;}
483 
484  // Initialize and access colour tag information.
485  void initColTag(int colTag = 0) {maxColTag = max( colTag,startColTag);}
486  int lastColTag() const {return maxColTag;}
487  int nextColTag() {return ++maxColTag;}
488 
489  // Access scale for which event as a whole is defined.
490  void scale( double scaleIn) {scaleSave = scaleIn;}
491  double scale() const {return scaleSave;}
492 
493  // Need a second scale if two hard interactions in event.
494  void scaleSecond( double scaleSecondIn) {scaleSecondSave = scaleSecondIn;}
495  double scaleSecond() const {return scaleSecondSave;}
496 
497  // Find complete list of daughters.
498  // Note: temporarily retained for CMS compatibility. Do not use!
499  vector<int> daughterList(int i) const {return entry[i].daughterList();}
500 
501  // Member functions for rotations and boosts of an event.
502  void rot(double theta, double phi)
503  {for (int i = 0; i < size(); ++i) entry[i].rot(theta, phi);}
504  void bst(double betaX, double betaY, double betaZ)
505  {for (int i = 0; i < size(); ++i) entry[i].bst(betaX, betaY, betaZ);}
506  void bst(double betaX, double betaY, double betaZ, double gamma)
507  {for (int i = 0; i < size(); ++i) entry[i].bst(betaX, betaY, betaZ,
508  gamma);}
509  void bst(const Vec4& vec)
510  {for (int i = 0; i < size(); ++i) entry[i].bst(vec);}
511  void rotbst(const RotBstMatrix& M)
512  {for (int i = 0; i < size(); ++i) entry[i].rotbst(M);}
513 
514  // Clear the list of junctions.
515  void clearJunctions() {junction.resize(0);}
516 
517  // Add a junction to the list, study it or extra input.
518  int appendJunction( int kind, int col0, int col1, int col2)
519  { junction.push_back( Junction( kind, col0, col1, col2) );
520  return junction.size() - 1;}
521  int appendJunction(Junction junctionIn) {junction.push_back(junctionIn);
522  return junction.size() - 1;}
523  int sizeJunction() const {return junction.size();}
524  bool remainsJunction(int i) const {return junction[i].remains();}
525  void remainsJunction(int i, bool remainsIn) {junction[i].remains(remainsIn);}
526  int kindJunction(int i) const {return junction[i].kind();}
527  int colJunction( int i, int j) const {return junction[i].col(j);}
528  void colJunction( int i, int j, int colIn) {junction[i].col(j, colIn);}
529  int endColJunction( int i, int j) const {return junction[i].endCol(j);}
530  void endColJunction( int i, int j, int endColIn)
531  {junction[i].endCol(j, endColIn);}
532  int statusJunction( int i, int j) const {return junction[i].status(j);}
533  void statusJunction( int i, int j, int statusIn)
534  {junction[i].status(j, statusIn);}
535  Junction& getJunction(int i) {return junction[i];}
536  const Junction& getJunction(int i) const {return junction[i];}
537  void eraseJunction(int i);
538 
539  // Save or restore the size of the junction list (throwing at the end).
540  void saveJunctionSize() {savedJunctionSize = junction.size();}
541  void restoreJunctionSize() {junction.resize(savedJunctionSize);}
542 
543  // List any junctions in the event; for debug mainly.
544  void listJunctions() const;
545 
546  // Save event record size at Parton Level, i.e. before hadronization.
547  void savePartonLevelSize() {savedPartonLevelSize = entry.size();}
548 
549  // Operator overloading allows to append one event to an existing one.
550  // Warning: particles should be OK, but some other information unreliable.
551  Event& operator+=(const Event& addEvent);
552 
553 private:
554 
555  // The Particle class needs to access particle data.
556  friend class Particle;
557 
558  // Constants: could only be changed in the code itself.
559  static const int IPERLINE;
560 
561  // Initialization data, normally only set once.
562  int startColTag;
563 
564  // The event: a vector containing all particles (entries).
565  // The explicit use of Pythia8:: qualifier patches a limitation in ROOT.
566  vector<Pythia8::Particle> entry;
567 
568  // The list of junctions.
569  // The explicit use of Pythia8:: qualifier patches a limitation in ROOT.
570  vector<Pythia8::Junction> junction;
571 
572  // The maximum colour tag of the event so far.
573  int maxColTag;
574 
575  // Saved entry and junction list sizes, for simple restoration.
576  int savedSize, savedJunctionSize, savedPartonLevelSize;
577 
578  // The scale of the event; linear quantity in GeV.
579  double scaleSave, scaleSecondSave;
580 
581  // Header specification in event listing (at most 40 characters wide).
582  string headerList;
583 
584  // Pointer to the particle data table.
585  // The //! below is ROOT notation that this member should not be saved.
586  ParticleData* particleDataPtr;
587 
588 };
589 
590 //==========================================================================
591 
592 } // end namespace Pythia8
593 
594 #endif // end Pythia8_Event_H
Definition: AgUStep.h:26