PYTHIA  8.312
Public Member Functions | Protected Attributes | Static Protected Attributes | List of all members
Particle Class Reference

#include <Event.h>

Inheritance diagram for Particle:
Py8Particle ColourParticle HelicityParticle

Public Member Functions

 Particle ()
 Constructors.
 
 Particle (int idIn, int statusIn=0, int mother1In=0, int mother2In=0, int daughter1In=0, int daughter2In=0, int colIn=0, int acolIn=0, double pxIn=0., double pyIn=0., double pzIn=0., double eIn=0., double mIn=0., double scaleIn=0., double polIn=9.)
 
 Particle (int idIn, int statusIn, int mother1In, int mother2In, int daughter1In, int daughter2In, int colIn, int acolIn, Vec4 pIn, double mIn=0., double scaleIn=0., double polIn=9.)
 
 Particle (const Particle &pt)
 
Particleoperator= (const Particle &pt)
 
virtual ~Particle ()
 Destructor.
 
void setEvtPtr (Event *evtPtrIn)
 Member functions to set the Event and ParticleDataEntry pointers.
 
void setPDEPtr (ParticleDataEntryPtr pdePtrIn=nullptr)
 Set pointer to the particle data species of the particle.
 
void id (int idIn)
 Member functions for input.
 
void status (int statusIn)
 
void statusPos ()
 
void statusNeg ()
 
void statusCode (int statusIn)
 
void mother1 (int mother1In)
 
void mother2 (int mother2In)
 
void mothers (int mother1In=0, int mother2In=0)
 
void daughter1 (int daughter1In)
 
void daughter2 (int daughter2In)
 
void daughters (int daughter1In=0, int daughter2In=0)
 
void col (int colIn)
 
void acol (int acolIn)
 
void cols (int colIn=0, int acolIn=0)
 
void p (Vec4 pIn)
 
void p (double pxIn, double pyIn, double pzIn, double eIn)
 
void px (double pxIn)
 
void py (double pyIn)
 
void pz (double pzIn)
 
void e (double eIn)
 
void m (double mIn)
 
void scale (double scaleIn)
 
void pol (double polIn)
 
void vProd (Vec4 vProdIn)
 
void vProd (double xProdIn, double yProdIn, double zProdIn, double tProdIn)
 
void xProd (double xProdIn)
 
void yProd (double yProdIn)
 
void zProd (double zProdIn)
 
void tProd (double tProdIn)
 
void vProdAdd (Vec4 vProdIn)
 
void tau (double tauIn)
 
int id () const
 Member functions for output.
 
int status () const
 
int mother1 () const
 
int mother2 () const
 
int daughter1 () const
 
int daughter2 () const
 
int col () const
 
int acol () const
 
Vec4 p () const
 
double px () const
 
double py () const
 
double pz () const
 
double e () const
 
double m () const
 
double scale () const
 
double pol () const
 
bool hasVertex () const
 
Vec4 vProd () const
 
double xProd () const
 
double yProd () const
 
double zProd () const
 
double tProd () const
 
double tau () const
 
int idAbs () const
 Member functions for output; derived int and bool quantities.
 
int statusAbs () const
 
bool isFinal () const
 
int intPol () const
 Find out if polarization is (close to) an integer.
 
bool isRescatteredIncoming () const
 
double m2 () const
 Member functions for output; derived double quantities.
 
double mCalc () const
 
double m2Calc () const
 
double eCalc () const
 
double pT () const
 
double pT2 () const
 
double mT () const
 
double mT2 () const
 
double pAbs () const
 
double pAbs2 () const
 
double eT () const
 
double eT2 () const
 
double theta () const
 
double phi () const
 
double thetaXZ () const
 
double pPos () const
 
double pNeg () const
 
double y () const
 Functions for rapidity and pseudorapidity.
 
double eta () const
 
double y (double mCut) const
 Rapidity with minimal transverse mass, and after rotation/boost.
 
double y (double mCut, RotBstMatrix &M) const
 
Vec4 vDec () const
 
double xDec () const
 
double yDec () const
 
double zDec () const
 
double tDec () const
 
virtual int index () const
 Methods that can refer back to the event the particle belongs to. More...
 
int iTopCopy () const
 Trace the first and last copy of one and the same particle.
 
int iBotCopy () const
 
int iTopCopyId (bool simplify=false) const
 
int iBotCopyId (bool simplify=false) const
 
vector< int > motherList () const
 Find complete list of mothers. More...
 
vector< int > daughterList () const
 Find complete list of daughters. More...
 
vector< int > daughterListRecursive () const
 
vector< int > sisterList (bool traceTopBot=false) const
 
bool isAncestor (int iAncestor) const
 
int statusHepMC () const
 Convert internal Pythia status codes to the HepMC status conventions. More...
 
bool isFinalPartonLevel () const
 Check if particle belonged to the final state at the Parton Level.
 
bool undoDecay ()
 
int colHV () const
 Get and set Hidden Valley colours, referring back to the event.
 
int acolHV () const
 
void colHV (int colHVin)
 
void acolHV (int acolHVin)
 
void colsHV (int colHVin, int acolHVin)
 
string name () const
 Further output, based on a pointer to a ParticleDataEntry object.
 
string nameWithStatus (int maxLen=20) const
 Particle name, with status but imposed maximum length -> may truncate. More...
 
int spinType () const
 
int chargeType () const
 
double charge () const
 
bool isCharged () const
 
bool isNeutral () const
 
int colType () const
 
double m0 () const
 
double mWidth () const
 
double mMin () const
 
double mMax () const
 
double mSel () const
 
double constituentMass () const
 
double tau0 () const
 
bool mayDecay () const
 
bool canDecay () const
 
bool doExternalDecay () const
 
bool isResonance () const
 
bool isVisible () const
 
bool isLepton () const
 
bool isQuark () const
 
bool isGluon () const
 
bool isDiquark () const
 
bool isParton () const
 
bool isHadron () const
 
bool isExotic () const
 
ParticleDataEntryparticleDataEntry () const
 
void rescale3 (double fac)
 Member functions that perform operations.
 
void rescale4 (double fac)
 
void rescale5 (double fac)
 
void rot (double thetaIn, double phiIn)
 
void bst (double betaX, double betaY, double betaZ)
 
void bst (double betaX, double betaY, double betaZ, double gamma)
 
void bst (const Vec4 &pBst)
 
void bst (const Vec4 &pBst, double mBst)
 
void bstback (const Vec4 &pBst)
 
void bstback (const Vec4 &pBst, double mBst)
 
void rotbst (const RotBstMatrix &M, bool boostVertex=true)
 
void offsetHistory (int minMother, int addMother, int minDaughter, int addDaughter)
 Add offsets to mother and daughter pointers (must be non-negative).
 
void offsetCol (int addCol)
 Add offsets to colour and anticolour (must be positive).
 

Protected Attributes

int idSave
 Properties of the current particle.
 
int statusSave
 
int mother1Save
 
int mother2Save
 
int daughter1Save
 
int daughter2Save
 
int colSave
 
int acolSave
 
Vec4 pSave
 
double mSave
 
double scaleSave
 
double polSave
 
bool hasVertexSave
 
Vec4 vProdSave
 
double tauSave
 
ParticleDataEntryPtr pdePtr
 
EventevtPtr
 ! More...
 

Static Protected Attributes

static const double TINY = 1e-20
 Constants: could only be changed in the code itself. More...
 

Detailed Description

Particle class. This class holds info on a particle in general.

Member Function Documentation

vector< int > daughterList ( ) const

Find complete list of daughters.

Vector of all the daughters; created empty. Done if no event pointer.

Simple cases: no or one daughter.

A range of daughters.

Two separated daughters.

Special case for two incoming beams: attach further initiators and remnants that have beam as mother.

Done.

vector< int > daughterListRecursive ( ) const

Find complete list of daughters recursively, i.e. including subsequent generations. Is intended specifically for resonance decays.

Vector of all the daughters; created empty. Done if no event pointer.

Find first generation of daughters.

Recursively add daughters of unstable particles.

Done.

int iBotCopyId ( bool  simplify = false) const

Check that particle belongs to event record. Initialize.

Simple solution when only first and last daughter are studied.

Else full solution where all daughters are studied.

int index ( ) const
virtual

Methods that can refer back to the event the particle belongs to.

Method to find the index of the particle in the event record.

Reimplemented in HelicityParticle, and Py8Particle.

bool isAncestor ( int  iAncestor) const

Check whether a given particle is an arbitrarily-steps-removed mother to another. For the parton -> hadron transition, only first-rank hadrons are associated with the respective end quark.

Begin loop to trace upwards from the daughter.

If positive match then done.

If out of range then failed to find match.

If unique mother then keep on moving up the chain.

If many mothers, except hadronization, then fail tracing.

For hadronization step, fail if not first rank, else move up.

Fail for ministring -> one hadron and for junctions.

End of loop. Should never reach beyond here.

int iTopCopyId ( bool  simplify = false) const

Trace the first and last copy of one and the same particle, also through shower branchings, making use of flavour matches. Stops tracing when this gives ambiguities.

Check that particle belongs to event record. Initialize.

Simple solution when only first and last mother are studied.

Else full solution where all mothers are studied.

vector< int > motherList ( ) const

Find complete list of mothers.

Vector of all the mothers; created empty. Done if no event pointer.

Special cases in the beginning, where the meaning of zero is unclear.

One mother or a carbon copy.

A range of mothers from string fragmentation.

Two separate mothers.

Done.

string nameWithStatus ( int  maxLen = 20) const

Particle name, with status but imposed maximum length -> may truncate.

Remove from end, excluding closing bracket and charge.

vector< int > sisterList ( bool  traceTopBot = false) const

Find complete list of sisters. Optionally traces up with iTopCopy and down with iBotCopy to give sisters at same level of evolution.

Vector of all the sisters; created empty.

Find all daughters of the mother.

Copy all daughters, excepting the input particle itself.

Done.

int statusHepMC ( ) const

Convert internal Pythia status codes to the HepMC status conventions.

Positive codes are final particles. Status -12 are beam particles.

Hadrons, muons, taus that decay normally are status 2.

Particle should not decay into itself (e.g. Bose-Einstein).

Other acceptable negative codes as their positive counterpart.

Unacceptable codes as 0.

bool undoDecay ( )

Recursively remove the decay products of a particle, update it to be undecayed, and update all mother/daughter indices to be correct. Warning: assumes that decay chains are nicely ordered.

Do not remove daughters of a parton, i.e. entry carrying colour.

Find range of daughters to remove.

Refuse if any of the daughters have other mothers.

Initialize range arrays for daughters and granddaughters.

Begin recursive search through all decay chains.

Find new daughter range, if present.

Check if the range duplicates or contradicts existing ones.

Add new range where relevant. Keep ranges ordered.

End of recursive search all decay chains.

Join adjacent ranges to reduce number of erase steps.

Iterate over relevant ranges, from bottom up.

Remove daughters in each range and update mother and daughter indices.

Update mother that has been undecayed.

Done.

Member Data Documentation

Event* evtPtr
protected

!

Pointer to the whole event record to which the particle belongs (if any). As above it should not be saved.

ParticleDataEntryPtr pdePtr
protected

Pointer to properties of the particle species. Should no be saved in a persistent copy of the event record. The ///! below is ROOT notation that this member should not be saved. Event::restorePtrs() can be called to restore the missing information.

const double TINY = 1e-20
staticprotected

Constants: could only be changed in the code itself.

Small number to avoid division by zero.

Particle class. This class holds info on a particle in general. Constants: could be changed here if desired, but normally should not. These are of technical nature, as described for each.


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