]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGJE/UserTasks/AliAnalysisTaskIDFFTCF.h
Merge branch 'master' of https://git.cern.ch/reps/AliRoot
[u/mrichter/AliRoot.git] / PWGJE / UserTasks / AliAnalysisTaskIDFFTCF.h
1 // *************************************************************************
2 // * Task for Fragmentation Function Analysis in PWG4 Jet Task Force Train *
3 // *************************************************************************
4
5 #ifndef ALIANALYSISTASKIDFFTCFN_H
6 #define ALIANALYSISTASKIDFFTCFN_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 AliAODJets;
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 AliAODTrack;
26 class AliAODMCParticle;
27
28 #include "AliAnalysisTaskSE.h"
29 #include "TAxis.h"
30 #include "THnSparse.h"
31 #include <TTreeStream.h>
32   
33 class AliAnalysisTaskIDFFTCF : public AliAnalysisTaskSE {
34
35  public:
36   
37   //----------------------------------------
38   class AliFragFuncHistos : public TObject
39   {
40     
41     public:
42     
43     AliFragFuncHistos(const char* name = "FFhistos", 
44                       Int_t nJetPt = 0, Float_t jetPtMin = 0, Float_t jetPtMax = 0,
45                       Int_t nPt = 0, Float_t ptMin = 0, Float_t ptMax = 0,
46                       Int_t nXi = 0, Float_t xiMin = 0, Float_t xiMax = 0,
47                       Int_t nZ  = 0, Float_t zMin  = 0, Float_t zMax  = 0);
48     
49     AliFragFuncHistos(const AliFragFuncHistos& copy);
50     AliFragFuncHistos& operator=(const AliFragFuncHistos &o);
51     virtual ~AliFragFuncHistos();
52     
53     virtual void DefineHistos();
54     virtual void FillFF(Float_t trackPt, Float_t jetPt,Bool_t incrementJetPt, Float_t norm = 0, Bool_t scaleStrangeness = kFALSE, Float_t scaleFacStrangeness = 1.);
55
56     virtual void AddToOutput(TList* list) const;
57
58   private:
59
60     Int_t   fNBinsJetPt; // FF histos bins
61     Float_t fJetPtMin;   // FF histos limits
62     Float_t fJetPtMax;   // FF histos limits
63     Int_t   fNBinsPt;    // FF histos bins
64     Float_t fPtMin;      // FF histos limits
65     Float_t fPtMax;      // FF histos limits
66     Int_t   fNBinsXi;    // FF histos bins
67     Float_t fXiMin;      // FF histos limits
68     Float_t fXiMax;      // FF histos limits
69     Int_t   fNBinsZ;     // FF histos bins
70     Float_t fZMin;       // FF histos limits
71     Float_t fZMax;       // FF histos limits
72   
73     TH2F*   fh2TrackPt;   //! FF: track transverse momentum 
74     TH2F*   fh2Xi;        //! FF: xi 
75     TH2F*   fh2Z;         //! FF: z  
76     TH1F*   fh1JetPt;     //! jet pt 
77
78     TString fNameFF;      // histo names prefix
79     
80     ClassDef(AliFragFuncHistos, 1);
81   };
82   
83   //----------------------------------------
84   class AliFragFuncQAJetHistos : public TObject
85   {
86
87   public:
88  
89     AliFragFuncQAJetHistos(const char* name = "QAJethistos",
90                 Int_t nPt  = 0, Float_t ptMin  = 0, Float_t ptMax  = 0,
91                 Int_t nEta = 0, Float_t etaMin = 0, Float_t etaMax = 0,
92                 Int_t nPhi = 0, Float_t phiMin = 0, Float_t phiMax = 0);
93       
94     AliFragFuncQAJetHistos(const AliFragFuncQAJetHistos& copy);
95     AliFragFuncQAJetHistos& operator=(const AliFragFuncQAJetHistos &o);
96     virtual ~AliFragFuncQAJetHistos();
97     virtual void DefineHistos();
98     virtual void FillJetQA(Float_t eta, Float_t phi, Float_t pt);
99     virtual void AddToOutput(TList* list) const;
100
101   private:
102     
103     Int_t   fNBinsPt;    // jet QA histos bins
104     Float_t fPtMin;      // jet QA histos limits
105     Float_t fPtMax;      // jet QA histos limits
106     Int_t   fNBinsEta;   // jet QA histos bins
107     Float_t fEtaMin;     // jet QA histos limits
108     Float_t fEtaMax;     // jet QA histos limits
109     Int_t   fNBinsPhi;   // jet QA histos bins
110     Float_t fPhiMin;     // jet QA histos limits
111     Float_t fPhiMax;     // jet QA histos limits
112     
113     TH2F*   fh2EtaPhi;   //! jet phi vs eta 
114     TH1F*   fh1Pt;       //! jet transverse momentum 
115     TString fNameQAJ;    // histo names prefix
116     
117     ClassDef(AliFragFuncQAJetHistos, 1);
118   };
119   
120   //----------------------------------------
121   class AliFragFuncQATrackHistos : public TObject
122   {
123
124   public:
125
126     AliFragFuncQATrackHistos(const char* name = "QATrackhistos", 
127                   Int_t nPt  = 0, Float_t ptMin  = 0, Float_t ptMax  = 0,
128                   Int_t nEta = 0, Float_t etaMin = 0, Float_t etaMax = 0,
129                   Int_t nPhi = 0, Float_t phiMin = 0, Float_t phiMax = 0, 
130                   Float_t ptThresh = 0);
131     
132     AliFragFuncQATrackHistos(const AliFragFuncQATrackHistos& copy);
133     AliFragFuncQATrackHistos& operator=(const AliFragFuncQATrackHistos &o);
134     virtual ~AliFragFuncQATrackHistos();
135     virtual void DefineHistos();
136     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.);
137     virtual void AddToOutput(TList* list) const;
138
139   private:
140     
141     Int_t   fNBinsPt;    // track QA histos bins in pt
142     Float_t fPtMin;      // track QA histos limits in pt
143     Float_t fPtMax;      // track QA histos limits in pt
144     Int_t   fNBinsEta;   // track QA histos bins in eta
145     Float_t fEtaMin;     // track QA histos limits in eta
146     Float_t fEtaMax;     // track QA histos limits in eta
147     Int_t   fNBinsPhi;   // track QA histos bins in phi
148     Float_t fPhiMin;     // track QA histos limits in phi
149     Float_t fPhiMax;     // track QA histos limits in phi
150
151     Float_t fHighPtThreshold; //  high pt track phi vs eta distribution
152
153     TH2F*   fh2EtaPhi;        //! track phi vs eta 
154     TH1F*   fh1Pt;            //! track transverse momentum 
155     TH2F*   fh2HighPtEtaPhi;  //! phi vs eta for high pt (>fgHighPtThreshold) tracks
156     TH2F*   fh2PhiPt;         //! track phi vs pt
157
158     TString fNameQAT;         // histo names prefix
159     
160     ClassDef(AliFragFuncQATrackHistos, 2);
161   };
162   
163   enum TPCCUTMODE{
164     kPIDNone = 0, 
165     kPIDN,
166     kMIGeo
167   };
168   static Bool_t fkDump;  //=1: enable debug streamer; =0 : not.
169
170   AliAnalysisTaskIDFFTCF(); 
171   AliAnalysisTaskIDFFTCF(const char *name);
172   AliAnalysisTaskIDFFTCF(const  AliAnalysisTaskIDFFTCF &copy);
173   AliAnalysisTaskIDFFTCF& operator=(const  AliAnalysisTaskIDFFTCF &o);
174   virtual ~AliAnalysisTaskIDFFTCF();
175   
176   virtual void   UserCreateOutputObjects();
177   virtual void   Init();
178   virtual void   LocalInit() {Init();}
179
180   virtual void   UserExec(Option_t *option);
181   virtual void   Terminate(Option_t* );
182   virtual Bool_t Notify();
183
184   virtual void   SetNonStdFile(char* c){fNonStdFile = c;} 
185
186   virtual void   SetTrackTypeGen(Int_t i){fTrackTypeGen = i;}
187   virtual void   SetJetTypeGen(Int_t i){fJetTypeGen = i;}
188   virtual void   SetJetTypeRecEff(Int_t i){fJetTypeRecEff = i;}
189
190   virtual void   SetBranchGenJets(const char* c){fBranchGenJets = c;}
191   virtual void   SetBranchRecJets(const char* c){fBranchRecJets = c;}
192
193   virtual void   SetTrackCuts(Float_t trackPt = 0.15, Float_t trackEtaMin = -0.9, Float_t trackEtaMax = 0.9, 
194                               Float_t trackPhiMin = 0., Float_t trackPhiMax = 2*TMath::Pi())
195   {fTrackPtCut = trackPt; fTrackEtaMin = trackEtaMin; fTrackEtaMax = trackEtaMax; 
196     fTrackPhiMin = trackPhiMin; fTrackPhiMax = trackPhiMax;}
197
198   virtual void   UseAODInputJets(Bool_t b) {fUseAODInputJets = b;}  
199   virtual void   SetFilterMask(UInt_t i) {fFilterMask = i;}
200   virtual void   UsePhysicsSelection(Bool_t b) {fUsePhysicsSelection = b;}
201   virtual void   SetEventSelectionMask(UInt_t i){fEvtSelectionMask = i;}
202   virtual void   SetEventClass(Int_t i){fEventClass = i;}
203   virtual void   SetMaxVertexZ(Float_t z){fMaxVertexZ = z;}
204   virtual void   UseLeadingJet(Bool_t b){fLeadingJets = b;}
205
206   virtual void   SetJetCuts(Float_t jetPt = 5., Float_t jetEtaMin = -0.5, Float_t jetEtaMax = 0.5, 
207                             Float_t jetPhiMin = 0., Float_t jetPhiMax = 2*TMath::Pi())
208   {fJetPtCut = jetPt; fJetEtaMin = jetEtaMin; fJetEtaMax = jetEtaMax; 
209     fJetPhiMin = jetPhiMin; fJetPhiMax = jetPhiMax;}
210
211   virtual void   SetFFRadius(Float_t r = 0.4) { fFFRadius = r; }
212   virtual void   SetFFMinLTrackPt(Float_t pt = -1) { fFFMinLTrackPt = pt; }
213   virtual void   SetFFMaxTrackPt(Float_t pt = -1) { fFFMaxTrackPt = pt; }
214   virtual void   SetFFMinNTracks(Int_t nTracks = 0) { fFFMinnTracks = nTracks; }
215   virtual void   SetQAMode(Int_t qa = 3)      {fQAMode = qa;}
216   virtual void   SetFFMode(Int_t ff = 1)      {fFFMode = ff;}
217   virtual void   SetEffMode(Int_t eff = 1)    {fEffMode = eff;}
218
219   static  void   SetProperties(TH1* h,const char* x, const char* y);
220   static  void   SetProperties(TH1* h,const char* x, const char* y,const char* z);
221   static  void   SetProperties(THnSparse* h, Int_t dim, const char** labels);
222
223   void SetTPCCutMode(Int_t mode){ fTPCCutMode = mode; }
224   Int_t GetTPCCutMode(){return fTPCCutMode; }
225
226   void SetTOFCutMode(Int_t mode){ fTOFCutMode = mode; }
227   Int_t GetTOFCutMode(){return fTOFCutMode; }
228
229   void   SetHighPtThreshold(Float_t pt = 5.) { fQATrackHighPtThreshold = pt; }
230
231   void   SetFFHistoBins(Int_t nJetPt = 245, Float_t jetPtMin = 5, Float_t jetPtMax = 250, 
232                         Int_t nPt = 200, Float_t ptMin = 0., Float_t ptMax = 200., 
233                         Int_t nXi = 70, Float_t xiMin = 0., Float_t xiMax = 7.,
234                         Int_t nZ = 22,  Float_t zMin = 0.,  Float_t zMax = 1.1)
235   { fFFNBinsJetPt = nJetPt; fFFJetPtMin = jetPtMin; fFFJetPtMax = jetPtMax; 
236     fFFNBinsPt = nPt; fFFPtMin = ptMin; fFFPtMax = ptMax;
237     fFFNBinsXi = nXi; fFFXiMin = xiMin; fFFXiMax = xiMax;
238     fFFNBinsZ  = nZ;  fFFZMin  = zMin;  fFFZMax  = zMax; }
239
240   void  SetQAJetHistoBins(Int_t nPt = 300, Float_t ptMin = 0., Float_t ptMax = 300.,
241                           Int_t nEta = 20, Float_t etaMin = -1.0, Float_t etaMax = 1.0,
242                           Int_t nPhi = 60, Float_t phiMin = 0., Float_t phiMax = 2*TMath::Pi())
243     { fQAJetNBinsPt = nPt; fQAJetPtMin = ptMin; fQAJetPtMax = ptMax;
244       fQAJetNBinsEta = nEta; fQAJetEtaMin = etaMin; fQAJetEtaMax = etaMax;
245       fQAJetNBinsPhi = nPhi; fQAJetPhiMin = phiMin; fQAJetPhiMax = phiMax; }
246   
247   void  SetQATrackHistoBins(Int_t nPt = 200, Float_t ptMin = 0., Float_t ptMax = 200.,
248                             Int_t nEta = 20, Float_t etaMin = -1.0, Float_t etaMax = 1.0,
249                             Int_t nPhi = 60, Float_t phiMin = 0., Float_t phiMax = 2*TMath::Pi())
250   { fQATrackNBinsPt = nPt; fQATrackPtMin = ptMin; fQATrackPtMax = ptMax;
251     fQATrackNBinsEta = nEta; fQATrackEtaMin = etaMin; fQATrackEtaMax = etaMax;
252     fQATrackNBinsPhi = nPhi; fQATrackPhiMin = phiMin; fQATrackPhiMax = phiMax; }
253   
254
255   Float_t  GetFFRadius() const { return fFFRadius; }
256   Float_t  GetFFMinLTrackPt() const { return fFFMinLTrackPt; }
257   Float_t  GetFFMaxTrackPt() const { return fFFMaxTrackPt; }
258   Float_t  GetFFMinNTracks() const { return fFFMinnTracks; }
259
260   void     GetJetTracksTrackrefs(TList* l, const AliAODJet* j, const Double_t minPtL, Double_t maxPt, Bool_t& isBadPt);
261   void     GetJetTracksPointing(TList* in, TList* out, const AliAODJet* j, Double_t r, Double_t& sumPt, Double_t minPtL, Double_t maxPt, Bool_t& isBadPt);  
262
263   void     AssociateGenRec(TList* tracksAODMCCharged,TList* tracksRec, TArrayI& indexAODTr,TArrayI& indexMCTr,TArrayS& isRefGen,TH2F* fh2PtRecVsGen);
264
265   void     FillSingleTrackHistosRecGen(AliFragFuncQATrackHistos* trackQAGen, AliFragFuncQATrackHistos* trackQARec, TList* tracksGen, 
266                                        const TArrayI& indexAODTr, const TArrayS& isRefGen, Int_t pdg = 0, Bool_t scaleGFL = kFALSE);
267
268
269   void     FillJetTrackHistosRec(AliFragFuncHistos* histRec,  AliAODJet* jet, 
270                                  TList* jetTrackList, const TList* tracksGen, const TList* tracksRec, const TArrayI& indexAODTr,
271                                  const TArrayS& isRefGen, TList* jetTrackListTR = 0, Int_t pdg = 0, Bool_t scaleGFL = kFALSE);
272
273
274   Float_t  CalcJetArea(Float_t etaJet, Float_t rc) const;
275  
276   void     BookQAHistos(TList* list = 0, AliFragFuncQATrackHistos** rec = 0, TString strTitRec = "", AliFragFuncQATrackHistos** gen = 0, TString strTitGen = "",
277                         AliFragFuncQATrackHistos** sec = 0, TString strTitSec = "");
278
279   void     BookFFHistos(TList* list, AliFragFuncHistos** rec = 0, TString strTitRec = "", AliFragFuncHistos** gen = 0, TString strTitGen = "",
280                         AliFragFuncHistos** sec = 0, TString strTitSec = "");
281
282   Double_t  TrackingPtGeantFlukaCorrectionPrMinus(Double_t pTmc);
283   Double_t  TrackingPtGeantFlukaCorrectionKaMinus(Double_t pTmc);
284     
285
286   // Consts
287   enum {kTrackUndef=0, kTrackAOD, kTrackAODQualityCuts, kTrackAODCuts,  
288         kTrackKineAll, kTrackKineCharged, kTrackKineChargedAcceptance, 
289         kTrackAODMCAll, kTrackAODMCCharged, kTrackAODMCChargedAcceptance, kTrackAODMCChargedSec, kTrackAOCMCChargedPrimAcceptance};
290   enum {kJetsUndef=0, kJetsRec, kJetsRecAcceptance, kJetsGen, kJetsGenAcceptance, kJetsKine, kJetsKineAcceptance};
291
292  
293  protected:
294   
295   Int_t   GetListOfTracks(TList* list, Int_t type);
296   Int_t   GetListOfJets(TList* list, Int_t type);
297
298   AliESDEvent* fESD;      // ESD event
299   AliAODEvent* fAOD;      // AOD event
300   AliAODEvent* fAODJets;  // AOD event with jet branch (case we have AOD both in input and output)
301   AliAODExtension  *fAODExtension; //! where we take the jets from can be input or output AOD
302   TString       fNonStdFile; // name of delta aod file to catch the extension
303  
304  
305   TString fBranchRecJets;         // branch name for reconstructed jets
306   TString fBranchGenJets;         // branch name for generated jets
307
308   Int_t   fTrackTypeGen;        // type of generated tracks
309   Int_t   fJetTypeGen;          // type of generated jets
310
311   Int_t   fJetTypeRecEff;       // type of jets used for filling reconstruction efficiency histos
312
313   Bool_t  fUseAODInputJets;     // take jets from in/output - only relevant if AOD event both in input AND output and we want to use output
314   UInt_t  fFilterMask;          // filter bit for selected tracks
315   Bool_t  fUsePhysicsSelection; // switch for event selection
316   UInt_t  fEvtSelectionMask;    // trigger class selection
317   Int_t   fEventClass;          // centrality class selection
318   Float_t fMaxVertexZ;          // maximum abs(z) position of primiary vertex [cm]
319   Bool_t  fLeadingJets;         // leading/all jets
320
321
322   Int_t fTPCCutMode;      //mode for cutting TPC for good dE/dx
323   Int_t fTOFCutMode;      //mode for cutting TOF
324   TTreeStream * fStream; //debug streamer
325   TTree * fTree;         //tree of streamer
326
327   // track cuts
328   Float_t fTrackPtCut;    // track transverse momentum cut
329   Float_t fTrackEtaMin;   // track eta cut
330   Float_t fTrackEtaMax;   // track eta cut
331   Float_t fTrackPhiMin;   // track phi cut
332   Float_t fTrackPhiMax;   // track phi cut
333   
334
335   // jet cuts
336   Float_t fJetPtCut;      // jet transverse momentum cut
337   Float_t fJetEtaMin;     // jet eta cut
338   Float_t fJetEtaMax;     // jet eta cut
339   Float_t fJetPhiMin;     // jet phi cut
340   Float_t fJetPhiMax;     // jet phi cut
341
342   Float_t fFFRadius;        // if radius > 0 construct FF from tracks within cone around jet axis, otherwise use trackRefs  
343   Float_t fFFMinLTrackPt;   // reject jets with leading track with pt smaller than this value
344   Float_t fFFMaxTrackPt;    // reject jets containing any track with pt larger than this value
345   Int_t   fFFMinnTracks;    // reject jets with less tracks than this value
346   Int_t   fQAMode;          // QA mode: 0x00=0 none, 0x01=1 track qa, 0x10=2 track qa, 0x11=3 both
347   Int_t   fFFMode;          // fragmentation function mode
348   Int_t   fEffMode;         // efficiency mode
349
350   Float_t fAvgTrials;       // average number of trials per event
351   
352   TList* fTracksRecCuts;           //! reconstructed tracks after cuts
353   TList* fTracksGen;               //! generated tracks 
354   TList* fTracksAODMCCharged;      //! AOD MC tracks 
355   TList* fTracksAODMCChargedSec;   //! AOD MC tracks - secondaries 
356   TList* fTracksRecQualityCuts;    //! reconstructed tracks after quality cuts, no acceptance/pt cut
357
358   TList* fJetsRec;        //! jets from reconstructed tracks
359   TList* fJetsRecCuts;    //! jets from reonstructed tracks after jet cuts 
360   TList* fJetsGen;        //! jets from generated tracks
361   TList* fJetsRecEff;     //! jets used for reconstruction efficiency histos 
362
363    
364   AliFragFuncQATrackHistos* fQATrackHistosRecCuts;  //! track QA: reconstructed tracks after cuts
365   AliFragFuncQATrackHistos* fQATrackHistosGen;      //! track QA: generated tracks
366   
367   AliFragFuncQAJetHistos*  fQAJetHistosRec;             //! jet QA: jets from reconstructed tracks
368   AliFragFuncQAJetHistos*  fQAJetHistosRecCuts;         //! jet QA: jets from reconstructed tracks after jet cuts 
369   AliFragFuncQAJetHistos*  fQAJetHistosRecCutsLeading;  //! jet QA: leading jet from reconstructed tracks after jet cuts 
370   AliFragFuncQAJetHistos*  fQAJetHistosGen;             //! jet QA: jets from generated tracks  
371   AliFragFuncQAJetHistos*  fQAJetHistosGenLeading;      //! jet QA: leading jet from generated tracks  
372   AliFragFuncQAJetHistos*  fQAJetHistosRecEffLeading;   //! jet QA: leading jet used for reconstruction efficiency histos  
373   
374
375   AliFragFuncHistos*  fFFHistosRecCutsInc;       //! inclusive FF (all jets) 
376   AliFragFuncHistos*  fFFHistosRecCutsIncPi;     //! inclusive FF (all jets) 
377   AliFragFuncHistos*  fFFHistosRecCutsIncPro;    //! inclusive FF (all jets) 
378   AliFragFuncHistos*  fFFHistosRecCutsIncK;      //! inclusive FF (all jets) 
379   AliFragFuncHistos*  fFFHistosRecCutsIncEl;     //! inclusive FF (all jets) 
380   AliFragFuncHistos*  fFFHistosRecCutsIncMu;     //! inclusive FF (all jets) 
381
382   AliFragFuncHistos*  fFFHistosRecLeadingTrack; //! FF reconstructed tracks after cuts: leading track pt / jet pt (all jets)
383
384   AliFragFuncHistos*  fFFHistosGenInc;       //! inclusive FF (all jets) 
385   AliFragFuncHistos*  fFFHistosGenIncPi;     //! inclusive FF (all jets) 
386   AliFragFuncHistos*  fFFHistosGenIncPro;    //! inclusive FF (all jets) 
387   AliFragFuncHistos*  fFFHistosGenIncK;      //! inclusive FF (all jets) 
388   AliFragFuncHistos*  fFFHistosGenIncEl;     //! inclusive FF (all jets) 
389   AliFragFuncHistos*  fFFHistosGenIncMu;     //! inclusive FF (all jets) 
390   AliFragFuncHistos*  fFFHistosGenLeadingTrack; //! FF reconstructed tracks after cuts: leading track pt / jet pt (all jets)
391
392   Float_t  fQATrackHighPtThreshold;       // track QA high transverse momentum threshold
393   
394   THnSparseD * fTHnIDFF;                //! tracks in jets
395   THnSparseD * fTHnIncl;                //! inclusive tracks
396
397   // histogram bins  
398
399   Int_t   fFFNBinsJetPt;    // FF histos bins
400   Float_t fFFJetPtMin;      // FF histos limits
401   Float_t fFFJetPtMax;      // FF histos limits
402
403   Int_t   fFFNBinsPt;       // FF histos bins
404   Float_t fFFPtMin;         // FF histos limits
405   Float_t fFFPtMax;         // FF histos limits
406
407   Int_t   fFFNBinsXi;       // FF histos bins
408   Float_t fFFXiMin;         // FF histos limits
409   Float_t fFFXiMax;         // FF histos limits
410
411   Int_t   fFFNBinsZ;        // FF histos bins
412   Float_t fFFZMin;          // FF histos limits
413   Float_t fFFZMax;          // FF histos limits
414
415   Int_t   fQAJetNBinsPt;    // jet QA histos bins
416   Float_t fQAJetPtMin;      // jet QA histos limits
417   Float_t fQAJetPtMax;      // jet QA histos limits
418   
419   Int_t   fQAJetNBinsEta;   // jet QA histos bins
420   Float_t fQAJetEtaMin;     // jet QA histos limits
421   Float_t fQAJetEtaMax;     // jet QA histos limits
422   
423   Int_t   fQAJetNBinsPhi;   // jet QA histos bins
424   Float_t fQAJetPhiMin;     // jet QA histos limits
425   Float_t fQAJetPhiMax;     // jet QA histos limits
426   
427   Int_t   fQATrackNBinsPt;  // track QA histos bins
428   Float_t fQATrackPtMin;    // track QA histos limits
429   Float_t fQATrackPtMax;    // track QA histos limits
430   
431   Int_t   fQATrackNBinsEta; // track QA histos bins
432   Float_t fQATrackEtaMin;   // track QA histos limits
433   Float_t fQATrackEtaMax;   // track QA histos limits
434   
435   Int_t   fQATrackNBinsPhi; // track QA histos bins
436   Float_t fQATrackPhiMin;   // track QA histos limits
437   Float_t fQATrackPhiMax;   // track QA histos limits
438   
439   // Histograms
440   TList *fCommonHistList;         // List of common histos
441   
442   TH1F  *fh1EvtSelection;         //! event cuts 
443   TH1F  *fh1VertexNContributors;  //! NContributors to prim vertex
444   TH1F  *fh1VertexZ;              //! prim vertex z distribution
445   TH1F  *fh1EvtMult;              //! number of reconstructed tracks after cuts 
446   TH1F  *fh1EvtCent;              //! centrality percentile 
447
448   TProfile* fh1Xsec;              //! pythia cross section and trials
449   TH1F*     fh1Trials;            //! sum of trials
450   TH1F*     fh1PtHard;            //! pt hard of the event
451   TH1F*     fh1PtHardTrials;      //! pt hard of the event
452
453   TH1F  *fh1nRecJetsCuts;         //! number of jets from reconstructed tracks per event 
454   TH1F  *fh1nGenJets;             //! number of jets from generated tracks per event
455   TH1F  *fh1nRecEffJets;          //! number of jets for reconstruction eff per event
456
457   TH2F  *fh2PtRecVsGenPrim;       //! association rec/gen MC: rec vs gen pt, primaries 
458   TH2F  *fh2PtRecVsGenSec;        //! association rec/gen MC: rec vs gen pt, secondaries 
459
460   // tracking efficiency / secondaries
461   
462   AliFragFuncQATrackHistos* fQATrackHistosRecEffGen;      //! tracking efficiency: generated primaries 
463   AliFragFuncQATrackHistos* fQATrackHistosRecEffRec;      //! tracking efficiency: reconstructed primaries
464   AliFragFuncQATrackHistos* fQATrackHistosSecRec;         //! reconstructed secondaries
465
466   AliFragFuncQATrackHistos* fQATrackHistosRecEffGenPi;     //! tracking efficiency: generated primaries 
467   AliFragFuncQATrackHistos* fQATrackHistosRecEffGenPro;    //! tracking efficiency: generated primaries 
468   AliFragFuncQATrackHistos* fQATrackHistosRecEffGenK;      //! tracking efficiency: generated primaries 
469   AliFragFuncQATrackHistos* fQATrackHistosRecEffGenEl;     //! tracking efficiency: generated primaries 
470   AliFragFuncQATrackHistos* fQATrackHistosRecEffGenMu;     //! tracking efficiency: generated primaries 
471
472   AliFragFuncQATrackHistos* fQATrackHistosRecEffRecPi;       //! tracking efficiency: generated primaries 
473   AliFragFuncQATrackHistos* fQATrackHistosRecEffRecPro;      //! tracking efficiency: generated primaries 
474   AliFragFuncQATrackHistos* fQATrackHistosRecEffRecK;        //! tracking efficiency: generated primaries 
475   AliFragFuncQATrackHistos* fQATrackHistosRecEffRecEl;       //! tracking efficiency: generated primaries 
476   AliFragFuncQATrackHistos* fQATrackHistosRecEffRecMu;       //! tracking efficiency: generated primaries 
477   AliFragFuncQATrackHistos* fQATrackHistosRecEffRecProGFL;   //! tracking efficiency: generated primaries 
478   AliFragFuncQATrackHistos* fQATrackHistosRecEffRecKGFL;     //! tracking efficiency: generated primaries 
479
480   AliFragFuncQATrackHistos* fQATrackHistosSecRecPi;       //! tracking efficiency: generated primaries 
481   AliFragFuncQATrackHistos* fQATrackHistosSecRecPro;      //! tracking efficiency: generated primaries 
482   AliFragFuncQATrackHistos* fQATrackHistosSecRecK;        //! tracking efficiency: generated primaries 
483   AliFragFuncQATrackHistos* fQATrackHistosSecRecEl;       //! tracking efficiency: generated primaries 
484   AliFragFuncQATrackHistos* fQATrackHistosSecRecMu;       //! tracking efficiency: generated primaries 
485   AliFragFuncQATrackHistos* fQATrackHistosSecRecProGFL;   //! tracking efficiency: generated primaries 
486   AliFragFuncQATrackHistos* fQATrackHistosSecRecKGFL;     //! tracking efficiency: generated primaries 
487
488
489
490   AliFragFuncHistos*  fFFHistosRecEffRec;                 //! tracking efficiency: FF reconstructed primaries
491   AliFragFuncHistos*  fFFHistosSecRec;                    //! secondary contamination: FF reconstructed secondaries 
492
493   AliFragFuncHistos*  fFFHistosRecEffRecPi;               //! tracking efficiency: FF reconstructed primaries
494   AliFragFuncHistos*  fFFHistosRecEffRecPro;              //! tracking efficiency: FF reconstructed primaries
495   AliFragFuncHistos*  fFFHistosRecEffRecK;                //! tracking efficiency: FF reconstructed primaries
496   AliFragFuncHistos*  fFFHistosRecEffRecEl;               //! tracking efficiency: FF reconstructed primaries
497   AliFragFuncHistos*  fFFHistosRecEffRecMu;               //! tracking efficiency: FF reconstructed primaries
498   AliFragFuncHistos*  fFFHistosRecEffRecProGFL;           //! tracking efficiency: FF reconstructed primaries
499   AliFragFuncHistos*  fFFHistosRecEffRecKGFL;             //! tracking efficiency: FF reconstructed primaries
500
501
502   AliFragFuncHistos*  fFFHistosSecRecPi;                  //! secondary contamination: FF reconstructed secondaries 
503   AliFragFuncHistos*  fFFHistosSecRecPro;                 //! secondary contamination: FF reconstructed secondaries 
504   AliFragFuncHistos*  fFFHistosSecRecK;                   //! secondary contamination: FF reconstructed secondaries 
505   AliFragFuncHistos*  fFFHistosSecRecEl;                  //! secondary contamination: FF reconstructed secondaries 
506   AliFragFuncHistos*  fFFHistosSecRecMu;                  //! secondary contamination: FF reconstructed secondaries 
507   AliFragFuncHistos*  fFFHistosSecRecProGFL;              //! secondary contamination: FF reconstructed secondaries 
508   AliFragFuncHistos*  fFFHistosSecRecKGFL;                //! secondary contamination: FF reconstructed secondaries 
509
510
511
512   TRandom3*                   fRandom;          // TRandom3 for background estimation 
513
514   ClassDef(AliAnalysisTaskIDFFTCF, 1);
515 };
516
517 #endif