]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGCF/FEMTOSCOPY/ESE/AliAnalysisTaskFemtoESE.h
Merge branch 'feature-movesplit'
[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, Int_t z, 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 SetMinQinv(Double_t min){fQinvMin = min;}
68   void SetMaxDCA(Double_t maxxy, Double_t maxz){fMaxDcaXY = maxxy; fMaxDcaZ = maxz;}
69   void SetPtCuts(Double_t min, Double_t max){fPtMin = min; fPtMax = max;}
70   void SetMaxEta(Double_t eta){fEtaMax = eta;}
71   void SetShareFraction(Double_t val) {fShareFraction = val;}
72   Double_t GetShareFraction() {return fShareFraction;}
73   void SetShareQuality(Double_t val) {fShareQuality = val;}
74   Double_t GetShareQuality() {return fShareQuality;}
75   void SetQPercCuts(Double_t min, Double_t max) {fMinQPerc = min; fMaxQPerc = max;}
76   Double_t GetMinQPerc(){return fMinQPerc;}
77   Double_t GetMaxQPerc(){return fMaxQPerc;}
78   void SetQPercDetector(Int_t det){fQPercDet = det;};
79   void SetEPDetector(Int_t det){fEPDet = det;};
80   void SetNMixingTracks(Int_t n){fMixingTracks = n;};
81   void SetQBinning(Int_t n, Double_t q){qbins = n; qlimit = q;};
82
83   void SetKtBins(Int_t n, Double_t* bins);
84   void SetEPBins(Int_t n);
85   void SetCentBins(Int_t n, Double_t* bins);
86   void SetVzBins(Int_t n, Double_t* bins);
87
88   //Double_t GetQPercLHC11h(Double_t qvec);
89
90   Double_t GetCentralityWeight(Double_t cent);
91
92  private:
93   Double_t GetQinv(Double_t[], Double_t[]);
94   void GetQosl(Double_t[], Double_t[], Double_t&, Double_t&, Double_t&);
95   Bool_t TrackCut(AliAODTrack* ftrack);
96   Bool_t EventCut(/*AliAODEvent* fevent*/);
97   Bool_t PairCut(AliFemtoESEBasicParticle* ftrack1, AliFemtoESEBasicParticle* ftrack2, Bool_t mix);
98   Double_t DeltaPhiStar(AliAODTrack* ftrack1, AliAODTrack* ftrack2, Double_t r);
99   TObjArray* CloneAndReduceTrackList(TObjArray* tracks, Double_t psi);
100   Double_t GetDeltaPhiEP(Double_t px1, Double_t py1, Double_t px2, Double_t py2, Double_t psi);
101   Bool_t FindBin(Double_t kt, Double_t phi, Double_t cent, Int_t& a, Int_t& b, Int_t&c);
102
103   AliAODEvent            *fAOD; //!    // AOD object
104   TList                  *fOutputList; //! Compact Output list
105   //AliPIDResponse         *fPIDResponse; //! PID response object; equivalent to AliAODpidUtil
106   AliHelperPID*     fHelperPID;      // points to class for PID
107   AliEventPoolManager**     fPoolMgr;         //![2] event pool manager
108   AliSpectraAODEventCuts* fEventCuts;
109   AliSpectraAODTrackCuts* fTrackCuts;
110
111   Int_t          fFilterBit;         // track selection cuts
112   UInt_t         fSelectBit;            // Select events according to AliAnalysisTaskJetServices bit maps 
113   Bool_t         bIsLHC10h;
114   Int_t fEventCounter;
115   Int_t fMixingTracks;
116   Double_t fBfield;
117   Double_t fMinSepPairEta;
118   Double_t fMinSepPairPhi;
119   Double_t fQinvMin;
120   Double_t fMaxDcaXY;
121   Double_t fMaxDcaZ;
122   Double_t fPtMin;
123   Double_t fPtMax;
124   Double_t fEtaMax;
125   Double_t fShareQuality;
126   Double_t fShareFraction;
127
128   Int_t nCountSamePairs;
129   Int_t nCountMixedPairs;
130   Int_t nCountTracks;
131
132   Double_t fMinQPerc;
133   Double_t fMaxQPerc;
134
135   Double_t qlimit;
136   Int_t qbins;
137
138   Int_t fQPercDet; // detector used for q-vector (0-V0A, 1-V0C)
139   Int_t fEPDet; // detector used for event plane (0-V0A, 1-V0C)
140
141   // binning for histograms
142   Int_t nKtBins;
143   Int_t nKtBins1;
144   Double_t* ktBins; //[nKtBins1]
145   Int_t nEPBins;
146   Int_t nEPBins1;
147   Double_t* epBins; //[nEPBins1]
148   Int_t nEPBinsMix;
149   Int_t nEPBinsMix1;
150   Double_t* epBinsMix; //[nEPBinsMix1]
151   Int_t nCentBins;
152   Int_t nCentBins1;
153   Double_t* centBins; //[nCentBins1]
154   Int_t nVzBins;
155   Int_t nVzBins1;
156   Double_t* vzBins; //[nVzBins1]
157
158   Double_t vertex[3];
159
160   TH3F***** hq;
161   TH3F***** hqmix;
162   TH3F***** hqinv;
163
164   Int_t nqPercBinsLHC11h;
165   Double_t* qPercBinsLHC11h; //[nqPercBinsLHC11h]
166
167   //TStopwatch *stopwatch;
168
169   // histograms
170   TH1D *hpx;
171   TH1D *hpy;
172   TH1D *hpz;
173   TH1D *hpt;
174   TH1D *hE;
175   TH2D *hphieta;
176   TH2D *hphieta_pid;
177   TH1D *hpt_pid;
178   TH2D *hvzcent;
179   TH1D *hcent;
180   TH1D* hcentUnweighted;
181   TH2D *hcentn;
182   TH3D *hphistaretapair10;
183   TH3D *hphistaretapair16;
184   TH3D *hphistaretapair10a;
185   TH3D *hphistaretapair16a;
186   TH3D *hphistaretapair10b;
187   TH3D *hphistaretapair16b;
188   TH3D *hphietapair;
189   TH3D *hphietapair2;
190   TH1D *hpidid;
191   TH2D *hkt;
192   TH1D *hktcheck;
193   TH3D *hkt3;
194   TH2D* hdcaxy;
195   TH1D* hdcaz;
196   TH1D* hsharequal;
197   TH1D* hsharefrac;
198   TH1D* hsharequalmix;
199   TH1D* hsharefracmix;
200   TH2D *hPsiTPC;
201   TH2D *hPsiV0A;
202   TH2D *hPsiV0C;
203   TH2D* hShiftTPC;
204   TH2D* hShiftV0A;
205   TH2D* hShiftV0C;
206   TH2D *hPsiMix;
207   TH2D *hCheckEPA;
208   TH2D *hCheckEPC;
209   TH2D* hCheckEPmix;
210   TH3D* hAvDphi;
211   TH3D* hNpairs;
212   TH1D* hPairDphi;
213   TH1D* hPairDphiMix;
214   TH2D *hcentq;
215   TH2D* hMixedDistTracks;
216   TH2D* hMixedDistEvents;
217   TH2D *hQvecV0A;
218   TH2D *hQvecV0C;
219   TH1D *hresV0ATPC;
220   TH1D *hresV0CTPC;
221   TH1D *hresV0AV0C;
222   TH3F* hqinvcheck;
223
224   TH1F* hktbins;
225   TH1F* hcentbins;
226   TH1F* hepbins;
227
228   ClassDef(AliAnalysisTaskFemtoESE, 1);
229 };
230
231
232 class AliFemtoESEBasicParticle : public AliAODTrack
233 {
234  public:
235  AliFemtoESEBasicParticle(Double_t En, Double_t px, Double_t py, Double_t pz, Short_t charge, Double_t phi, Double_t eta)
236    : fE(En), fPx(px), fPy(py), fPz(pz), fCharge(charge), fPhi(phi), fEta(eta), fPsiEP(0.0)
237   {
238   }
239   ~AliFemtoESEBasicParticle() {}
240
241   // kinematics
242   virtual Double_t Px() const { return fPx;}
243   virtual Double_t Py() const { return fPy; }
244   virtual Double_t Pz() const { return fPz; }
245   virtual Double_t Pt() const { return sqrt(fPx*fPx+fPy*fPy); }
246   virtual Double_t P()  const { return sqrt(fPx*fPx+fPy*fPy+fPz*fPz);; }
247   virtual Bool_t   PxPyPz(Double_t[3]) const { AliFatal("Not implemented"); return 0; }
248   
249   virtual Double_t Xv() const { AliFatal("Not implemented"); return 0; }
250   virtual Double_t Yv() const { AliFatal("Not implemented"); return 0; }
251   virtual Double_t Zv() const { AliFatal("Not implemented"); return 0; }
252   virtual Bool_t   XvYvZv(Double_t[3]) const { AliFatal("Not implemented"); return 0; }
253   
254   virtual Double_t OneOverPt()  const { AliFatal("Not implemented"); return 0; }
255   virtual Double_t Phi()        const { return fPhi; }
256   virtual Double_t Theta()      const { AliFatal("Not implemented"); return 0; }
257   
258   virtual Double_t E()          const { return fE; }
259   virtual Double_t M()          const { AliFatal("Not implemented"); return 0; }
260   
261   virtual Double_t Eta()        const { return fEta; }
262   virtual Double_t Y()          const { AliFatal("Not implemented"); return 0; }
263   
264   virtual Short_t Charge()      const { return fCharge; }
265   virtual Int_t   GetLabel()    const { AliFatal("Not implemented"); return 0; }
266   // PID
267   virtual Int_t   PdgCode()     const { AliFatal("Not implemented"); return 0; }      
268   virtual const Double_t *PID() const { AliFatal("Not implemented"); return 0; }
269
270   void SetPsiEP(Double_t psi) {fPsiEP = psi;}
271   Double_t GetPsiEP() {return fPsiEP;}
272
273   /*void SetTPCClusterMap(TBits b){fClusterMap = TBits(b);};
274   void SetTPCSharedMap(TBits b){fSharedMap = TBits(b);};
275   TBits GetTPCClusterMap(){return fClusterMap;};
276   TBits GetTPCSharedMap(){return fSharedMap;};*/
277   
278 private:
279   Double_t fE;
280   Double_t fPx;
281   Double_t fPy;
282   Double_t fPz;
283   Short_t fCharge;
284   Double_t fPhi;
285   Double_t fEta;
286   Double_t fPsiEP;
287   //TBits fClusterMap;
288   //TBits fSharedMap;
289
290   
291   ClassDef( AliFemtoESEBasicParticle, 1);
292 };
293
294 #endif