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