]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGJE/UserTasks/AliAnalysisTaskJetProtonCorr.h
Updates from Filip
[u/mrichter/AliRoot.git] / PWGJE / UserTasks / AliAnalysisTaskJetProtonCorr.h
1 #include "TH1.h"
2 #include "TH2.h"
3 #include "TH3.h"
4
5 #include "AliLog.h"
6
7 #include "AliAnalysisTaskSE.h"
8
9 #define ID(x) x, #x
10 #define LAB(x) x + 1, #x
11
12 class TList;
13 class TClonesArray;
14 class AliVParticle;
15 class AliVTrack;
16 class AliPIDResponse;
17 class AliEventPoolManager;
18 class AliEventPool;
19 class AliESDtrackCuts;
20
21 class AliAnalysisTaskJetProtonCorr :
22   public AliAnalysisTaskSE
23 {
24 public:
25   AliAnalysisTaskJetProtonCorr(const char *name = "jets_trg_trd");
26   ~AliAnalysisTaskJetProtonCorr();
27
28   // analysis operations
29   virtual void   UserCreateOutputObjects();
30   virtual Bool_t Notify();
31   virtual void   UserExec(Option_t *option);
32   virtual void   Terminate(const Option_t *option);
33
34   // task configuration
35   void SetJetBranchName(const char* const branchName) { strncpy(fJetBranchName, branchName, fgkStringLength-1); }
36   const char* GetJetBranchName() const { return fJetBranchName; }
37
38   void SetPtThrPart(Float_t minPt, Float_t maxPt) { fTrgPartPtMin = minPt; fTrgPartPtMax = maxPt; }
39   Float_t GetPtMinPart() const { return fTrgPartPtMin; }
40   Float_t GetPtMaxPart() const { return fTrgPartPtMax; }
41   void SetPtThrJet(Float_t minPt, Float_t maxPt) { fTrgJetPtMin = minPt; fTrgJetPtMax = maxPt; }
42   Float_t GetPtMinJet() const { return fTrgJetPtMin; }
43   Float_t GetPtMaxJet() const { return fTrgJetPtMax; }
44   void SetPtThrAss(Float_t minPt, Float_t maxPt) { fAssPartPtMin = minPt; fAssPartPtMax = maxPt; }
45   Float_t GetPtMinAss() const { return fAssPartPtMin; }
46   Float_t GetPtMaxAss() const { return fAssPartPtMax; }
47
48   void PrintTask(Option_t *option, Int_t indent) const;
49
50   static Double_t TOFsignal(Double_t *x, Double_t *par)
51   {
52     Double_t norm = par[0];
53     Double_t mean = par[1];
54     Double_t sigma = par[2];
55     Double_t tail = par[3];
56
57     if (x[0] <= (tail + mean))
58       return norm * TMath::Gaus(x[0], mean, sigma);
59     else
60       return norm * TMath::Gaus(tail + mean, mean, sigma) * TMath::Exp(-tail * (x[0] - tail - mean) / (sigma * sigma));
61   }
62
63   // histograms
64   enum Hist_t {
65       kHistStat = 0,
66       kHistCentrality,
67       kHistCentralityUsed,
68       kHistCentralityCheck,
69       kHistCentralityCheckUsed,
70       kHistSignalTPC,
71       kHistSignalTOF,
72       kHistBetaTOF,
73       kHistDeltaTPC,
74       kHistDeltaTPCSemi,
75       kHistDeltaTOF,
76       kHistDeltaTOFSemi,
77       kHistExpSigmaTOFe,
78       kHistExpSigmaTOFmu,
79       kHistExpSigmaTOFpi,
80       kHistExpSigmaTOFk,
81       kHistExpSigmaTOFp,
82       kHistExpSigmaTOFd,
83       kHistExpSigmaTOFeSemi,
84       kHistExpSigmaTOFmuSemi,
85       kHistExpSigmaTOFpiSemi,
86       kHistExpSigmaTOFkSemi,
87       kHistExpSigmaTOFpSemi,
88       kHistExpSigmaTOFdSemi,
89       kHistCmpSigmaTOFe,
90       kHistCmpSigmaTOFmu,
91       kHistCmpSigmaTOFpi,
92       kHistCmpSigmaTOFk,
93       kHistCmpSigmaTOFp,
94       kHistCmpSigmaTOFd,
95       kHistCmpSigmaTOFeSemi,
96       kHistCmpSigmaTOFmuSemi,
97       kHistCmpSigmaTOFpiSemi,
98       kHistCmpSigmaTOFkSemi,
99       kHistCmpSigmaTOFpSemi,
100       kHistCmpSigmaTOFdSemi,
101       kHistNsigmaTPCe,
102       kHistNsigmaTPCmu,
103       kHistNsigmaTPCpi,
104       kHistNsigmaTPCk,
105       kHistNsigmaTPCp,
106       kHistNsigmaTPCd,
107       kHistNsigmaTPCe_e,
108       kHistNsigmaTOFe,
109       kHistNsigmaTOFmu,
110       kHistNsigmaTOFpi,
111       kHistNsigmaTOFk,
112       kHistNsigmaTOFp,
113       kHistNsigmaTOFd,
114       kHistNsigmaTOFmismatch,
115       kHistNsigmaTOFmismatch2,
116       kHistDeltaTOFe,
117       kHistDeltaTOFmu,
118       kHistDeltaTOFpi,
119       kHistDeltaTOFk,
120       kHistDeltaTOFp,
121       kHistDeltaTOFd,
122       kHistNsigmaTPCeSemi,
123       kHistNsigmaTPCmuSemi,
124       kHistNsigmaTPCpiSemi,
125       kHistNsigmaTPCkSemi,
126       kHistNsigmaTPCpSemi,
127       kHistNsigmaTPCdSemi,
128       kHistNsigmaTPCe_eSemi,
129       kHistNsigmaTOFeSemi,
130       kHistNsigmaTOFmuSemi,
131       kHistNsigmaTOFpiSemi,
132       kHistNsigmaTOFkSemi,
133       kHistNsigmaTOFpSemi,
134       kHistNsigmaTOFdSemi,
135       kHistNsigmaTOFmismatchSemi,
136       kHistNsigmaTOFmismatch2Semi,
137       kHistDeltaTOFeSemi,
138       kHistDeltaTOFmuSemi,
139       kHistDeltaTOFpiSemi,
140       kHistDeltaTOFkSemi,
141       kHistDeltaTOFpSemi,
142       kHistDeltaTOFdSemi,
143       kHistNsigmaTPCTOF,
144       kHistNsigmaTPCTOFPt,
145       kHistNsigmaTPCTOFUsed,
146       kHistNsigmaTPCTOFUsedCentral,
147       kHistNsigmaTPCTOFUsedSemiCentral,
148       kHistNsigmaTPCTOFUsedPt,
149       kHistNsigmaTPCTOFUsedPtCentral,
150       kHistNsigmaTPCTOFUsedPtSemiCentral,
151       kHistEvPlane,
152       kHistEvPlaneUsed,
153       kHistEvPlaneCheck,
154       kHistEvPlaneCheckUsed,
155       kHistEvPlaneCorr,
156       kHistJetPtCentral,
157       kHistJetPtSemi,
158       kHistEtaPhiTrgHad,
159       kHistEtaPhiTrgJet,
160       kHistEtaPhiAssHad,
161       kHistEtaPhiAssProt,
162       kHistLast
163   };
164
165   // statistics
166   enum Stat_t {
167       kStatSeen = 1,
168       kStatTrg,
169       kStatCent,
170       kStatEvPlane,
171       kStatPID,
172       kStatUsed,
173       kStatEvCuts,
174       kStatCentral,
175       kStatSemiCentral,
176       kStatLast
177   };
178
179   // trigger conditions
180   enum Trigger_t {
181       kTriggerMB = 0,
182       kTriggerInt,
183       kTriggerLast
184   };
185
186   // classification
187   enum CorrType_t {
188     kCorrHadHad = 0,
189     kCorrHadProt,
190     kCorrJetHad,
191     kCorrJetProt,
192     kCorrLast
193   };
194
195   enum Class_t {
196     kClCentral = 0,
197     kClSemiCentral,
198     kClDijet,
199     kClLast
200   };
201
202   enum Trg_t {
203     kTrgHad = 0,
204     kTrgJet,
205     kTrgLast
206   };
207
208   enum Ass_t {
209     kAssHad,
210     kAssProt,
211     kAssLast
212   };
213
214   enum Ev_t {
215     kEvSame = 0,
216     kEvMix,
217     kEvLast
218   };
219
220   class AliHistCorr : public TNamed {
221   public:
222     AliHistCorr(TString name, TList *outputList = 0x0);
223     ~AliHistCorr();
224
225     void Trigger(Float_t weight = 1.) { fHistStat->Fill(1., weight); }
226     void Fill(AliVParticle *trgPart, AliVParticle *assPart, Float_t weight = 1.);
227
228   protected:
229     TList *fOutputList;
230
231     TH1F *fHistStat;
232
233     TH1F *fHistCorrPhi;
234     TH2F *fHistCorrPhi2;
235     TH2F *fHistCorrEtaPhi;
236
237     AliHistCorr(const AliHistCorr &rhs);
238     AliHistCorr& operator=(const AliHistCorr &rhs);
239
240     ClassDef(AliHistCorr, 1);
241   };
242
243   AliHistCorr*& GetHistCorr(CorrType_t corr, Class_t cl, Ev_t ev) { return fHistCorr[kEvLast*(kClLast*corr + cl) + ev]; }
244   AliEventPoolManager*& GetPoolMgr(Trg_t trg, Ass_t ass) { return fPoolMgr[kTrgLast * ass + trg]; }
245   AliEventPool*& GetPool(Class_t cls, Trg_t trg, Ass_t ass) { return fPool[kClLast * (kTrgLast * ass + trg) + cls]; }
246
247 protected:
248   AliMCEvent  *fMCEvent; //!
249   AliESDEvent *fESDEvent; //!
250   AliAODEvent *fAODEvent; //!
251
252   UInt_t fTriggerMask;          //! internal representation of trigger conditions
253   UInt_t fClassMask;            //! internal representation of event classes
254   Float_t fCentrality; //!
255   Float_t fCentralityCheck; //!
256   Float_t fZvtx; //!
257   AliPIDResponse *fPIDResponse; //!
258   Float_t fEventPlane; //!
259   Float_t fEventPlaneCheck; //!
260   TObjArray *fPrimTrackArray; //!
261   TClonesArray *fJetArray; //!
262
263   AliEventPoolManager *fPoolMgr[kTrgLast * kAssLast]; //!
264   AliEventPool *fPool[kClLast * kTrgLast * kAssLast]; //!
265
266   AliHistCorr **fHistCorr; //! [kCorrLast*kEvLast*kClLast]; //!
267
268   Bool_t DetectTriggers();
269   void   MarkTrigger(Trigger_t trg) { fTriggerMask |= (1 << trg); }
270   Bool_t IsTrigger(Trigger_t trg) const { return (fTriggerMask & (1 << trg)); }
271
272   Bool_t DetectClasses();
273   void   MarkClass(Class_t cl) { fClassMask |= (1 << cl); }
274   Bool_t IsClass(Class_t cl) const { return (fClassMask & (1 << cl)); }
275
276   Bool_t PrepareEvent();
277   Bool_t CleanUpEvent();
278
279   Float_t GetCentrality() const { return fCentrality; }
280   Float_t GetEventPlane() const { return fEventPlane; }
281   AliPIDResponse* GetPID() const { return fPIDResponse; }
282   Bool_t IsCentral() { return ((fCentrality >= 0.) && (fCentrality <= 10.)); }
283   Bool_t IsSemiCentral() { return ((fCentrality >= 30.) && (fCentrality <= 50.)); }
284
285   AliVTrack* GetLeadingTrack(AliAODJet *jet) const;
286
287   Float_t GetDPhiStar(Float_t phi1, Float_t pt1, Float_t charge1,
288                       Float_t phi2, Float_t pt2, Float_t charge2,
289                       Float_t radius, Float_t bSign);
290
291   Bool_t AcceptTrigger(AliVTrack *trg);
292   Bool_t AcceptTrigger(AliAODJet *trg);
293   Bool_t AcceptAssoc(AliVTrack *trk);
294   Bool_t IsProton(AliVTrack *trk);
295   Bool_t AcceptAngleToEvPlane(Float_t phi, Float_t psi);
296   Bool_t AcceptTwoTracks(AliVParticle *trgPart, AliVParticle *assPart);
297
298   TObjArray* CloneTracks(TObjArray *tracks) const;
299
300   Bool_t Correlate(CorrType_t corr, Class_t cl, Ev_t ev,
301                    TCollection *trgArray, TCollection *assArray, Float_t weight = 1.);
302
303   Bool_t Correlate(Trg_t trg, Ass_t ass, Class_t cl, Ev_t ev,
304                    TCollection *trgArray, TCollection *assArray, Float_t weight = 1.);
305
306   // output objects
307   TList *fOutputList;           //! list of output objects
308
309   // histogram management
310   TH1  *fHist[kHistLast];       //! pointers to histogram
311   const char *fShortTaskId;     //! short identifier for the task
312
313   TH1*&  GetHistogram(Hist_t hist, const Int_t idx = 0) { return fHist[hist + idx]; }
314
315   TH1*   AddHistogram(Hist_t hist, const char *hid, TString title,
316                       Int_t xbins, Float_t xmin, Float_t xmax, Int_t binType = 1);
317   TH2*   AddHistogram(Hist_t hist, const char *hid, TString title,
318                       Int_t xbins, Float_t xmin, Float_t xmax,
319                       Int_t ybins, Float_t ymin, Float_t ymax, Int_t binType = 1);
320   TH3*   AddHistogram(Hist_t hist, const char *hid, TString title,
321                       Int_t xbins, Float_t xmin, Float_t xmax,
322                       Int_t ybins, Float_t ymin, Float_t ymax,
323                       Int_t zbins, Float_t zmin, Float_t zmax, Int_t binType = 1);
324
325   void    FillH1(Hist_t hist, Float_t x, Float_t weight = 1., Int_t idx = 0)
326   { GetHistogram(hist, idx)->Fill(x, weight); }
327   void    FillH2(Hist_t hist, Float_t x, Float_t y, Float_t weight = 1., Int_t idx = 0)
328   { ((TH2*) GetHistogram(hist, idx))->Fill(x, y, weight); }
329   void    FillH3(Hist_t hist, Float_t x, Float_t y, Float_t z, Float_t weight = 1., Int_t idx = 0)
330   { ((TH3*) GetHistogram(hist, idx))->Fill(x, y, z, weight); }
331
332   const char* fkCorrTypeName[kCorrLast]; //!
333   const char* fkClassName[kClLast]; //!
334   const char* fkEvName[kEvLast]; //!
335
336   // task configuration
337   static const Int_t fgkStringLength = 100; // max length for the jet branch name
338   char fJetBranchName[fgkStringLength];     // jet branch name
339
340   Bool_t fUseStandardCuts;
341
342   AliESDtrackCuts *fCutsPrim;   // track cuts for primary particles
343   Float_t fCutsTwoTrackEff;
344
345   Float_t fTrgPartPtMin;
346   Float_t fTrgPartPtMax;
347   Float_t fTrgJetPtMin;
348   Float_t fTrgJetPtMax;
349   Float_t fTrgJetLeadTrkPtMin;
350   Float_t fAssPartPtMin;
351   Float_t fAssPartPtMax;
352   Float_t fTrgAngleToEvPlane;
353
354   // not implemented
355   AliAnalysisTaskJetProtonCorr(const AliAnalysisTaskJetProtonCorr &rhs);
356   AliAnalysisTaskJetProtonCorr& operator=(const AliAnalysisTaskJetProtonCorr &rhs);
357
358   ClassDef(AliAnalysisTaskJetProtonCorr, 1);
359 };