]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGCF/FEMTOSCOPY/ESE/AliAnalysisTaskFemtoESE.h
Merge branch 'master' of https://git.cern.ch/reps/AliRoot
[u/mrichter/AliRoot.git] / PWGCF / FEMTOSCOPY / ESE / AliAnalysisTaskFemtoESE.h
1 #ifndef AliAnalysisTaskFemtoESE_cxx
2 #define AliAnalysisTaskFemtoESE_cxx
3
4
5 class TH1F;
6 class TH2F;
7 class TH3F;
8 class TH1D;
9 class TH2D;
10 class TH3D;
11
12 class TProfile;
13 class TProfile2D;
14
15 class AliESDEvent;
16 class AliAODEvent;
17 class AliESDtrackCuts;
18 class AliESDpid;
19 class AliHelperPID;
20 class AliEventPoolManager;
21 class AliSpectraAODEventCuts;
22 class AliSpectraAODTrackCuts;
23 class AliAODTrack;
24 class AliEventPool;
25
26 //class TStopwatch;
27
28 class AliFemtoESEBasicParticle;
29
30 #include "AliAnalysisTask.h"
31 #include "AliAnalysisTaskSE.h"
32 #include "AliESDpid.h"
33 #include "AliAODPid.h"
34
35 #include "AliAODEvent.h"
36 #include "AliAnalysisManager.h"
37 #include "AliAODInputHandler.h"
38 #include "AliSpectraAODEventCuts.h"
39 #include "AliSpectraAODTrackCuts.h"
40
41 class AliAnalysisTaskFemtoESE : public AliAnalysisTaskSE {
42  public:
43
44   AliAnalysisTaskFemtoESE();
45   AliAnalysisTaskFemtoESE(const char* name);
46   virtual ~AliAnalysisTaskFemtoESE();
47   AliAnalysisTaskFemtoESE(const AliAnalysisTaskFemtoESE &/*obj*/); 
48   AliAnalysisTaskFemtoESE &operator=(const AliAnalysisTaskFemtoESE &/*obj*/);
49
50   virtual void   UserCreateOutputObjects();
51   virtual void   UserExec(Option_t *option);
52   void TrackLoop(TObjArray *tracks, AliEventPool *pool, Double_t psiEP, Float_t centralityPercentile);
53   virtual void   Terminate(Option_t *);
54
55   AliHelperPID* GetHelperPID() { return fHelperPID; }
56   void SetHelperPID(AliHelperPID* pid){ fHelperPID = pid; }
57   AliSpectraAODEventCuts* GetEventCuts() {return fEventCuts;}
58   void SetEventCuts(AliSpectraAODEventCuts* cuts) {fEventCuts = cuts;}
59   AliSpectraAODTrackCuts* GetTrackCuts() {return fTrackCuts;}
60   void SetTrackCuts(AliSpectraAODTrackCuts* cuts) {fTrackCuts = cuts;}
61
62   Int_t GetTrackFilterBit(){return fFilterBit;}
63   void SetTrackFilterBit(Int_t bit){fFilterBit=bit;}
64   void SetEventSelectionBit( UInt_t val ) {fSelectBit = val;}
65   void SetIsLHC10h(Bool_t val) {bIsLHC10h = val;}
66   void SetMinSepPair(Double_t eta,Double_t phi){fMinSepPairEta = eta; fMinSepPairPhi = phi;}
67   void SetShareFraction(Double_t val) {fShareFraction = val;}
68   Double_t GetShareFraction() {return fShareFraction;}
69   void SetShareQuality(Double_t val) {fShareQuality = val;}
70   Double_t GetShareQuality() {return fShareQuality;}
71   void SetQPercCuts(Double_t min, Double_t max) {fMinQPerc = min; fMaxQPerc = max;}
72   Double_t GetMinQPerc(){return fMinQPerc;}
73   Double_t GetMaxQPerc(){return fMaxQPerc;}
74   void SetQPercDetector(Int_t det){fQPercDet = det;};
75   void SetEPDetector(Int_t det){fEPDet = det;};
76   void SetNMixingTracks(Int_t n){fMixingTracks = n;};
77
78   void SetKtBins(Int_t n, Double_t* bins);
79   void SetEPBins(Int_t n, Double_t min, Double_t max);
80   void SetCentBins(Int_t n, Double_t* bins);
81   void SetVzBins(Int_t n, Double_t* bins);
82
83   Double_t GetQPercLHC11h(Double_t qvec);
84
85  private:
86   Double_t GetQinv2(Double_t[], Double_t[]);
87   void GetQosl(Double_t[], Double_t[], Double_t&, Double_t&, Double_t&);
88   Bool_t TrackCut(AliAODTrack* ftrack);
89   Bool_t EventCut(/*AliAODEvent* fevent*/);
90   Bool_t PairCut(AliFemtoESEBasicParticle* ftrack1, AliFemtoESEBasicParticle* ftrack2, Bool_t mix);
91   Double_t DeltaPhiStar(AliAODTrack* ftrack1, AliAODTrack* ftrack2, Double_t r);
92   TObjArray* CloneAndReduceTrackList(TObjArray* tracks, Double_t psi);
93   Double_t GetDeltaPhiEP(Double_t px1, Double_t py1, Double_t px2, Double_t py2, Double_t psi);
94   Bool_t FindBin(Double_t kt, Double_t phi, Double_t cent, Int_t& a, Int_t& b, Int_t&c);
95
96   AliAODEvent            *fAOD; //!    // AOD object
97   TList                  *fOutputList; //! Compact Output list
98   //AliPIDResponse         *fPIDResponse; //! PID response object; equivalent to AliAODpidUtil
99   AliHelperPID*     fHelperPID;      // points to class for PID
100   AliEventPoolManager**     fPoolMgr;         //[2] event pool manager
101   AliSpectraAODEventCuts* fEventCuts;
102   AliSpectraAODTrackCuts* fTrackCuts;
103
104   Int_t          fFilterBit;         // track selection cuts
105   UInt_t         fSelectBit;            // Select events according to AliAnalysisTaskJetServices bit maps 
106   Bool_t         bIsLHC10h;
107   Int_t fEventCounter;
108   Int_t fMixingTracks;
109   Double_t fBfield;
110   Double_t fMinSepPairEta;
111   Double_t fMinSepPairPhi;
112   Double_t fShareQuality;
113   Double_t fShareFraction;
114
115   Int_t nCountSamePairs;
116   Int_t nCountMixedPairs;
117   Int_t nCountTracks;
118
119   Double_t fMinQPerc;
120   Double_t fMaxQPerc;
121
122   Int_t fQPercDet; // detector used for q-vector (0-V0A, 1-V0C)
123   Int_t fEPDet; // detector used for event plane (0-V0A, 1-V0C)
124
125   Double_t fPsiEPmix;
126   Double_t fPsiEPmixtemp;
127
128   // binning for histograms
129   Int_t nKtBins;
130   Int_t nKtBins1;
131   Double_t* ktBins; //[nKtBins1]
132   Int_t nEPBins;
133   Int_t nEPBins1;
134   Double_t* epBins; //[nEPBins1]
135   Int_t nCentBins;
136   Int_t nCentBins1;
137   Double_t* centBins; //[nCentBins1]
138   Int_t nVzBins;
139   Int_t nVzBins1;
140   Double_t* vzBins; //[nVzBins1]
141
142   Double_t vertex[3];
143
144   TH3F**** hq;
145   TH3F**** hqmix;
146
147   Int_t nqPercBinsLHC11h;
148   Double_t* qPercBinsLHC11h; //[nqPercBinsLHC11h]
149
150   //TStopwatch *stopwatch;
151
152   // histograms
153   TH1D *hpx;
154   TH1D *hpy;
155   TH1D *hpz;
156   TH1D *hpt;
157   TH1D *hE;
158   TH2D *hphieta;
159   TH2D *hphieta_pid;
160   TH1D *hpt_pid;
161   TH2D *hvzcent;
162   TH1D *hcent;
163   TH2D *hcentn;
164   TH3D *hphistaretapair10;
165   TH3D *hphistaretapair16;
166   TH3D *hphistaretapair10a;
167   TH3D *hphistaretapair16a;
168   TH3D *hphistaretapair10b;
169   TH3D *hphistaretapair16b;
170   TH3D *hphietapair;
171   TH3D *hphietapair2;
172   TH1D *hpidid;
173   TH1D *hkt;
174   TH1D *hktcheck;
175   TH3D *hkt3;
176   TH2D* hdcaxy;
177   TH1D* hdcaz;
178   TH1D* hsharequal;
179   TH1D* hsharefrac;
180   TH1D* hsharequalmix;
181   TH1D* hsharefracmix;
182   TH1D *hPsiTPC;
183   TH1D *hPsiV0A;
184   TH1D *hPsiV0C;
185   TH1D *hCheckEPA;
186   TH1D *hCheckEPC;
187   TH2D* hCheckEPmix;
188   TH2D *hcentq;
189   TH2D* hMixedDist;
190   TH2D *hQvecV0A;
191   TH2D *hQvecV0C;
192   TH1D *hresV0ATPC;
193   TH1D *hresV0CTPC;
194   TH1D *hresV0AV0C;
195   TH1F* hktbins;
196   TH1F* hcentbins;
197   TH1F* hepbins;
198
199   ClassDef(AliAnalysisTaskFemtoESE, 1);
200 };
201
202
203 class AliFemtoESEBasicParticle : public AliAODTrack
204 {
205  public:
206  AliFemtoESEBasicParticle(Double_t En, Double_t px, Double_t py, Double_t pz, Short_t charge, Double_t phi, Double_t eta)
207    : fE(En), fPx(px), fPy(py), fPz(pz), fCharge(charge), fPhi(phi), fEta(eta), fPsiEP(0.0)
208   {
209   }
210   ~AliFemtoESEBasicParticle() {}
211
212   // kinematics
213   virtual Double_t Px() const { return fPx;}
214   virtual Double_t Py() const { return fPy; }
215   virtual Double_t Pz() const { return fPz; }
216   virtual Double_t Pt() const { return sqrt(fPx*fPx+fPy*fPy); }
217   virtual Double_t P()  const { return sqrt(fPx*fPx+fPy*fPy+fPz*fPz);; }
218   virtual Bool_t   PxPyPz(Double_t[3]) const { AliFatal("Not implemented"); return 0; }
219   
220   virtual Double_t Xv() const { AliFatal("Not implemented"); return 0; }
221   virtual Double_t Yv() const { AliFatal("Not implemented"); return 0; }
222   virtual Double_t Zv() const { AliFatal("Not implemented"); return 0; }
223   virtual Bool_t   XvYvZv(Double_t[3]) const { AliFatal("Not implemented"); return 0; }
224   
225   virtual Double_t OneOverPt()  const { AliFatal("Not implemented"); return 0; }
226   virtual Double_t Phi()        const { return fPhi; }
227   virtual Double_t Theta()      const { AliFatal("Not implemented"); return 0; }
228   
229   virtual Double_t E()          const { return fE; }
230   virtual Double_t M()          const { AliFatal("Not implemented"); return 0; }
231   
232   virtual Double_t Eta()        const { return fEta; }
233   virtual Double_t Y()          const { AliFatal("Not implemented"); return 0; }
234   
235   virtual Short_t Charge()      const { return fCharge; }
236   virtual Int_t   GetLabel()    const { AliFatal("Not implemented"); return 0; }
237   // PID
238   virtual Int_t   PdgCode()     const { AliFatal("Not implemented"); return 0; }      
239   virtual const Double_t *PID() const { AliFatal("Not implemented"); return 0; }
240
241   void SetPsiEP(Double_t psi) {fPsiEP = psi;}
242   Double_t GetPsiEP() {return fPsiEP;}
243
244   /*void SetTPCClusterMap(TBits b){fClusterMap = TBits(b);};
245   void SetTPCSharedMap(TBits b){fSharedMap = TBits(b);};
246   TBits GetTPCClusterMap(){return fClusterMap;};
247   TBits GetTPCSharedMap(){return fSharedMap;};*/
248   
249 private:
250   Double_t fE;
251   Double_t fPx;
252   Double_t fPy;
253   Double_t fPz;
254   Short_t fCharge;
255   Double_t fPhi;
256   Double_t fEta;
257   Double_t fPsiEP;
258   //TBits fClusterMap;
259   //TBits fSharedMap;
260
261   
262   ClassDef( AliFemtoESEBasicParticle, 1);
263 };
264
265 #endif