]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG4/AliAODParticleCorrelation.h
New more general analysis implemention for particle identification and correlation...
[u/mrichter/AliRoot.git] / PWG4 / AliAODParticleCorrelation.h
1 #ifndef AliAODParticleCorrelation_H
2 #define AliAODParticleCorrelation_H
3 /* Copyright(c) 1998-2007, ALICE Experiment at CERN, All rights reserved. *
4  * See cxx source for full Copyright notice                               */
5
6 /* $Id: AliAODParticleCorrelation.h 20609 2007-08-31 11:17:02Z morsch $ */
7
8 //-------------------------------------------------------------------------
9 //     Copy of AOD photon class
10 //     Author: Yves Schutz, CERN, Gustavo Conesa, INFN
11 //-------------------------------------------------------------------------
12
13 #include "TLorentzVector.h"
14 #include "AliVParticle.h"
15 //#include "AliAODVertex.h"
16 #include "AliAODJet.h"
17 #include "TString.h"
18
19 class AliAODParticleCorrelation : public AliVParticle {
20
21  public:
22     AliAODParticleCorrelation();
23     AliAODParticleCorrelation(Double_t px, Double_t py, Double_t pz, Double_t e);
24     AliAODParticleCorrelation(TLorentzVector & p);  
25     virtual ~AliAODParticleCorrelation();
26     AliAODParticleCorrelation(const AliAODParticleCorrelation& photon); 
27     AliAODParticleCorrelation& operator=(const AliAODParticleCorrelation& photon);
28
29     // AliVParticle methods
30     virtual Double_t Px()         const { return fMomentum->Px();      }
31     virtual Double_t Py()         const { return fMomentum->Py();      }
32     virtual Double_t Pz()         const { return fMomentum->Pz();      }
33     virtual Double_t Pt()         const { return fMomentum->Pt();      }
34     virtual Double_t P()          const { return fMomentum->P();       }
35     virtual Bool_t   PxPyPz(Double_t p[3]) const { p[0] = Px(); p[1] = Py(); p[2] = Pz(); return kTRUE; }
36     virtual Double_t OneOverPt()  const { return 1. / fMomentum->Pt(); }
37     virtual Double_t Phi()        const;
38     virtual Double_t Theta()      const { return fMomentum->Theta();   }
39     virtual Double_t E()          const { return fMomentum->E();       }
40     virtual Double_t M()          const { return fMomentum->M();       }
41     virtual Double_t Eta()        const { return fMomentum->Eta();     }
42     virtual Double_t Y()          const { return fMomentum->Rapidity();}
43     virtual Double_t Xv()         const {return -999.;} // put reasonable values here
44     virtual Double_t Yv()         const {return -999.;} //
45     virtual Double_t Zv()         const {return -999.;} //
46     virtual Bool_t   XvYvZv(Double_t x[3]) const { x[0] = Xv(); x[1] = Yv(); x[2] = Zv(); return kTRUE; }  
47     virtual void     Print(Option_t* /*option*/) const;
48
49     //
50     //Dummy
51     virtual Short_t Charge()      const { return 0;}
52     virtual const Double_t* PID() const { return NULL;}
53     //
54
55     virtual Float_t GetR()   const {return fR ; }
56     virtual Int_t GetPdg()   const {return fPdg ; }
57     virtual Int_t GetTag()   const {return fTag ; }
58     virtual Int_t GetLabel()   const {return fLabel ; }
59     virtual TString GetDetector()   const {return fDetector ; }
60
61     virtual void SetR(Float_t r)   {fR = r ; }
62     virtual void SetPdg(Int_t pdg)   {fPdg = pdg ; }
63     virtual void SetTag(Int_t tag)   {fTag = tag ; }
64     virtual void SetLabel(Int_t l)   {fLabel = l ; }
65     virtual void SetDetector(TString d)   {fDetector = d ; }
66
67     virtual TRefArray* GetRefTracks()           const { return  fRefTracks;}
68     virtual void     AddTrack(TObject *tr) {fRefTracks->Add(tr);}
69     TObject* GetTrack(Int_t i) {return fRefTracks->At(i);}
70
71     virtual TRefArray* GetRefClusters()           const { return  fRefClusters;}
72     virtual void     AddCluster(TObject *tr) {fRefClusters->Add(tr);}
73     TObject* GetCluster(Int_t i) {return fRefClusters->At(i);}
74
75     virtual TRefArray* GetRefIsolationConeTracks()           const { return  fRefIsolationConeTracks;}
76     virtual void     AddIsolationConeTrack(TObject *tr) {fRefIsolationConeTracks->Add(tr);}
77     TObject* GetIsolationConeTrack(Int_t i) {return fRefIsolationConeTracks->At(i);}
78
79     virtual TRefArray* GetRefIsolationConeClusters()           const { return  fRefIsolationConeClusters;}
80     virtual void     AddIsolationConeCluster(TObject *tr) {fRefIsolationConeClusters->Add(tr);}
81     TObject* GetIsolationConeCluster(Int_t i) {return fRefIsolationConeClusters->At(i);}
82
83     virtual TRefArray* GetRefBackgroundTracks()           const { return  fRefBackgroundTracks;}
84     virtual void     AddBackgroundTrack(TObject *tr) {fRefBackgroundTracks->Add(tr);}
85     TObject* GetBackgroundTrack(Int_t i) {return fRefBackgroundTracks->At(i);}
86
87     virtual TRefArray* GetRefBackgroundClusters()           const { return  fRefBackgroundClusters;}
88     virtual void     AddBackgroundCluster(TObject *tr) {fRefBackgroundClusters->Add(tr);}
89     TObject* GetBackgroundCluster(Int_t i) {return fRefBackgroundClusters->At(i);}
90
91     virtual void SetLeadingDetector(TString d)   {fLeadingDetector = d ; }
92     virtual TString GetLeadingDetector()   const {return fLeadingDetector ; }
93     
94     virtual TLorentzVector * GetLeading()           const { return  fLeading;}
95     virtual void  SetLeading(TLorentzVector *lead) {fLeading = lead;}
96
97     virtual TLorentzVector * GetCorrelatedJet()           const { return  fCorrJet;}
98     virtual void  SetCorrelatedJet(TLorentzVector *jet) {fCorrJet = jet;}
99
100     void SetRefJet(AliAODJet* jet)  { fRefJet = jet;}
101     //AliAODJet* GetJet() {return ((AliAODJet*) fRefJet);}
102     TRef GetRefJet() {return fRefJet;}
103
104  private:
105     TLorentzVector* fMomentum;           // Photon 4-momentum vector
106     Int_t          fPdg; // id of particle
107     Int_t          fTag; // tag of particle (decay, fragment, prompt photon)
108     Int_t          fLabel; // MC label
109     TString          fDetector; // Detector where particle was measured.
110     Float_t        fR ; // Isolation cone size
111     TRefArray*     fRefTracks;  // array of references to the tracks belonging to the jet / all selected hadrons  
112     TRefArray*     fRefClusters; // array of references to the clusters belonging to the jet / all selected hadrons  
113     
114     TRefArray*     fRefIsolationConeTracks;  // array of references to the tracks belonging to the cone around direct particle candidate  
115     TRefArray*     fRefIsolationConeClusters; // array of references to the clusters belonging to the  cone around direct particle candidate  
116
117     TRefArray*     fRefBackgroundTracks;  // array of references to the tracks for background stimation
118     TRefArray*     fRefBackgroundClusters; // array of references to the clusters for background stimation 
119
120     TString     fLeadingDetector; // Detector where leading particle was measured.
121     TLorentzVector* fLeading;     // Leading Particle 4-momentum vector
122     
123     TLorentzVector* fCorrJet;     // Jet  4-momentum vector
124     
125     TRef        fRefJet; // Rerence to jet found with JETAN and correlated with particle
126
127   ClassDef(AliAODParticleCorrelation,1);
128 };
129
130 inline Double_t AliAODParticleCorrelation::Phi() const
131 {
132     // Return phi
133     Double_t phi = fMomentum->Phi();
134     if (phi < 0.) phi += 2. * TMath::Pi();
135     return phi;
136 }
137
138 #endif