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) 2014 Torbjorn Sjostrand.
3 // PYTHIA is licenced under the GNU GPL version 2, 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  bool isRescatteredIncoming() const {return
152  (statusSave == -34 || statusSave == -45 ||
153  statusSave == -46 || statusSave == -54);}
154 
155  // Member functions for output; derived double quantities.
156  double m2() const {return (mSave >= 0.) ? mSave*mSave
157  : -mSave*mSave;}
158  double mCalc() const {return pSave.mCalc();}
159  double m2Calc() const {return pSave.m2Calc();}
160  double eCalc() const {return sqrt(abs(m2() + pSave.pAbs2()));}
161  double pT() const {return pSave.pT();}
162  double pT2() const {return pSave.pT2();}
163  double mT() const {double temp = m2() + pSave.pT2();
164  return (temp >= 0.) ? sqrt(temp) : -sqrt(-temp);}
165  double mT2() const {return m2() + pSave.pT2();}
166  double pAbs() const {return pSave.pAbs();}
167  double pAbs2() const {return pSave.pAbs2();}
168  double eT() const {return pSave.eT();}
169  double eT2() const {return pSave.eT2();}
170  double theta() const {return pSave.theta();}
171  double phi() const {return pSave.phi();}
172  double thetaXZ() const {return pSave.thetaXZ();}
173  double pPos() const {return pSave.pPos();}
174  double pNeg() const {return pSave.pNeg();}
175  double y() const;
176  double eta() const;
177  Vec4 vDec() const {return (tauSave > 0. && mSave > 0.)
178  ? vProdSave + tauSave * pSave / mSave : vProdSave;}
179  double xDec() const {return (tauSave > 0. && mSave > 0.)
180  ? vProdSave.px() + tauSave * pSave.px() / mSave : vProdSave.px();}
181  double yDec() const {return (tauSave > 0. && mSave > 0.)
182  ? vProdSave.py() + tauSave * pSave.py() / mSave : vProdSave.py();}
183  double zDec() const {return (tauSave > 0. && mSave > 0.)
184  ? vProdSave.pz() + tauSave * pSave.pz() / mSave : vProdSave.pz();}
185  double tDec() const {return (tauSave > 0. && mSave > 0.)
186  ? vProdSave.e() + tauSave * pSave.e() / mSave : vProdSave.e();}
187 
188  // Methods that can refer back to the event the particle belongs to.
189  int index() const;
190  int statusHepMC() const;
191  int iTopCopy() const;
192  int iBotCopy() const;
193  int iTopCopyId() const;
194  int iBotCopyId() const;
195  vector<int> motherList() const;
196  vector<int> daughterList() const;
197  vector<int> sisterList(bool traceTopBot = false) const;
198  bool isAncestor(int iAncestor) const;
199  bool undoDecay();
200 
201  // Further output, based on a pointer to a ParticleDataEntry object.
202  string name() const {
203  return (pdePtr != 0) ? pdePtr->name(idSave) : " ";}
204  string nameWithStatus(int maxLen = 20) const;
205  int spinType() const {
206  return (pdePtr != 0) ? pdePtr->spinType() : 0;}
207  int chargeType() const {
208  return (pdePtr != 0) ? pdePtr->chargeType(idSave) : 0;}
209  double charge() const {
210  return (pdePtr != 0) ? pdePtr->charge(idSave) : 0;}
211  bool isCharged() const {
212  return (pdePtr != 0) ? (pdePtr->chargeType(idSave) != 0) : false;}
213  bool isNeutral() const {
214  return (pdePtr != 0) ? (pdePtr->chargeType(idSave) == 0) : false;}
215  int colType() const {
216  return (pdePtr != 0) ? pdePtr->colType(idSave) : 0;}
217  double m0() const {
218  return (pdePtr != 0) ? pdePtr->m0() : 0.;}
219  double mWidth() const {
220  return (pdePtr != 0) ? pdePtr->mWidth() : 0.;}
221  double mMin() const {
222  return (pdePtr != 0) ? pdePtr->mMin() : 0.;}
223  double mMax() const {
224  return (pdePtr != 0) ? pdePtr->mMax() : 0.;}
225  double mSel() const {
226  return (pdePtr != 0) ? pdePtr->mSel() : 0.;}
227  double constituentMass() const {
228  return (pdePtr != 0) ? pdePtr->constituentMass() : 0.;}
229  double tau0() const {
230  return (pdePtr != 0) ? pdePtr->tau0() : 0.;}
231  bool mayDecay() const {
232  return (pdePtr != 0) ? pdePtr->mayDecay() : false;}
233  bool canDecay() const {
234  return (pdePtr != 0) ? pdePtr->canDecay() : false;}
235  bool doExternalDecay() const {
236  return (pdePtr != 0) ? pdePtr->doExternalDecay() : false;}
237  bool isResonance() const {
238  return (pdePtr != 0) ? pdePtr->isResonance() : false;}
239  bool isVisible() const {
240  return (pdePtr != 0) ? pdePtr->isVisible() : false;}
241  bool isLepton() const {
242  return (pdePtr != 0) ? pdePtr->isLepton() : false;}
243  bool isQuark() const {
244  return (pdePtr != 0) ? pdePtr->isQuark() : false;}
245  bool isGluon() const {
246  return (pdePtr != 0) ? pdePtr->isGluon() : false;}
247  bool isDiquark() const {
248  return (pdePtr != 0) ? pdePtr->isDiquark() : false;}
249  bool isParton() const {
250  return (pdePtr != 0) ? pdePtr->isParton() : false;}
251  bool isHadron() const {
252  return (pdePtr != 0) ? pdePtr->isHadron() : false;}
253  ParticleDataEntry& particleDataEntry() const {return *pdePtr;}
254 
255  // Member functions that perform operations.
256  void rescale3(double fac) {pSave.rescale3(fac);}
257  void rescale4(double fac) {pSave.rescale4(fac);}
258  void rescale5(double fac) {pSave.rescale4(fac); mSave *= fac;}
259  void rot(double thetaIn, double phiIn) {pSave.rot(thetaIn, phiIn);
260  if (hasVertexSave) vProdSave.rot(thetaIn, phiIn);}
261  void bst(double betaX, double betaY, double betaZ) {
262  pSave.bst(betaX, betaY, betaZ);
263  if (hasVertexSave) vProdSave.bst(betaX, betaY, betaZ);}
264  void bst(double betaX, double betaY, double betaZ, double gamma) {
265  pSave.bst(betaX, betaY, betaZ, gamma);
266  if (hasVertexSave) vProdSave.bst(betaX, betaY, betaZ, gamma);}
267  void bst(const Vec4& pBst) {pSave.bst(pBst);
268  if (hasVertexSave) vProdSave.bst(pBst);}
269  void bst(const Vec4& pBst, double mBst) {pSave.bst(pBst, mBst);
270  if (hasVertexSave) vProdSave.bst(pBst, mBst);}
271  void bstback(const Vec4& pBst) {pSave.bstback(pBst);
272  if (hasVertexSave) vProdSave.bstback(pBst);}
273  void bstback(const Vec4& pBst, double mBst) {pSave.bstback(pBst, mBst);
274  if (hasVertexSave) vProdSave.bstback(pBst, mBst);}
275  void rotbst(const RotBstMatrix& M) {pSave.rotbst(M);
276  if (hasVertexSave) vProdSave.rotbst(M);}
277  void offsetHistory( int minMother, int addMother, int minDaughter,
278  int addDaughter);
279  void offsetCol( int addCol);
280 
281 private:
282 
283  // Constants: could only be changed in the code itself.
284  static const double TINY;
285 
286  // Properties of the current particle.
287  int idSave, statusSave, mother1Save, mother2Save, daughter1Save,
288  daughter2Save, colSave, acolSave;
289  Vec4 pSave;
290  double mSave, scaleSave, polSave;
291  bool hasVertexSave;
292  Vec4 vProdSave;
293  double tauSave;
294 
295  // Pointer to properties of the particle species.
296  // Should no be saved in a persistent copy of the event record.
297  // The //! below is ROOT notation that this member should not be saved.
298  // Event::restorePtrs() can be called to restore the missing information.
299  ParticleDataEntry* pdePtr;
300 
301  // Pointer to the whole event record to which the particle belongs (if any).
302  // As above it should not be saved.
303  Event* evtPtr;
304 
305 };
306 
307 // Invariant mass of a pair and its square.
308 // (Not part of class proper, but tightly linked.)
309 double m(const Particle&, const Particle&);
310 double m2(const Particle&, const Particle&);
311 
312 //==========================================================================
313 
314 // The junction class stores what kind of junction it is, the colour indices
315 // of the legs at the junction and as far out as legs have been traced,
316 // and the status codes assigned for fragmentation of each leg.
317 
318 class Junction {
319 
320 public:
321 
322  // Constructors.
323  Junction() : remainsSave(true), kindSave(0) {
324  for (int j = 0; j < 3; ++j) {
325  colSave[j] = 0; endColSave[j] = 0; statusSave[j] = 0; } }
326  Junction( int kindIn, int col0In, int col1In, int col2In)
327  : remainsSave(true), kindSave(kindIn) {colSave[0] = col0In;
328  colSave[1] = col1In; colSave[2] = col2In;
329  for (int j = 0; j < 3; ++j) {
330  endColSave[j] = colSave[j]; statusSave[j] = 0; } }
331  Junction(const Junction& ju) : remainsSave(ju.remainsSave),
332  kindSave(ju.kindSave) { for (int j = 0; j < 3; ++j) {
333  colSave[j] = ju.colSave[j]; endColSave[j] = ju.endColSave[j];
334  statusSave[j] = ju.statusSave[j]; } }
335  Junction& operator=(const Junction& ju) {if (this != &ju) {
336  remainsSave = ju.remainsSave; kindSave = ju.kindSave;
337  for (int j = 0; j < 3; ++j) { colSave[j] = ju.colSave[j];
338  endColSave[j] = ju.endColSave[j]; statusSave[j] = ju.statusSave[j]; } }
339  return *this; }
340 
341  // Set values.
342  void remains(bool remainsIn) {remainsSave = remainsIn;}
343  void col(int j, int colIn) {colSave[j] = colIn; endColSave[j] = colIn;}
344  void cols(int j, int colIn, int endColIn) {colSave[j] = colIn;
345  endColSave[j] = endColIn;}
346  void endCol(int j, int endColIn) {endColSave[j] = endColIn;}
347  void status(int j, int statusIn) {statusSave[j] = statusIn;}
348 
349  // Read out value.
350  bool remains() const {return remainsSave;}
351  int kind() const {return kindSave;}
352  int col(int j) const {return colSave[j];}
353  int endCol(int j) const {return endColSave[j];}
354  int status(int j) const {return statusSave[j];}
355 
356 private:
357 
358  // Kind, positions of the three ends and their status codes.
359  bool remainsSave;
360  int kindSave, colSave[3], endColSave[3], statusSave[3];
361 
362 };
363 
364 //==========================================================================
365 
366 // The Event class holds all info on the generated event.
367 
368 class Event {
369 
370 public:
371 
372  // Constructors.
373  Event(int capacity = 100) : startColTag(100), maxColTag(100),
374  savedSize(0), savedJunctionSize(0), scaleSave(0.), scaleSecondSave(0.),
375  headerList("----------------------------------------"),
376  particleDataPtr(0) { entry.reserve(capacity); }
377  Event& operator=(const Event& oldEvent);
378  Event(const Event& oldEvent) {*this = oldEvent;}
379 
380  // Initialize header for event listing, particle data table, and colour.
381  void init( string headerIn = "", ParticleData* particleDataPtrIn = 0,
382  int startColTagIn = 100) {
383  headerList.replace(0, headerIn.length() + 2, headerIn + " ");
384  particleDataPtr = particleDataPtrIn; startColTag = startColTagIn;}
385 
386  // Clear event record.
387  void clear() {entry.resize(0); maxColTag = startColTag; scaleSave = 0.;
388  scaleSecondSave = 0.; clearJunctions();}
389 
390  // Clear event record, and set first particle empty.
391  void reset() {clear(); append(90, -11, 0, 0, 0., 0., 0., 0., 0.);}
392 
393  // Overload index operator to access element of event record.
394  Particle& operator[](int i) {return entry[i];}
395  const Particle& operator[](int i) const {return entry[i];}
396 
397  // Implement standard references to elements in the particle array.
398  Particle& front() {return entry.front();}
399  Particle& at(int i) {return entry.at(i);}
400  Particle& back() {return entry.back();}
401 
402  // Event record size.
403  int size() const {return entry.size();}
404 
405  // Put a new particle at the end of the event record; return index.
406  int append(Particle entryIn) {
407  entry.push_back(entryIn); setEvtPtr();
408  if (entryIn.col() > maxColTag) maxColTag = entryIn.col();
409  if (entryIn.acol() > maxColTag) maxColTag = entryIn.acol();
410  return entry.size() - 1;
411  }
412  int append(int id, int status, int mother1, int mother2, int daughter1,
413  int daughter2, int col, int acol, double px, double py, double pz,
414  double e, double m = 0., double scaleIn = 0., double polIn = 9.) {
415  entry.push_back( Particle(id, status, mother1, mother2, daughter1,
416  daughter2, col, acol, px, py, pz, e, m, scaleIn, polIn) ); setEvtPtr();
417  if (col > maxColTag) maxColTag = col;
418  if (acol > maxColTag) maxColTag = acol;
419  return entry.size() - 1;
420  }
421  int append(int id, int status, int mother1, int mother2, int daughter1,
422  int daughter2, int col, int acol, Vec4 p, double m = 0.,
423  double scaleIn = 0., double polIn = 9.) {
424  entry.push_back( Particle(id, status, mother1, mother2, daughter1,
425  daughter2, col, acol, p, m, scaleIn, polIn) ); setEvtPtr();
426  if (col > maxColTag) maxColTag = col;
427  if (acol > maxColTag) maxColTag = acol;
428  return entry.size() - 1;
429  }
430 
431  // Brief versions of append: no mothers and no daughters.
432  int append(int id, int status, int col, int acol, double px, double py,
433  double pz, double e, double m = 0., double scaleIn = 0.,
434  double polIn = 9.) { entry.push_back( Particle(id, status, 0, 0, 0, 0,
435  col, acol, px, py, pz, e, m, scaleIn, polIn) ); setEvtPtr();
436  if (col > maxColTag) maxColTag = col;
437  if (acol > maxColTag) maxColTag = acol;
438  return entry.size() - 1;
439  }
440  int append(int id, int status, int col, int acol, Vec4 p, double m = 0.,
441  double scaleIn = 0., double polIn = 9.) {entry.push_back( Particle(id,
442  status, 0, 0, 0, 0, col, acol, p, m, scaleIn, polIn) ); setEvtPtr();
443  if (col > maxColTag) maxColTag = col;
444  if (acol > maxColTag) maxColTag = acol;
445  return entry.size() - 1;
446  }
447 
448  // Set pointer to the event for a particle, by default latest one.
449  void setEvtPtr(int iSet = -1) {if (iSet < 0) iSet = entry.size() - 1;
450  entry[iSet].setEvtPtr( this);}
451 
452  // Add a copy of an existing particle at the end of the event record.
453  int copy(int iCopy, int newStatus = 0);
454 
455  // List the particles in an event.
456  void list() const;
457  void list(ostream& os) const;
458  void list(bool showScaleAndVertex, bool showMothersAndDaughters = false)
459  const;
460  void list(bool showScaleAndVertex, bool showMothersAndDaughters,
461  ostream& os) const;
462 
463  // Remove last n entries.
464  void popBack(int nRemove = 1) { if (nRemove ==1) entry.pop_back();
465  else {int newSize = max( 0, size() - nRemove);
466  entry.resize(newSize);} }
467 
468  // Remove entries from iFirst to iLast, including endpoints.
469  void remove(int iFirst, int iLast) {
470  if (iFirst < 0 || iLast >= int(entry.size()) || iLast < iFirst) return;
471  entry.erase( entry.begin() + iFirst, entry.begin() + iLast + 1);
472  }
473 
474  // Undo the decay of a single particle (where daughters well-defined).
475  bool undoDecay(int i);
476 
477  // Restore all ParticleDataEntry* pointers in the Particle vector.
478  // Useful when a persistent copy of the event record is read back in.
479  void restorePtrs() { for (int i = 0; i < size(); ++i) setEvtPtr(i); }
480 
481  // Save or restore the size of the event record (throwing at the end).
482  void saveSize() {savedSize = entry.size();}
483  void restoreSize() {entry.resize(savedSize);}
484  int savedSizeValue() {return savedSize;}
485 
486  // Initialize and access colour tag information.
487  void initColTag(int colTag = 0) {maxColTag = max( colTag,startColTag);}
488  int lastColTag() const {return maxColTag;}
489  int nextColTag() {return ++maxColTag;}
490 
491  // Access scale for which event as a whole is defined.
492  void scale( double scaleIn) {scaleSave = scaleIn;}
493  double scale() const {return scaleSave;}
494 
495  // Need a second scale if two hard interactions in event.
496  void scaleSecond( double scaleSecondIn) {scaleSecondSave = scaleSecondIn;}
497  double scaleSecond() const {return scaleSecondSave;}
498 
499  // Find complete list of daughters and mothers.
500  vector<int> motherList(int i) const;
501  vector<int> daughterList(int i) const;
502 
503  // Convert to HepMC status code conventions.
504  int statusHepMC(int i) const;
505 
506  // Trace the first and last copy of one and the same particle.
507  int iTopCopy(int i) const;
508  int iBotCopy(int i) const;
509 
510  // Trace the first and last copy of a particle, using flavour match.
511  int iTopCopyId(int i) const;
512  int iBotCopyId(int i) const;
513 
514  // Find list of sisters, also tracking up and down identical copies.
515  vector<int> sisterList(int i) const;
516  vector<int> sisterListTopBot(int i, bool widenSearch = true) const;
517 
518  // Check whether two particles have a direct mother-daughter relation.
519  bool isAncestor(int i, int iAncestor) const;
520 
521  // Member functions for rotations and boosts of an event.
522  void rot(double theta, double phi)
523  {for (int i = 0; i < size(); ++i) entry[i].rot(theta, phi);}
524  void bst(double betaX, double betaY, double betaZ)
525  {for (int i = 0; i < size(); ++i) entry[i].bst(betaX, betaY, betaZ);}
526  void bst(double betaX, double betaY, double betaZ, double gamma)
527  {for (int i = 0; i < size(); ++i) entry[i].bst(betaX, betaY, betaZ,
528  gamma);}
529  void bst(const Vec4& vec)
530  {for (int i = 0; i < size(); ++i) entry[i].bst(vec);}
531  void rotbst(const RotBstMatrix& M)
532  {for (int i = 0; i < size(); ++i) entry[i].rotbst(M);}
533 
534  // Clear the list of junctions.
535  void clearJunctions() {junction.resize(0);}
536 
537  // Add a junction to the list, study it or extra input.
538  int appendJunction( int kind, int col0, int col1, int col2)
539  { junction.push_back( Junction( kind, col0, col1, col2) );
540  return junction.size() - 1;}
541  int appendJunction(Junction junctionIn) {junction.push_back(junctionIn);
542  return junction.size() - 1;}
543  int sizeJunction() const {return junction.size();}
544  bool remainsJunction(int i) const {return junction[i].remains();}
545  void remainsJunction(int i, bool remainsIn) {junction[i].remains(remainsIn);}
546  int kindJunction(int i) const {return junction[i].kind();}
547  int colJunction( int i, int j) const {return junction[i].col(j);}
548  void colJunction( int i, int j, int colIn) {junction[i].col(j, colIn);}
549  int endColJunction( int i, int j) const {return junction[i].endCol(j);}
550  void endColJunction( int i, int j, int endColIn)
551  {junction[i].endCol(j, endColIn);}
552  int statusJunction( int i, int j) const {return junction[i].status(j);}
553  void statusJunction( int i, int j, int statusIn)
554  {junction[i].status(j, statusIn);}
555  Junction& getJunction(int i) {return junction[i];}
556  const Junction& getJunction(int i) const {return junction[i];}
557  void eraseJunction(int i);
558 
559  // Save or restore the size of the junction list (throwing at the end).
560  void saveJunctionSize() {savedJunctionSize = junction.size();}
561  void restoreJunctionSize() {junction.resize(savedJunctionSize);}
562 
563  // List any junctions in the event; for debug mainly.
564  void listJunctions(ostream& os = cout) const;
565 
566  // Operator overloading allows to append one event to an existing one.
567  // Warning: particles should be OK, but some other information unreliable.
568  Event& operator+=(const Event& addEvent);
569 
570 private:
571 
572  // The Particle class needs to access particle data.
573  friend class Particle;
574 
575  // Constants: could only be changed in the code itself.
576  static const int IPERLINE;
577 
578  // Initialization data, normally only set once.
579  int startColTag;
580 
581  // The event: a vector containing all particles (entries).
582  // The explicit use of Pythia8:: qualifier patches a limitation in ROOT.
583  vector<Pythia8::Particle> entry;
584 
585  // The list of junctions.
586  // The explicit use of Pythia8:: qualifier patches a limitation in ROOT.
587  vector<Pythia8::Junction> junction;
588 
589  // The maximum colour tag of the event so far.
590  int maxColTag;
591 
592  // Saved entry and junction list sizes, for simple restoration.
593  int savedSize, savedJunctionSize;
594 
595  // The scale of the event; linear quantity in GeV.
596  double scaleSave, scaleSecondSave;
597 
598  // Header specification in event listing (at most 40 characters wide).
599  string headerList;
600 
601  // Pointer to the particle data table.
602  // The //! below is ROOT notation that this member should not be saved.
603  ParticleData* particleDataPtr;
604 
605 };
606 
607 //==========================================================================
608 
609 } // end namespace Pythia8
610 
611 #endif // end Pythia8_Event_H
Definition: AgUStep.h:26