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