00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016 #ifndef ROOT_TParticle
00017 #define ROOT_TParticle
00018
00019 #ifndef ROOT_TNamed
00020 #include "TNamed.h"
00021 #endif
00022 #ifndef ROOT_TAttLine
00023 #include "TAttLine.h"
00024 #endif
00025 #ifndef ROOT_TAtt3D
00026 #include "TAtt3D.h"
00027 #endif
00028 #ifndef ROOT_TLorentzVector
00029 #include "TLorentzVector.h"
00030 #endif
00031
00032 class TParticlePDG;
00033
00034 class TParticle : public TObject, public TAttLine, public TAtt3D {
00035
00036
00037 protected:
00038
00039 Int_t fPdgCode;
00040 Int_t fStatusCode;
00041 Int_t fMother[2];
00042 Int_t fDaughter[2];
00043 Float_t fWeight;
00044
00045 Double_t fCalcMass;
00046
00047 Double_t fPx;
00048 Double_t fPy;
00049 Double_t fPz;
00050 Double_t fE;
00051
00052 Double_t fVx;
00053 Double_t fVy;
00054 Double_t fVz;
00055 Double_t fVt;
00056
00057 Double_t fPolarTheta;
00058 Double_t fPolarPhi;
00059
00060 TParticlePDG* fParticlePDG;
00061
00062
00063
00064
00065
00066
00067 public:
00068
00069 TParticle();
00070
00071 TParticle(Int_t pdg, Int_t status,
00072 Int_t mother1, Int_t mother2,
00073 Int_t daughter1, Int_t daughter2,
00074 Double_t px, Double_t py, Double_t pz, Double_t etot,
00075 Double_t vx, Double_t vy, Double_t vz, Double_t time);
00076
00077 TParticle(Int_t pdg, Int_t status,
00078 Int_t mother1, Int_t mother2,
00079 Int_t daughter1, Int_t daughter2,
00080 const TLorentzVector &p,
00081 const TLorentzVector &v);
00082
00083 TParticle(const TParticle &part);
00084
00085 virtual ~TParticle();
00086
00087 TParticle& operator=(const TParticle&);
00088
00089
00090
00091
00092 Int_t GetStatusCode () const { return fStatusCode; }
00093 Int_t GetPdgCode () const { return fPdgCode; }
00094 Int_t GetFirstMother () const { return fMother[0]; }
00095 Int_t GetMother (Int_t i) const { return fMother[i]; }
00096 Int_t GetSecondMother () const { return fMother[1]; }
00097 Bool_t IsPrimary () const { return fMother[0]>-1 ? kFALSE : kTRUE; }
00098 Int_t GetFirstDaughter() const { return fDaughter[0]; }
00099 Int_t GetDaughter (Int_t i) const { return fDaughter[i]; }
00100 Int_t GetLastDaughter () const { return fDaughter[1]; }
00101 Double_t GetCalcMass () const { return fCalcMass; }
00102 Double_t GetMass ();
00103 Int_t GetNDaughters () const { return fDaughter[1]>0 ? fDaughter[1]-fDaughter[0]+1 : 0;}
00104 Float_t GetWeight () const { return fWeight; }
00105 void GetPolarisation(TVector3 &v);
00106 TParticlePDG* GetPDG (Int_t mode = 0);
00107 Int_t Beauty ();
00108 Int_t Charm ();
00109 Int_t Strangeness ();
00110 void Momentum(TLorentzVector &v) const { v.SetPxPyPzE(fPx,fPy,fPz,fE); }
00111 void ProductionVertex(TLorentzVector &v) const { v.SetXYZT(fVx,fVy,fVz,fVt); }
00112
00113
00114
00115 Double_t Vx () const { return fVx; }
00116 Double_t Vy () const { return fVy; }
00117 Double_t Vz () const { return fVz; }
00118 Double_t T () const { return fVt; }
00119 Double_t R () const { return TMath::Sqrt(fVx*fVx+fVy*fVy); }
00120 Double_t Rho () const { return TMath::Sqrt(fVx*fVx+fVy*fVy+fVz*fVz); }
00121 Double_t Px () const { return fPx; }
00122 Double_t Py () const { return fPy; }
00123 Double_t Pz () const { return fPz; }
00124 Double_t P () const { return TMath::Sqrt(fPx*fPx+fPy*fPy+fPz*fPz); }
00125 Double_t Pt () const { return TMath::Sqrt(fPx*fPx+fPy*fPy); }
00126 Double_t Energy () const { return fE; }
00127 Double_t Eta () const
00128 {
00129 Double_t pmom = P();
00130 if (pmom != TMath::Abs(fPz)) return 0.5*TMath::Log((pmom+fPz)/(pmom-fPz));
00131 else return 1.e30;
00132 }
00133 Double_t Y () const
00134 {
00135 if (fE != TMath::Abs(fPz)) return 0.5*TMath::Log((fE+fPz)/(fE-fPz));
00136 else return 1.e30;
00137 }
00138
00139 Double_t Phi () const { return TMath::Pi()+TMath::ATan2(-fPy,-fPx); }
00140 Double_t Theta () const { return (fPz==0)?TMath::PiOver2():TMath::ACos(fPz/P()); }
00141
00142
00143
00144 void SetFirstMother (int code) { fMother[0] = code ; }
00145 void SetMother (int i, int code) { fMother[i] = code ; }
00146 void SetLastMother (int code) { fMother[1] = code ; }
00147 void SetFirstDaughter(int code) { fDaughter[0] = code ; }
00148 void SetDaughter(int i, int code) { fDaughter[i] = code ; }
00149 void SetLastDaughter(int code) { fDaughter[1] = code ; }
00150 void SetCalcMass(Double_t mass) { fCalcMass=mass;}
00151 void SetPdgCode(Int_t pdg);
00152 void SetPolarisation(Double_t polx, Double_t poly, Double_t polz);
00153 void SetPolarisation(const TVector3& v) {SetPolarisation(v.X(), v.Y(), v.Z());}
00154 void SetStatusCode(int status) {fStatusCode = status;}
00155 void SetWeight(Float_t weight = 1) {fWeight = weight; }
00156 void SetMomentum(Double_t px, Double_t py, Double_t pz, Double_t e) {fPx=px; fPy=py; fPz=pz; fE=e;}
00157 void SetMomentum(const TLorentzVector& p) {SetMomentum(p.Px(),p.Py(),p.Pz(),p.Energy());}
00158 void SetProductionVertex(Double_t vx, Double_t vy, Double_t vz, Double_t t) {fVx=vx; fVy=vy; fVz=vz; fVt=t;}
00159 void SetProductionVertex(const TLorentzVector& v) {SetProductionVertex(v.X(),v.Y(),v.Z(),v.T());}
00160
00161
00162
00163 virtual void Paint(Option_t *option = "");
00164 virtual void Print(Option_t *option = "") const;
00165 virtual void Sizeof3D() const;
00166 virtual Int_t DistancetoPrimitive(Int_t px, Int_t py);
00167 virtual void ExecuteEvent(Int_t event, Int_t px, Int_t py);
00168 virtual const char *GetName() const;
00169 virtual const char *GetTitle() const;
00170
00171
00172 ClassDef(TParticle,2)
00173 };
00174
00175 #endif