]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HBTAN/AliHBTParticle.h
72a55f713f578fec29ed815ce1b4b942832ba59a
[u/mrichter/AliRoot.git] / HBTAN / AliHBTParticle.h
1 #ifndef ALIHBTPARTICLE
2 #define ALIHBTPARTICLE
3 //___________________________________________________________
4 /////////////////////////////////////////////////////////////
5 //
6 // class AliHBTParticle
7 //
8 // Ali HBT Particle: simplified class TParticle
9 // Simplified in order to minimize the size of object
10 //  - we want to keep a lot of such a objects in memory
11 // Additionaly adjusted for HBT Analysies purposes
12 //
13 /////////////////////////////////////////////////////////////
14
15 #include <TObject.h>
16 #include <TLorentzVector.h>
17 #include <TMath.h>
18 #include <TDatabasePDG.h>
19
20 #include "AliConst.h"
21
22 class TParticle;
23 //class AliHBTTrackPoints;
24
25 class AliHBTParticle : public TObject
26 {
27
28 public:
29                                 // ****** constructors and destructor
30   AliHBTParticle();
31   AliHBTParticle(const AliHBTParticle& in); 
32  
33   AliHBTParticle(Int_t pdg, Int_t idx, Double_t px, Double_t py, Double_t pz, Double_t etot,
34                  Double_t vx, Double_t vy, Double_t vz, Double_t time);
35
36   AliHBTParticle(Int_t pdg, Float_t prob, Int_t idx, Double_t px, Double_t py, Double_t pz, Double_t etot,
37                  Double_t vx, Double_t vy, Double_t vz, Double_t time);
38
39   AliHBTParticle(const TParticle& p,Int_t idx);
40
41   virtual ~AliHBTParticle();
42   
43   AliHBTParticle& operator=(const AliHBTParticle& in); 
44   
45   void           SetPIDprobability(Int_t pdg, Float_t prob = 1.0);
46   Float_t        GetPIDprobability(Int_t pdg);
47   
48   Int_t          GetPdgCode      () const { return (fPids)?fPids[fPdgIdx]:0;}
49   Int_t          GetPid          () const { return GetPdgCode();}
50   Float_t        GetPidProb      () const { return (fPidProb)?fPidProb[fPdgIdx]:0;}
51   
52   Int_t          GetUID          () const { return fIdxInEvent;}
53   Int_t          GetNumberOfPids () const { return fNPids;}
54   Int_t          GetNthPid         (Int_t idx) const;
55   Float_t        GetNthPidProb     (Int_t idx) const;
56       
57   void           SetPdgCode(Int_t pdg, Float_t prob = 1.0);
58   Double_t       GetCalcMass     () const { return fCalcMass; }
59   Double_t       GetMass         ()       { return (GetPDG())?GetPDG()->Mass():-1.;}
60
61   TParticlePDG*  GetPDG          (){return TDatabasePDG::Instance()->GetParticle(GetPdgCode());}
62
63   Int_t          Beauty          ()  { return GetPDG()->Beauty(); }
64   Int_t          Charm           ()  { return GetPDG()->Charm(); }
65   Int_t          Strangeness     ()  { return GetPDG()->Strangeness();}
66   void ProductionVertex(TLorentzVector &v) { v.SetXYZT(fVx,fVy,fVz,fVt);}
67
68
69   Double_t         Vx    () const { return fVx;}
70   Double_t         Vy    () const { return fVy;}
71   Double_t         Vz    () const { return fVz;}
72   Double_t         T     () const { return fVt;}
73
74   Double_t         Px    () const { return fPx; } //X coordinate of the momentum
75   Double_t         Py    () const { return fPy; } //Y coordinate of the momentum
76   Double_t         Pz    () const { return fPz; } //Z coordinate of the momentum
77   Double_t         P     () const                 //momentum
78     { return TMath::Sqrt(fPx*fPx+fPy*fPy+fPz*fPz); }
79   
80   void Momentum(TLorentzVector &v) { v.SetPxPyPzE(fPx,fPy,fPz,fE);}
81     
82   Double_t         Pt    () const  //transverse momentum
83     { return TMath::Sqrt(fPx*fPx+fPy*fPy); }
84   Double_t         Energy() const { return fE; }
85   
86                                    //Pseudo Rapidity
87   Double_t         Eta   () const { if (P() != fPz) return 0.5*TMath::Log((P()+fPz)/(P()-fPz)); 
88                                     else return 1.e30;}
89
90                                    //Rapidity
91   Double_t         Y     () const { if (fE  != fPz) return 0.5*TMath::Log((fE+fPz)/(fE-fPz));
92                                     else return 1.e30;}
93
94   Double_t         Phi   () const { return TMath::Pi()+TMath::ATan2(-fPy,-fPx); }
95
96   Double_t         Theta () const { return (fPz==0)?TMath::PiOver2():TMath::ACos(fPz/P()); }
97
98   // setters
99
100   void           SetMomentum(Double_t px, Double_t py, Double_t pz, Double_t e)
101                              {fPx=px; fPy=py; fPz=pz; fE=e;}
102   void           SetMomentum(const TLorentzVector& p)
103                              {SetMomentum(p.Px(),p.Py(),p.Pz(),p.Energy());}
104
105   void           SetProductionVertex(Double_t vx, Double_t vy, Double_t vz, Double_t t)
106                              {fVx=vx; fVy=vy; fVz=vz; fVt=t;}
107   void           SetProductionVertex(const TLorentzVector& v)
108                              {SetProductionVertex(v.X(),v.Y(),v.Z(),v.T());}
109   void           SetCalcMass(Double_t mass) {fCalcMass = mass;}
110   
111   void           SetUID(Int_t id){fIdxInEvent = id;}
112   
113   const Char_t*  GetName() const; 
114   void           Print() const;
115
116   static void    SetDebug(Int_t dbg=1){fgDebug=dbg;}
117   static Int_t   GetDebug(){return fgDebug;}
118   static Int_t   fgDebug; //debug printout level
119   
120 protected:
121   Int_t          GetPidSlot(Int_t pdg) const;//returns position of the given PID in fPids (and fPidProb) array.
122
123 private:
124   Char_t         fPdgIdx;               // index of PDG code of the particle in fPids
125   Int_t          fIdxInEvent;           // index of a particle: the same particle can appear in the event
126                                         //  many times with different pid's. Idx allows to check that they are the same particles
127   Int_t          fNPids;                // number of non-zero proboble Pids
128   Int_t         *fPids;                 // [fNPids] Array with PIDs
129   Float_t       *fPidProb;              // [fNPids] PIDs probabilities
130   Double_t       fCalcMass;             // Calculated mass
131
132   Double_t       fPx;                   // x component of momentum
133   Double_t       fPy;                   // y component of momentum
134   Double_t       fPz;                   // z component of momentum
135   Double_t       fE;                    // Energy
136
137   Double_t       fVx;                   // x of production vertex
138   Double_t       fVy;                   // y of production vertex
139   Double_t       fVz;                   // z of production vertex
140   Double_t       fVt;                   // t of production vertex
141
142   ClassDef(AliHBTParticle,2)  // TParticle vertex particle information
143 };
144
145 #endif