]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGJE/UserTasks/AliAnalysisTaskIDFragmentationFunction.h
- Added histo for jet pT correlation with pp V0M mult percentile estimator - Also...
[u/mrichter/AliRoot.git] / PWGJE / UserTasks / AliAnalysisTaskIDFragmentationFunction.h
1 // *************************************************************************
2 // * Task for Fragmentation Function Analysis in PWG4 Jet Task Force Train *
3 // *************************************************************************
4
5 #ifndef ALIANALYSISTASKIDFRAGMENTATIONFUNCTION_H
6 #define ALIANALYSISTASKIDFRAGMENTATIONFUNCTION_H
7
8 /* Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
9  * See cxx source for full Copyright notice                               */
10
11 /* $Id$ */
12
13 class AliESDEvent;
14 class AliAODEvent;
15 class AliAODJet;
16 class AliAODExtension;
17 class TList;
18 class TH1F;
19 class TH2F;
20 class TH3F;
21 class TProfile;
22 class THnSparse; 
23 class TRandom3;
24 class TArrayS;
25 class AliAnalysisUtils;
26 class AliAODTrack;
27 class AliAODMCParticle;
28
29 #include "AliAnalysisTaskSE.h"
30 #include "AliPID.h"
31   
32 class AliAnalysisTaskIDFragmentationFunction : public AliAnalysisTaskSE {
33
34  public:
35   
36   //----------------------------------------
37   class AliFragFuncHistos : public TObject
38   {
39     
40     public:
41     
42     AliFragFuncHistos(const char* name = "FFhistos", 
43                       Int_t nJetPt = 0, Float_t jetPtMin = 0, Float_t jetPtMax = 0,
44                       Int_t nPt = 0, Float_t ptMin = 0, Float_t ptMax = 0,
45                       Int_t nXi = 0, Float_t xiMin = 0, Float_t xiMax = 0,
46                       Int_t nZ  = 0, Float_t zMin  = 0, Float_t zMax  = 0);
47     
48     AliFragFuncHistos(const AliFragFuncHistos& copy);
49     AliFragFuncHistos& operator=(const AliFragFuncHistos &o);
50     virtual ~AliFragFuncHistos();
51     
52     virtual void DefineHistos();
53     virtual void FillFF(Float_t trackPt, Float_t jetPt,Bool_t incrementJetPt, Float_t norm = 0, Bool_t scaleStrangeness = kFALSE, Float_t scaleFacStrangeness = 1.);
54
55     virtual void AddToOutput(TList* list) const;
56
57   private:
58
59     Int_t   fNBinsJetPt; // FF histos bins
60     Float_t fJetPtMin;   // FF histos limits
61     Float_t fJetPtMax;   // FF histos limits
62     Int_t   fNBinsPt;    // FF histos bins
63     Float_t fPtMin;      // FF histos limits
64     Float_t fPtMax;      // FF histos limits
65     Int_t   fNBinsXi;    // FF histos bins
66     Float_t fXiMin;      // FF histos limits
67     Float_t fXiMax;      // FF histos limits
68     Int_t   fNBinsZ;     // FF histos bins
69     Float_t fZMin;       // FF histos limits
70     Float_t fZMax;       // FF histos limits
71   
72     TH2F*   fh2TrackPt;   //! FF: track transverse momentum 
73     TH2F*   fh2Xi;        //! FF: xi 
74     TH2F*   fh2Z;         //! FF: z  
75     TH1F*   fh1JetPt;     //! jet pt 
76
77     TString fNameFF;      // histo names prefix
78     
79     ClassDef(AliFragFuncHistos, 1);
80   };
81   
82   //----------------------------------------
83   class AliFragFuncQAJetHistos : public TObject
84   {
85
86   public:
87  
88     AliFragFuncQAJetHistos(const char* name = "QAJethistos",
89                 Int_t nPt  = 0, Float_t ptMin  = 0, Float_t ptMax  = 0,
90                 Int_t nEta = 0, Float_t etaMin = 0, Float_t etaMax = 0,
91                 Int_t nPhi = 0, Float_t phiMin = 0, Float_t phiMax = 0);
92       
93     AliFragFuncQAJetHistos(const AliFragFuncQAJetHistos& copy);
94     AliFragFuncQAJetHistos& operator=(const AliFragFuncQAJetHistos &o);
95     virtual ~AliFragFuncQAJetHistos();
96     virtual void DefineHistos();
97     virtual void FillJetQA(Float_t eta, Float_t phi, Float_t pt);
98     virtual void AddToOutput(TList* list) const;
99
100   private:
101     
102     Int_t   fNBinsPt;    // jet QA histos bins
103     Float_t fPtMin;      // jet QA histos limits
104     Float_t fPtMax;      // jet QA histos limits
105     Int_t   fNBinsEta;   // jet QA histos bins
106     Float_t fEtaMin;     // jet QA histos limits
107     Float_t fEtaMax;     // jet QA histos limits
108     Int_t   fNBinsPhi;   // jet QA histos bins
109     Float_t fPhiMin;     // jet QA histos limits
110     Float_t fPhiMax;     // jet QA histos limits
111     
112     TH2F*   fh2EtaPhi;   //! jet phi vs eta 
113     TH1F*   fh1Pt;       //! jet transverse momentum 
114     TString fNameQAJ;    // histo names prefix
115     
116     ClassDef(AliFragFuncQAJetHistos, 1);
117   };
118   
119   //----------------------------------------
120   class AliFragFuncQATrackHistos : public TObject
121   {
122
123   public:
124
125     AliFragFuncQATrackHistos(const char* name = "QATrackhistos", 
126                   Int_t nPt  = 0, Float_t ptMin  = 0, Float_t ptMax  = 0,
127                   Int_t nEta = 0, Float_t etaMin = 0, Float_t etaMax = 0,
128                   Int_t nPhi = 0, Float_t phiMin = 0, Float_t phiMax = 0, 
129                   Float_t ptThresh = 0);
130     
131     AliFragFuncQATrackHistos(const AliFragFuncQATrackHistos& copy);
132     AliFragFuncQATrackHistos& operator=(const AliFragFuncQATrackHistos &o);
133     virtual ~AliFragFuncQATrackHistos();
134     virtual void DefineHistos();
135     virtual void FillTrackQA(Float_t eta, Float_t phi, Float_t pt, Bool_t weightPt = kFALSE, Float_t norm = 0., Bool_t scaleStrangeness = kFALSE, Float_t scaleFacStrangeness = 1.);
136     virtual void AddToOutput(TList* list) const;
137
138   private:
139     
140     Int_t   fNBinsPt;    // track QA histos bins in pt
141     Float_t fPtMin;      // track QA histos limits in pt
142     Float_t fPtMax;      // track QA histos limits in pt
143     Int_t   fNBinsEta;   // track QA histos bins in eta
144     Float_t fEtaMin;     // track QA histos limits in eta
145     Float_t fEtaMax;     // track QA histos limits in eta
146     Int_t   fNBinsPhi;   // track QA histos bins in phi
147     Float_t fPhiMin;     // track QA histos limits in phi
148     Float_t fPhiMax;     // track QA histos limits in phi
149
150     Float_t fHighPtThreshold; //  high pt track phi vs eta distribution
151
152     TH2F*   fh2EtaPhi;        //! track phi vs eta 
153     TH1F*   fh1Pt;            //! track transverse momentum 
154     TH2F*   fh2HighPtEtaPhi;  //! phi vs eta for high pt (>fgHighPtThreshold) tracks
155     TH2F*   fh2PhiPt;         //! track phi vs pt
156
157     TString fNameQAT;         // histo names prefix
158     
159     ClassDef(AliFragFuncQATrackHistos, 2);
160   };
161   
162
163   AliAnalysisTaskIDFragmentationFunction(); 
164   AliAnalysisTaskIDFragmentationFunction(const char *name);
165   AliAnalysisTaskIDFragmentationFunction(const  AliAnalysisTaskIDFragmentationFunction &copy);
166   AliAnalysisTaskIDFragmentationFunction& operator=(const  AliAnalysisTaskIDFragmentationFunction &o);
167   virtual ~AliAnalysisTaskIDFragmentationFunction();
168   
169   virtual void   UserCreateOutputObjects();
170   virtual void   Init();
171   virtual void   LocalInit() {Init();}
172
173   virtual void   UserExec(Option_t *option);
174   virtual void   Terminate(Option_t* );
175   virtual Bool_t Notify();
176
177   virtual void   SetNonStdFile(char* c){fNonStdFile = c;} 
178
179   virtual void   SetTrackTypeGen(Int_t i){fTrackTypeGen = i;}
180   virtual void   SetJetTypeGen(Int_t i){fJetTypeGen = i;}
181   virtual void   SetJetTypeRecEff(Int_t i){fJetTypeRecEff = i;}
182
183   virtual void   SetBranchRecBackClusters(const char* c){fBranchRecBckgClusters = c;}
184   virtual void   SetBranchGenJets(const char* c){fBranchGenJets = c;}
185   virtual void   SetBranchRecJets(const char* c){fBranchRecJets = c;}
186   virtual void   SetBranchEmbeddedJets(const char* c){fBranchEmbeddedJets = c;}
187
188   virtual void   SetTrackCuts(Float_t trackPt = 0.15, Float_t trackEtaMin = -0.9, Float_t trackEtaMax = 0.9, 
189                               Float_t trackPhiMin = 0., Float_t trackPhiMax = 2*TMath::Pi())
190   {fTrackPtCut = trackPt; fTrackEtaMin = trackEtaMin; fTrackEtaMax = trackEtaMax; 
191     fTrackPhiMin = trackPhiMin; fTrackPhiMax = trackPhiMax;}
192
193   virtual void   UseExtraTracks()        { fUseExtraTracks =  1;}
194   virtual void   UseExtraonlyTracks()    { fUseExtraTracks = -1;}
195
196   virtual void   UseExtraTracksBgr()     { fUseExtraTracksBgr =  1;}
197   virtual void   UseExtraonlyTracksBgr() { fUseExtraTracksBgr = -1;}
198
199   virtual void   SetCutFractionPtEmbedded(Float_t cut = 0) { fCutFractionPtEmbedded = cut; }
200   virtual void   SetUseEmbeddedJetAxis(Bool_t b = kTRUE)   { fUseEmbeddedJetAxis = b; }
201   virtual void   SetUseEmbeddedJetPt(Bool_t  b = kTRUE)    { fUseEmbeddedJetPt   = b; }
202
203   virtual void   UseAODInputJets(Bool_t b) {fUseAODInputJets = b;}  
204   virtual void   SetFilterMask(UInt_t i) {fFilterMask = i;}
205   virtual void   UsePhysicsSelection(Bool_t b) {fUsePhysicsSelection = b;}
206   virtual void   SetEventSelectionMask(UInt_t i){fEvtSelectionMask = i;}
207   virtual void   SetEventClass(Int_t i){fEventClass = i;}
208   virtual void   SetMaxVertexZ(Float_t z){fMaxVertexZ = z;}
209   virtual void   SetJetCuts(Float_t jetPt = 5., Float_t jetEtaMin = -0.5, Float_t jetEtaMax = 0.5, 
210                             Float_t jetPhiMin = 0., Float_t jetPhiMax = 2*TMath::Pi())
211   {fJetPtCut = jetPt; fJetEtaMin = jetEtaMin; fJetEtaMax = jetEtaMax; 
212     fJetPhiMin = jetPhiMin; fJetPhiMax = jetPhiMax;}
213
214   virtual void   SetFFRadius(Float_t r = 0.4) { fFFRadius = r; }
215   virtual void   SetFFMinLTrackPt(Float_t pt = -1) { fFFMinLTrackPt = pt; }
216   virtual void   SetFFMaxTrackPt(Float_t pt = -1) { fFFMaxTrackPt = pt; }
217   virtual void   SetFFMinNTracks(Int_t nTracks = 0) { fFFMinnTracks = nTracks; }
218   virtual void   SetFFBckgRadius(Float_t r = 0.7) { fFFBckgRadius = r; }
219   virtual void   SetBckgMode(Bool_t bg = 1) { fBckgMode = bg; }
220   virtual void   SetBckgType(Int_t bg0 = 0, Int_t bg1 = 0,Int_t bg2 = 0, Int_t bg3 = 0, Int_t bg4 = 0) 
221   { fBckgType[0] = bg0; fBckgType[1] = bg1; fBckgType[2] = bg2; fBckgType[3] = bg3; fBckgType[4] = bg4; }
222   virtual void   SetQAMode(Int_t qa = 3)      {fQAMode = qa;}
223   virtual void   SetFFMode(Int_t ff = 1)      {fFFMode = ff;}
224   virtual void   SetIDFFMode(Int_t idff = 0)      {fIDFFMode = idff;}
225   virtual void   SetEffMode(Int_t eff = 1)    {fEffMode = eff;}
226   virtual void   SetJSMode(Int_t js = 1)      {fJSMode = js;}
227
228   static  void   SetProperties(TH1* h,const char* x, const char* y);
229   static  void   SetProperties(TH1* h,const char* x, const char* y,const char* z);
230   static  void   SetProperties(THnSparse* h,Int_t dim, const char** labels);
231
232   void   SetHighPtThreshold(Float_t pt = 5.) { fQATrackHighPtThreshold = pt; }
233
234   void   SetFFHistoBins(Int_t nJetPt = 245, Float_t jetPtMin = 5, Float_t jetPtMax = 250, 
235                         Int_t nPt = 200, Float_t ptMin = 0., Float_t ptMax = 200., 
236                         Int_t nXi = 70, Float_t xiMin = 0., Float_t xiMax = 7.,
237                         Int_t nZ = 22,  Float_t zMin = 0.,  Float_t zMax = 1.1)
238   { fFFNBinsJetPt = nJetPt; fFFJetPtMin = jetPtMin; fFFJetPtMax = jetPtMax; 
239     fFFNBinsPt = nPt; fFFPtMin = ptMin; fFFPtMax = ptMax;
240     fFFNBinsXi = nXi; fFFXiMin = xiMin; fFFXiMax = xiMax;
241     fFFNBinsZ  = nZ;  fFFZMin  = zMin;  fFFZMax  = zMax; }
242
243   void  SetQAJetHistoBins(Int_t nPt = 300, Float_t ptMin = 0., Float_t ptMax = 300.,
244                           Int_t nEta = 20, Float_t etaMin = -1.0, Float_t etaMax = 1.0,
245                           Int_t nPhi = 60, Float_t phiMin = 0., Float_t phiMax = 2*TMath::Pi())
246     { fQAJetNBinsPt = nPt; fQAJetPtMin = ptMin; fQAJetPtMax = ptMax;
247       fQAJetNBinsEta = nEta; fQAJetEtaMin = etaMin; fQAJetEtaMax = etaMax;
248       fQAJetNBinsPhi = nPhi; fQAJetPhiMin = phiMin; fQAJetPhiMax = phiMax; }
249   
250   void  SetQATrackHistoBins(Int_t nPt = 200, Float_t ptMin = 0., Float_t ptMax = 200.,
251                             Int_t nEta = 20, Float_t etaMin = -1.0, Float_t etaMax = 1.0,
252                             Int_t nPhi = 60, Float_t phiMin = 0., Float_t phiMax = 2*TMath::Pi())
253   { fQATrackNBinsPt = nPt; fQATrackPtMin = ptMin; fQATrackPtMax = ptMax;
254     fQATrackNBinsEta = nEta; fQATrackEtaMin = etaMin; fQATrackEtaMax = etaMax;
255     fQATrackNBinsPhi = nPhi; fQATrackPhiMin = phiMin; fQATrackPhiMax = phiMax; }
256   
257
258
259   Float_t  GetFFRadius() const { return fFFRadius; }
260   Float_t  GetFFMinLTrackPt() const { return fFFMinLTrackPt; }
261   Float_t  GetFFMaxTrackPt() const { return fFFMaxTrackPt; }
262   Float_t  GetFFMinNTracks() const { return fFFMinnTracks; }
263   Float_t  GetFFBckgRadius() const { return fFFBckgRadius; }
264   void     GetJetTracksTrackrefs(TList* l, const AliAODJet* j, Double_t minPtL, Double_t maxPt, Bool_t& isBadPt);
265   void     GetJetTracksPointing(TList* in, TList* out, const AliAODJet* j, Double_t r, Double_t& sumPt, Double_t minPtL, Double_t maxPt,
266                                 Bool_t& isBadPt);  
267   void     GetTracksOutOfNJets(Int_t nCases, TList* in, TList* out, const TList* jets, Double_t& pt);
268   void     GetTracksOutOfNJetsStat(Int_t nCases, TList* in, TList* out, const TList* jets, Double_t& pt, Double_t &normFactor);
269   void     GetTracksTiltedwrpJetAxis(Float_t alpha, TList* inputlist, TList* outputlist, const AliAODJet* jet, Double_t radius, Double_t& sumPt);
270   void     GetTracksTiltedwrpJetAxisWindow(Float_t alpha, TList* inputlist, TList* outputlist, const AliAODJet* jet, Double_t radius, Double_t& sumPt, Double_t &normFactor);
271
272   void     AssociateGenRec(TList* tracksAODMCCharged,TList* tracksRec, TArrayI& indexAODTr,TArrayI& indexMCTr,TArrayS& isRefGen,TH2F* fh2PtRecVsGen);
273
274   void     FillSingleTrackHistosRecGen(AliFragFuncQATrackHistos* trackQAGen, AliFragFuncQATrackHistos* trackQARec, TList* tracksGen, 
275                                        const TArrayI& indexAODTr, const TArrayS& isRefGen, Bool_t scaleStrangeness = kFALSE);
276
277
278   void     FillJetTrackHistosRec(AliFragFuncHistos* histRec,  AliAODJet* jet, 
279                                  TList* jetTrackList, const TList* tracksGen, const TList* tracksRec, const TArrayI& indexAODTr,
280                                  const TArrayS& isRefGen, TList* jetTrackListTR = 0, Bool_t scaleStrangeness = kFALSE,
281                                  Bool_t fillJS = kFALSE, TProfile* hProNtracksLeadingJet = 0, TProfile** hProDelRPtSum = 0, TProfile* hProDelR80pcPt = 0);
282
283
284   Float_t  CalcJetArea(Float_t etaJet, Float_t rc) const;
285   void     GetClusterTracksOutOf1Jet(AliAODJet* jet, TList* outputlist, Double_t &normFactor);
286   void     GetClusterTracksMedian(TList* outputlist, Double_t &normFactor);
287
288   void     FillBckgHistos(Int_t type, TList* inputtracklist, TList* inputjetlist, AliAODJet* jet, 
289                           AliFragFuncHistos* ffbckghistocuts,AliFragFuncQATrackHistos* qabckghistos,TH1F* fh1Mult = 0); 
290  
291   Double_t GetMCStrangenessFactor(Double_t pt) const;
292   Double_t GetMCStrangenessFactorCMS(AliAODMCParticle* daughter) const;
293   
294   Bool_t IsSecondaryWithStrangeMotherMC(AliAODMCParticle* part);
295
296   void FillJetShape(AliAODJet* jet, TList* list,  TProfile* hProNtracksLeadingJet, TProfile** hProDelRPtSum, TProfile* hProDelR80pcPt=0, Double_t dPhiUE=0, Double_t normUE = 0, Bool_t scaleStrangeness = kFALSE);
297
298   const TString* GetNamesOfInclusivePIDtasks() const { return fNameInclusivePIDtask; };
299   void SetNamesOfInclusivePIDtasks(Int_t numNames, const TString* names);
300   
301   const TString* GetNamesOfJetPIDtasks() const { return fNameJetPIDtask; };
302   void SetNamesOfJetPIDtasks(Int_t numNames, const TString* names);
303         
304         
305         Bool_t GetIsPP() const { return fIsPP; };
306         void SetIsPP(Bool_t flag) { fIsPP = flag; };
307   
308   Bool_t GetOnlyLeadingJets() const { return fOnlyLeadingJets; }
309   void SetOnlyLeadingJets(Bool_t onlyLeadingJets) { fOnlyLeadingJets = onlyLeadingJets; }
310   
311   // Consts
312   enum {kTrackUndef=0, kTrackAOD, kTrackAODQualityCuts, kTrackAODCuts, 
313         kTrackAODExtra, kTrackAODExtraonly, kTrackAODExtraCuts, kTrackAODExtraonlyCuts, 
314         kTrackKineAll, kTrackKineCharged, kTrackKineChargedAcceptance, 
315         kTrackAODMCAll, kTrackAODMCCharged, kTrackAODMCChargedAcceptance, kTrackAODMCChargedSecS, kTrackAODMCChargedSecNS, kTrackAOCMCChargedPrimAcceptance};
316   enum {kJetsUndef=0, kJetsRec, kJetsRecAcceptance, kJetsGen, kJetsGenAcceptance, kJetsKine, kJetsKineAcceptance,kJetsEmbedded};
317   enum {kBckgNone=0, kBckgPerp, kBckgOutLJ, kBckgOut2J, kBckgClusters, kBckgClustersOutLeading, kBckgOut3J, kBckgOutAJ, kBckgOutLJStat, 
318         kBckgOut2JStat, kBckgOut3JStat, kBckgOutAJStat,  kBckgASide, kBckgASideWindow, kBckgPerpWindow, kBckgPerp2, kBckgPerp2Area};
319
320  
321  protected:
322   
323   Int_t   GetListOfTracks(TList* list, Int_t type);
324   Int_t   GetListOfJets(TList* list, Int_t type);
325   Int_t   GetListOfBckgJets(TList *list, Int_t type);
326
327   AliESDEvent* fESD;      //! ESD event
328   AliAODEvent* fAOD;      //! AOD event
329   AliAODEvent* fAODJets;  //! AOD event with jet branch (case we have AOD both in input and output)
330   AliAODExtension  *fAODExtension; //! where we take the jets from can be input or output AOD
331   TString       fNonStdFile; // name of delta aod file to catch the extension
332  
333  
334   TString fBranchRecJets;         // branch name for reconstructed jets
335   TString fBranchRecBckgClusters; // branch name for reconstructed background clusters 
336   TString fBranchGenJets;         // branch name for generated jets
337   TString fBranchEmbeddedJets;    // branch name for embedded jets
338
339   Int_t   fTrackTypeGen;        // type of generated tracks
340   Int_t   fJetTypeGen;          // type of generated jets
341
342   Int_t   fJetTypeRecEff;       // type of jets used for filling reconstruction efficiency histos
343
344   Bool_t  fUseAODInputJets;     // take jets from in/output - only relevant if AOD event both in input AND output and we want to use output
345   UInt_t  fFilterMask;          // filter bit for selected tracks
346   Bool_t  fUsePhysicsSelection; // switch for event selection
347   UInt_t  fEvtSelectionMask;    // trigger class selection
348   Int_t   fEventClass;          // centrality class selection
349   Float_t fMaxVertexZ;          // maximum abs(z) position of primiary vertex [cm]
350
351   // track cuts
352   Float_t fTrackPtCut;    // track transverse momentum cut
353   Float_t fTrackEtaMin;   // track eta cut
354   Float_t fTrackEtaMax;   // track eta cut
355   Float_t fTrackPhiMin;   // track phi cut
356   Float_t fTrackPhiMax;   // track phi cut
357   
358   Int_t   fUseExtraTracks;         // +/- 1: embedded extra/extra only tracks, default: 0 (ignore extra tracks)
359   Int_t   fUseExtraTracksBgr;      // +/- 1: background: use embedded extra/extra only tracks, default: 0 (ignore extra tracks)
360   Float_t fCutFractionPtEmbedded;  // cut on ratio of embedded pt found in jet
361   Bool_t  fUseEmbeddedJetAxis;     // use axis of embedded jet for FF
362   Bool_t  fUseEmbeddedJetPt;       // use axis of embedded jet for FF
363
364   // jet cuts
365   Float_t fJetPtCut;      // jet transverse momentum cut
366   Float_t fJetEtaMin;     // jet eta cut
367   Float_t fJetEtaMax;     // jet eta cut
368   Float_t fJetPhiMin;     // jet phi cut
369   Float_t fJetPhiMax;     // jet phi cut
370
371   Float_t fFFRadius;        // if radius > 0 construct FF from tracks within cone around jet axis, otherwise use trackRefs  
372   Float_t fFFMinLTrackPt;   // reject jets with leading track with pt smaller than this value
373   Float_t fFFMaxTrackPt;    // reject jets containing any track with pt larger than this value
374   Int_t   fFFMinnTracks;    // reject jets with less tracks than this value
375   Float_t fFFBckgRadius;    // compute background outside cone of this radius around jet axes
376   Bool_t  fBckgMode;        // Set background subtraction mode
377   Int_t   fBckgType[5];     // Set background subtraction mode
378   Int_t   fQAMode;          // QA mode: 0x00=0 none, 0x01=1 track qa, 0x10=2 track qa, 0x11=3 both
379   Int_t   fFFMode;          // fragmentation function mode
380   Int_t   fIDFFMode;        // identified fragmentation function mode (implicitely enables also normal fragmentation function mode (but not effeciency, background etc.)!)
381   Int_t   fEffMode;         // efficiency mode
382   Int_t   fJSMode;          // jet shape mode
383
384   Float_t fAvgTrials;       // average number of trials per event
385   
386   TList* fTracksRecCuts;           //! reconstructed tracks after cuts
387   TList* fTracksRecCutsEfficiency; //! reconstructed tracks after cuts for efficiency
388   TList* fTracksGen;               //! generated tracks 
389   TList* fTracksAODMCCharged;      //! AOD MC tracks 
390   TList* fTracksAODMCChargedSecNS; //! AOD MC tracks - secondaries (non-strangeness) 
391   TList* fTracksAODMCChargedSecS;  //! AOD MC tracks - secondaries (from strangeness)
392   TList* fTracksRecQualityCuts;    //! reconstructed tracks after quality cuts, no acceptance/pt cut
393
394   
395   TList* fJetsRec;        //! jets from reconstructed tracks
396   TList* fJetsRecCuts;    //! jets from reonstructed tracks after jet cuts 
397   TList* fJetsGen;        //! jets from generated tracks
398   TList* fJetsRecEff;     //! jets used for reconstruction efficiency histos 
399   TList* fJetsEmbedded;   //! jets used for embedding
400
401   TList* fBckgJetsRec;      //! jets from reconstructed tracks
402   TList* fBckgJetsRecCuts;  //! jets from reonstructed tracks after jet cuts
403   TList* fBckgJetsGen;      //! jets from generated tracks
404  
405   
406   AliFragFuncQATrackHistos* fQATrackHistosRecCuts;  //! track QA: reconstructed tracks after cuts
407   AliFragFuncQATrackHistos* fQATrackHistosGen;      //! track QA: generated tracks
408   
409   AliFragFuncQAJetHistos*  fQAJetHistosRec;             //! jet QA: jets from reconstructed tracks
410   AliFragFuncQAJetHistos*  fQAJetHistosRecCuts;         //! jet QA: jets from reconstructed tracks after jet cuts 
411   AliFragFuncQAJetHistos*  fQAJetHistosRecCutsLeading;  //! jet QA: leading jet from reconstructed tracks after jet cuts 
412   AliFragFuncQAJetHistos*  fQAJetHistosGen;             //! jet QA: jets from generated tracks  
413   AliFragFuncQAJetHistos*  fQAJetHistosGenLeading;      //! jet QA: leading jet from generated tracks  
414   AliFragFuncQAJetHistos*  fQAJetHistosRecEffLeading;   //! jet QA: leading jet used for reconstruction efficiency histos  
415   
416
417   AliFragFuncHistos*  fFFHistosRecCuts;         //! FF reconstructed tracks after cuts (leading jet) 
418   AliFragFuncHistos*  fFFHistosRecCutsInc;      //! inclusive FF (all jets) 
419   AliFragFuncHistos*  fFFHistosRecLeadingTrack; //! FF reconstructed tracks after cuts: leading track pt / jet pt (all jets)
420
421   AliFragFuncHistos*  fFFHistosGen;             //! FF generated tracks after cuts 
422   AliFragFuncHistos*  fFFHistosGenInc;          //! inclusive FF (all jets) 
423   AliFragFuncHistos*  fFFHistosGenLeadingTrack; //! FF reconstructed tracks after cuts: leading track pt / jet pt (all jets)
424
425   Float_t  fQATrackHighPtThreshold;       // track QA high transverse momentum threshold
426   
427   // histogram bins  
428
429   Int_t   fFFNBinsJetPt;    // FF histos bins
430   Float_t fFFJetPtMin;      // FF histos limits
431   Float_t fFFJetPtMax;      // FF histos limits
432
433   Int_t   fFFNBinsPt;       // FF histos bins
434   Float_t fFFPtMin;         // FF histos limits
435   Float_t fFFPtMax;         // FF histos limits
436
437   Int_t   fFFNBinsXi;       // FF histos bins
438   Float_t fFFXiMin;         // FF histos limits
439   Float_t fFFXiMax;         // FF histos limits
440
441   Int_t   fFFNBinsZ;        // FF histos bins
442   Float_t fFFZMin;          // FF histos limits
443   Float_t fFFZMax;          // FF histos limits
444
445   Int_t   fQAJetNBinsPt;    // jet QA histos bins
446   Float_t fQAJetPtMin;      // jet QA histos limits
447   Float_t fQAJetPtMax;      // jet QA histos limits
448   
449   Int_t   fQAJetNBinsEta;   // jet QA histos bins
450   Float_t fQAJetEtaMin;     // jet QA histos limits
451   Float_t fQAJetEtaMax;     // jet QA histos limits
452   
453   Int_t   fQAJetNBinsPhi;   // jet QA histos bins
454   Float_t fQAJetPhiMin;     // jet QA histos limits
455   Float_t fQAJetPhiMax;     // jet QA histos limits
456   
457   Int_t   fQATrackNBinsPt;  // track QA histos bins
458   Float_t fQATrackPtMin;    // track QA histos limits
459   Float_t fQATrackPtMax;    // track QA histos limits
460   
461   Int_t   fQATrackNBinsEta; // track QA histos bins
462   Float_t fQATrackEtaMin;   // track QA histos limits
463   Float_t fQATrackEtaMax;   // track QA histos limits
464   
465   Int_t   fQATrackNBinsPhi; // track QA histos bins
466   Float_t fQATrackPhiMin;   // track QA histos limits
467   Float_t fQATrackPhiMax;   // track QA histos limits
468   
469   // Histograms
470   TList *fCommonHistList;         // List of common histos
471   
472   TH1F  *fh1EvtSelection;         //! event cuts 
473   TH1F  *fh1VtxSelection;         //! type of accepted vertices
474   TH1F  *fh1VertexNContributors;  //! NContributors to prim vertex
475   TH1F  *fh1VertexZ;              //! prim vertex z distribution
476   TH1F  *fh1EvtMult;              //! number of reconstructed tracks after cuts 
477   TH1F  *fh1EvtCent;              //! centrality percentile 
478
479   TProfile* fh1Xsec;              //! pythia cross section and trials
480   TH1F*     fh1Trials;            //! sum of trials
481   TH1F*     fh1PtHard;            //! pt hard of the event
482   TH1F*     fh1PtHardTrials;      //! pt hard of the event
483
484   TH1F  *fh1nRecJetsCuts;         //! number of jets from reconstructed tracks per event 
485   TH1F  *fh1nGenJets;             //! number of jets from generated tracks per event
486   TH1F  *fh1nRecEffJets;          //! number of jets for reconstruction eff per event
487   TH1F  *fh1nEmbeddedJets;        //! number of embedded jets per event
488
489   TH1F  *fh1nRecBckgJetsCuts;     //! number of jets from reconstructed tracks per event
490   TH1F  *fh1nGenBckgJets;         //! number of jets from generated tracks per event
491   TH2F  *fh2PtRecVsGenPrim;       //! association rec/gen MC: rec vs gen pt, primaries 
492   TH2F  *fh2PtRecVsGenSec;        //! association rec/gen MC: rec vs gen pt, secondaries 
493   
494   TH2F  *fhDCA_XY;                //! DCA XY for all rec. particles
495   TH2F  *fhDCA_Z;                 //! DCA Z for all rec. particles
496   
497   TH2F  *fhJetPtRefMultEta5;      //! Jet pT vs. reference multiplicity (|eta|<0.5)
498   TH2F  *fhJetPtRefMultEta8;      //! Jet pT vs. reference multiplicity (|eta|<0.8)
499   TH2F  *fhJetPtMultPercent;      //! Jet pT vs. multiplicity percentile (usually V0M)
500
501   TH2F  *fhDCA_XY_prim_MCID[AliPID::kSPECIES];   //! DCA XY for all rec. prim. particles sorted by MC-ID
502   TH2F  *fhDCA_Z_prim_MCID[AliPID::kSPECIES];    //! DCA Z for all rec. prim. particles sorted by MC-ID
503  
504   TH2F  *fhDCA_XY_sec_MCID[AliPID::kSPECIES];    //! DCA XY for all rec. sec. particles sorted by MC-ID
505   TH2F  *fhDCA_Z_sec_MCID[AliPID::kSPECIES];     //! DCA Z for all rec. sec. particles sorted by MC-ID
506  
507   // tracking efficiency / secondaries
508   
509   AliFragFuncQATrackHistos* fQATrackHistosRecEffGen;      //! tracking efficiency: generated primaries 
510   AliFragFuncQATrackHistos* fQATrackHistosRecEffRec;      //! tracking efficiency: reconstructed primaries
511   AliFragFuncQATrackHistos* fQATrackHistosSecRecNS;       //! reconstructed secondaries (non-strangeness)
512   AliFragFuncQATrackHistos* fQATrackHistosSecRecS;        //! reconstructed secondaries (strange mothers)
513   AliFragFuncQATrackHistos* fQATrackHistosSecRecSsc;      //! reconstructed secondaries (strange mothers) - scale factor
514
515   AliFragFuncHistos*  fFFHistosRecEffRec;                 //! tracking efficiency: FF reconstructed primaries
516   AliFragFuncHistos*  fFFHistosSecRecNS;                  //! secondary contamination: FF reconstructed secondaries (non-strangeness)
517   AliFragFuncHistos*  fFFHistosSecRecS;                   //! secondary contamination: FF reconstructed secondaries (strange mothers)
518   AliFragFuncHistos*  fFFHistosSecRecSsc;                 //! secondary contamination: FF reconstructed secondaries (strange mothers) - scale factor
519
520
521   // Background
522   TH1F  *fh1BckgMult0; //! background multiplicity
523   TH1F  *fh1BckgMult1; //! background multiplicity
524   TH1F  *fh1BckgMult2; //! background multiplicity
525   TH1F  *fh1BckgMult3; //! background multiplicity
526   TH1F  *fh1BckgMult4; //! background multiplicity
527
528   // embedding
529   TH1F* fh1FractionPtEmbedded;         //! ratio embedded pt in rec jet to embedded jet pt 
530   TH1F* fh1IndexEmbedded;              //! index embedded jet matching to leading rec jet 
531   TH2F* fh2DeltaPtVsJetPtEmbedded;     //! delta pt rec - embedded jet
532   TH2F* fh2DeltaPtVsRecJetPtEmbedded;  //! delta pt rec - embedded jet
533   TH1F* fh1DeltaREmbedded;             //! delta R  rec - embedded jet
534
535
536   AliFragFuncQATrackHistos* fQABckgHisto0RecCuts;  //! track QA: reconstructed tracks after cuts
537   AliFragFuncQATrackHistos* fQABckgHisto0Gen;      //! track QA: generated tracks
538   AliFragFuncQATrackHistos* fQABckgHisto1RecCuts;  //! track QA: reconstructed tracks after cuts
539   AliFragFuncQATrackHistos* fQABckgHisto1Gen;      //! track QA: generated tracks
540   AliFragFuncQATrackHistos* fQABckgHisto2RecCuts;  //! track QA: reconstructed tracks after cuts
541   AliFragFuncQATrackHistos* fQABckgHisto2Gen;      //! track QA: generated tracks
542   AliFragFuncQATrackHistos* fQABckgHisto3RecCuts;  //! track QA: reconstructed tracks after cuts
543   AliFragFuncQATrackHistos* fQABckgHisto3Gen;      //! track QA: generated tracks
544   AliFragFuncQATrackHistos* fQABckgHisto4RecCuts;  //! track QA: reconstructed tracks after cuts
545   AliFragFuncQATrackHistos* fQABckgHisto4Gen;      //! track QA: generated tracks
546   
547   AliFragFuncHistos*  fFFBckgHisto0RecCuts;       //! Bckg (outside leading jet or 2 jets or more) FF reconstructed tracks after cuts 
548   AliFragFuncHistos*  fFFBckgHisto0Gen;           //! Bckg (outside leading jet or 2 jets or more) FF generated tracks after cuts 
549   AliFragFuncHistos*  fFFBckgHisto1RecCuts;       //! Bckg (outside leading jet or 2 jets or more) FF reconstructed tracks after cuts 
550   AliFragFuncHistos*  fFFBckgHisto1Gen;           //! Bckg (outside leading jet or 2 jets or more) FF generated tracks after cuts 
551   AliFragFuncHistos*  fFFBckgHisto2RecCuts;       //! Bckg (outside leading jet or 2 jets or more) FF reconstructed tracks after cuts 
552   AliFragFuncHistos*  fFFBckgHisto2Gen;           //! Bckg (outside leading jet or 2 jets or more) FF generated tracks after cuts 
553   AliFragFuncHistos*  fFFBckgHisto3RecCuts;       //! Bckg (outside leading jet or 3 jets or more) FF reconstructed tracks after cuts 
554   AliFragFuncHistos*  fFFBckgHisto3Gen;           //! Bckg (outside leading jet or 3 jets or more) FF generated tracks after cuts 
555   AliFragFuncHistos*  fFFBckgHisto4RecCuts;       //! Bckg (outside leading jet or 4 jets or more) FF reconstructed tracks after cuts 
556   AliFragFuncHistos*  fFFBckgHisto4Gen;           //! Bckg (outside leading jet or 4 jets or more) FF generated tracks after cuts 
557
558   AliFragFuncHistos*  fFFBckgHisto0RecEffRec;     //! Bckg (outside leading jet or 2 jets or more) FF reconstructed primaries after cuts 
559   AliFragFuncHistos*  fFFBckgHisto0SecRecNS;      //! secondary contamination: FF reconstructed secondaries (non-strangeness)
560   AliFragFuncHistos*  fFFBckgHisto0SecRecS;       //! secondary contamination: FF reconstructed secondaries (strange mothers)
561   AliFragFuncHistos*  fFFBckgHisto0SecRecSsc;     //! secondary contamination: FF reconstructed secondaries (strange mothers) - scale factor
562
563   TProfile* fProNtracksLeadingJet;          //! jet shape 
564   TProfile* fProDelR80pcPt;                 //! jet shape 
565   TProfile* fProDelRPtSum[5];               //! jet shape 
566
567   TProfile* fProNtracksLeadingJetGen;       //! jet shape 
568   TProfile* fProDelR80pcPtGen;              //! jet shape 
569   TProfile* fProDelRPtSumGen[5];            //! jet shape 
570
571   TProfile* fProNtracksLeadingJetBgrPerp2;  //! jet shape 
572   TProfile* fProDelRPtSumBgrPerp2[5];       //! jet shape 
573
574   TProfile* fProNtracksLeadingJetRecPrim;   //! jet shape 
575   TProfile* fProDelR80pcPtRecPrim;          //! jet shape 
576   TProfile* fProDelRPtSumRecPrim[5];        //! jet shape 
577
578   TProfile* fProNtracksLeadingJetRecSecNS;  //! jet shape 
579   TProfile* fProDelRPtSumRecSecNS[5];       //! jet shape 
580
581   TProfile* fProNtracksLeadingJetRecSecS;   //! jet shape 
582   TProfile* fProDelRPtSumRecSecS[5];        //! jet shape 
583
584   TProfile* fProNtracksLeadingJetRecSecSsc; //! jet shape 
585   TProfile* fProDelRPtSumRecSecSsc[5];      //! jet shape 
586   
587
588   TRandom3* fRandom;                        //! TRandom3 for background estimation 
589   
590   Bool_t fOnlyLeadingJets;                  // Flag indicating whether some histos are filled with leading jets only or all jets
591   
592   AliAnalysisUtils *fAnaUtils;              //! Object to use analysis utils like pile-up rejection
593   
594   // PID framework
595   Int_t fNumInclusivePIDtasks;              // Number of inclusive PID tasks used 
596   Int_t fNumJetPIDtasks;                    // Number of jet PID tasks used
597   
598   TString* fNameInclusivePIDtask;           //[fNumInclusivePIDtasks] Names of the tasks for inclusive PID spectra
599   TString* fNameJetPIDtask;                 //[fNumJetPIDtasks] Names of the tasks for jet PID spectra
600   
601   AliAnalysisTaskPID** fInclusivePIDtask;   //! Pointer to tasks for inclusive PID spectra
602   AliAnalysisTaskPID** fJetPIDtask;         //! Pointer to tasks for jet PID spectra
603   
604   Bool_t fUseInclusivePIDtask;              // Process inclusive PID spectra?
605   Bool_t fUseJetPIDtask;                    // Process jet PID spectra?
606   
607   Bool_t fIsPP;                             // Is pp collision system? -> If yes, centrality will be set to -1
608   
609   AliFragFuncHistos* fIDFFHistosRecCuts[AliPID::kSPECIES];    //! Identified FF reconstructed tracks after cuts 
610   AliFragFuncHistos* fIDFFHistosGen[AliPID::kSPECIES];    //! Identified FF generated tracks after cuts 
611
612   ClassDef(AliAnalysisTaskIDFragmentationFunction, 20);
613 };
614
615
616 inline void AliAnalysisTaskIDFragmentationFunction::SetNamesOfInclusivePIDtasks(Int_t numNames, const TString* names)
617 {
618   delete [] fNameInclusivePIDtask;
619   fNameInclusivePIDtask = 0x0;
620   
621   if (!names || numNames < 0) {
622     fNumInclusivePIDtasks = 0;
623     return;
624   }
625   
626   fNumInclusivePIDtasks = numNames;
627   
628   if (numNames > 0) {
629     fNameInclusivePIDtask = new TString[numNames];
630     
631     for (Int_t i = 0; i < numNames; i++) {
632       fNameInclusivePIDtask[i] = names[i];
633     }
634   }  
635 }
636
637 inline void AliAnalysisTaskIDFragmentationFunction::SetNamesOfJetPIDtasks(Int_t numNames, const TString* names)
638 {
639   delete [] fNameJetPIDtask;
640   fNameJetPIDtask = 0x0;
641   
642   if (!names || numNames < 0) {
643     fNumJetPIDtasks = 0;
644     return;
645   }
646   
647   fNumJetPIDtasks = numNames;
648   
649   if (numNames > 0) {
650     fNameJetPIDtask = new TString[numNames];
651     
652     for (Int_t i = 0; i < numNames; i++) {
653       fNameJetPIDtask[i] = names[i];
654     }
655   }  
656 }
657
658 #endif