]> 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
25 #include "AliAnalysisTask.h"
26 #include "AliAnalysisTaskSE.h"
27 #include "AliESDpid.h"
28 #include "AliAODPid.h"
29
30 #include "AliAODEvent.h"
31 #include "AliAnalysisManager.h"
32 #include "AliAODInputHandler.h"
33 #include "AliSpectraAODEventCuts.h"
34 #include "AliSpectraAODTrackCuts.h"
35
36 class AliAnalysisTaskFemtoESE : public AliAnalysisTaskSE {
37  public:
38
39   AliAnalysisTaskFemtoESE();
40   AliAnalysisTaskFemtoESE(const char* name);
41   virtual ~AliAnalysisTaskFemtoESE();
42   AliAnalysisTaskFemtoESE(const AliAnalysisTaskFemtoESE &/*obj*/); 
43   AliAnalysisTaskFemtoESE &operator=(const AliAnalysisTaskFemtoESE &/*obj*/);
44
45   virtual void   UserCreateOutputObjects();
46   virtual void   UserExec(Option_t *option);
47   virtual void   Terminate(Option_t *);
48
49   AliHelperPID* GetHelperPID() { return fHelperPID; }
50   void SetHelperPID(AliHelperPID* pid){ fHelperPID = pid; }
51   AliSpectraAODEventCuts* GetEventCuts() {return fEventCuts;}
52   void SetEventCuts(AliSpectraAODEventCuts* cuts) {fEventCuts = cuts;}
53   AliSpectraAODTrackCuts* GetTrackCuts() {return fTrackCuts;}
54   void SetTrackCuts(AliSpectraAODTrackCuts* cuts) {fTrackCuts = cuts;}
55
56   Int_t GetTrackFilterBit(){return fFilterBit;}
57   void SetTrackFilterBit(Int_t bit){fFilterBit=bit;}
58   void SetEventSelectionBit( UInt_t val ) {fSelectBit = val;}
59   void SetMinSepPair(Double_t eta,Double_t phi){fMinSepPairEta = eta; fMinSepPairPhi = phi;}
60   void SetShareFraction(Double_t val) {fShareFraction = val;}
61   Double_t GetShareFraction() {return fShareFraction;}
62   void SetShareQuality(Double_t val) {fShareQuality = val;}
63   Double_t GetShareQuality() {return fShareQuality;}
64   void SetQPercCuts(Double_t min, Double_t max) {fMinQPerc = min; fMaxQPerc = max;}
65   Double_t GetMinQPerc(){return fMinQPerc;}
66   Double_t GetMaxQPerc(){return fMaxQPerc;}
67   void SetQPercDetector(Int_t det){fQPercDet = det;};
68   void SetEPDetector(Int_t det){fEPDet = det;};
69
70
71   void SetKtBins(Int_t n, Double_t* bins);
72   void SetEPBins(Int_t n, Double_t min, Double_t max);
73   void SetCentBins(Int_t n, Double_t* bins);
74   void SetVzBins(Int_t n, Double_t* bins);
75
76  private:
77   Double_t GetQinv(Double_t[], Double_t[]);
78   void GetQosl(Double_t[], Double_t[], Double_t&, Double_t&, Double_t&);
79   Bool_t TrackCut(AliAODTrack* ftrack);
80   Bool_t EventCut(/*AliAODEvent* fevent*/);
81   Bool_t PairCut(AliAODTrack* ftrack1, AliAODTrack* ftrack2, Bool_t mix);
82   Double_t DeltaPhiStar(AliAODTrack* ftrack1, AliAODTrack* ftrack2, Double_t r);
83   TObjArray* CloneAndReduceTrackList(TObjArray* tracks, Double_t psi);
84   Double_t GetDeltaPhiEP(Double_t px1, Double_t py1, Double_t px2, Double_t py2, Double_t psi);
85   Bool_t FindBin(Double_t kt, Double_t phi, Double_t cent, Double_t vz, Int_t& a, Int_t& b, Int_t&c, Int_t& d);
86
87   AliAODEvent            *fAOD; //!    // AOD object
88   TList                  *fOutputList; //! Compact Output list
89   //AliPIDResponse         *fPIDResponse; //! PID response object; equivalent to AliAODpidUtil
90   AliHelperPID*     fHelperPID;      // points to class for PID
91   AliEventPoolManager*     fPoolMgr;         //! event pool manager
92   AliSpectraAODEventCuts* fEventCuts;
93   AliSpectraAODTrackCuts* fTrackCuts;
94
95   Int_t          fFilterBit;         // track selection cuts
96   UInt_t                fSelectBit;            // Select events according to AliAnalysisTaskJetServices bit maps 
97   Int_t fEventCounter;
98   Int_t fMixingTracks;
99   Double_t fBfield;
100   Double_t fMinSepPairEta;
101   Double_t fMinSepPairPhi;
102   Double_t fShareQuality;
103   Double_t fShareFraction;
104
105   Int_t nCountSamePairs;
106   Int_t nCountMixedPairs;
107   Int_t nCountTracks;
108
109   Double_t fMinQPerc;
110   Double_t fMaxQPerc;
111
112   Int_t fQPercDet; // detector used for q-vector (0-V0A, 1-V0C)
113   Int_t fEPDet; // detector used for event plane (0-V0A, 1-V0C)
114
115   Double_t fPsiEPmix;
116   Double_t fPsiEPmixtemp;
117
118   // binning for histograms
119   Int_t nKtBins;
120   Int_t nKtBins1;
121   Double_t* ktBins; //[nKtBins1]
122   Int_t nEPBins;
123   Int_t nEPBins1;
124   Double_t* epBins; //[nEPBins1]
125   Int_t nCentBins;
126   Int_t nCentBins1;
127   Double_t* centBins; //[nCentBins1]
128   Int_t nVzBins;
129   Int_t nVzBins1;
130   Double_t* vzBins; //[nVzBins1]
131
132   Double_t vertex[3];
133
134   TH3F***** hq;
135   TH3F***** hqmix;
136
137   ClassDef(AliAnalysisTaskFemtoESE, 1);
138 };
139
140
141 class AliFemtoESEBasicParticle : public AliAODTrack
142 {
143  public:
144  AliFemtoESEBasicParticle(Double_t En, Double_t px, Double_t py, Double_t pz, Short_t charge, Double_t phi, Double_t eta)
145    : fE(En), fPx(px), fPy(py), fPz(pz), fCharge(charge), fPhi(phi), fEta(eta), fPsiEP(0.0)
146   {
147   }
148   ~AliFemtoESEBasicParticle() {}
149
150   // kinematics
151   virtual Double_t Px() const { return fPx;}
152   virtual Double_t Py() const { return fPy; }
153   virtual Double_t Pz() const { return fPz; }
154   virtual Double_t Pt() const { return sqrt(fPx*fPx+fPy*fPy); }
155   virtual Double_t P()  const { return sqrt(fPx*fPx+fPy*fPy+fPz*fPz);; }
156   virtual Bool_t   PxPyPz(Double_t[3]) const { AliFatal("Not implemented"); return 0; }
157   
158   virtual Double_t Xv() const { AliFatal("Not implemented"); return 0; }
159   virtual Double_t Yv() const { AliFatal("Not implemented"); return 0; }
160   virtual Double_t Zv() const { AliFatal("Not implemented"); return 0; }
161   virtual Bool_t   XvYvZv(Double_t[3]) const { AliFatal("Not implemented"); return 0; }
162   
163   virtual Double_t OneOverPt()  const { AliFatal("Not implemented"); return 0; }
164   virtual Double_t Phi()        const { return fPhi; }
165   virtual Double_t Theta()      const { AliFatal("Not implemented"); return 0; }
166   
167   virtual Double_t E()          const { return fE; }
168   virtual Double_t M()          const { AliFatal("Not implemented"); return 0; }
169   
170   virtual Double_t Eta()        const { return fEta; }
171   virtual Double_t Y()          const { AliFatal("Not implemented"); return 0; }
172   
173   virtual Short_t Charge()      const { return fCharge; }
174   virtual Int_t   GetLabel()    const { AliFatal("Not implemented"); return 0; }
175   // PID
176   virtual Int_t   PdgCode()     const { AliFatal("Not implemented"); return 0; }      
177   virtual const Double_t *PID() const { AliFatal("Not implemented"); return 0; }
178
179   void SetPsiEP(Double_t psi) {fPsiEP = psi;}
180   Double_t GetPsiEP() {return fPsiEP;}
181
182   /*void SetTPCClusterMap(TBits b){fClusterMap = TBits(b);};
183   void SetTPCSharedMap(TBits b){fSharedMap = TBits(b);};
184   TBits GetTPCClusterMap(){return fClusterMap;};
185   TBits GetTPCSharedMap(){return fSharedMap;};*/
186   
187 private:
188   Double_t fE;
189   Double_t fPx;
190   Double_t fPy;
191   Double_t fPz;
192   Short_t fCharge;
193   Double_t fPhi;
194   Double_t fEta;
195   Double_t fPsiEP;
196   //TBits fClusterMap;
197   //TBits fSharedMap;
198
199   
200   ClassDef( AliFemtoESEBasicParticle, 1);
201 };
202
203 #endif