PYTHIA  8.312
Public Member Functions | Friends | List of all members
Event Class Reference

The Event class holds all info on the generated event. More...

#include <Event.h>

Public Member Functions

 Event (int capacity=100)
 Constructors.
 
Eventoperator= (const Event &oldEvent)
 Copy all information from one event record to another. More...
 
 Event (const Event &oldEvent)
 
void init (string headerIn="", ParticleData *particleDataPtrIn=0, int startColTagIn=100)
 Initialize header for event listing, particle data table, and colour.
 
void clear ()
 Clear event record.
 
void free ()
 
void reset ()
 Clear event record, and set first particle empty.
 
Particleoperator[] (int i)
 Overload index operator to access element of event record.
 
const Particleoperator[] (int i) const
 
Particlefront ()
 Implement standard references to elements in the particle array.
 
Particleat (int i)
 
Particleback ()
 
const Particlefront () const
 
const Particleat (int i) const
 
const Particleback () const
 
vector< Pythia8::Particle >::iterator begin ()
 Implement iterators for the particle array.
 
vector< Pythia8::Particle >::iterator end ()
 
vector< Pythia8::Particle >::const_iterator begin () const
 
vector< Pythia8::Particle >::const_iterator end () const
 
int size () const
 Event record size.
 
int append (Particle entryIn)
 Put a new particle at the end of the event record; return index.
 
int append (int id, int status, int mother1, int mother2, int daughter1, int daughter2, int col, int acol, double px, double py, double pz, double e, double m=0., double scaleIn=0., double polIn=9.)
 
int append (int id, int status, int mother1, int mother2, int daughter1, int daughter2, int col, int acol, Vec4 p, double m=0., double scaleIn=0., double polIn=9.)
 
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.
 
int append (int id, int status, int col, int acol, Vec4 p, double m=0., double scaleIn=0., double polIn=9.)
 
void setEvtPtr (int iSet=-1)
 Set pointer to the event for a particle, by default latest one.
 
int copy (int iCopy, int newStatus=0)
 Add a copy of an existing particle at the end of the event record. More...
 
void list (bool showScaleAndVertex=false, bool showMothersAndDaughters=false, int precision=3) const
 List the particles in an event. More...
 
void popBack (int nRemove=1)
 Remove last n entries.
 
void remove (int iFirst, int iLast, bool shiftHistory=true)
 
void restorePtrs ()
 
void saveSize ()
 Save or restore the size of the event record (throwing at the end).
 
void restoreSize ()
 
int savedSizeValue ()
 
void initColTag (int colTag=0)
 Initialize and access colour tag information.
 
int lastColTag () const
 
int nextColTag ()
 
void scale (double scaleIn)
 Access scale for which event as a whole is defined.
 
double scale () const
 
void scaleSecond (double scaleSecondIn)
 Need a second scale if two hard interactions in event.
 
double scaleSecond () const
 
vector< int > daughterList (int i) const
 
int nFinal (bool chargedOnly=false) const
 Return number of final-state particles, optionally charged only.
 
double dyAbs (int i1, int i2) const
 Find separation in y, eta, phi or R between two particles.
 
double detaAbs (int i1, int i2) const
 
double dphiAbs (int i1, int i2) const
 
double RRapPhi (int i1, int i2) const
 
double REtaPhi (int i1, int i2) const
 
void rot (double theta, double phi)
 Member functions for rotations and boosts of an event.
 
void bst (double betaX, double betaY, double betaZ)
 
void bst (double betaX, double betaY, double betaZ, double gamma)
 
void bst (const Vec4 &vec)
 
void rotbst (const RotBstMatrix &M, bool boostVertices=true)
 
void clearJunctions ()
 Clear the list of junctions.
 
int appendJunction (int kind, int col0, int col1, int col2)
 Add a junction to the list, study it or extra input.
 
int appendJunction (Junction junctionIn)
 
int sizeJunction () const
 
bool remainsJunction (int i) const
 
void remainsJunction (int i, bool remainsIn)
 
int kindJunction (int i) const
 
int colJunction (int i, int j) const
 
void colJunction (int i, int j, int colIn)
 
int endColJunction (int i, int j) const
 
void endColJunction (int i, int j, int endColIn)
 
int statusJunction (int i, int j) const
 
void statusJunction (int i, int j, int statusIn)
 
JunctiongetJunction (int i)
 
const JunctiongetJunction (int i) const
 
void eraseJunction (int i)
 Erase junction stored in specified slot and move up the ones under.
 
void saveJunctionSize ()
 Save or restore the size of the junction list (throwing at the end).
 
void restoreJunctionSize ()
 
void listJunctions () const
 List any junctions in the event; for debug mainly. More...
 
bool hasHVcols () const
 Tell whether event has Hidden Valley colours stored.
 
void listHVcols () const
 List any Hidden Valley colours. Find largest HV colour. More...
 
int maxHVcols () const
 Find largest Hidden Valley (anti)colour in an event.
 
void saveHVcolsSize ()
 Save or restore the size of the HV-coloured list (throwing at the end).
 
void restoreHVcolsSize ()
 
void savePartonLevelSize ()
 Save event record size at Parton Level, i.e. before hadronization.
 
void addStringBreaks (StringBreaks &sbIn)
 Add information about string breaks for a single string.
 
void clearStringBreaks ()
 Clear string breaks information.
 
const vector< StringBreaks > & getStringBreaks () const
 Get information about all string breaks in the event.
 
Eventoperator+= (const Event &addEvent)
 Operator overloading allows to append one event to an existing one. More...
 
const vector< Particle > * particles () const
 Direct access to the particles via constant pointer.
 

Friends

class Particle
 The Particle class needs to access particle data.
 

Detailed Description

The Event class holds all info on the generated event.

Member Function Documentation

int copy ( int  iCopy,
int  newStatus = 0 
)

Add a copy of an existing particle at the end of the event record.

Add a copy of an existing particle at the end of the event record; return index. Three cases, depending on sign of new status code: Positive: copy is viewed as daughter, status of original is negated. Negative: copy is viewed as mother, status of original is unchanged. Zero: the new is a perfect carbon copy (maybe to be changed later).

Protect against attempt to copy negative entries (e.g., junction codes) or entries beyond end of record.

Simple carbon copy.

Set up to make new daughter of old.

Set up to make new mother of old.

Done.

vector<int> daughterList ( int  i) const
inline

Find complete list of daughters. Note: temporarily retained for CMS compatibility. Do not use!

void list ( bool  showScaleAndVertex = false,
bool  showMothersAndDaughters = false,
int  precision = 3 
) const

List the particles in an event.

Print an event.

Header.

Precision. At high energy switch to scientific format for momenta.

Listing of complete event.

Basic line for a particle, always printed.

Optional extra line for scale value, polarization and production vertex.

Optional extra line, giving a complete list of mothers and daughters.

Extra blank separation line when each particle spans more than one line.

Statistics on momentum and charge.

Line with sum charge, momentum, energy and invariant mass.

Listing finished.

void listHVcols ( ) const

List any Hidden Valley colours. Find largest HV colour.

Print the Hidden Valley colours in an event.

void listJunctions ( ) const

List any junctions in the event; for debug mainly.

Print the junctions in an event.

Header.

Loop through junctions in event and list them.

Alternative if no junctions. Listing finished.

Event & operator+= ( const Event addEvent)

Operator overloading allows to append one event to an existing one.

Operator overloading allows to append one event to an existing one. Warning: particles should be OK, but some other information unreliable.

Find offsets. One less since won't copy line 0.

Add energy to zeroth line and calculate new invariant mass.

Read out particles from line 1 (not 0) onwards.

Add offset to nonzero mother, daughter and colour indices.

Append particle to summed event.

Read out junctions one by one.

Add colour offsets to all three legs.

Append junction to summed event.

Append Hidden Valley colour information.

Set header that indicates character as sum of events.

Done.

Event & operator= ( const Event oldEvent)

Copy all information from one event record to another.

Do not copy if same.

Reset all current info in the event.

Copy particle data table; needed for individual particles.

Copy all the particles one by one.

Copy all the junctions one by one.

Copy the Hidden Valley colour information.

Copy all other values.

Done.

void remove ( int  iFirst,
int  iLast,
bool  shiftHistory = true 
)

Remove entries from iFirst to iLast, including endpoints, and fix history. (To the extent possible; history pointers in removed range are zeroed.)

Check that removal range is sensible.

Remove the entries.

Loop over remaining particles; read out mothers and daughters.

Shift mother and daughter indices according to removed number. Set zero if in removed range.

Set the new values.

void restorePtrs ( )
inline

Restore all ParticleDataEntryPtr pointers in the Particle vector. Useful when a persistent copy of the event record is read back in.


The documentation for this class was generated from the following files: