]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGJE/UserTasks/AliAnalysisTaskJetProtonCorr.h
- Don't change binning of multiplicity, if not stored (in order to maintain jet code...
[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       kHistNsigmaTPCTOFUsedCentralMCe,
211       kHistNsigmaTPCTOFUsedCentralMCmu,
212       kHistNsigmaTPCTOFUsedCentralMCpi,
213       kHistNsigmaTPCTOFUsedCentralMCk,
214       kHistNsigmaTPCTOFUsedCentralMCp,
215       kHistNsigmaTPCTOFUsedCentralMCd,
216
217       kHistNsigmaTPCTOFUsedSemiCentralMCe,
218       kHistNsigmaTPCTOFUsedSemiCentralMCmu,
219       kHistNsigmaTPCTOFUsedSemiCentralMCpi,
220       kHistNsigmaTPCTOFUsedSemiCentralMCk,
221       kHistNsigmaTPCTOFUsedSemiCentralMCp,
222       kHistNsigmaTPCTOFUsedSemiCentralMCd,
223
224       kHistEvPlane,
225       kHistEvPlaneUsed,
226       kHistEvPlaneCheck,
227       kHistEvPlaneCheckUsed,
228       kHistEvPlane3,
229       kHistEvPlaneCorr,
230       kHistEvPlaneCorrNoTrgJets,
231       kHistEvPlaneCorrNoTrgJetsTrgd,
232       kHistJetPtCentral,
233       kHistJetPtSemi,
234       kHistEtaPhiTrgHad,
235       kHistEtaPhiTrgJet,
236       kHistEtaPhiAssHad,
237       kHistEtaPhiAssProt,
238       kHistPhiTrgJetEvPlane,
239       kHistPhiTrgHadEvPlane,
240       kHistPhiRndTrgJetEvPlane,
241       kHistPhiRndTrgHadEvPlane,
242       kHistPhiAssHadEvPlane,
243       kHistPhiAssHadVsEvPlane,
244       kHistPhiAssProtEvPlane,
245       kHistPhiTrgJetEvPlane3,
246       kHistPhiTrgHadEvPlane3,
247       kHistPhiAssHadEvPlane3,
248       kHistPhiAssProtEvPlane3,
249       kHistLast
250   };
251
252   // statistics
253   enum Stat_t {
254       kStatSeen = 1,
255       kStatTrg,
256       kStatVtx,
257       kStatCent,
258       kStatEvPlane,
259       kStatPID,
260       kStatUsed,
261       kStatEvCuts,
262       kStatCentral,
263       kStatSemiCentral,
264       kStatLast
265   };
266
267   // trigger conditions
268   enum Trigger_t {
269       kTriggerMB = 0,
270       kTriggerInt,
271       kTriggerLast
272   };
273
274   // classification
275   enum CorrType_t {
276     kCorrHadHad = 0,
277     kCorrHadProt,
278     kCorrJetHad,
279     kCorrJetProt,
280
281     kCorrRndHadHad,
282     kCorrRndHadProt,
283     kCorrRndJetHad,
284     kCorrRndJetProt,
285
286     kCorrRndHadExcHad,
287     kCorrRndHadExcProt,
288     kCorrRndJetExcHad,
289     kCorrRndJetExcProt,
290
291     kCorrLast
292   };
293
294   enum Class_t {
295     kClCentral = 0,
296     kClSemiCentral,
297     // kClDijet,
298     kClLast
299   };
300
301   enum Trg_t {
302     kTrgHad = 0,
303     kTrgJet,
304     kTrgHadRnd,
305     kTrgJetRnd,
306     kTrgLast
307   };
308
309   enum Ass_t {
310     kAssHad = 0,
311     kAssProt,
312     kAssHadJetExc,
313     kAssProtJetExc,
314     kAssHadHadExc,
315     kAssProtHadExc,
316     kAssLast
317   };
318
319   enum Ev_t {
320     kEvSame = 0,
321     kEvMix,
322     kEvLast
323   };
324
325   class AliPartCorr : public AliVParticle {
326   public:
327     AliPartCorr(Float_t eta = 0., Float_t phi = 0., Float_t pt = 0., Float_t charge = 0) :
328       fPt(pt), fEta(eta), fPhi(phi), fCharge(charge) {}
329     AliPartCorr(const AliVParticle &rhs) :
330       fPt(rhs.Pt()), fEta(rhs.Eta()), fPhi(rhs.Phi()), fCharge(rhs.Charge()) {}
331     virtual ~AliPartCorr() {}
332     
333     // kinematics
334     virtual Double_t Px() const { AliFatal("not implemented"); return 0; }
335     virtual Double_t Py() const { AliFatal("not implemented"); return 0; }
336     virtual Double_t Pz() const { AliFatal("not implemented"); return 0; }
337     virtual Double_t Pt() const { return fPt; }
338     virtual Double_t P()  const { AliFatal("not implemented"); return 0; }
339     virtual Bool_t   PxPyPz(Double_t[3]) const { AliFatal("not implemented"); return 0; }
340
341     virtual Double_t Xv() const { AliFatal("not implemented"); return 0; }
342     virtual Double_t Yv() const { AliFatal("not implemented"); return 0; }
343     virtual Double_t Zv() const { AliFatal("not implemented"); return 0; }
344     virtual Bool_t   XvYvZv(Double_t[3]) const { AliFatal("not implemented"); return 0; }
345
346     virtual Double_t OneOverPt()  const { AliFatal("not implemented"); return 0; }
347     virtual Double_t Phi()        const { return fPhi; }
348     virtual Double_t Theta()      const { AliFatal("not implemented"); return 0; }
349
350     virtual Double_t E()          const { AliFatal("not implemented"); return 0; }
351     virtual Double_t M()          const { AliFatal("not implemented"); return 0; }
352
353     virtual Double_t Eta()        const { return fEta; }
354     virtual Double_t Y()          const { AliFatal("not implemented"); return 0; }
355
356     virtual Short_t Charge()      const { return fCharge; }
357     virtual Int_t   GetLabel()    const { AliFatal("not implemented"); return 0; }
358
359     virtual Int_t   PdgCode()     const { AliFatal("not implemented"); return 0; }
360     virtual const Double_t *PID() const { AliFatal("not implemented"); return 0; }
361
362   protected:
363     Float_t fPt;
364     Float_t fEta;
365     Float_t fPhi;
366     Short_t fCharge;
367
368     ClassDef(AliPartCorr, 1);
369   };
370
371   class AliHistCorr : public TNamed {
372   public:
373     AliHistCorr(TString name, TList *outputList = 0x0);
374     ~AliHistCorr();
375
376     void Trigger(Float_t phi, Float_t eta, Float_t weight = 1.) {
377       fHistStat->Fill(1., weight);
378       if (fHistCorrTrgEtaPhi)
379         fHistCorrTrgEtaPhi->Fill(phi, eta, weight);
380     }
381     void Ass(Float_t phi, Float_t eta, Float_t weight = 1.) {
382       if (fHistCorrAssEtaPhi)
383         fHistCorrAssEtaPhi->Fill(phi, eta, weight);
384     }
385     void Fill(AliVParticle *trgPart, AliVParticle *assPart, Float_t weight = 1.);
386     void Fill(TLorentzVector *trgPart, AliVParticle *assPart, Float_t weight = 1.);
387     void Fill(TLorentzVector *trgPart, TLorentzVector *assPart, Float_t weight = 1.);
388
389   protected:
390     TList *fOutputList;
391
392     TH1F *fHistStat;
393
394     TH1F *fHistCorrPhi;
395     TH2F *fHistCorrPhi2;
396     TH2F *fHistCorrEtaPhi;
397     TH2F *fHistCorrAvgEtaPhi;
398     TH2F *fHistCorrTrgEtaPhi;
399     TH2F *fHistCorrAssEtaPhi;
400
401     const Float_t fHistDphiLo;
402     const Int_t   fHistDphiNbins;
403     const Int_t   fHistDetaNbins;
404
405     AliHistCorr(const AliHistCorr &rhs);
406     AliHistCorr& operator=(const AliHistCorr &rhs);
407
408     ClassDef(AliHistCorr, 1);
409   };
410
411   AliHistCorr*& GetHistCorr(CorrType_t corr, Class_t cl, Ev_t ev) { return fHistCorr[kEvLast*(kClLast*corr + cl) + ev]; }
412   AliEventPoolManager*& GetPoolMgr(Ass_t ass) { return fPoolMgr[ass]; }
413   AliEventPool*& GetPool(Ass_t ass) { return fPool[ass]; }
414
415 protected:
416   AliMCEvent  *fMCEvent; //!
417   AliESDEvent *fESDEvent; //!
418   AliAODEvent *fAODEvent; //!
419
420   Int_t fRunNumber; //! current run number
421   AliOADBContainer *fOADBContainerTOF; //! container for OADB entry with TOF parameters
422   AliTOFPIDParams *fParamsTOF; //! TOF parametrization
423
424   AliEventplane *fEventplane; //! pointer to the event plane
425
426   UInt_t fTriggerMask;          //! internal representation of trigger conditions
427   UInt_t fClassMask;            //! internal representation of event classes
428   Float_t fCentrality; //!
429   Float_t fCentralityCheck; //!
430   Float_t fZvtx; //!
431   AliPIDResponse *fPIDResponse; //!
432   Float_t fEventPlaneAngle; //!
433   Float_t fEventPlaneAngleCheck; //!
434   Float_t fEventPlaneAngle3; //!
435
436   TObjArray *fPrimTrackArrayAss; //!
437   TObjArray *fPrimTrackArrayTrg; //!
438   TClonesArray *fPrimConstrainedTrackArray; //!
439   TClonesArray *fJetArray; //!
440
441   AliEventPoolManager *fPoolMgr[kAssProt + 1]; //!
442   AliEventPool *fPool[kAssProt + 1]; //!
443
444   AliHistCorr **fHistCorr; //! [kCorrLast*kEvLast*kClLast]; //!
445
446   Int_t fErrorMsg; //! remaining number of error messages to be printed
447
448   Bool_t DetectTriggers();
449   void   MarkTrigger(Trigger_t trg) { fTriggerMask |= (1 << trg); }
450   Bool_t IsTrigger(Trigger_t trg) const { return (fTriggerMask & (1 << trg)); }
451
452   Bool_t DetectClasses();
453   void   MarkClass(Class_t cl) { fClassMask |= (1 << cl); }
454   Bool_t IsClass(Class_t cl) const { return (fClassMask & (1 << cl)); }
455
456   Bool_t PrepareEvent();
457   Bool_t CleanUpEvent();
458
459   Float_t GetCentrality() const { return fCentrality; }
460   Float_t GetEventPlaneAngle() const { return fEventPlaneAngle; }
461   AliPIDResponse* GetPID() const { return fPIDResponse; }
462   Bool_t IsCentral() const { return ((fCentrality >= 0.) && (fCentrality <= 10.)); }
463   Bool_t IsSemiCentral() const { return ((fCentrality >= 30.) && (fCentrality <= 50.)); }
464
465   AliVTrack* GetLeadingTrack(const AliAODJet *jet) const;
466
467   Float_t GetDPhiStar(Float_t phi1, Float_t pt1, Float_t charge1,
468                       Float_t phi2, Float_t pt2, Float_t charge2,
469                       Float_t radius, Float_t bSign) const;
470
471   Bool_t AcceptTrigger(const AliVTrack *trg);
472   Bool_t AcceptTrigger(const AliAODJet *trg);
473   Bool_t AcceptAssoc(const AliVTrack *trk) const;
474   Bool_t IsProton(const AliVTrack *trk) const;
475   Bool_t AcceptAngleToEvPlane(Float_t phi, Float_t psi) const;
476   Bool_t AcceptTwoTracks(const AliVParticle *trgPart, const AliVParticle *assPart) const;
477
478   TObjArray* CloneTracks(TObjArray *tracks) const;
479
480   Bool_t GenerateRandom(TCollection *trgJetArray, TCollection *trgHadArray,
481                         TCollection *assHadJetArray, TCollection *assProtJetArray,
482                         TCollection *assHadHadArray, TCollection *assProtHadArray,
483                         Float_t pFraction = .5);
484
485   Bool_t Correlate(CorrType_t corr, Class_t cl, Ev_t ev,
486                    TCollection *trgArray, TCollection *assArray, Float_t weight = 1.);
487
488   Bool_t Correlate(Trg_t trg, Ass_t ass, Class_t cl, Ev_t ev,
489                    TCollection *trgArray, TCollection *assArray, Float_t weight = 1.);
490
491   // output objects
492   TList *fOutputList;           //! list of output objects
493
494   // histogram management
495   TH1  *fHist[kHistLast];       //! pointers to histogram
496   const char *fShortTaskId;     //! short identifier for the task
497
498   TH1*&  GetHistogram(Hist_t hist, const Int_t idx = 0) { return fHist[hist + idx]; }
499
500   TH1*   AddHistogram(Hist_t hist, const char *hid, TString title,
501                       Int_t xbins, Float_t xmin, Float_t xmax, Int_t binType = 1);
502   TH2*   AddHistogram(Hist_t hist, const char *hid, TString title,
503                       Int_t xbins, Float_t xmin, Float_t xmax,
504                       Int_t ybins, Float_t ymin, Float_t ymax, Int_t binType = 1);
505   TH3*   AddHistogram(Hist_t hist, const char *hid, TString title,
506                       Int_t xbins, Float_t xmin, Float_t xmax,
507                       Int_t ybins, Float_t ymin, Float_t ymax,
508                       Int_t zbins, Float_t zmin, Float_t zmax, Int_t binType = 1);
509
510   void    FillH1(Hist_t hist, Float_t x, Float_t weight = 1., Int_t idx = 0)
511   { GetHistogram(hist, idx)->Fill(x, weight); }
512   void    FillH2(Hist_t hist, Float_t x, Float_t y, Float_t weight = 1., Int_t idx = 0)
513   { ((TH2*) GetHistogram(hist, idx))->Fill(x, y, weight); }
514   void    FillH3(Hist_t hist, Float_t x, Float_t y, Float_t z, Float_t weight = 1., Int_t idx = 0)
515   { ((TH3*) GetHistogram(hist, idx))->Fill(x, y, z, weight); }
516
517   const char* fkCorrTypeName[kCorrLast]; //!
518   const char* fkClassName[kClLast]; //!
519   const char* fkEvName[kEvLast]; //!
520
521   // task configuration
522   static const Int_t fgkStringLength = 100; // max length for the jet branch name
523   char fJetBranchName[fgkStringLength];     // jet branch name
524
525   const Bool_t fUseStandardCuts;
526   Bool_t fUseEvplaneV0;
527
528   AliESDtrackCuts *fCutsPrimTrg;        // track cuts for primary particles (trigger)
529   AliESDtrackCuts *fCutsPrimTrgConstrain;       // track cuts for primary particles (trigger)
530   AliESDtrackCuts *fCutsPrimAss;        // track cuts for primary particles (associate)
531   Float_t fCutsTwoTrackEff;
532
533   UInt_t  fAssFilterMask;
534   Bool_t  fRequirePID;
535   Float_t fTrgJetEtaMax;
536   Float_t fHadEtaMax;
537   Float_t fTrgPartPtMin;
538   Float_t fTrgPartPtMax;
539   Float_t fTrgJetPtMin;
540   Float_t fTrgJetPtMax;
541   Float_t fTrgJetLeadTrkPtMin;
542   Float_t fTrgJetLeadTrkPtMax;
543   Float_t fTrgJetAreaMin;
544   Float_t fAssPartPtMin;
545   Float_t fAssPartPtMax;
546   Float_t fTrgAngleToEvPlane;
547
548   Float_t fToyMeanNoPart;
549   Float_t fToyRadius;
550   Float_t fToySmearPhi;
551
552   TF1 *fTrgJetPhiModCent;
553   TF1 *fTrgJetPhiModSemi;
554   TF1 *fTrgHadPhiModCent;
555   TF1 *fTrgHadPhiModSemi;
556
557   // not implemented
558   AliAnalysisTaskJetProtonCorr(const AliAnalysisTaskJetProtonCorr &rhs);
559   AliAnalysisTaskJetProtonCorr& operator=(const AliAnalysisTaskJetProtonCorr &rhs);
560
561   ClassDef(AliAnalysisTaskJetProtonCorr, 1);
562 };