]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGJE/UserTasks/AliAnalysisTaskJetProtonCorr.h
add correation with V0mult to Rho task
[u/mrichter/AliRoot.git] / PWGJE / UserTasks / AliAnalysisTaskJetProtonCorr.h
1 #include "TH1.h"
2 #include "TH2.h"
3 #include "TH3.h"
4 #include "TF1.h"
5
6 #include "AliLog.h"
7
8 #include "AliVParticle.h"
9 #include "AliAnalysisTaskSE.h"
10 #include "AliESDtrackCuts.h"
11
12 #define ID(x) x, #x
13 #define LAB(x) x + 1, #x
14
15 class TList;
16 class TClonesArray;
17 class AliOADBContainer;
18 class AliTOFPIDParams;
19 class AliVTrack;
20 class AliPIDResponse;
21 class AliEventPoolManager;
22 class AliEventPool;
23 class AliEventplane;
24 class TSpline;
25
26 class AliAnalysisTaskJetProtonCorr :
27   public AliAnalysisTaskSE
28 {
29 public:
30   AliAnalysisTaskJetProtonCorr(const char *name = "jets_trg_trd");
31   ~AliAnalysisTaskJetProtonCorr();
32
33   // analysis operations
34   virtual void   UserCreateOutputObjects();
35   virtual Bool_t Notify();
36   virtual void   UserExec(Option_t *option);
37   virtual void   Terminate(const Option_t *option);
38
39   void SetParamsTOF();
40
41   // task configuration
42   void SetJetBranchName(const char* branchName) { strncpy(fJetBranchName, branchName, fgkStringLength-1); }
43   const char* GetJetBranchName() const { return fJetBranchName; }
44
45   void SetPtThrPart(Float_t minPt, Float_t maxPt) { fTrgPartPtMin = minPt; fTrgPartPtMax = maxPt; }
46   Float_t GetPtMinPart() const { return fTrgPartPtMin; }
47   Float_t GetPtMaxPart() const { return fTrgPartPtMax; }
48   void SetPtThrJet(Float_t minPt, Float_t maxPt) { fTrgJetPtMin = minPt; fTrgJetPtMax = maxPt; }
49   Float_t GetPtMinJet() const { return fTrgJetPtMin; }
50   Float_t GetPtMaxJet() const { return fTrgJetPtMax; }
51   void SetPtThrAss(Float_t minPt, Float_t maxPt) { fAssPartPtMin = minPt; fAssPartPtMax = maxPt; }
52   Float_t GetPtMinAss() const { return fAssPartPtMin; }
53   Float_t GetPtMaxAss() const { return fAssPartPtMax; }
54
55   void SetTwoTrackCut(Float_t cut) { fCutsTwoTrackEff = cut; }
56   Float_t GetTwoTrackCut() const { return fCutsTwoTrackEff; }
57
58   void SetTrackCutsAss(const AliESDtrackCuts &cuts) { *fCutsPrimAss = cuts; }
59   void SetTrackCutsTrg(const AliESDtrackCuts &cuts) { *fCutsPrimTrg = cuts; }
60   void SetTrackCutsTrgConstrain(const AliESDtrackCuts &cuts) { *fCutsPrimTrgConstrain = cuts; }
61
62   void SetTrgJetEtaMax(Float_t etamax) { fTrgJetEtaMax = etamax; }
63   Float_t GetTrgJetEtaMax() const { return fTrgJetEtaMax; }
64   void SetHadEtaMax(Float_t etamax) { fHadEtaMax = etamax; }
65   Float_t GetHadEtaMax() const { return fHadEtaMax; }
66
67   void SetUseEvplaneV0(Bool_t usev0 = kTRUE) { fUseEvplaneV0 = usev0; }
68   Bool_t GetUseEvplaneV0() const { return fUseEvplaneV0; }
69
70   void SetJetV2(Float_t v2Cent, Float_t v2Semi) {
71     fTrgJetPhiModCent->SetParameter(0, v2Cent);
72     fTrgJetPhiModSemi->SetParameter(0, v2Semi);
73   }
74   void SetHadV2(Float_t v2Cent, Float_t v2Semi) {
75     fTrgHadPhiModCent->SetParameter(0, v2Cent);
76     fTrgHadPhiModSemi->SetParameter(0, v2Semi);
77   }
78
79   void SetLeadingTrackPtMin(Float_t pt) { fTrgJetLeadTrkPtMin = pt; }
80   void SetLeadingTrackPtMax(Float_t pt) { fTrgJetLeadTrkPtMax = pt; }
81
82   void SetJetAreaMin(Float_t area) { fTrgJetAreaMin = area; }
83   Float_t GetJetAreaMin() const { return fTrgJetAreaMin; }
84
85   void SetRequirePID(Bool_t req = kTRUE) { fRequirePID = req; }
86   Bool_t GetRequirePID() const { return fRequirePID; }
87
88   void SetFilterMask(Int_t mask) { fAssFilterMask = mask; }
89   Int_t GetFilterMask() const { return fAssFilterMask; }
90
91   void SetErrorCount(Int_t cnt) { fErrorMsg = cnt; }
92   Int_t GetErrorCount() const { return fErrorMsg; }
93
94   void SetTrgAngleToEvPlane(Float_t angle) { fTrgAngleToEvPlane = angle; }
95   Float_t GetTrgAngleToEvPlane() const { return fTrgAngleToEvPlane; }
96
97   void SetToyMeanNoPart(Float_t mean) { fToyMeanNoPart = mean; }
98   Float_t GetToyMeanNoPart() const { return fToyMeanNoPart; }
99   void SetToyRadius(Float_t radius) { fToyRadius = radius; }
100   Float_t GetToyRadius() const { return fToyRadius; }
101   void SetToySmearPhi(Float_t sigma) { fToySmearPhi = sigma; }
102   Float_t GetToySmearPhi() const { return fToySmearPhi; }
103
104   void SetEventPlaneResSpline(TSpline *spline) { fSplineEventPlaneRes = spline; }
105   const TSpline *GetEventPlaneResSpline() const { return fSplineEventPlaneRes; }
106
107   void PrintTask(Option_t *option, Int_t indent) const;
108
109   static Double_t TOFsignal(Double_t *x, Double_t *par)
110   {
111     Double_t norm = par[0];
112     Double_t mean = par[1];
113     Double_t sigma = par[2];
114     Double_t tail = par[3];
115
116     if (x[0] <= (tail + mean))
117       return norm * TMath::Gaus(x[0], mean, sigma);
118     else
119       return norm * TMath::Gaus(tail + mean, mean, sigma) * TMath::Exp(-tail * (x[0] - tail - mean) / (sigma * sigma));
120   }
121
122   // histograms
123   enum Hist_t {
124       kHistStat = 0,
125       kHistVertexNctb,
126       kHistVertexZ,
127       kHistCentrality,
128       kHistCentralityUsed,
129       kHistCentralityCheck,
130       kHistCentralityCheckUsed,
131       kHistCentralityVsMult,
132       kHistSignalTPC,
133       kHistSignalTOF,
134       kHistBetaTOF,
135       kHistDeltaTPC,
136       kHistDeltaTPCSemi,
137       kHistDeltaTOF,
138       kHistDeltaTOFSemi,
139       // kHistExpSigmaTOFe,
140       // kHistExpSigmaTOFmu,
141       // kHistExpSigmaTOFpi,
142       // kHistExpSigmaTOFk,
143       // kHistExpSigmaTOFp,
144       // kHistExpSigmaTOFd,
145       // kHistExpSigmaTOFeSemi,
146       // kHistExpSigmaTOFmuSemi,
147       // kHistExpSigmaTOFpiSemi,
148       // kHistExpSigmaTOFkSemi,
149       // kHistExpSigmaTOFpSemi,
150       // kHistExpSigmaTOFdSemi,
151       // kHistCmpSigmaTOFe,
152       // kHistCmpSigmaTOFmu,
153       // kHistCmpSigmaTOFpi,
154       // kHistCmpSigmaTOFk,
155       // kHistCmpSigmaTOFp,
156       // kHistCmpSigmaTOFd,
157       // kHistCmpSigmaTOFeSemi,
158       // kHistCmpSigmaTOFmuSemi,
159       // kHistCmpSigmaTOFpiSemi,
160       // kHistCmpSigmaTOFkSemi,
161       // kHistCmpSigmaTOFpSemi,
162       // kHistCmpSigmaTOFdSemi,
163       kHistNsigmaTPCe,
164       kHistNsigmaTPCmu,
165       kHistNsigmaTPCpi,
166       kHistNsigmaTPCk,
167       kHistNsigmaTPCp,
168       kHistNsigmaTPCd,
169       kHistNsigmaTPCe_e,
170       kHistNsigmaTOFe,
171       kHistNsigmaTOFmu,
172       kHistNsigmaTOFpi,
173       kHistNsigmaTOFk,
174       kHistNsigmaTOFp,
175       kHistNsigmaTOFd,
176       kHistNsigmaTOFmismatch,
177       kHistDeltaTOFe,
178       kHistDeltaTOFmu,
179       kHistDeltaTOFpi,
180       kHistDeltaTOFk,
181       kHistDeltaTOFp,
182       kHistDeltaTOFd,
183
184       kHistNsigmaTPCeSemi,
185       kHistNsigmaTPCmuSemi,
186       kHistNsigmaTPCpiSemi,
187       kHistNsigmaTPCkSemi,
188       kHistNsigmaTPCpSemi,
189       kHistNsigmaTPCdSemi,
190       kHistNsigmaTPCe_eSemi,
191       kHistNsigmaTOFeSemi,
192       kHistNsigmaTOFmuSemi,
193       kHistNsigmaTOFpiSemi,
194       kHistNsigmaTOFkSemi,
195       kHistNsigmaTOFpSemi,
196       kHistNsigmaTOFdSemi,
197       kHistNsigmaTOFmismatchSemi,
198       kHistDeltaTOFeSemi,
199       kHistDeltaTOFmuSemi,
200       kHistDeltaTOFpiSemi,
201       kHistDeltaTOFkSemi,
202       kHistDeltaTOFpSemi,
203       kHistDeltaTOFdSemi,
204
205       kHistNsigmaTPCTOF,
206       kHistNsigmaTPCTOFPt,
207       kHistNsigmaTPCTOFUsed,
208       kHistNsigmaTPCTOFUsedCentral,
209       kHistNsigmaTPCTOFUsedSemiCentral,
210       kHistNsigmaTPCTOFUsedPt,
211       kHistNsigmaTPCTOFUsedPtCentral,
212       kHistNsigmaTPCTOFUsedPtSemiCentral,
213
214       kHistNsigmaTPCTOFUsedCentralMCe,
215       kHistNsigmaTPCTOFUsedCentralMCmu,
216       kHistNsigmaTPCTOFUsedCentralMCpi,
217       kHistNsigmaTPCTOFUsedCentralMCk,
218       kHistNsigmaTPCTOFUsedCentralMCp,
219       kHistNsigmaTPCTOFUsedCentralMCd,
220
221       kHistNsigmaTPCTOFUsedSemiCentralMCe,
222       kHistNsigmaTPCTOFUsedSemiCentralMCmu,
223       kHistNsigmaTPCTOFUsedSemiCentralMCpi,
224       kHistNsigmaTPCTOFUsedSemiCentralMCk,
225       kHistNsigmaTPCTOFUsedSemiCentralMCp,
226       kHistNsigmaTPCTOFUsedSemiCentralMCd,
227
228       kHistNevMix,
229
230       kHistEvPlane,
231       kHistEvPlaneRes,
232       kHistEvPlaneUsed,
233       kHistEvPlaneCheck,
234       kHistEvPlaneCheckUsed,
235       kHistEvPlane3,
236       kHistEvPlaneCorr,
237       kHistEvPlaneCross,
238       kHistEvPlaneCorrNoTrgJets,
239       kHistEvPlaneCorrNoTrgJetsTrgd,
240       kHistJetPtCentral,
241       kHistJetPtSemi,
242       kHistEtaPhiTrgHad,
243       kHistEtaPhiTrgJet,
244       kHistEtaPhiAssHad,
245       kHistEtaPhiAssProt,
246       kHistPhiTrgJetEvPlane,
247       kHistPhiTrgHadEvPlane,
248       kHistPhiRndTrgJetEvPlane,
249       kHistPhiRndTrgHadEvPlane,
250       kHistPhiAssHadEvPlane,
251       kHistPhiAssHadVsEvPlane,
252       kHistPhiAssProtEvPlane,
253       kHistPhiTrgJetEvPlane3,
254       kHistPhiTrgHadEvPlane3,
255       kHistPhiAssHadEvPlane3,
256       kHistPhiAssProtEvPlane3,
257       kHistLast
258   };
259
260   // statistics
261   enum Stat_t {
262       kStatSeen = 1,
263       kStatTrg,
264       kStatVtx,
265       kStatCent,
266       kStatEvPlane,
267       kStatPID,
268       kStatUsed,
269       kStatEvCuts,
270       kStatCentral,
271       kStatSemiCentral,
272       kStatLast
273   };
274
275   // trigger conditions
276   enum Trigger_t {
277       kTriggerMB = 0,
278       kTriggerInt,
279       kTriggerLast
280   };
281
282   // classification
283   enum CorrType_t {
284     kCorrHadHad = 0,
285     kCorrHadProt,
286     kCorrJetHad,
287     kCorrJetProt,
288
289     kCorrRndHadHad,
290     kCorrRndHadProt,
291     kCorrRndJetHad,
292     kCorrRndJetProt,
293
294     kCorrRndHadExcHad,
295     kCorrRndHadExcProt,
296     kCorrRndJetExcHad,
297     kCorrRndJetExcProt,
298
299     kCorrLast
300   };
301
302   enum Class_t {
303     kClCentral = 0,
304     kClSemiCentral,
305     // kClDijet,
306     kClLast
307   };
308
309   enum Trg_t {
310     kTrgHad = 0,
311     kTrgJet,
312     kTrgHadRnd,
313     kTrgJetRnd,
314     kTrgLast
315   };
316
317   enum Ass_t {
318     kAssHad = 0,
319     kAssProt,
320     kAssHadJetExc,
321     kAssProtJetExc,
322     kAssHadHadExc,
323     kAssProtHadExc,
324     kAssLast
325   };
326
327   enum Ev_t {
328     kEvSame = 0,
329     kEvMix,
330     kEvLast
331   };
332
333   class AliPartCorr : public AliVParticle {
334   public:
335     AliPartCorr(Float_t eta = 0., Float_t phi = 0., Float_t pt = 0., Float_t charge = 0) :
336       fPt(pt), fEta(eta), fPhi(phi), fCharge(charge) {}
337     AliPartCorr(const AliVParticle &rhs) :
338       fPt(rhs.Pt()), fEta(rhs.Eta()), fPhi(rhs.Phi()), fCharge(rhs.Charge()) {}
339     virtual ~AliPartCorr() {}
340     
341     // kinematics
342     virtual Double_t Px() const { AliFatal("not implemented"); return 0; }
343     virtual Double_t Py() const { AliFatal("not implemented"); return 0; }
344     virtual Double_t Pz() const { AliFatal("not implemented"); return 0; }
345     virtual Double_t Pt() const { return fPt; }
346     virtual Double_t P()  const { AliFatal("not implemented"); return 0; }
347     virtual Bool_t   PxPyPz(Double_t[3]) const { AliFatal("not implemented"); return 0; }
348
349     virtual Double_t Xv() const { AliFatal("not implemented"); return 0; }
350     virtual Double_t Yv() const { AliFatal("not implemented"); return 0; }
351     virtual Double_t Zv() const { AliFatal("not implemented"); return 0; }
352     virtual Bool_t   XvYvZv(Double_t[3]) const { AliFatal("not implemented"); return 0; }
353
354     virtual Double_t OneOverPt()  const { AliFatal("not implemented"); return 0; }
355     virtual Double_t Phi()        const { return fPhi; }
356     virtual Double_t Theta()      const { AliFatal("not implemented"); return 0; }
357
358     virtual Double_t E()          const { AliFatal("not implemented"); return 0; }
359     virtual Double_t M()          const { AliFatal("not implemented"); return 0; }
360
361     virtual Double_t Eta()        const { return fEta; }
362     virtual Double_t Y()          const { AliFatal("not implemented"); return 0; }
363
364     virtual Short_t Charge()      const { return fCharge; }
365     virtual Int_t   GetLabel()    const { AliFatal("not implemented"); return 0; }
366
367     virtual Int_t   PdgCode()     const { AliFatal("not implemented"); return 0; }
368     virtual const Double_t *PID() const { AliFatal("not implemented"); return 0; }
369
370   protected:
371     Float_t fPt;
372     Float_t fEta;
373     Float_t fPhi;
374     Short_t fCharge;
375
376     ClassDef(AliPartCorr, 1);
377   };
378
379   class AliHistCorr : public TNamed {
380   public:
381     AliHistCorr(TString name, TList *outputList = 0x0);
382     ~AliHistCorr();
383
384     void Trigger(Float_t phi, Float_t eta, Float_t qpt, Float_t weight) {
385       fHistStat->Fill(1., weight);
386       if (fHistCorrTrgEtaPhi)
387         fHistCorrTrgEtaPhi->Fill(phi, eta, weight);
388       if (fHistCorrTrgEtaPhiQpt)
389         fHistCorrTrgEtaPhiQpt->Fill(phi, eta, qpt, weight);
390     }
391     void Ass(Float_t phi, Float_t eta, Float_t qpt, Float_t weight) {
392       if (fHistCorrAssEtaPhi)
393         fHistCorrAssEtaPhi->Fill(phi, eta, weight);
394       if (fHistCorrAssEtaPhiQpt)
395         fHistCorrAssEtaPhiQpt->Fill(phi, eta, qpt, weight);
396     }
397     void Fill(AliVParticle *trgPart, AliVParticle *assPart, Float_t weight = 1.);
398     void Fill(TLorentzVector *trgPart, AliVParticle *assPart, Float_t weight = 1.);
399     void Fill(TLorentzVector *trgPart, TLorentzVector *assPart, Float_t weight = 1.);
400
401   protected:
402     TList *fOutputList;
403
404     TH1F *fHistStat;
405
406     TH1F *fHistCorrPhi;
407     TH2F *fHistCorrPhi2;
408     TH2F *fHistCorrEtaPhi;
409     TH2F *fHistCorrAvgEtaPhi;
410     TH2F *fHistCorrTrgEtaPhi;
411     TH2F *fHistCorrAssEtaPhi;
412     TH3F *fHistCorrTrgEtaPhiQpt;
413     TH3F *fHistCorrAssEtaPhiQpt;
414
415     const Float_t fHistDphiLo;
416     const Int_t   fHistDphiNbins;
417     const Int_t   fHistDetaNbins;
418
419     AliHistCorr(const AliHistCorr &rhs);
420     AliHistCorr& operator=(const AliHistCorr &rhs);
421
422     ClassDef(AliHistCorr, 1);
423   };
424
425   AliHistCorr*& GetHistCorr(CorrType_t corr, Class_t cl, Ev_t ev) { return fHistCorr[kEvLast*(kClLast*corr + cl) + ev]; }
426   AliEventPoolManager*& GetPoolMgr(Ass_t ass) { return fPoolMgr[ass]; }
427   AliEventPool*& GetPool(Ass_t ass) { return fPool[ass]; }
428
429 protected:
430   AliMCEvent  *fMCEvent; //!
431   AliESDEvent *fESDEvent; //!
432   AliAODEvent *fAODEvent; //!
433
434   Int_t fRunNumber; //! current run number
435   AliOADBContainer *fOADBContainerTOF; //! container for OADB entry with TOF parameters
436   AliTOFPIDParams *fParamsTOF; //! TOF parametrization
437
438   AliEventplane *fEventplane; //! pointer to the event plane
439
440   UInt_t fTriggerMask;          //! internal representation of trigger conditions
441   UInt_t fClassMask;            //! internal representation of event classes
442   Float_t fCentrality; //!
443   Float_t fCentralityCheck; //!
444   Float_t fZvtx; //!
445   AliPIDResponse *fPIDResponse; //!
446   Float_t fEventPlaneAngle; //!
447   Float_t fEventPlaneRes; //!
448   Float_t fEventPlaneAngleCheck; //!
449   Float_t fEventPlaneAngle3; //!
450
451   TObjArray *fPrimTrackArrayAss; //!
452   TObjArray *fPrimTrackArrayTrg; //!
453   TClonesArray *fPrimConstrainedTrackArray; //!
454   TClonesArray *fJetArray; //!
455
456   AliEventPoolManager *fPoolMgr[kAssProt + 1]; //!
457   AliEventPool *fPool[kAssProt + 1]; //!
458
459   AliHistCorr **fHistCorr; //! [kCorrLast*kEvLast*kClLast]; //!
460
461   Int_t fErrorMsg; //! remaining number of error messages to be printed
462
463   Bool_t DetectTriggers();
464   void   MarkTrigger(Trigger_t trg) { fTriggerMask |= (1 << trg); }
465   Bool_t IsTrigger(Trigger_t trg) const { return (fTriggerMask & (1 << trg)); }
466
467   Bool_t DetectClasses();
468   void   MarkClass(Class_t cl) { fClassMask |= (1 << cl); }
469   Bool_t IsClass(Class_t cl) const { return (fClassMask & (1 << cl)); }
470
471   Bool_t PrepareEvent();
472   Bool_t CleanUpEvent();
473
474   Float_t GetCentrality() const { return fCentrality; }
475   Float_t GetEventPlaneAngle() const { return fEventPlaneAngle; }
476   AliPIDResponse* GetPID() const { return fPIDResponse; }
477   Bool_t IsCentral() const { return ((fCentrality >= 0.) && (fCentrality <= 10.)); }
478   Bool_t IsSemiCentral() const { return ((fCentrality >= 30.) && (fCentrality <= 50.)); }
479
480   AliVTrack* GetLeadingTrack(const AliAODJet *jet) const;
481
482   Float_t GetDPhiStar(Float_t phi1, Float_t pt1, Float_t charge1,
483                       Float_t phi2, Float_t pt2, Float_t charge2,
484                       Float_t radius, Float_t bSign) const;
485
486   Bool_t AcceptTrigger(AliVTrack *trg);
487   Bool_t AcceptTrigger(AliAODJet *trg);
488   Bool_t AcceptAssoc(const AliVTrack *trk) const;
489   Bool_t IsProton(const AliVTrack *trk) const;
490   Bool_t AcceptAngleToEvPlane(Float_t phi, Float_t psi) const;
491   Bool_t AcceptTwoTracks(const AliVParticle *trgPart, const AliVParticle *assPart) const;
492
493   Float_t GetPhiRel2(AliVParticle *part) const;
494
495   TObjArray* CloneTracks(TObjArray *tracks) const;
496
497   Bool_t GenerateRandom(TCollection *trgJetArray, TCollection *trgHadArray,
498                         TCollection *assHadJetArray, TCollection *assProtJetArray,
499                         TCollection *assHadHadArray, TCollection *assProtHadArray,
500                         Float_t pFraction = .5);
501
502   Bool_t Correlate(CorrType_t corr, Class_t cl, Ev_t ev,
503                    TCollection *trgArray, TCollection *assArray, Float_t weight = 1.);
504
505   Bool_t Correlate(Trg_t trg, Ass_t ass, Class_t cl, Ev_t ev,
506                    TCollection *trgArray, TCollection *assArray, Float_t weight = 1.);
507
508   // output objects
509   TList *fOutputList;           //! list of output objects
510
511   // histogram management
512   TH1  *fHist[kHistLast];       //! pointers to histogram
513   const char *fShortTaskId;     //! short identifier for the task
514
515   TH1*&  GetHistogram(Hist_t hist, const Int_t idx = 0) { return fHist[hist + idx]; }
516
517   TH1*   AddHistogram(Hist_t hist, const char *hid, TString title,
518                       Int_t xbins, Float_t xmin, Float_t xmax, Int_t binType = 1);
519   TH2*   AddHistogram(Hist_t hist, const char *hid, TString title,
520                       Int_t xbins, Float_t xmin, Float_t xmax,
521                       Int_t ybins, Float_t ymin, Float_t ymax, Int_t binType = 1);
522   TH3*   AddHistogram(Hist_t hist, const char *hid, TString title,
523                       Int_t xbins, Float_t xmin, Float_t xmax,
524                       Int_t ybins, Float_t ymin, Float_t ymax,
525                       Int_t zbins, Float_t zmin, Float_t zmax, Int_t binType = 1);
526
527   void    FillH1(Hist_t hist, Float_t x, Float_t weight = 1., Int_t idx = 0)
528   { GetHistogram(hist, idx)->Fill(x, weight); }
529   void    FillH2(Hist_t hist, Float_t x, Float_t y, Float_t weight = 1., Int_t idx = 0)
530   { ((TH2*) GetHistogram(hist, idx))->Fill(x, y, weight); }
531   void    FillH3(Hist_t hist, Float_t x, Float_t y, Float_t z, Float_t weight = 1., Int_t idx = 0)
532   { ((TH3*) GetHistogram(hist, idx))->Fill(x, y, z, weight); }
533
534   const char* fkCorrTypeName[kCorrLast]; //!
535   const char* fkClassName[kClLast]; //!
536   const char* fkEvName[kEvLast]; //!
537
538   // task configuration
539   static const Int_t fgkStringLength = 100; // max length for the jet branch name
540   char fJetBranchName[fgkStringLength];     // jet branch name
541
542   const Bool_t fUseStandardCuts;
543   Bool_t fUseEvplaneV0;
544
545   AliESDtrackCuts *fCutsPrimTrg;        // track cuts for primary particles (trigger)
546   AliESDtrackCuts *fCutsPrimTrgConstrain;       // track cuts for primary particles (trigger)
547   AliESDtrackCuts *fCutsPrimAss;        // track cuts for primary particles (associate)
548   Float_t fCutsTwoTrackEff;
549
550   UInt_t  fAssFilterMask;
551   Bool_t  fRequirePID;
552   Float_t fTrgJetEtaMax;
553   Float_t fHadEtaMax;
554   Float_t fTrgPartPtMin;
555   Float_t fTrgPartPtMax;
556   Float_t fTrgJetPtMin;
557   Float_t fTrgJetPtMax;
558   Float_t fTrgJetLeadTrkPtMin;
559   Float_t fTrgJetLeadTrkPtMax;
560   Float_t fTrgJetAreaMin;
561   Float_t fAssPartPtMin;
562   Float_t fAssPartPtMax;
563   Float_t fTrgAngleToEvPlane;
564
565   Float_t fToyMeanNoPart;
566   Float_t fToyRadius;
567   Float_t fToySmearPhi;
568
569   TF1 *fTrgJetPhiModCent;
570   TF1 *fTrgJetPhiModSemi;
571   TF1 *fTrgHadPhiModCent;
572   TF1 *fTrgHadPhiModSemi;
573
574   Float_t fTrgJetV2Cent;
575   Float_t fTrgJetV2Semi;
576   Float_t fTrgHadV2Cent;
577   Float_t fTrgHadV2Semi;
578
579   TSpline *fSplineEventPlaneRes;
580
581   // not implemented
582   AliAnalysisTaskJetProtonCorr(const AliAnalysisTaskJetProtonCorr &rhs);
583   AliAnalysisTaskJetProtonCorr& operator=(const AliAnalysisTaskJetProtonCorr &rhs);
584
585   ClassDef(AliAnalysisTaskJetProtonCorr, 1);
586 };