Logo ROOT   6.30.04
Reference Guide
 All Namespaces Files Pages
TParticle.h
Go to the documentation of this file.
1 // @(#)root/eg:$Id$
2 // Author: Rene Brun , Federico Carminati 26/04/99
3 
4 /*************************************************************************
5  * Copyright (C) 1995-2000, Rene Brun and Fons Rademakers. *
6  * All rights reserved. *
7  * *
8  * For the licensing terms see $ROOTSYS/LICENSE. *
9  * For the list of contributors see $ROOTSYS/README/CREDITS. *
10  *************************************************************************/
11 //////////////////////////////////////////////////////////////////////////
12 // //
13 // TParticle: defines equivalent of HEPEVT particle //
14 //////////////////////////////////////////////////////////////////////////
15 
16 #ifndef ROOT_TParticle
17 #define ROOT_TParticle
18 
19 #include "TNamed.h"
20 #include "TAttLine.h"
21 #include "TAtt3D.h"
22 #include "TLorentzVector.h"
23 
24 class TParticlePDG;
25 
26 class TParticle : public TObject, public TAttLine, public TAtt3D {
27 
28 
29 protected:
30 
31  Int_t fPdgCode; // PDG code of the particle
32  Int_t fStatusCode; // generation status code
33  Int_t fMother[2]; // Indices of the mother particles
34  Int_t fDaughter[2]; // Indices of the daughter particles
35  Float_t fWeight; // particle weight
36 
37  Double_t fCalcMass; // Calculated mass
38 
39  Double_t fPx; // x component of momentum
40  Double_t fPy; // y component of momentum
41  Double_t fPz; // z component of momentum
42  Double_t fE; // Energy
43 
44  Double_t fVx; // x of production vertex
45  Double_t fVy; // y of production vertex
46  Double_t fVz; // z of production vertex
47  Double_t fVt; // t of production vertex
48 
49  Double_t fPolarTheta; // Polar angle of polarisation
50  Double_t fPolarPhi; // azymutal angle of polarisation
51 
52  mutable TParticlePDG* fParticlePDG; //! reference to the particle record in PDG database
53  //----------------------------------------------------------------------------
54  // functions
55  //----------------------------------------------------------------------------
56 public:
57  // ****** constructors and destructor
58  TParticle();
59 
60  TParticle(Int_t pdg, Int_t status,
61  Int_t mother1, Int_t mother2,
62  Int_t daughter1, Int_t daughter2,
63  Double_t px, Double_t py, Double_t pz, Double_t etot,
64  Double_t vx, Double_t vy, Double_t vz, Double_t time);
65 
66  TParticle(Int_t pdg, Int_t status,
67  Int_t mother1, Int_t mother2,
68  Int_t daughter1, Int_t daughter2,
69  const TLorentzVector &p,
70  const TLorentzVector &v);
71 
72  TParticle(const TParticle &part);
73 
74  virtual ~TParticle();
75 
76  TParticle& operator=(const TParticle&);
77 
78 // virtual TString* Name () const { return fName.Data(); }
79 // virtual char* GetName() const { return fName.Data(); }
80 
81  Double_t Ek () const { return fE-fCalcMass; }
82  Int_t GetStatusCode () const { return fStatusCode; }
83  Int_t GetPdgCode () const { return fPdgCode; }
84  Int_t GetFirstMother () const { return fMother[0]; }
85  Int_t GetMother (Int_t i) const { return fMother[i]; }
86  Int_t GetSecondMother () const { return fMother[1]; }
87  Bool_t IsPrimary () const { return fMother[0]>-1 ? kFALSE : kTRUE; } //Is this particle primary one?
88  Int_t GetFirstDaughter() const { return fDaughter[0]; }
89  Int_t GetDaughter (Int_t i) const { return fDaughter[i]; }
90  Int_t GetLastDaughter () const { return fDaughter[1]; }
91  Double_t GetCalcMass () const { return fCalcMass; }
92  Double_t GetMass () const;
93  Int_t GetNDaughters () const { return fDaughter[1]>0 ? fDaughter[1]-fDaughter[0]+1 : 0;}
94  Float_t GetWeight () const { return fWeight; }
95  void GetPolarisation(TVector3 &v) const;
96  TParticlePDG* GetPDG (Int_t mode = 0) const;
97  Int_t Beauty () const;
98  Int_t Charm () const;
99  Int_t Strangeness () const;
100  Double_t PhiX () const { return TMath::Pi()+TMath::ATan2(-fPz,-fPy); } // note that PhiX() returns an angle between 0 and 2pi
101  Double_t PhiY () const { return TMath::Pi()+TMath::ATan2(-fPx,-fPz); } // note that PhiY() returns an angle between 0 and 2pi
102  Double_t PhiZ () const { return TMath::Pi()+TMath::ATan2(-fPy,-fPx); } // note that PhiZ() returns an angle between 0 and 2pi
103  Double_t ThetaX () const { return (fPx==0)?TMath::PiOver2():TMath::ACos(fPx/P()); }
104  Double_t ThetaY () const { return (fPy==0)?TMath::PiOver2():TMath::ACos(fPy/P()); }
105  Double_t ThetaZ () const { return (fPz==0)?TMath::PiOver2():TMath::ACos(fPz/P()); }
106  Double_t GetPolarTheta () const { return fPolarTheta; }
107  Double_t GetPolarPhi () const { return fPolarPhi; }
108  void GetPolarisation (Double_t &theta, Double_t &phi) const { theta = fPolarTheta; phi = fPolarPhi; }
109  void SetPolarTheta (Double_t theta) { fPolarTheta = theta; }
110  void SetPolarPhi (Double_t phi) { fPolarPhi = phi; }
111  void SetPolarisation (Double_t theta, Double_t phi) { fPolarTheta = theta; fPolarPhi = phi; }
112  void Momentum(TLorentzVector &v) const { v.SetPxPyPzE(fPx,fPy,fPz,fE); }
113  void ProductionVertex(TLorentzVector &v) const { v.SetXYZT(fVx,fVy,fVz,fVt); }
114 
115  Double_t Theta(const TParticle &p) // Returns the angle between momenta of particles
116  {
117  Double_t v = P()*p.P();
118  if (v == 0) v = 1; else v = (fPx*p.Px()+fPy*p.Py()+fPz*p.Pz())/v;
119  if (v > 1) v = 1; else if (v < -1) v = -1; // just a precaution
120  return TMath::ACos(v);
121  }
122 
123 // ****** redefine several most oftenly used methods of LORENTZ_VECTOR
124 
125  Double_t Vx () const { return fVx; }
126  Double_t Vy () const { return fVy; }
127  Double_t Vz () const { return fVz; }
128  Double_t T () const { return fVt; }
129  Double_t R () const { return TMath::Sqrt(fVx*fVx+fVy*fVy); } //Radius of production vertex in cylindrical system
130  Double_t Rho () const { return TMath::Sqrt(fVx*fVx+fVy*fVy+fVz*fVz); } //Radius of production vertex in spherical system
131  Double_t Px () const { return fPx; }
132  Double_t Py () const { return fPy; }
133  Double_t Pz () const { return fPz; }
134  Double_t P () const { return TMath::Sqrt(fPx*fPx+fPy*fPy+fPz*fPz); }
135  Double_t Pt () const { return TMath::Sqrt(fPx*fPx+fPy*fPy); }
136  Double_t Energy () const { return fE; }
137  Double_t Eta () const
138  {
139  Double_t pmom = P();
140  if (pmom != TMath::Abs(fPz)) return 0.5*TMath::Log((pmom+fPz)/(pmom-fPz));
141  else return 1.e30;
142  }
143  Double_t Y () const
144  {
145  if (fE != TMath::Abs(fPz)) return 0.5*TMath::Log((fE+fPz)/(fE-fPz));
146  else return 1.e30;
147  }
148 
149  Double_t Phi () const { return TMath::Pi()+TMath::ATan2(-fPy,-fPx); } // note that Phi() returns an angle between 0 and 2pi
150  Double_t Theta () const { return (fPz==0)?TMath::PiOver2():TMath::ACos(fPz/P()); }
151 
152  // setters
153 
154  void SetFirstMother (int code) { fMother[0] = code ; }
155  void SetMother (int i, int code) { fMother[i] = code ; }
156  void SetLastMother (int code) { fMother[1] = code ; }
157  void SetFirstDaughter(int code) { fDaughter[0] = code ; }
158  void SetDaughter(int i, int code) { fDaughter[i] = code ; }
159  void SetLastDaughter(int code) { fDaughter[1] = code ; }
160  void SetCalcMass(Double_t mass) { fCalcMass=mass;}
161  void SetPdgCode(Int_t pdg);
162  void SetPolarisation(Double_t polx, Double_t poly, Double_t polz);
163  void SetPolarisation(const TVector3& v) {SetPolarisation(v.X(), v.Y(), v.Z());}
164  void SetStatusCode(int status) {fStatusCode = status;}
165  void SetWeight(Float_t weight = 1) {fWeight = weight; }
166  void SetMomentum(Double_t px, Double_t py, Double_t pz, Double_t e) {fPx=px; fPy=py; fPz=pz; fE=e;}
167  void SetMomentum(const TLorentzVector& p) {SetMomentum(p.Px(),p.Py(),p.Pz(),p.Energy());}
168  void SetProductionVertex(Double_t vx, Double_t vy, Double_t vz, Double_t t) {fVx=vx; fVy=vy; fVz=vz; fVt=t;}
169  void SetProductionVertex(const TLorentzVector& v) {SetProductionVertex(v.X(),v.Y(),v.Z(),v.T());}
170 
171  // ****** overloaded functions of TObject
172 
173  virtual void Paint(Option_t *option = "");
174  virtual void Print(Option_t *option = "") const;
175  virtual void Sizeof3D() const;
176  virtual Int_t DistancetoPrimitive(Int_t px, Int_t py);
177  virtual void ExecuteEvent(Int_t event, Int_t px, Int_t py);
178  virtual const char *GetName() const;
179  virtual const char *GetTitle() const;
180 
181 
182  ClassDef(TParticle,2) // TParticle vertex particle information
183 };
184 
185 #endif