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