11 #ifndef Pythia8_Event_H 12 #define Pythia8_Event_H 14 #include "Pythia8/Basics.h" 15 #include "Pythia8/ParticleData.h" 16 #include "Pythia8/PythiaStdlib.h" 23 class ParticleDataEntry;
24 class ResonanceWidths;
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.),
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.)),
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.),
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),
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;
85 void setPDEPtr(ParticleDataEntryPtr pdePtrIn =
nullptr);
89 void status(
int statusIn) {statusSave = statusIn;}
90 void statusPos() {statusSave = abs(statusSave);}
91 void statusNeg() {statusSave = -abs(statusSave);}
92 void statusCode(
int statusIn) {statusSave =
93 (statusSave > 0) ? abs(statusIn) : -abs(statusIn);}
94 void mother1(
int mother1In) {mother1Save = mother1In;}
95 void mother2(
int mother2In) {mother2Save = mother2In;}
96 void mothers(
int mother1In = 0,
int mother2In = 0)
97 {mother1Save = mother1In; mother2Save = mother2In;}
98 void daughter1(
int daughter1In) {daughter1Save = daughter1In;}
99 void daughter2(
int daughter2In) {daughter2Save = daughter2In;}
100 void daughters(
int daughter1In = 0,
int daughter2In = 0)
101 {daughter1Save = daughter1In; daughter2Save = daughter2In;}
102 void col(
int colIn) {colSave = colIn;}
103 void acol(
int acolIn) {acolSave = acolIn;}
104 void cols(
int colIn = 0,
int acolIn = 0) {colSave = colIn;
106 void p(
Vec4 pIn) {pSave = pIn;}
107 void p(
double pxIn,
double pyIn,
double pzIn,
double eIn)
108 {pSave.p(pxIn, pyIn, pzIn, eIn);}
109 void px(
double pxIn) {pSave.px(pxIn);}
110 void py(
double pyIn) {pSave.py(pyIn);}
111 void pz(
double pzIn) {pSave.pz(pzIn);}
112 void e(
double eIn) {pSave.e(eIn);}
113 void m(
double mIn) {mSave = mIn;}
114 void scale(
double scaleIn) {scaleSave = scaleIn;}
115 void pol(
double polIn) {polSave = polIn;}
116 void vProd(
Vec4 vProdIn) {vProdSave = vProdIn; hasVertexSave =
true;}
117 void vProd(
double xProdIn,
double yProdIn,
double zProdIn,
double tProdIn)
118 {vProdSave.p(xProdIn, yProdIn, zProdIn, tProdIn); hasVertexSave =
true;}
119 void xProd(
double xProdIn) {vProdSave.px(xProdIn); hasVertexSave =
true;}
120 void yProd(
double yProdIn) {vProdSave.py(yProdIn); hasVertexSave =
true;}
121 void zProd(
double zProdIn) {vProdSave.pz(zProdIn); hasVertexSave =
true;}
122 void tProd(
double tProdIn) {vProdSave.e(tProdIn); hasVertexSave =
true;}
123 void vProdAdd(
Vec4 vProdIn) {vProdSave += vProdIn; hasVertexSave =
true;}
124 void tau(
double tauIn) {tauSave = tauIn;}
128 int status()
const {
return statusSave;}
129 int mother1()
const {
return mother1Save;}
130 int mother2()
const {
return mother2Save;}
131 int daughter1()
const {
return daughter1Save;}
132 int daughter2()
const {
return daughter2Save;}
133 int col()
const {
return colSave;}
134 int acol()
const {
return acolSave;}
135 Vec4 p()
const {
return pSave;}
136 double px()
const {
return pSave.px();}
137 double py()
const {
return pSave.py();}
138 double pz()
const {
return pSave.pz();}
139 double e()
const {
return pSave.e();}
140 double m()
const {
return mSave;}
141 double scale()
const {
return scaleSave;}
142 double pol()
const {
return polSave;}
143 bool hasVertex()
const {
return hasVertexSave;}
144 Vec4 vProd()
const {
return vProdSave;}
145 double xProd()
const {
return vProdSave.px();}
146 double yProd()
const {
return vProdSave.py();}
147 double zProd()
const {
return vProdSave.pz();}
148 double tProd()
const {
return vProdSave.e();}
149 double tau()
const {
return tauSave;}
153 int statusAbs()
const {
return abs(statusSave);}
154 bool isFinal()
const {
return (statusSave > 0);}
156 bool isRescatteredIncoming()
const {
return 157 (statusSave == -34 || statusSave == -45 ||
158 statusSave == -46 || statusSave == -54);}
161 double m2()
const {
return (mSave >= 0.) ? mSave*mSave
163 double mCalc()
const {
return pSave.mCalc();}
164 double m2Calc()
const {
return pSave.m2Calc();}
165 double eCalc()
const {
return sqrt(abs(
m2() + pSave.pAbs2()));}
166 double pT()
const {
return pSave.pT();}
167 double pT2()
const {
return pSave.pT2();}
168 double mT()
const {
double temp =
m2() + pSave.pT2();
169 return (temp >= 0.) ? sqrt(temp) : -sqrt(-temp);}
170 double mT2()
const {
return m2() + pSave.pT2();}
171 double pAbs()
const {
return pSave.pAbs();}
172 double pAbs2()
const {
return pSave.pAbs2();}
173 double eT()
const {
return pSave.eT();}
174 double eT2()
const {
return pSave.eT2();}
175 double theta()
const {
return pSave.theta();}
176 double phi()
const {
return pSave.phi();}
177 double thetaXZ()
const {
return pSave.thetaXZ();}
178 double pPos()
const {
return pSave.pPos();}
179 double pNeg()
const {
return pSave.pNeg();}
182 double y(
double mCut)
const;
184 Vec4 vDec()
const {
return (tauSave > 0. && mSave > 0.)
185 ? vProdSave + tauSave * pSave / mSave : vProdSave;}
186 double xDec()
const {
return (tauSave > 0. && mSave > 0.)
187 ? vProdSave.px() + tauSave * pSave.px() / mSave : vProdSave.px();}
188 double yDec()
const {
return (tauSave > 0. && mSave > 0.)
189 ? vProdSave.py() + tauSave * pSave.py() / mSave : vProdSave.py();}
190 double zDec()
const {
return (tauSave > 0. && mSave > 0.)
191 ? vProdSave.pz() + tauSave * pSave.pz() / mSave : vProdSave.pz();}
192 double tDec()
const {
return (tauSave > 0. && mSave > 0.)
193 ? vProdSave.e() + tauSave * pSave.e() / mSave : vProdSave.e();}
196 virtual int index()
const;
198 int iBotCopy()
const;
204 vector<int>
sisterList(
bool traceTopBot =
false)
const;
213 void colHV(
int colHVin);
214 void acolHV(
int acolHVin);
215 void colsHV(
int colHVin,
int acolHVin);
221 int spinType()
const {
223 int chargeType()
const {
225 double charge()
const {
227 bool isCharged()
const {
229 bool isNeutral()
const {
231 int colType()
const {
235 double mWidth()
const {
237 double mMin()
const {
239 double mMax()
const {
241 double mSel()
const {
243 double constituentMass()
const {
245 double tau0()
const {
247 bool mayDecay()
const {
249 bool canDecay()
const {
251 bool doExternalDecay()
const {
252 return (
pdePtr != 0) ?
pdePtr->doExternalDecay() :
false;}
253 bool isResonance()
const {
255 bool isVisible()
const {
257 bool isLepton()
const {
259 bool isQuark()
const {
261 bool isGluon()
const {
263 bool isDiquark()
const {
265 bool isParton()
const {
267 bool isHadron()
const {
269 bool isExotic()
const {
275 void rescale4(
double fac) {pSave.rescale4(fac);}
276 void rescale5(
double fac) {pSave.rescale4(fac); mSave *= fac;}
277 void rot(
double thetaIn,
double phiIn) {pSave.
rot(thetaIn, phiIn);
278 if (hasVertexSave) vProdSave.
rot(thetaIn, phiIn);}
279 void bst(
double betaX,
double betaY,
double betaZ) {
280 pSave.
bst(betaX, betaY, betaZ);
281 if (hasVertexSave) vProdSave.
bst(betaX, betaY, betaZ);}
282 void bst(
double betaX,
double betaY,
double betaZ,
double gamma) {
283 pSave.
bst(betaX, betaY, betaZ, gamma);
284 if (hasVertexSave) vProdSave.
bst(betaX, betaY, betaZ, gamma);}
285 void bst(
const Vec4& pBst) {pSave.
bst(pBst);
286 if (hasVertexSave) vProdSave.
bst(pBst);}
287 void bst(
const Vec4& pBst,
double mBst) {pSave.
bst(pBst, mBst);
288 if (hasVertexSave) vProdSave.
bst(pBst, mBst);}
289 void bstback(
const Vec4& pBst) {pSave.
bstback(pBst);
290 if (hasVertexSave) vProdSave.
bstback(pBst);}
291 void bstback(
const Vec4& pBst,
double mBst) {pSave.
bstback(pBst, mBst);
292 if (hasVertexSave) vProdSave.
bstback(pBst, mBst);}
294 if (hasVertexSave && boostVertex) vProdSave.
rotbst(M);}
295 void offsetHistory(
int minMother,
int addMother,
int minDaughter,
305 int idSave, statusSave, mother1Save, mother2Save, daughter1Save,
306 daughter2Save, colSave, acolSave;
308 double mSave, scaleSave, polSave;
343 Junction() : remainsSave(true), kindSave(0), colSave(), endColSave(),
346 Junction(
int kindIn,
int col0In,
int col1In,
int col2In)
347 : remainsSave(
true), kindSave(kindIn), colSave(), endColSave(),
349 colSave[0] = col0In; colSave[1] = col1In; colSave[2] = col2In;
350 for (
int j = 0; j < 3; ++j) {
351 endColSave[j] = colSave[j]; } }
353 kindSave(ju.kindSave), colSave(), endColSave(), statusSave() {
354 for (
int j = 0; j < 3; ++j) {
355 colSave[j] = ju.colSave[j]; endColSave[j] = ju.endColSave[j];
356 statusSave[j] = ju.statusSave[j]; } }
358 remainsSave = ju.remainsSave; kindSave = ju.kindSave;
359 for (
int j = 0; j < 3; ++j) { colSave[j] = ju.colSave[j];
360 endColSave[j] = ju.endColSave[j]; statusSave[j] = ju.statusSave[j]; } }
364 void remains(
bool remainsIn) {remainsSave = remainsIn;}
365 void col(
int j,
int colIn) {colSave[j] = colIn; endColSave[j] = colIn;}
366 void cols(
int j,
int colIn,
int endColIn) {colSave[j] = colIn;
367 endColSave[j] = endColIn;}
368 void endCol(
int j,
int endColIn) {endColSave[j] = endColIn;}
369 void status(
int j,
int statusIn) {statusSave[j] = statusIn;}
373 int kind()
const {
return kindSave;}
374 int col(
int j)
const {
return colSave[j];}
375 int endCol(
int j)
const {
return endColSave[j];}
376 int status(
int j)
const {
return statusSave[j];}
382 int kindSave, colSave[3], endColSave[3], statusSave[3];
396 HVcols(
int iHVin,
int colHVin,
int acolHVin) : iHV(iHVin),
397 colHV(colHVin), acolHV(acolHVin) {}
413 Event(
int capacity = 100) : startColTag(100), maxColTag(100),
414 savedSize(0), savedJunctionSize(0), savedHVcolsSize(0),
415 savedPartonLevelSize(0), scaleSave(0.), scaleSecondSave(0.),
416 headerList(
"----------------------------------------"),
417 particleDataPtr(0) { entry.reserve(capacity); }
419 Event(
const Event& oldEvent) {*
this = oldEvent;}
423 int startColTagIn = 100) {
424 headerList.replace(0, headerIn.length() + 2, headerIn +
" ");
425 particleDataPtr = particleDataPtrIn; startColTag = startColTagIn;}
428 void clear() {entry.resize(0); maxColTag = startColTag;
429 savedPartonLevelSize = 0; scaleSave = 0.; scaleSecondSave = 0.;
430 clearJunctions(); clearHV();}
431 void free() {vector<Particle>().swap(entry); maxColTag = startColTag;
432 savedPartonLevelSize = 0; scaleSave = 0.; scaleSecondSave = 0.;
433 clearJunctions(); clearHV();}
436 void reset() {clear(); append(90, -11, 0, 0, 0., 0., 0., 0., 0.);}
440 const Particle& operator[](
int i)
const {
return entry.at(i);}
444 Particle& at(
int i) {
return entry.at(i);}
445 Particle& back() {
return entry.back();}
446 const Particle& front()
const {
return entry.front();}
447 const Particle& at(
int i)
const {
return entry.at(i);}
448 const Particle& back()
const {
return entry.back();}
451 vector<Pythia8::Particle>::iterator
begin() {
return entry.begin(); }
452 vector<Pythia8::Particle>::iterator end() {
return entry.end(); }
453 vector<Pythia8::Particle>::const_iterator begin()
const 454 {
return entry.begin(); }
455 vector<Pythia8::Particle>::const_iterator end()
const 456 {
return entry.end(); }
459 int size()
const {
return entry.size();}
464 if (entryIn.col() > maxColTag) maxColTag = entryIn.col();
465 if (entryIn.acol() > maxColTag) maxColTag = entryIn.acol();
466 return entry.size() - 1;
468 int append(
int id,
int status,
int mother1,
int mother2,
int daughter1,
469 int daughter2,
int col,
int acol,
double px,
double py,
double pz,
470 double e,
double m = 0.,
double scaleIn = 0.,
double polIn = 9.) {
471 entry.push_back(
Particle(
id, status, mother1, mother2, daughter1,
472 daughter2, col, acol, px, py, pz, e,
m, scaleIn, polIn) );
setEvtPtr();
473 if (col > maxColTag) maxColTag = col;
474 if (acol > maxColTag) maxColTag = acol;
475 return entry.size() - 1;
477 int append(
int id,
int status,
int mother1,
int mother2,
int daughter1,
478 int daughter2,
int col,
int acol,
Vec4 p,
double m = 0.,
479 double scaleIn = 0.,
double polIn = 9.) {
480 entry.push_back(
Particle(
id, status, mother1, mother2, daughter1,
481 daughter2, col, acol, p,
m, scaleIn, polIn) );
setEvtPtr();
482 if (col > maxColTag) maxColTag = col;
483 if (acol > maxColTag) maxColTag = acol;
484 return entry.size() - 1;
488 int append(
int id,
int status,
int col,
int acol,
double px,
double py,
489 double pz,
double e,
double m = 0.,
double scaleIn = 0.,
490 double polIn = 9.) { entry.push_back(
Particle(
id, status, 0, 0, 0, 0,
491 col, acol, px, py, pz, e,
m, scaleIn, polIn) );
setEvtPtr();
492 if (col > maxColTag) maxColTag = col;
493 if (acol > maxColTag) maxColTag = acol;
494 return entry.size() - 1;
496 int append(
int id,
int status,
int col,
int acol,
Vec4 p,
double m = 0.,
497 double scaleIn = 0.,
double polIn = 9.) {entry.push_back(
Particle(
id,
498 status, 0, 0, 0, 0, col, acol, p,
m, scaleIn, polIn) );
setEvtPtr();
499 if (col > maxColTag) maxColTag = col;
500 if (acol > maxColTag) maxColTag = acol;
501 return entry.size() - 1;
505 void setEvtPtr(
int iSet = -1) {
if (iSet < 0) iSet = entry.size() - 1;
506 entry[iSet].setEvtPtr(
this);}
509 int copy(
int iCopy,
int newStatus = 0);
512 void list(
bool showScaleAndVertex =
false,
513 bool showMothersAndDaughters =
false,
int precision = 3)
const;
516 void popBack(
int nRemove = 1) {
if (nRemove ==1) entry.pop_back();
517 else {
int newSize = max( 0, size() - nRemove);
518 entry.resize(newSize);} }
522 void remove(
int iFirst,
int iLast,
bool shiftHistory =
true);
530 void restoreSize() {entry.resize(savedSize);}
531 int savedSizeValue() {
return savedSize;}
534 void initColTag(
int colTag = 0) {maxColTag = max( colTag,startColTag);}
535 int lastColTag()
const {
return maxColTag;}
536 int nextColTag() {
return ++maxColTag;}
539 void scale(
double scaleIn) {scaleSave = scaleIn;}
540 double scale()
const {
return scaleSave;}
543 void scaleSecond(
double scaleSecondIn) {scaleSecondSave = scaleSecondIn;}
544 double scaleSecond()
const {
return scaleSecondSave;}
548 vector<int>
daughterList(
int i)
const {
return entry[i].daughterList();}
551 int nFinal(
bool chargedOnly =
false)
const {
553 for (
int i = 0; i < size(); ++i)
554 if (entry[i].isFinal() && (!chargedOnly || entry[i].isCharged()))
559 double dyAbs(
int i1,
int i2)
const {
560 return abs( entry[i1].
y() - entry[i2].
y() ); }
561 double detaAbs(
int i1,
int i2)
const {
562 return abs( entry[i1].eta() - entry[i2].eta() ); }
563 double dphiAbs(
int i1,
int i2)
const {
564 double dPhiTmp = abs( entry[i1].
phi() - entry[i2].
phi() );
565 if (dPhiTmp > M_PI) dPhiTmp = 2. * M_PI - dPhiTmp;
567 double RRapPhi(
int i1,
int i2)
const {
568 return sqrt(
pow2(dyAbs(i1, i2)) +
pow2(dphiAbs(i1, i2)) ); }
569 double REtaPhi(
int i1,
int i2)
const {
570 return sqrt(
pow2(detaAbs(i1, i2)) +
pow2(dphiAbs(i1, i2)) ); }
573 void rot(
double theta,
double phi)
574 {
for (
int i = 0; i < size(); ++i) entry[i].rot(theta, phi);}
575 void bst(
double betaX,
double betaY,
double betaZ)
576 {
for (
int i = 0; i < size(); ++i) entry[i].bst(betaX, betaY, betaZ);}
577 void bst(
double betaX,
double betaY,
double betaZ,
double gamma)
578 {
for (
int i = 0; i < size(); ++i) entry[i].bst(betaX, betaY, betaZ,
580 void bst(
const Vec4& vec)
581 {
for (
int i = 0; i < size(); ++i) entry[i].bst(vec);}
582 void rotbst(
const RotBstMatrix& M,
bool boostVertices =
true)
583 {
for (
int i = 0; i < size(); ++i) entry[i].rotbst(M, boostVertices);}
590 { junction.push_back(
Junction( kind, col0, col1, col2) );
591 return junction.size() - 1;}
592 int appendJunction(
Junction junctionIn) {junction.push_back(junctionIn);
593 return junction.size() - 1;}
594 int sizeJunction()
const {
return junction.size();}
595 bool remainsJunction(
int i)
const {
return junction[i].remains();}
596 void remainsJunction(
int i,
bool remainsIn) {junction[i].remains(remainsIn);}
597 int kindJunction(
int i)
const {
return junction[i].kind();}
598 int colJunction(
int i,
int j)
const {
return junction[i].col(j);}
599 void colJunction(
int i,
int j,
int colIn) {junction[i].col(j, colIn);}
600 int endColJunction(
int i,
int j)
const {
return junction[i].endCol(j);}
601 void endColJunction(
int i,
int j,
int endColIn)
602 {junction[i].endCol(j, endColIn);}
603 int statusJunction(
int i,
int j)
const {
return junction[i].status(j);}
604 void statusJunction(
int i,
int j,
int statusIn)
605 {junction[i].status(j, statusIn);}
606 Junction& getJunction(
int i) {
return junction[i];}
607 const Junction& getJunction(
int i)
const {
return junction[i];}
608 void eraseJunction(
int i);
612 void restoreJunctionSize() {junction.resize(savedJunctionSize);}
615 void listJunctions()
const;
619 for (
const HVcols& col: hvCols) {
if (at(col.iHV).isFinal())
return true;}
623 void listHVcols()
const;
624 int maxHVcols()
const;
628 void restoreHVcolsSize() {hvCols.resize(savedHVcolsSize);}
638 const vector<Particle>*
particles()
const {
return &entry;}
646 static const int IPERLINE;
653 vector<Pythia8::Particle> entry;
657 vector<Pythia8::Junction> junction;
661 vector<Pythia8::HVcols> hvCols;
664 bool findIndexHV(
int iIn) {
if (iIn > 0 && iIn == iEventHV)
return true;
665 for (
int i = 0; i < int(hvCols.size()); ++i)
if (hvCols[i].iHV == iIn)
666 {iEventHV = iIn; iIndexHV = i;
return true; }
668 int iEventHV, iIndexHV;
671 void clearHV() {hvCols.resize(0); iEventHV = -1; iIndexHV = -1;}
677 int savedSize, savedJunctionSize, savedHVcolsSize, savedPartonLevelSize;
680 double scaleSave, scaleSecondSave;
virtual ~Particle()
Destructor.
Definition: Event.h:81
ParticleDataEntryPtr pdePtr
Definition: Event.h:317
constexpr double pow2(const double &x)
Powers of small integers - for balance speed/code clarity.
Definition: PythiaStdlib.h:182
const vector< Particle > * particles() const
Direct access to the particles via constant pointer.
Definition: Event.h:638
int append(int id, int status, int col, int acol, double px, double py, double pz, double e, double m=0., double scaleIn=0., double polIn=9.)
Brief versions of append: no mothers and no daughters.
Definition: Event.h:488
This class holds info on a single particle species.
Definition: ParticleData.h:125
void clearJunctions()
Clear the list of junctions.
Definition: Event.h:586
bool isAncestor(int iAncestor) const
Definition: Event.cc:332
void bstback(const Vec4 &pIn)
Boost given by a Vec4 p; boost in opposite direction.
Definition: Basics.cc:509
void rescale3(double fac)
Member functions that perform operations.
Definition: Basics.h:95
The HVcols class stores Hidden Valley colours for HV-coloured particles.
Definition: Event.h:390
The Event class holds all info on the generated event.
Definition: Event.h:408
bool isFinalPartonLevel() const
Check if particle belonged to the final state at the Parton Level.
Definition: Event.cc:411
double m2() const
Member functions for output; derived double quantities.
Definition: Event.h:161
double RRapPhi(const Vec4 &v1, const Vec4 &v2)
R is distance in cylindrical (y/eta, phi) coordinates.
Definition: Basics.cc:764
int iBotCopyId(bool simplify=false) const
Definition: Event.cc:153
void setEvtPtr(int iSet=-1)
Set pointer to the event for a particle, by default latest one.
Definition: Event.h:505
HVcols()
Default constructor and constructor with standard input.
Definition: Event.h:395
vector< int > sisterList(bool traceTopBot=false) const
Definition: Event.cc:301
void id(int idIn)
Member functions for input.
Definition: Event.h:88
void clear()
Clear event record.
Definition: Event.h:428
void rot(double theta, double phi)
Member functions for rotations and boosts of an event.
Definition: Event.h:573
Junction()
Constructors.
Definition: Event.h:343
void rotbst(const RotBstMatrix &M)
Arbitrary combination of rotations and boosts defined by 4 * 4 matrix.
Definition: Basics.cc:551
int colHV() const
Get and set Hidden Valley colours, referring back to the event.
Definition: Event.cc:522
void scale(double scaleIn)
Access scale for which event as a whole is defined.
Definition: Event.h:539
void remains(bool remainsIn)
Set values.
Definition: Event.h:364
static const double TINY
Constants: could only be changed in the code itself.
Definition: Event.h:302
void saveSize()
Save or restore the size of the event record (throwing at the end).
Definition: Event.h:529
vector< int > daughterListRecursive() const
Definition: Event.cc:270
void reset()
Clear event record, and set first particle empty.
Definition: Event.h:436
void rot(double thetaIn, double phiIn)
Rotation (simple).
Definition: Basics.cc:368
int append(Particle entryIn)
Put a new particle at the end of the event record; return index.
Definition: Event.h:462
int statusHepMC() const
Convert internal Pythia status codes to the HepMC status conventions.
Definition: Event.cc:383
void saveHVcolsSize()
Save or restore the size of the HV-coloured list (throwing at the end).
Definition: Event.h:627
vector< int > daughterList() const
Find complete list of daughters.
Definition: Event.cc:225
void initColTag(int colTag=0)
Initialize and access colour tag information.
Definition: Event.h:534
double REtaPhi(const Vec4 &v1, const Vec4 &v2)
Distance in cylindrical (eta, phi) coordinates.
Definition: Basics.cc:775
Particle & operator[](int i)
Overload index operator to access element of event record.
Definition: Event.h:439
bool undoDecay()
Definition: Event.cc:426
int iHV
Index of HV-particle and its HV-colour and HV-anticolour.
Definition: Event.h:400
double y() const
Functions for rapidity and pseudorapidity.
Definition: Event.cc:54
int id() const
Member functions for output.
Definition: Event.h:127
Particle()
Constructors.
Definition: Event.h:37
double dyAbs(int i1, int i2) const
Find separation in y, eta, phi or R between two particles.
Definition: Event.h:559
int idSave
Properties of the current particle.
Definition: Event.h:305
int iTopCopy() const
Trace the first and last copy of one and the same particle.
Definition: Event.cc:95
int appendJunction(int kind, int col0, int col1, int col2)
Add a junction to the list, study it or extra input.
Definition: Event.h:589
void bst(double betaX, double betaY, double betaZ)
Boost (simple).
Definition: Basics.cc:434
void popBack(int nRemove=1)
Remove last n entries.
Definition: Event.h:516
int size() const
Event record size.
Definition: Event.h:459
bool hasHVcols() const
Tell whether event has Hidden Valley colours stored.
Definition: Event.h:618
vector< int > motherList() const
Find complete list of mothers.
Definition: Event.cc:189
Event(int capacity=100)
Constructors.
Definition: Event.h:413
int idAbs() const
Member functions for output; derived int and bool quantities.
Definition: Event.h:152
Particle & front()
Implement standard references to elements in the particle array.
Definition: Event.h:443
vector< Pythia8::Particle >::iterator begin()
Implement iterators for the particle array.
Definition: Event.h:451
void setPDEPtr(ParticleDataEntryPtr pdePtrIn=nullptr)
Set pointer to the particle data species of the particle.
Definition: Event.cc:30
double phi(const Vec4 &v1, const Vec4 &v2)
phi is azimuthal angle between v1 and v2 around z axis.
Definition: Basics.cc:704
void scaleSecond(double scaleSecondIn)
Need a second scale if two hard interactions in event.
Definition: Event.h:543
void rescale3(double fac)
Member functions that perform operations.
Definition: Event.h:274
int nFinal(bool chargedOnly=false) const
Return number of final-state particles, optionally charged only.
Definition: Event.h:551
string name() const
Further output, based on a pointer to a ParticleDataEntry object.
Definition: Event.h:218
void setEvtPtr(Event *evtPtrIn)
Member functions to set the Event and ParticleDataEntry pointers.
Definition: Event.h:84
Header for classes to set beam momentum and interaction vertex spread.
Definition: Analysis.h:20
double m(const Vec4 &v1)
Invariant mass and its square.
Definition: Basics.cc:588
vector< int > daughterList(int i) const
Definition: Event.h:548
void saveJunctionSize()
Save or restore the size of the junction list (throwing at the end).
Definition: Event.h:611
void init(string headerIn="", ParticleData *particleDataPtrIn=0, int startColTagIn=100)
Initialize header for event listing, particle data table, and colour.
Definition: Event.h:422
void savePartonLevelSize()
Save event record size at Parton Level, i.e. before hadronization.
Definition: Event.h:631
void restorePtrs()
Definition: Event.h:526
int intPol() const
Find out if polarization is (close to) an integer.
Definition: Event.cc:40
void offsetHistory(int minMother, int addMother, int minDaughter, int addDaughter)
Add offsets to mother and daughter pointers (must be non-negative).
Definition: Event.cc:575
bool remains() const
Read out value.
Definition: Event.h:372
void offsetCol(int addCol)
Add offsets to colour and anticolour (must be positive).
Definition: Event.cc:590
int iTopCopyId(bool simplify=false) const
Definition: Event.cc:121
This class holds a map of all ParticleDataEntries.
Definition: ParticleData.h:422
string nameWithStatus(int maxLen=20) const
Particle name, with status but imposed maximum length -> may truncate.
Definition: Event.cc:558
Event * evtPtr
!
Definition: Event.h:321
virtual int index() const
Methods that can refer back to the event the particle belongs to.
Definition: Event.cc:87