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