A++ » INC » HepParticle

class HepParticle: public TObject


 Abstract HEP particle class

 Note for developers: Some often used values like Pt, Eta etc. are
 stored as transient data members for fast access. Their computation
 from the persistent data members is either done in the normal ctor
 or in SetPtEtaPhiE(), or triggered by the fCompute flag in case of
 using the default ctor.

 Also included here is a function allowing to print TLorentzVector
 Cartesian/(pt,eta,phi) coordinates.

 !!! IMPORTANT: Do not acces the transient data members directly,  !!!
 !!! but use ALWAYS their corresponding access function like Pt(), !!!
 !!! Eta()etc. to ensure that the variables are computed properly. !!!

 !!! IMPORTANT: Do not use Pt(), Pmag(), Thethe(), Phi(), Eta() in !!!
 !!! TTree::Draw() or TTree::Scan() since these functions make use !!!
 !!! of a cache which will NOT be updated when jumping to the next !!!
 !!! event inside TTree::Draw() or TTree::Scan().                  !!!


 Author: Oliver Maria Kind <mailto: kind@mail.desy.de>
 Update: $Id: HepParticle.cxx,v 1.28 2014/02/18 12:15:18 kind Exp $
 Copyright: 2008 (C) Oliver Maria Kind

Function Members (Methods)

public:
virtual~HepParticle()
voidTObject::AbstractMethod(const char* method) const
virtual voidTObject::AppendPad(Option_t* option = "")
virtual voidTObject::Browse(TBrowser* b)
Float_tCharge() const
static TClass*Class()
virtual const char*TObject::ClassName() const
virtual voidClear(Option_t* option = "")
virtual TObject*TObject::Clone(const char* newname = "") const
virtual Int_tCompare(const TObject* obj) const
virtual voidTObject::Copy(TObject& object) const
Float_tCosTheta()
virtual voidTObject::Delete(Option_t* option = "")MENU
Float_tDeltaEta(HepParticle* prt)
Float_tDeltaEta(HepJet* jet)
Float_tDeltaPhi(HepParticle* prt) const
Float_tDeltaPhi(HepJet* jet) const
Float_tDeltaPtFrac(HepParticle* prt)
Float_tDeltaPtFrac(HepJet* jet)
Float_tDeltaR(HepParticle* prt)
Float_tDeltaR(HepJet* jet)
virtual Int_tTObject::DistancetoPrimitive(Int_t px, Int_t py)
virtual voidTObject::Draw(Option_t* option = "")
virtual voidTObject::DrawClass() constMENU
virtual TObject*TObject::DrawClone(Option_t* option = "") constMENU
virtual voidTObject::Dump() constMENU
Float_tE() const
virtual voidTObject::Error(const char* method, const char* msgfmt) const
Float_tEt() const
Float_tEta()
virtual voidTObject::Execute(const char* method, const char* params, Int_t* error = 0)
virtual voidTObject::Execute(TMethod* method, TObjArray* params, Int_t* error = 0)
virtual voidTObject::ExecuteEvent(Int_t event, Int_t px, Int_t py)
virtual voidTObject::Fatal(const char* method, const char* msgfmt) const
virtual TObject*TObject::FindObject(const char* name) const
virtual TObject*TObject::FindObject(const TObject* obj) const
virtual voidGetCovMatrixPtEtaPhi(TMatrixD& CovMatPtEtaPhi) const
virtual Option_t*TObject::GetDrawOption() const
static Long_tTObject::GetDtorOnly()
virtual const char*TObject::GetIconName() const
Int_tGetId() const
virtual const char*TObject::GetName() const
virtual char*TObject::GetObjectInfo(Int_t px, Int_t py) const
static Bool_tTObject::GetObjectStat()
virtual Option_t*TObject::GetOption() const
Int_tGetPdgCode() const
const char*GetPdgName() const
virtual const char*TObject::GetTitle() const
virtual UInt_tTObject::GetUniqueID() const
virtual Bool_tTObject::HandleTimer(TTimer* timer)
virtual ULong_tTObject::Hash() const
Bool_tHasInvalidPdgCode() const
HepParticle()
HepParticle(const HepParticle&)
HepParticle(Int_t Id, Float_t Px, Float_t Py, Float_t Pz, Float_t E, Int_t PdgCode)
virtual voidTObject::Info(const char* method, const char* msgfmt) const
virtual Bool_tTObject::InheritsFrom(const char* classname) const
virtual Bool_tTObject::InheritsFrom(const TClass* cl) const
virtual voidTObject::Inspect() constMENU
Float_tInvDeltaPtFrac(HepParticle* prt)
voidTObject::InvertBit(UInt_t f)
virtual TClass*IsA() const
Bool_tIsBeauty() const
Bool_tIsBeautyBar() const
Bool_tIsBeautyPrime() const
Bool_tIsBeautyPrimeBar() const
Bool_tIsBeautyPrimeQuark() const
Bool_tIsBeautyQuark() const
Bool_tIsChargedLepton() const
Bool_tIsCharm() const
Bool_tIsCharmBar() const
Bool_tIsCharmQuark() const
Bool_tIsDown() const
Bool_tIsDownBar() const
Bool_tIsDownQuark() const
Bool_tIsElectron() const
Bool_tIsElectronNeutrino() const
Bool_tIsEMinus() const
Bool_tIsEPlus() const
virtual Bool_tTObject::IsEqual(const TObject* obj) const
virtual Bool_tTObject::IsFolder() const
Bool_tIsKaon() const
Bool_tIsKMinus() const
Bool_tIsKPlus() const
Bool_tIsLightQuark() const
Bool_tIsMuMinus() const
Bool_tIsMuon() const
Bool_tIsMuonNeutrino() const
Bool_tIsMuPlus() const
Bool_tIsNegative() const
Bool_tIsNeutral() const
Bool_tIsNeutrino() const
Bool_tTObject::IsOnHeap() const
Bool_tIsPhoton() const
Bool_tIsPiMinus() const
Bool_tIsPion() const
Bool_tIsPiPlus() const
Bool_tIsPositive() const
Bool_tIsProton() const
Bool_tIsQuark() const
virtual Bool_tIsSortable() const
Bool_tIsStrange() const
Bool_tIsStrangeBar() const
Bool_tIsStrangeQuark() const
Bool_tIsStringOrCluster() const
Bool_tIsTau() const
Bool_tIsTauMinus() const
Bool_tIsTauNeutrino() const
Bool_tIsTauPlus() const
Bool_tIsTop() const
Bool_tIsTopBar() const
Bool_tIsTopPrime() const
Bool_tIsTopPrimeBar() const
Bool_tIsTopPrimeQuark() const
Bool_tIsTopQuark() const
Bool_tIsUp() const
Bool_tIsUpBar() const
Bool_tIsUpQuark() const
Bool_tIsWBoson() const
Bool_tIsWMinus() const
Bool_tIsWPlus() const
Bool_tIsZ0Boson() const
Bool_tTObject::IsZombie() const
virtual voidTObject::ls(Option_t* option = "") const
Double_tM(Option_t* option = "PDG") const
Double_tMass(Option_t* option = "PDG") const
voidTObject::MayNotUse(const char* method) const
Float_tMperp() const
Float_tMperp2() const
Float_tMt() const
Float_tMt2() const
virtual Bool_tTObject::Notify()
voidTObject::Obsolete(const char* method, const char* asOfVers, const char* removedFromVers) const
voidTObject::operator delete(void* ptr)
voidTObject::operator delete(void* ptr, void* vp)
voidTObject::operator delete[](void* ptr)
voidTObject::operator delete[](void* ptr, void* vp)
void*TObject::operator new(size_t sz)
void*TObject::operator new(size_t sz, void* vp)
void*TObject::operator new[](size_t sz)
void*TObject::operator new[](size_t sz, void* vp)
HepParticle&operator=(const HepParticle&)
Bool_toperator==(HepParticle Partner)
const TLorentzVector&P() const
TVector3P3() const
virtual voidTObject::Paint(Option_t* option = "")
Float_tPhi()
Float_tPmag()
virtual voidTObject::Pop()
virtual voidPrint(Option_t* option = "")
static voidPrintFooter()
static voidPrintHeader()
static voidPrintLorentzVector(TLorentzVector vec)
Float_tPt()
Float_tPx() const
Float_tPy() const
Float_tPz() const
virtual Int_tTObject::Read(const char* name)
virtual voidTObject::RecursiveRemove(TObject* obj)
voidTObject::ResetBit(UInt_t f)
virtual voidTObject::SaveAs(const char* filename = "", Option_t* option = "") constMENU
virtual voidTObject::SavePrimitive(ostream& out, Option_t* option = "")
voidTObject::SetBit(UInt_t f)
voidTObject::SetBit(UInt_t f, Bool_t set)
virtual voidTObject::SetDrawOption(Option_t* option = "")MENU
static voidTObject::SetDtorOnly(void* obj)
voidSetId(Int_t Id)
static voidTObject::SetObjectStat(Bool_t stat)
voidSetP(const TLorentzVector& p)
voidSetPdgCode(Int_t pdg)
voidSetPtEtaPhiE(Float_t Pt, Float_t Eta, Float_t Phi, Float_t E)
voidSetPtEtaPhiM(Float_t Pt, Float_t Eta, Float_t Phi, Float_t M)
voidSetPxPyPzE(Float_t Px, Float_t Py, Float_t Pz, Float_t E)
virtual voidTObject::SetUniqueID(UInt_t uid)
virtual voidShowMembers(TMemberInspector& insp) const
virtual voidStreamer(TBuffer&)
voidStreamerNVirtual(TBuffer& ClassDef_StreamerNVirtual_b)
virtual voidTObject::SysError(const char* method, const char* msgfmt) const
Bool_tTObject::TestBit(UInt_t f) const
Int_tTObject::TestBits(UInt_t f) const
Float_tTheta()
virtual voidTObject::UseCurrentStyle()
virtual voidTObject::Warning(const char* method, const char* msgfmt) const
virtual Int_tTObject::Write(const char* name = 0, Int_t option = 0, Int_t bufsize = 0)
virtual Int_tTObject::Write(const char* name = 0, Int_t option = 0, Int_t bufsize = 0) const
protected:
virtual voidTObject::DoError(int level, const char* location, const char* fmt, va_list va) const
voidTObject::MakeZombie()
private:
voidComputeTransientVars()

Data Members

public:
static TObject::(anonymous)TObject::kBitMask
static TObject::EStatusBitsTObject::kCanDelete
static TObject::EStatusBitsTObject::kCannotPick
static TObject::EStatusBitsTObject::kHasUUID
static TObject::EStatusBitsTObject::kInvalidObject
static TObject::(anonymous)TObject::kIsOnHeap
static TObject::EStatusBitsTObject::kIsReferenced
static TObject::EStatusBitsTObject::kMustCleanup
static TObject::EStatusBitsTObject::kNoContextMenu
static TObject::(anonymous)TObject::kNotDeleted
static TObject::EStatusBitsTObject::kObjInCanvas
static TObject::(anonymous)TObject::kOverwrite
static TObject::(anonymous)TObject::kSingleKey
static TObject::(anonymous)TObject::kWriteDelete
static TObject::(anonymous)TObject::kZombie
protected:
static HepDatabasePDG*fDbasePDG! PDG database
Int_tfIdId number (for convenience only)
TLorentzVectorfP4-momentum (Px, Py, Pz, E) in GeV
Int_tfPdgCodePDG code
private:
Bool_tfCompute! Compute the transient variables
Float_tfCosTheta! Cosine of the polar angle
Float_tfEta! Pseudo-rapidity of the 3-momentum
Float_tfPhi! Azimuth of the 3-momentum (rad)
Float_tfPmag! Magnitude of the 3-mom. vector (GeV)
Float_tfPt! Transverse momentum (GeV)
Float_tfTheta! Polar angle of the 3-momentum (rad)

Class Charts

Inheritance Chart:
TObject
HepParticle
HepDecayParticle
Hep2BodyDecay
AtlD0DecayKPi
AtlDstarDecayDPi
AtlK0sDecayPiPi
AtlLambdaDecayPPi
AtlPhiDecayKK
AtlPhotonConv
HepK0sDecay
HepTopDecay
 [more...]

Function documentation

HepParticle()
 Default constructor

 The fCompute flag triggers the computation of the transient
 member variables the next time one of them is accessed. This
 mechanism is neccessary since the values of the persistent data
 members needed for the computation are filles AFTER the call to
 this ctor

HepParticle(Int_t Id, Float_t Px, Float_t Py, Float_t Pz, Float_t E, Int_t PdgCode)
 Normal constructor

~HepParticle()
 Default destructor

void Clear(Option_t* option = "")
 Clear object

Double_t Mass(Option_t* option = "PDG") const
 Particle mass

 Options available:
    "PDG" = mass as given in the PDG dbase (default)
    "REC" = 'reconstructed' invariant mass stored in
             the 4-momentum vector

Bool_t operator==(HepParticle Partner)
 Compare Particle with Partner w.r.t. PdgCode and 4-Momentum
 return kTRUE in case of Matching

 Comparison is made between Objects, not Pointers

Int_t Compare(const TObject* obj) const
 Compares transverse momentum of this particles with another.
 Necessary for sorting

Bool_t HasInvalidPdgCode() const
 Test validity of PDG code
 Compounds like eg. deuteron, tritium or alphas etc. are
 unknown to the PDG database

const char* GetPdgName() const
 Return particle's name

Float_t Charge() const
 Get charge in units of e
 (Note that in the PDG dbase the charges are stored in 1/3 of e)

Float_t DeltaPhi(HepParticle* prt) const
 Compute azimuth difference between this and the given particle
 in the range from -pi to +pi

Float_t DeltaPhi(HepJet* jet) const
 Compute azimuth difference between this particle and the given
 jet in the range from -pi to +pi

Float_t DeltaEta(HepParticle* prt)
 Compute (signed) pseudo-rapidity difference between this and
 the given particle

Float_t DeltaEta(HepJet* jet)
 Compute (signed) pseudo-rapidity difference between this
 particle and the given jet

Float_t DeltaR(HepParticle* prt)
 Compute distance in eta-phi plane between this and the given
 particle

Float_t DeltaR(HepJet* jet)
 Compute distance in eta-phi plane between this particle and the
 given jet

Float_t Mperp2() const
 Mt^2 = Et^2 - Pt^2

 Use this function for the Jacobian-peak method

Float_t DeltaPtFrac(HepParticle* prt)
 DeltaPtFrac = (Pt_1 - Pt_2)Pt_1

 Unsigned, for better checking systematic errors

Float_t DeltaPtFrac(HepJet* jet)
Float_t InvDeltaPtFrac(HepParticle* prt)
 An "Inverse" DeltaPtFrac is needed, due to Multiple Scattering in
 the Inner Detector
 InvDeltaPtFrac = ( 1/Pt_1 - 1/Pt_2)*Pt_1

 Unsigned again, for better checking systematic errors
 I dont know yet if I can write the same algorithm for the jets

void ComputeTransientVars()
 Compute transient variables like Pt, Theta, Eta etc. for fast
 access

void Print(Option_t* option = "")
 Print object information

 Options available:
   "nohead" - No header containing the variable names is
              displayed. Useful when printing a whole table
              for a list of particles

void PrintHeader()
 Print information header


 GINUS denotes (in this order):

 IsGenerator()
 IsGenInteracting()
 IsGenNonInteracting()
 IsGenSimulStable()
 IsGenStable()

void PrintFooter()
 Print footer

void GetCovMatrixPtEtaPhi(TMatrixD& CovMatPtEtaPhi) const
 Dummy for returning the covariance matrix in some derived
 classes (reconstructed electrons and muons, ...)

void PrintLorentzVector(TLorentzVector vec)
 Print the Cartesian and (pt, eta, phi) coordinates of a
 TLorentVector. Needed, since this functionality is not
 provided by RooT. Guards Eta() from Pt = 0.

HepParticle()
Int_t GetId() const
{ return fId; }
Int_t GetPdgCode() const
{ return fPdgCode; }
Bool_t IsPositive() const
{ return ( Charge() > 0. ) ? kTRUE : kFALSE; }
Bool_t IsNegative() const
{ return ( Charge() < 0. ) ? kTRUE : kFALSE; }
Bool_t IsNeutral() const
{ return ( Charge() == 0. ) ? kTRUE : kFALSE; }
Bool_t IsSortable() const
{ return kTRUE; }
Bool_t IsChargedLepton() const
 Particle is electron, mu or tau
Bool_t IsNeutrino() const
 Particle is neutrino or antineutrino
|| IsTauNeutrino()
 Particle is a neutrino
Bool_t IsPhoton() const
 Particle is Gamma
Bool_t IsElectron() const
 Particle is e+ or e-
Bool_t IsEMinus() const
 Particle is e-
Bool_t IsEPlus() const
 Particle is e+
Bool_t IsElectronNeutrino() const
 Particle is electron neutrino
Bool_t IsMuon() const
 Particle is mu+ or mu-
Bool_t IsMuPlus() const
 Particle is mu+-
Bool_t IsMuMinus() const
 Particle is mu-
Bool_t IsMuonNeutrino() const
 Particle is muon neutrino
Bool_t IsTau() const
 Particle is tau+ or tau-
Bool_t IsTauMinus() const
 Particle is tau-
Bool_t IsTauPlus() const
 Particle is tau+
Bool_t IsPion() const
 Particle is Pi+ or Pi-
Bool_t IsPiMinus() const
 Particle is Pi-
Bool_t IsPiPlus() const
 Particle is Pi+
Bool_t IsKaon() const
 Particle is K+ or K-
Bool_t IsProton() const
 Particle is p or antiproton
Bool_t IsKMinus() const
 Particle is K-
Bool_t IsKPlus() const
 Particle is Ki+
Bool_t IsStringOrCluster() const
 Particle is cluster (91) or string (92) on hadronisation level
Bool_t IsCharmQuark() const
 Particle is c or c_bar
Bool_t IsCharm() const
 Particle is c
Bool_t IsCharmBar() const
 Particle is c_bar
Bool_t IsBeautyQuark() const
 Particle is b or b_bar
Bool_t IsBeauty() const
 Particle is b
Bool_t IsBeautyBar() const
 Particle is b_bar
Bool_t IsDown() const
 Particle is d
Bool_t IsDownBar() const
 Particle is d_bar
Bool_t IsDownQuark() const
 Particle is d or d_bar quark
Bool_t IsUp() const
 Particle is u
Bool_t IsUpBar() const
 Particle is u_bar
Bool_t IsUpQuark() const
 Particle is u or u_bar quark
Bool_t IsStrange() const
 Particle is s
Bool_t IsStrangeBar() const
 Particle is s_bar
Bool_t IsStrangeQuark() const
 Particle is s or s_bar quark
Bool_t IsTopQuark() const
 Particle is t or t_bar
Bool_t IsTop() const
 Particle is t
Bool_t IsTopBar() const
 Particle is t_bar
Bool_t IsLightQuark() const
 Is any of u, u_bar, d, d_bar or s, s_bar
Bool_t IsQuark() const
 Is any of u, u_bar, d, d_bar, s, s_bar, c, c_bar, b, b_bar, t, t_bar
Bool_t IsBeautyPrimeQuark() const
 Particle is b' or b'_bar
Bool_t IsBeautyPrime() const
 Particle is b'
Bool_t IsBeautyPrimeBar() const
 Particle is b'_bar
Bool_t IsTopPrimeQuark() const
 Particle is t' or t'_bar
Bool_t IsTopPrime() const
 Particle is t'
Bool_t IsTopPrimeBar() const
 Particle is t'_bar
Bool_t IsWBoson() const
 Particle is W+ or W-
Bool_t IsWPlus() const
 Particle is W+
Bool_t IsWMinus() const
 Particle is W-
Bool_t IsZ0Boson() const
 Particle is Z0
Double_t M(Option_t* option = "PDG") const
 Alias for Mass()
Float_t Mt2() const
     Mt^2 = E^2 - Pz^2 = M^2 + Pt^2    (GeV^2)

 !!! Warning: Do not use this function for   !!!
 !!! the reconstruction of heavy particles   !!!
 !!! like the W-boson with the Jacobian-peak !!!
 !!! method. Use Mperp() instead             !!!

Float_t Mt() const
     Mt = Sqrt(E^2 - Pz^2 = M^2 + Pt^2)   (GeV)

 !!! Warning: Do not use this function for   !!!
 !!! the reconstruction of heavy particles   !!!
 !!! like the W-boson with the Jacobian-peak !!!
 !!! method. Use Mperp() instead             !!!

Float_t Mperp() const
 Mt = Sqrt(Et^2 - Pt^2)     (GeV)
 Use this function for the Jacobian-peak method

 Note that the mass is signed
const TLorentzVector& P() const
 Return 4-momentum vector (GeV)
TVector3 P3() const
 Return 3-momentum vector (GeV)
Float_t Pmag()
 Return magnitude of 3-momentum (GeV)

 !!! The return value is cached. Do not use this function !!!
 !!! inside TTree::Draw() or TTree::Scan().               !!!
 !!! Use eg "fElectrons.fP.Mag()" instead.                !!!

Float_t Pt()
 Return transvere momentum (GeV)

 !!! The return value is cached. Do not use this function !!!
 !!! inside TTree::Draw() or TTree::Scan().               !!!
 !!! Use eg "fElectrons.fP.Perp()" instead.               !!!

Float_t Px() const
 Return momentum x component (GeV)
Float_t Py() const
 Return momentum y component (GeV)
Float_t Pz() const
 Return momentum z component (GeV)
Float_t E() const
 Return energy (GeV)
Float_t Et() const
 Return transverse energy (GeV)
Float_t Eta()
 Return pseudo-rapidity

 !!! The return value is cached. Do not use this function !!!
 !!! inside TTree::Draw() or TTree::Scan().               !!!
 !!! Use eg "fElectrons.fP.Eta()" instead.                !!!

Float_t Phi()
 Return azimuth (rad)

 !!! The return value is cached. Do not use this function !!!
 !!! inside TTree::Draw() or TTree::Scan().               !!!
 !!! Use eg "fElectrons.fP.Phi()" instead.                !!!

Float_t Theta()
 Return polar angle (rad)

 !!! The return value is cached. Do not use this function !!!
 !!! inside TTree::Draw() or TTree::Scan().               !!!
 !!! Use eg "fElectrons.fP.Theta()" instead.              !!!

Float_t CosTheta()
 Return cosine of the polar angle

 !!! The return value is cached. Do not use this function !!!
 !!! inside TTree::Draw() or TTree::Scan().               !!!
 !!! Use eg "fElectrons.fP.CosTheta()" instead.           !!!

void SetId(Int_t Id)
{ fId = Id; }
void SetPdgCode(Int_t pdg)
{ fPdgCode = pdg; }
void SetPxPyPzE(Float_t Px, Float_t Py, Float_t Pz, Float_t E)
void SetPtEtaPhiE(Float_t Pt, Float_t Eta, Float_t Phi, Float_t E)
void SetPtEtaPhiM(Float_t Pt, Float_t Eta, Float_t Phi, Float_t M)
void SetP(const TLorentzVector& p)