]> git.uio.no Git - u/mrichter/AliRoot.git/blob - CORRFW/AliCFTrackQualityCuts.cxx
updated merging function (Michele)
[u/mrichter/AliRoot.git] / CORRFW / AliCFTrackQualityCuts.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
7  * Permission to use, copy, modify and distribute this software and its   *
8  * documentation strictly for non-commercial purposes is hereby granted   *
9  * without fee, provided that the above copyright notice appears in all   *
10  * copies and that both the copyright notice and this permission notice   *
11  * appear in the supporting documentation. The authors make no claims     *
12  * about the suitability of this software for any purpose. It is          *
13  * provided "as is" without express or implied warranty.                  *
14  **************************************************************************/
15
16 // The class AliCFTrackQualityCuts is designed to select reconstructed tracks
17 // of high quality and to provide corresponding QA histograms.
18 // This class inherits from the Analysis' Framework abstract base class
19 // AliAnalysisCuts and is a part of the Correction Framework.
20 // This class acts on single, reconstructed tracks, it is applicable on
21 // ESD and AOD data.
22 // It mainly consists of a IsSelected function that returns a boolean.
23 // This function checks whether the considered track passes a set of cuts:
24 // - number of clusters in the ITS
25 // - number of clusters in the TPC
26 // - number of clusters in the TRD
27 // - ratio of found / finable number of clusters in the TPC
28 // - number of tracklets in the TRD
29 // - number TRD tracklets used for pid
30 // - chi2 / cluster in the ITS
31 // - chi2 / cluster in the TPC
32 // - chi2 / tracklet in the TRD
33 // - number of clusters in the TPC used for dEdx calculation
34 // - successful TPC refit
35 // - successful ITS refit
36 // - covariance matrix diagonal elements
37 //
38 // The cut values for these cuts are set with the corresponding set functions.
39 // All cut classes provided by the correction framework are supposed to be
40 // added in the Analysis Framwork's class AliAnalysisFilter and applied by
41 // the filter via a loop.
42 //
43 // author: I. Kraus (Ingrid.Kraus@cern.ch)
44 // idea taken form
45 // AliESDtrackCuts writte by Jan Fiete Grosse-Oetringhaus and
46 // AliRsnDaughterCut class written by A. Pulvirenti.
47
48
49 #include <TCanvas.h>
50 #include <TDirectory.h>
51 #include <TH2.h>
52 #include <TBits.h>
53
54 #include <AliESDtrack.h>
55 #include <AliESDtrackCuts.h>
56 #include <AliLog.h>
57 #include "AliCFTrackQualityCuts.h"
58 #include "AliAODTrack.h"
59
60 ClassImp(AliCFTrackQualityCuts)
61
62 //__________________________________________________________________________________
63 AliCFTrackQualityCuts::AliCFTrackQualityCuts() :
64   AliCFCutBase(),
65   fMinNClusterTPC(0),
66   fMinNClusterITS(0),
67   fMinNClusterTRD(0),
68   fMinFoundClusterTPC(0),
69   fMinNTrackletTRD(0),
70   fMinNTrackletTRDpid(0),
71   fMaxChi2PerClusterTPC(0),
72   fMaxChi2PerClusterITS(0),
73   fMaxChi2PerTrackletTRD(0),
74   fMinNdEdxClusterTPC(0),
75   fCovariance11Max(0),
76   fCovariance22Max(0),
77   fCovariance33Max(0),
78   fCovariance44Max(0),
79   fCovariance55Max(0),
80   fStatus(0),
81   fhCutStatistics(0),
82   fhCutCorrelation(0),
83   fBitmap(0x0),
84   fTrackCuts(0x0),
85   fhNBinsClusterTPC(0),
86   fhNBinsClusterITS(0),
87   fhNBinsClusterTRD(0),
88   fhNBinsFoundClusterTPC(0),
89   fhNBinsTrackletTRD(0),
90   fhNBinsTrackletTRDpid(0),
91   fhNBinsChi2TPC(0),
92   fhNBinsChi2ITS(0),
93   fhNBinsChi2TRD(0),
94   fhNBinsdEdxClusterTPC(0),
95   fhNBinsCovariance11(0),
96   fhNBinsCovariance22(0),
97   fhNBinsCovariance33(0),
98   fhNBinsCovariance44(0),
99   fhNBinsCovariance55(0),
100   fhBinLimClusterTPC(0x0),
101   fhBinLimClusterITS(0x0),
102   fhBinLimClusterTRD(0x0),
103   fhBinLimFoundClusterTPC(0x0),
104   fhBinLimTrackletTRD(0x0),
105   fhBinLimTrackletTRDpid(0x0),
106   fhBinLimChi2TPC(0x0),
107   fhBinLimChi2ITS(0x0),
108   fhBinLimChi2TRD(0x0),
109   fhBinLimdEdxClusterTPC(0x0),
110   fhBinLimCovariance11(0x0),
111   fhBinLimCovariance22(0x0),
112   fhBinLimCovariance33(0x0),
113   fhBinLimCovariance44(0x0),
114   fhBinLimCovariance55(0x0)
115 {
116   //
117   // Default constructor
118   //
119   Initialise();
120 }
121 //__________________________________________________________________________________
122 AliCFTrackQualityCuts::AliCFTrackQualityCuts(Char_t* name, Char_t* title) :
123   AliCFCutBase(name,title),
124   fMinNClusterTPC(0),
125   fMinNClusterITS(0),
126   fMinNClusterTRD(0),
127   fMinFoundClusterTPC(0),
128   fMinNTrackletTRD(0),
129   fMinNTrackletTRDpid(0),
130   fMaxChi2PerClusterTPC(0),
131   fMaxChi2PerClusterITS(0),
132   fMaxChi2PerTrackletTRD(0),
133   fMinNdEdxClusterTPC(0),
134   fCovariance11Max(0),
135   fCovariance22Max(0),
136   fCovariance33Max(0),
137   fCovariance44Max(0),
138   fCovariance55Max(0),
139   fStatus(0),
140   fhCutStatistics(0),
141   fhCutCorrelation(0),
142   fBitmap(0x0),
143   fTrackCuts(0x0),
144   fhNBinsClusterTPC(0),
145   fhNBinsClusterITS(0),
146   fhNBinsClusterTRD(0),
147   fhNBinsFoundClusterTPC(0),
148   fhNBinsTrackletTRD(0),
149   fhNBinsTrackletTRDpid(0),
150   fhNBinsChi2TPC(0),
151   fhNBinsChi2ITS(0),
152   fhNBinsChi2TRD(0),
153   fhNBinsdEdxClusterTPC(0),
154   fhNBinsCovariance11(0),
155   fhNBinsCovariance22(0),
156   fhNBinsCovariance33(0),
157   fhNBinsCovariance44(0),
158   fhNBinsCovariance55(0),
159   fhBinLimClusterTPC(0x0),
160   fhBinLimClusterITS(0x0),
161   fhBinLimClusterTRD(0x0),
162   fhBinLimFoundClusterTPC(0x0),
163   fhBinLimTrackletTRD(0x0),
164   fhBinLimTrackletTRDpid(0x0),
165   fhBinLimChi2TPC(0x0),
166   fhBinLimChi2ITS(0x0),
167   fhBinLimChi2TRD(0x0),
168   fhBinLimdEdxClusterTPC(0x0),
169   fhBinLimCovariance11(0x0),
170   fhBinLimCovariance22(0x0),
171   fhBinLimCovariance33(0x0),
172   fhBinLimCovariance44(0x0),
173   fhBinLimCovariance55(0x0)
174 {
175   //
176   // Constructor
177   //
178   Initialise();
179 }
180 //__________________________________________________________________________________
181 AliCFTrackQualityCuts::AliCFTrackQualityCuts(const AliCFTrackQualityCuts& c) :
182   AliCFCutBase(c),
183   fMinNClusterTPC(c.fMinNClusterTPC),
184   fMinNClusterITS(c.fMinNClusterITS),
185   fMinNClusterTRD(c.fMinNClusterTRD),
186   fMinFoundClusterTPC(c.fMinFoundClusterTPC),
187   fMinNTrackletTRD(c.fMinNTrackletTRD),
188   fMinNTrackletTRDpid(c.fMinNTrackletTRDpid),
189   fMaxChi2PerClusterTPC(c.fMaxChi2PerClusterTPC),
190   fMaxChi2PerClusterITS(c.fMaxChi2PerClusterITS),
191   fMaxChi2PerTrackletTRD(c.fMaxChi2PerTrackletTRD),
192   fMinNdEdxClusterTPC(c.fMinNdEdxClusterTPC),
193   fCovariance11Max(c.fCovariance11Max),
194   fCovariance22Max(c.fCovariance22Max),
195   fCovariance33Max(c.fCovariance33Max),
196   fCovariance44Max(c.fCovariance44Max),
197   fCovariance55Max(c.fCovariance55Max),
198   fStatus(c.fStatus),
199   fhCutStatistics(c.fhCutStatistics),
200   fhCutCorrelation(c.fhCutCorrelation),
201   fBitmap(c.fBitmap),
202   fTrackCuts(c.fTrackCuts),
203   fhNBinsClusterTPC(c.fhNBinsClusterTPC),
204   fhNBinsClusterITS(c.fhNBinsClusterITS),
205   fhNBinsClusterTRD(c.fhNBinsClusterTRD),
206   fhNBinsFoundClusterTPC(c.fhNBinsFoundClusterTPC),
207   fhNBinsTrackletTRD(c.fhNBinsTrackletTRD),
208   fhNBinsTrackletTRDpid(c.fhNBinsTrackletTRDpid),
209   fhNBinsChi2TPC(c.fhNBinsChi2TPC),
210   fhNBinsChi2ITS(c.fhNBinsChi2ITS),
211   fhNBinsChi2TRD(c.fhNBinsChi2TRD),
212   fhNBinsdEdxClusterTPC(c.fhNBinsdEdxClusterTPC),
213   fhNBinsCovariance11(c.fhNBinsCovariance11),
214   fhNBinsCovariance22(c.fhNBinsCovariance22),
215   fhNBinsCovariance33(c.fhNBinsCovariance33),
216   fhNBinsCovariance44(c.fhNBinsCovariance44),
217   fhNBinsCovariance55(c.fhNBinsCovariance55),
218   fhBinLimClusterTPC(c.fhBinLimClusterTPC),
219   fhBinLimClusterITS(c.fhBinLimClusterITS),
220   fhBinLimClusterTRD(c.fhBinLimClusterTRD),
221   fhBinLimFoundClusterTPC(c.fhBinLimFoundClusterTPC),
222   fhBinLimTrackletTRD(c.fhBinLimTrackletTRD),
223   fhBinLimTrackletTRDpid(c.fhBinLimTrackletTRDpid),
224   fhBinLimChi2TPC(c.fhBinLimChi2TPC),
225   fhBinLimChi2ITS(c.fhBinLimChi2ITS),
226   fhBinLimChi2TRD(c.fhBinLimChi2TRD),
227   fhBinLimdEdxClusterTPC(c.fhBinLimdEdxClusterTPC),
228   fhBinLimCovariance11(c.fhBinLimCovariance11),
229   fhBinLimCovariance22(c.fhBinLimCovariance22),
230   fhBinLimCovariance33(c.fhBinLimCovariance33),
231   fhBinLimCovariance44(c.fhBinLimCovariance44),
232   fhBinLimCovariance55(c.fhBinLimCovariance55)
233 {
234   //
235   // copy constructor
236   //
237   ((AliCFTrackQualityCuts &) c).Copy(*this);
238 }
239 //__________________________________________________________________________________
240 AliCFTrackQualityCuts& AliCFTrackQualityCuts::operator=(const AliCFTrackQualityCuts& c)
241 {
242   //
243   // Assignment operator
244   //
245   if (this != &c) {
246     AliCFCutBase::operator=(c) ;
247     fMinNClusterTPC = c.fMinNClusterTPC ;
248     fMinNClusterITS = c.fMinNClusterITS ;
249     fMinNClusterTRD = c.fMinNClusterTRD ;
250     fMinFoundClusterTPC = c.fMinFoundClusterTPC ;
251     fMinNTrackletTRD = c.fMinNTrackletTRD ;
252     fMinNTrackletTRDpid = c.fMinNTrackletTRDpid ;
253     fMaxChi2PerClusterTPC = c.fMaxChi2PerClusterTPC ;
254     fMaxChi2PerClusterITS = c.fMaxChi2PerClusterITS ;
255     fMaxChi2PerTrackletTRD = c.fMaxChi2PerTrackletTRD ;
256     fMinNdEdxClusterTPC = c.fMinNdEdxClusterTPC;
257     fCovariance11Max = c.fCovariance11Max ;
258     fCovariance22Max = c.fCovariance22Max ;
259     fCovariance33Max = c.fCovariance33Max ;
260     fCovariance44Max = c.fCovariance44Max ;
261     fCovariance55Max = c.fCovariance55Max ;
262     fStatus = c.fStatus ;
263     fhCutStatistics = c.fhCutStatistics ;
264     fhCutCorrelation = c.fhCutCorrelation ;
265     fBitmap =  c.fBitmap ;
266     fTrackCuts = c.fTrackCuts ;
267     fhNBinsClusterTPC = c.fhNBinsClusterTPC ;
268     fhNBinsClusterITS = c.fhNBinsClusterITS ;
269     fhNBinsClusterTRD = c.fhNBinsClusterTRD ;
270     fhNBinsFoundClusterTPC = c.fhNBinsFoundClusterTPC ;
271     fhNBinsTrackletTRD = c.fhNBinsTrackletTRD ;
272     fhNBinsTrackletTRDpid = c.fhNBinsTrackletTRDpid ;
273     fhNBinsChi2TPC = c.fhNBinsChi2TPC ;
274     fhNBinsChi2ITS = c.fhNBinsChi2ITS ;
275     fhNBinsChi2TRD = c.fhNBinsChi2TRD ;
276     fhNBinsdEdxClusterTPC = c.fhNBinsdEdxClusterTPC ;
277     fhNBinsCovariance11 = c.fhNBinsCovariance11 ;
278     fhNBinsCovariance22 = c.fhNBinsCovariance22 ;
279     fhNBinsCovariance33 = c.fhNBinsCovariance33 ;
280     fhNBinsCovariance44 = c.fhNBinsCovariance44 ;
281     fhNBinsCovariance55 = c.fhNBinsCovariance55 ;
282     fhBinLimClusterTPC = c.fhBinLimClusterTPC ;
283     fhBinLimClusterITS = c.fhBinLimClusterITS ;
284     fhBinLimClusterTRD = c.fhBinLimClusterTRD ;
285     fhBinLimFoundClusterTPC = c.fhBinLimFoundClusterTPC ;
286     fhBinLimTrackletTRD = c.fhBinLimTrackletTRD ;
287     fhBinLimTrackletTRDpid = c.fhBinLimTrackletTRDpid ;
288     fhBinLimChi2TPC = c.fhBinLimChi2TPC ;
289     fhBinLimChi2ITS = c.fhBinLimChi2ITS ;
290     fhBinLimChi2TRD = c.fhBinLimChi2TRD ;
291     fhBinLimdEdxClusterTPC = c.fhBinLimdEdxClusterTPC ;
292     fhBinLimCovariance11 = c.fhBinLimCovariance11 ;
293     fhBinLimCovariance22 = c.fhBinLimCovariance22 ;
294     fhBinLimCovariance33 = c.fhBinLimCovariance33 ;
295     fhBinLimCovariance44 = c.fhBinLimCovariance44 ;
296     fhBinLimCovariance55 = c.fhBinLimCovariance55 ;
297
298     for (Int_t i=0; i<c.kNHist; i++){
299       for (Int_t j=0; j<c.kNStepQA; j++){
300         if(c.fhQA[i][j]) fhQA[i][j] = (TH1F*)c.fhQA[i][j]->Clone();
301       }
302     }
303     ((AliCFTrackQualityCuts &) c).Copy(*this);
304   }
305   return *this;
306 }
307 //__________________________________________________________________________________
308 AliCFTrackQualityCuts::~AliCFTrackQualityCuts()
309 {
310   //
311   // destructor
312   //
313   if (fhCutStatistics)  delete fhCutStatistics;
314   if (fhCutCorrelation) delete fhCutCorrelation;
315
316   for (Int_t i=0; i<kNHist; i++){
317     for (Int_t j=0; j<kNStepQA; j++){
318       if(fhQA[i][j]) delete fhQA[i][j];
319     }
320   }
321   if(fBitmap) delete fBitmap;
322   if(fTrackCuts) delete fTrackCuts;
323   if(fhBinLimClusterTPC) delete fhBinLimClusterTPC;
324   if(fhBinLimClusterITS) delete fhBinLimClusterITS;
325   if(fhBinLimClusterTRD) delete fhBinLimClusterTRD;
326   if(fhBinLimFoundClusterTPC) delete fhBinLimFoundClusterTPC;
327   if(fhBinLimTrackletTRD) delete fhBinLimTrackletTRD;
328   if(fhBinLimTrackletTRDpid) delete fhBinLimTrackletTRDpid;
329   if(fhBinLimChi2TPC) delete fhBinLimChi2TPC;
330   if(fhBinLimChi2ITS) delete fhBinLimChi2ITS;
331   if(fhBinLimChi2TRD) delete fhBinLimChi2TRD;
332   if(fhBinLimdEdxClusterTPC) delete fhBinLimdEdxClusterTPC;
333   if(fhBinLimCovariance11) delete fhBinLimCovariance11;
334   if(fhBinLimCovariance22) delete fhBinLimCovariance22;
335   if(fhBinLimCovariance33) delete fhBinLimCovariance33;
336   if(fhBinLimCovariance44) delete fhBinLimCovariance44;
337   if(fhBinLimCovariance55) delete fhBinLimCovariance55;
338 }
339 //__________________________________________________________________________________
340 void AliCFTrackQualityCuts::Initialise()
341 {
342   //
343   // sets everything to zero
344   //
345   fMinNClusterTPC = 0;
346   fMinNClusterITS = 0;
347   fMinNClusterTRD = 0;
348   fMinFoundClusterTPC = 0;
349   fMinNTrackletTRD = 0;
350   fMinNTrackletTRDpid = 0;
351   fMaxChi2PerClusterTPC = 0;
352   fMaxChi2PerClusterITS = 0;
353   fMaxChi2PerTrackletTRD = 0;
354   fMinNdEdxClusterTPC = 0;
355   fCovariance11Max = 0;
356   fCovariance22Max = 0;
357   fCovariance33Max = 0;
358   fCovariance44Max = 0;
359   fCovariance55Max = 0;
360   fStatus = 0;
361
362   SetMinNClusterTPC();
363   SetMinNClusterITS();
364   SetMinNClusterTRD();
365   SetMinFoundClusterTPC();
366   SetMinNTrackletTRD();
367   SetMinNTrackletTRDpid();
368   SetMaxChi2PerClusterTPC();
369   SetMaxChi2PerClusterITS();
370   SetMaxChi2PerTrackletTRD();
371   SetMinNdEdxClusterTPC();
372   SetMaxCovDiagonalElements();
373   SetStatus();
374
375   for (Int_t i=0; i<kNHist; i++){
376     for (Int_t j=0; j<kNStepQA; j++)
377       fhQA[i][j] = 0x0;
378   }
379   fhCutStatistics = 0;
380   fhCutCorrelation = 0;
381   fBitmap=new TBits(0);
382   fTrackCuts=new AliESDtrackCuts("aliESDtrackCuts","aliESDtrackCuts");
383
384   //set default bining for QA histograms
385   SetHistogramBins(kCutClusterTPC,165,-0.5,164.5);
386   SetHistogramBins(kCutClusterITS,8,-0.5,7.5);
387   SetHistogramBins(kCutClusterTRD,120,-0.5,119.5);
388   SetHistogramBins(kCutMinFoundClusterTPC,110,-0.05,1.05);
389   SetHistogramBins(kCutTrackletTRD,7,-0.5,6.5);
390   SetHistogramBins(kCutTrackletTRDpid,7,-0.5,6.5);
391   SetHistogramBins(kCutChi2TPC,500,0.,10.);
392   SetHistogramBins(kCutChi2ITS,500,0.,10.);
393   SetHistogramBins(kCutChi2TRD,500,0.,10.);
394   SetHistogramBins(kCutdEdxClusterTPC,165,-0.5,164.5);
395   SetHistogramBins(kCutCovElement11,200,0.,20.);
396   SetHistogramBins(kCutCovElement22,200,0.,20.);
397   SetHistogramBins(kCutCovElement33,100,0.,1.);
398   SetHistogramBins(kCutCovElement44,100,0.,5.);
399   SetHistogramBins(kCutCovElement55,100,0.,5.);
400 }
401 //__________________________________________________________________________________
402 void AliCFTrackQualityCuts::Copy(TObject &c) const
403 {
404   //
405   // Copy function
406   //
407   AliCFTrackQualityCuts& target = (AliCFTrackQualityCuts &) c;
408
409   target.Initialise();
410
411   if (fhCutStatistics)  target.fhCutStatistics = (TH1F*) fhCutStatistics->Clone();
412   if (fhCutCorrelation) target.fhCutCorrelation = (TH2F*) fhCutCorrelation->Clone();
413
414   for (Int_t i=0; i<kNHist; i++){
415     for (Int_t j=0; j<kNStepQA; j++){
416       if(fhQA[i][j]) target.fhQA[i][j] = (TH1F*)fhQA[i][j]->Clone();
417     }
418   }
419   TNamed::Copy(c);
420 }
421 //__________________________________________________________________________________
422 void AliCFTrackQualityCuts::SelectionBitMap(TObject* obj)
423 {
424   //
425   // test if the track passes the single cuts
426   // and store the information in a bitmap
427   //
428
429   // bitmap stores the decision of each single cut
430   for(Int_t i=0; i<kNCuts; i++)fBitmap->SetBitNumber(i,kFALSE);
431
432   if (!obj) return;
433   if (!obj->InheritsFrom("AliVParticle")) {
434     AliError("object must derived from AliVParticle !");
435     return;
436   }
437
438   Bool_t isESDTrack = strcmp(obj->ClassName(),"AliESDtrack") == 0 ? kTRUE : kFALSE ;
439   Bool_t isAODTrack = strcmp(obj->ClassName(),"AliAODTrack") == 0 ? kTRUE : kFALSE ;
440
441   if (!(isESDTrack || isAODTrack)) {
442     AliError("object must be an ESDtrack or an AODtrack !");
443     return;
444   }
445
446   AliESDtrack * esdTrack = 0x0 ;
447   AliAODTrack * aodTrack = 0x0 ; 
448   if (isESDTrack) esdTrack = dynamic_cast<AliESDtrack*>(obj);
449   if (isAODTrack) aodTrack = dynamic_cast<AliAODTrack*>(obj);
450
451   fTrackCuts->SetMinNClustersTPC(fMinNClusterTPC);
452   fTrackCuts->SetMinNClustersITS(fMinNClusterITS);
453   fTrackCuts->SetMaxChi2PerClusterTPC(fMaxChi2PerClusterTPC);
454   fTrackCuts->SetMaxChi2PerClusterITS(fMaxChi2PerClusterITS);
455   fTrackCuts->SetMaxCovDiagonalElements(fCovariance11Max,fCovariance22Max,fCovariance33Max,fCovariance44Max,fCovariance55Max);
456
457 // // // remove following 5 lines when AliESDtrackCuts is updated
458    Int_t    nClustersTPC = 0;
459    Int_t    nClustersITS = 0 ;
460    Float_t  chi2PerClusterTPC =  0 ;
461    Float_t  chi2PerClusterITS = 0 ;
462    Double_t extCov[15]={0,0,0,0,0,0,0,0,0,0,0,0,0,0,0};
463
464   // get cut quantities
465   Int_t   nClustersTRD = 0;
466   Int_t   nTrackletsTRD = 0;
467   Float_t chi2PerTrackletTRD = 0;
468   Float_t fractionFoundClustersTPC = 0;
469
470   if (isESDTrack) {
471     nClustersTRD = esdTrack->GetTRDncls();
472     nTrackletsTRD = esdTrack->GetTRDntracklets();
473     if (nTrackletsTRD != 0) chi2PerTrackletTRD = esdTrack->GetTRDchi2() / Float_t(nTrackletsTRD);
474 // // // include following line when AliESDtrackCuts is updated
475 //     if (esdTrack->GetTPCNclsF() != 0) fractionFoundClustersTPC = fTrackCuts->GetCutVariable(2) / float(esdTrack->GetTPCNclsF());
476
477 // // // remove following 6 lines when AliESDtrackCuts is updated
478      nClustersTPC = esdTrack->GetTPCclusters(0x0);
479      nClustersITS = esdTrack->GetITSclusters(0x0);
480      if (nClustersTPC != 0) chi2PerClusterTPC = esdTrack->GetTPCchi2() / Float_t(nClustersTPC);
481      if (nClustersITS != 0) chi2PerClusterITS = esdTrack->GetITSchi2() / Float_t(nClustersITS);
482      esdTrack->GetExternalCovariance(extCov);
483      if (esdTrack->GetTPCNclsF() != 0) fractionFoundClustersTPC = float(nClustersTPC) / float(esdTrack->GetTPCNclsF());
484   }
485
486   // fill the bitmap
487   Int_t iCutBit = 0;
488
489 // // // include following lines when AliESDtrackCuts is updated
490 //   fBitmap->SetBitNumber(iCutBit,fTrackCuts->GetCutDecision(2)); iCutBit++; // nClustersTPC
491 //   fBitmap->SetBitNumber(iCutBit,fTrackCuts->GetCutDecision(3)); iCutBit++; // nClustersITS
492 // // // remove following 6 lines when AliESDtrackCuts is updated
493    if (nClustersTPC >= fMinNClusterTPC)
494      fBitmap->SetBitNumber(iCutBit,kTRUE);
495    iCutBit++;
496    if (nClustersITS >= fMinNClusterITS)
497      fBitmap->SetBitNumber(iCutBit,kTRUE);
498    iCutBit++;
499
500   if (nClustersTRD >= fMinNClusterTRD)
501     fBitmap->SetBitNumber(iCutBit,kTRUE);
502   iCutBit++;
503 // // //   if ((fMinFoundClusterTPC <= 0) || (fTrackCuts->GetCutVariable(2) > 0 && (fractionFoundClustersTPC >= fMinFoundClusterTPC)))
504   if ((fMinFoundClusterTPC <= 0) || (nClustersTPC > 0 && (fractionFoundClustersTPC >= fMinFoundClusterTPC)))
505     fBitmap->SetBitNumber(iCutBit,kTRUE);
506   iCutBit++;
507   if (nTrackletsTRD >= fMinNTrackletTRD)
508     fBitmap->SetBitNumber(iCutBit,kTRUE);
509   iCutBit++;
510   if (!isESDTrack || esdTrack->GetTRDntrackletsPID() >= fMinNTrackletTRDpid)
511     fBitmap->SetBitNumber(iCutBit,kTRUE);
512   iCutBit++;
513
514 // // // include following lines when AliESDtrackCuts is updated
515 //   fBitmap->SetBitNumber(iCutBit,fTrackCuts->GetCutDecision(4)); iCutBit++; // chi2PerClusterTPC
516 //   fBitmap->SetBitNumber(iCutBit,fTrackCuts->GetCutDecision(5)); iCutBit++; // chi2PerClusterITS
517 // // // remove following 6 lines when AliESDtrackCuts is updated
518    if (chi2PerClusterTPC <= fMaxChi2PerClusterTPC)
519      fBitmap->SetBitNumber(iCutBit,kTRUE);
520    iCutBit++;
521    if (chi2PerClusterITS <= fMaxChi2PerClusterITS)
522      fBitmap->SetBitNumber(iCutBit,kTRUE);
523    iCutBit++;
524
525   if (chi2PerTrackletTRD <= fMaxChi2PerTrackletTRD)
526     fBitmap->SetBitNumber(iCutBit,kTRUE);
527   iCutBit++;
528   if (!isESDTrack || esdTrack->GetTPCsignalN() >= fMinNdEdxClusterTPC)
529     fBitmap->SetBitNumber(iCutBit,kTRUE);
530   iCutBit++;
531
532 // // // include following lines when AliESDtrackCuts is updated
533 //   fBitmap->SetBitNumber(iCutBit,fTrackCuts->GetCutDecision(6)); iCutBit++; // extCov[0]
534 //   fBitmap->SetBitNumber(iCutBit,fTrackCuts->GetCutDecision(7)); iCutBit++; // extCov[2]
535 //   fBitmap->SetBitNumber(iCutBit,fTrackCuts->GetCutDecision(8)); iCutBit++; // extCov[5]
536 //   fBitmap->SetBitNumber(iCutBit,fTrackCuts->GetCutDecision(9)); iCutBit++; // extCov[9]
537 //   fBitmap->SetBitNumber(iCutBit,fTrackCuts->GetCutDecision(10)); iCutBit++; // extCov[14]
538 // // // remove following lines when AliESDtrackCuts is updated
539    if (extCov[0]  <= fCovariance11Max)
540      fBitmap->SetBitNumber(iCutBit,kTRUE);
541    iCutBit++;
542    if (extCov[2]  <= fCovariance22Max)
543      fBitmap->SetBitNumber(iCutBit,kTRUE);
544    iCutBit++;
545    if (extCov[5]  <= fCovariance33Max)
546      fBitmap->SetBitNumber(iCutBit,kTRUE);
547    iCutBit++;
548    if (extCov[9]  <= fCovariance44Max)
549      fBitmap->SetBitNumber(iCutBit,kTRUE);
550    iCutBit++;
551    if (extCov[14] <= fCovariance55Max)
552      fBitmap->SetBitNumber(iCutBit,kTRUE);
553    iCutBit++;
554
555
556   if (isESDTrack) {
557     if ( (esdTrack->GetStatus() & fStatus) == fStatus ) fBitmap->SetBitNumber(iCutBit,kTRUE);
558   }
559   else {
560     if ( (aodTrack->GetStatus() & fStatus) == fStatus ) fBitmap->SetBitNumber(iCutBit,kTRUE);
561   }
562
563   return;
564 }
565 //__________________________________________________________________________________
566 Bool_t AliCFTrackQualityCuts::IsSelected(TObject* obj) {
567   //
568   // loops over decisions of single cuts and returns if the track is accepted
569   //
570   SelectionBitMap(obj);
571
572   if (fIsQAOn) FillHistograms(obj,0);
573   Bool_t isSelected = kTRUE;
574
575   for (UInt_t icut=0; icut<fBitmap->GetNbits();icut++) {
576     if(!fBitmap->TestBitNumber(icut)) {
577         isSelected = kFALSE;
578         break;
579     }
580   }
581   if (!isSelected) return kFALSE ;
582   if (fIsQAOn) FillHistograms(obj,1);
583   return kTRUE;
584 }
585 //__________________________________________________________________________________
586 void AliCFTrackQualityCuts::SetHistogramBins(Int_t index, Int_t nbins, Double_t *bins)
587 {
588   //
589   // variable bin size
590   //
591   if(!fIsQAOn) return;
592
593   if (index<0 || index>=kNHist) {
594     Error("SetHistogramBins","could not determine histogram from index %d",index);
595     return;
596   }
597
598   switch(index){
599   case kCutClusterTPC:
600     fhNBinsClusterTPC=nbins+1;
601     fhBinLimClusterTPC=new Double_t[nbins+1];
602     for(Int_t i=0;i<nbins+1;i++)fhBinLimClusterTPC[i]=bins[i];
603     break;
604
605   case kCutClusterITS:
606     fhNBinsClusterITS=nbins+1;
607     fhBinLimClusterITS=new Double_t[nbins+1];
608     for(Int_t i=0;i<nbins+1;i++)fhBinLimClusterITS[i]=bins[i];
609     break;
610
611   case kCutClusterTRD:
612     fhNBinsClusterTRD=nbins+1;
613     fhBinLimClusterTRD=new Double_t[nbins+1];
614     for(Int_t i=0;i<nbins+1;i++)fhBinLimClusterTRD[i]=bins[i];
615     break;
616
617   case kCutMinFoundClusterTPC:
618     fhNBinsFoundClusterTPC=nbins+1;
619     fhBinLimFoundClusterTPC=new Double_t[nbins+1];
620     for(Int_t i=0;i<nbins+1;i++)fhBinLimFoundClusterTPC[i]=bins[i];
621     break;
622
623   case kCutTrackletTRD:
624     fhNBinsTrackletTRD=nbins+1;
625     fhBinLimTrackletTRD=new Double_t[nbins+1];
626     for(Int_t i=0;i<nbins+1;i++)fhBinLimTrackletTRD[i]=bins[i];
627     break;
628
629   case kCutTrackletTRDpid:
630     fhNBinsTrackletTRDpid=nbins+1;
631     fhBinLimTrackletTRDpid=new Double_t[nbins+1];
632     for(Int_t i=0;i<nbins+1;i++)fhBinLimTrackletTRDpid[i]=bins[i];
633     break;
634
635   case kCutChi2TPC:
636     fhNBinsChi2TPC=nbins+1;
637     fhBinLimChi2TPC=new Double_t[nbins+1];
638     for(Int_t i=0;i<nbins+1;i++)fhBinLimChi2TPC[i]=bins[i];
639     break;
640
641   case kCutChi2ITS:
642     fhNBinsChi2ITS=nbins+1;
643     fhBinLimChi2ITS=new Double_t[nbins+1];
644     for(Int_t i=0;i<nbins+1;i++)fhBinLimChi2ITS[i]=bins[i];
645     break;
646
647   case kCutChi2TRD:
648     fhNBinsChi2TRD=nbins+1;
649     fhBinLimChi2TRD=new Double_t[nbins+1];
650     for(Int_t i=0;i<nbins+1;i++)fhBinLimChi2TRD[i]=bins[i];
651     break;
652
653   case kCutdEdxClusterTPC:
654     fhNBinsdEdxClusterTPC=nbins+1;
655     fhBinLimdEdxClusterTPC=new Double_t[nbins+1];
656     for(Int_t i=0;i<nbins+1;i++)fhBinLimdEdxClusterTPC[i]=bins[i];
657     break;
658
659   case kCutCovElement11:
660     fhNBinsCovariance11=nbins+1;
661     fhBinLimCovariance11=new Double_t[nbins+1];
662     for(Int_t i=0;i<nbins+1;i++)fhBinLimCovariance11[i]=bins[i];
663     break;
664
665   case kCutCovElement22:
666     fhNBinsCovariance22=nbins+1;
667     fhBinLimCovariance22=new Double_t[nbins+1];
668     for(Int_t i=0;i<nbins+1;i++)fhBinLimCovariance22[i]=bins[i];
669     break;
670
671   case kCutCovElement33:
672     fhNBinsCovariance33=nbins+1;
673     fhBinLimCovariance33=new Double_t[nbins+1];
674     for(Int_t i=0;i<nbins+1;i++)fhBinLimCovariance33[i]=bins[i];
675     break;
676
677   case kCutCovElement44:
678     fhNBinsCovariance44=nbins+1;
679     fhBinLimCovariance44=new Double_t[nbins+1];
680     for(Int_t i=0;i<nbins+1;i++)fhBinLimCovariance44[i]=bins[i];
681     break;
682
683   case kCutCovElement55:
684     fhNBinsCovariance55=nbins+1;
685     fhBinLimCovariance55=new Double_t[nbins+1];
686     for(Int_t i=0;i<nbins+1;i++)fhBinLimCovariance55[i]=bins[i];
687     break;
688  }
689 }
690 //__________________________________________________________________________________
691 void AliCFTrackQualityCuts::SetHistogramBins(Int_t index, Int_t nbins, Double_t xmin, Double_t xmax)
692 {
693   //
694   // fixed bin size
695   //
696   if (index<0 || index>=kNHist) {
697     Error("SetHistogramBins","could not determine histogram from index %d",index);
698     return;
699   }
700
701   switch(index){
702   case kCutClusterTPC:
703     fhNBinsClusterTPC=nbins+1;
704     fhBinLimClusterTPC=new Double_t[nbins+1];
705     for(Int_t i=0;i<nbins+1;i++)fhBinLimClusterTPC[i]=xmin+i*(xmax-xmin)/Double_t(nbins);
706     break;
707
708   case kCutClusterITS:
709     fhNBinsClusterITS=nbins+1;
710     fhBinLimClusterITS=new Double_t[nbins+1];
711     for(Int_t i=0;i<nbins+1;i++)fhBinLimClusterITS[i]=xmin+i*(xmax-xmin)/Double_t(nbins);
712     break;
713
714   case kCutClusterTRD:
715     fhNBinsClusterTRD=nbins+1;
716     fhBinLimClusterTRD=new Double_t[nbins+1];
717     for(Int_t i=0;i<nbins+1;i++)fhBinLimClusterTRD[i]=xmin+i*(xmax-xmin)/Double_t(nbins);
718     break;
719
720   case kCutMinFoundClusterTPC:
721     fhNBinsFoundClusterTPC=nbins+1;
722     fhBinLimFoundClusterTPC=new Double_t[nbins+1];
723     for(Int_t i=0;i<nbins+1;i++)fhBinLimFoundClusterTPC[i]=xmin+i*(xmax-xmin)/Double_t(nbins);
724     break;
725
726   case kCutTrackletTRD:
727     fhNBinsTrackletTRD=nbins+1;
728     fhBinLimTrackletTRD=new Double_t[nbins+1];
729     for(Int_t i=0;i<nbins+1;i++)fhBinLimTrackletTRD[i]=xmin+i*(xmax-xmin)/Double_t(nbins);
730     break;
731
732   case kCutTrackletTRDpid:
733     fhNBinsTrackletTRDpid=nbins+1;
734     fhBinLimTrackletTRDpid=new Double_t[nbins+1];
735     for(Int_t i=0;i<nbins+1;i++)fhBinLimTrackletTRDpid[i]=xmin+i*(xmax-xmin)/Double_t(nbins);
736     break;
737
738   case kCutChi2TPC:
739     fhNBinsChi2TPC=nbins+1;
740     fhBinLimChi2TPC=new Double_t[nbins+1];
741     for(Int_t i=0;i<nbins+1;i++)fhBinLimChi2TPC[i]=xmin+i*(xmax-xmin)/Double_t(nbins);
742     break;
743
744   case kCutChi2ITS:
745     fhNBinsChi2ITS=nbins+1;
746     fhBinLimChi2ITS=new Double_t[nbins+1];
747     for(Int_t i=0;i<nbins+1;i++)fhBinLimChi2ITS[i]=xmin+i*(xmax-xmin)/Double_t(nbins);
748     break;
749
750   case kCutChi2TRD:
751     fhNBinsChi2TRD=nbins+1;
752     fhBinLimChi2TRD=new Double_t[nbins+1];
753     for(Int_t i=0;i<nbins+1;i++)fhBinLimChi2TRD[i]=xmin+i*(xmax-xmin)/Double_t(nbins);
754     break;
755
756   case kCutdEdxClusterTPC:
757     fhNBinsdEdxClusterTPC=nbins+1;
758     fhBinLimdEdxClusterTPC=new Double_t[nbins+1];
759     for(Int_t i=0;i<nbins+1;i++)fhBinLimdEdxClusterTPC[i]=xmin+i*(xmax-xmin)/Double_t(nbins);
760     break;
761
762   case kCutCovElement11:
763     fhNBinsCovariance11=nbins+1;
764     fhBinLimCovariance11=new Double_t[nbins+1];
765     for(Int_t i=0;i<nbins+1;i++)fhBinLimCovariance11[i]=xmin+i*(xmax-xmin)/Double_t(nbins);
766     break;
767
768   case kCutCovElement22:
769     fhNBinsCovariance22=nbins+1;
770     fhBinLimCovariance22=new Double_t[nbins+1];
771     for(Int_t i=0;i<nbins+1;i++)fhBinLimCovariance22[i]=xmin+i*(xmax-xmin)/Double_t(nbins);
772     break;
773
774   case kCutCovElement33:
775     fhNBinsCovariance33=nbins+1;
776     fhBinLimCovariance33=new Double_t[nbins+1];
777     for(Int_t i=0;i<nbins+1;i++)fhBinLimCovariance33[i]=xmin+i*(xmax-xmin)/Double_t(nbins);
778     break;
779
780   case kCutCovElement44:
781     fhNBinsCovariance44=nbins+1;
782     fhBinLimCovariance44=new Double_t[nbins+1];
783     for(Int_t i=0;i<nbins+1;i++)fhBinLimCovariance44[i]=xmin+i*(xmax-xmin)/Double_t(nbins);
784     break;
785
786   case kCutCovElement55:
787     fhNBinsCovariance55=nbins+1;
788     fhBinLimCovariance55=new Double_t[nbins+1];
789     for(Int_t i=0;i<nbins+1;i++)fhBinLimCovariance55[i]=xmin+i*(xmax-xmin)/Double_t(nbins);
790     break;
791  }
792 }
793 //__________________________________________________________________________________
794  void AliCFTrackQualityCuts::DefineHistograms() {
795   //
796   // histograms for cut variables, cut statistics and cut correlations
797   //
798   Int_t color = 2;
799
800   // book cut statistics and cut correlation histograms
801   fhCutStatistics = new TH1F(Form("%s_cut_statistics",GetName()), Form("%s cut statistics",GetName()), kNCuts,0.5,kNCuts+0.5);
802   fhCutStatistics->SetLineWidth(2);
803   fhCutStatistics->GetXaxis()->SetBinLabel(1, "nClustersTPC");
804   fhCutStatistics->GetXaxis()->SetBinLabel(2, "nClustersITS");
805   fhCutStatistics->GetXaxis()->SetBinLabel(3, "nClustersTRD");
806   fhCutStatistics->GetXaxis()->SetBinLabel(4, "fractionClustersTPC");
807   fhCutStatistics->GetXaxis()->SetBinLabel(5, "ntrackletsTRD");
808   fhCutStatistics->GetXaxis()->SetBinLabel(6, "ntrackletsTRDpid");
809   fhCutStatistics->GetXaxis()->SetBinLabel(7, "chi2PerClusterTPC");
810   fhCutStatistics->GetXaxis()->SetBinLabel(8, "chi2PerClusterITS");
811   fhCutStatistics->GetXaxis()->SetBinLabel(9, "chi2PerTrackletTRD");
812   fhCutStatistics->GetXaxis()->SetBinLabel(10, "ndEdxClusterTPC");
813   fhCutStatistics->GetXaxis()->SetBinLabel(11, "covElement11");
814   fhCutStatistics->GetXaxis()->SetBinLabel(12, "covElement22");
815   fhCutStatistics->GetXaxis()->SetBinLabel(13, "covElement33");
816   fhCutStatistics->GetXaxis()->SetBinLabel(14, "covElement44");
817   fhCutStatistics->GetXaxis()->SetBinLabel(15, "covElement55");
818   fhCutStatistics->GetXaxis()->SetBinLabel(16, "status");
819
820   fhCutCorrelation = new TH2F(Form("%s_cut_correlation",GetName()), Form("%s cut  correlation",GetName()), kNCuts,0.5,kNCuts+0.5,kNCuts,0.5,kNCuts+0.5);
821   fhCutCorrelation->SetLineWidth(2);
822   fhCutCorrelation->GetXaxis()->SetBinLabel(1, "nClustersTPC");
823   fhCutCorrelation->GetXaxis()->SetBinLabel(2, "nClustersITS");
824   fhCutCorrelation->GetXaxis()->SetBinLabel(3, "nClustersTRD");
825   fhCutCorrelation->GetXaxis()->SetBinLabel(4, "fractionClustersTPC");
826   fhCutCorrelation->GetXaxis()->SetBinLabel(5, "ntrackletsTRD");
827   fhCutCorrelation->GetXaxis()->SetBinLabel(6, "ntrackletsTRDpid");
828   fhCutCorrelation->GetXaxis()->SetBinLabel(7, "chi2PerClusterTPC");
829   fhCutCorrelation->GetXaxis()->SetBinLabel(8, "chi2PerClusterITS");
830   fhCutCorrelation->GetXaxis()->SetBinLabel(9, "chi2PerTrackletTRD");
831   fhCutCorrelation->GetXaxis()->SetBinLabel(10, "ndEdxClustersTPC");
832   fhCutCorrelation->GetXaxis()->SetBinLabel(11, "covElement11");
833   fhCutCorrelation->GetXaxis()->SetBinLabel(12, "covElement22");
834   fhCutCorrelation->GetXaxis()->SetBinLabel(13, "covElement33");
835   fhCutCorrelation->GetXaxis()->SetBinLabel(14, "covElement44");
836   fhCutCorrelation->GetXaxis()->SetBinLabel(15, "covElement55");
837   fhCutCorrelation->GetXaxis()->SetBinLabel(16, "status");
838
839   fhCutCorrelation->GetYaxis()->SetBinLabel(1, "nClustersTPC");
840   fhCutCorrelation->GetYaxis()->SetBinLabel(2, "nClustersITS");
841   fhCutCorrelation->GetYaxis()->SetBinLabel(3, "nClustersTRD");
842   fhCutCorrelation->GetYaxis()->SetBinLabel(4, "fractionClustersTPC");
843   fhCutCorrelation->GetYaxis()->SetBinLabel(5, "ntrackletsTRD");
844   fhCutCorrelation->GetYaxis()->SetBinLabel(6, "ntrackletsTRDpid");
845   fhCutCorrelation->GetYaxis()->SetBinLabel(7, "chi2PerClusterTPC");
846   fhCutCorrelation->GetYaxis()->SetBinLabel(8, "chi2PerClusterITS");
847   fhCutCorrelation->GetYaxis()->SetBinLabel(9, "chi2PerTrackletTRD");
848   fhCutCorrelation->GetYaxis()->SetBinLabel(10, "ndEdxClustersTPC");
849   fhCutCorrelation->GetYaxis()->SetBinLabel(11, "covElement11");
850   fhCutCorrelation->GetYaxis()->SetBinLabel(12, "covElement22");
851   fhCutCorrelation->GetYaxis()->SetBinLabel(13, "covElement33");
852   fhCutCorrelation->GetYaxis()->SetBinLabel(14, "covElement44");
853   fhCutCorrelation->GetYaxis()->SetBinLabel(15, "covElement55");
854   fhCutCorrelation->GetYaxis()->SetBinLabel(16, "status");
855
856
857   // book QA histograms
858   Char_t str[256];
859   for (Int_t i=0; i<kNStepQA; i++) {
860     if (i==0) sprintf(str," ");
861     else sprintf(str,"_cut");
862
863     fhQA[kCutClusterTPC][i]     = new TH1F(Form("%s_nClustersTPC%s",GetName(),str)     ,"",fhNBinsClusterTPC-1,fhBinLimClusterTPC);
864     fhQA[kCutClusterITS][i]     = new TH1F(Form("%s_nClustersITS%s",GetName(),str)     ,"",fhNBinsClusterITS-1,fhBinLimClusterITS);
865     fhQA[kCutClusterTRD][i]     = new TH1F(Form("%s_nClustersTRD%s",GetName(),str)     ,"",fhNBinsClusterTRD-1,fhBinLimClusterTRD);
866     fhQA[kCutMinFoundClusterTPC][i]     = new TH1F(Form("%s_fractionClustersTPC%s",GetName(),str)     ,"",fhNBinsFoundClusterTPC-1,fhBinLimFoundClusterTPC);
867     fhQA[kCutTrackletTRD][i]    = new TH1F(Form("%s_ntrackletsTRD%s",GetName(),str)     ,"",fhNBinsTrackletTRD-1,fhBinLimTrackletTRD);
868     fhQA[kCutTrackletTRDpid][i] = new TH1F(Form("%s_ntrackletsTRDpid%s",GetName(),str)     ,"",fhNBinsTrackletTRDpid-1,fhBinLimTrackletTRDpid);
869     fhQA[kCutChi2TPC][i]        = new TH1F(Form("%s_chi2PerClusterTPC%s",GetName(),str),"",fhNBinsChi2TPC-1,fhBinLimChi2TPC);
870     fhQA[kCutChi2ITS][i]        = new TH1F(Form("%s_chi2PerClusterITS%s",GetName(),str),"",fhNBinsChi2ITS-1,fhBinLimChi2ITS);
871     fhQA[kCutChi2TRD][i]        = new TH1F(Form("%s_chi2PerTrackletTRD%s",GetName(),str),"",fhNBinsChi2TRD-1,fhBinLimChi2TRD);
872     fhQA[kCutdEdxClusterTPC][i] = new TH1F(Form("%s_ndEdxClustersTPC%s",GetName(),str)     ,"",fhNBinsdEdxClusterTPC-1,fhBinLimdEdxClusterTPC);
873     fhQA[kCutCovElement11][i]   = new TH1F(Form("%s_covMatrixDiagonal11%s",GetName(),str),"",fhNBinsCovariance11-1,fhBinLimCovariance11);
874     fhQA[kCutCovElement22][i]   = new TH1F(Form("%s_covMatrixDiagonal22%s",GetName(),str),"",fhNBinsCovariance22-1,fhBinLimCovariance22);
875     fhQA[kCutCovElement33][i]   = new TH1F(Form("%s_covMatrixDiagonal33%s",GetName(),str),"",fhNBinsCovariance33-1,fhBinLimCovariance33);
876     fhQA[kCutCovElement44][i]   = new TH1F(Form("%s_covMatrixDiagonal44%s",GetName(),str),"",fhNBinsCovariance44-1,fhBinLimCovariance44);
877     fhQA[kCutCovElement55][i]   = new TH1F(Form("%s_covMatrixDiagonal55%s",GetName(),str),"",fhNBinsCovariance55-1,fhBinLimCovariance55);
878
879     fhQA[kCutClusterTPC][i]     ->SetXTitle("n TPC clusters");
880     fhQA[kCutClusterITS][i]     ->SetXTitle("n ITS clusters");
881     fhQA[kCutClusterTRD][i]     ->SetXTitle("n TRD clusters");
882     fhQA[kCutMinFoundClusterTPC][i]->SetXTitle("fraction TPC clusters");
883     fhQA[kCutTrackletTRD][i] ->SetXTitle("n tracklets TRD");
884     fhQA[kCutTrackletTRDpid][i]->SetXTitle("n tracklets TRD pid");
885     fhQA[kCutChi2TPC][i]        ->SetXTitle("#chi^{2} per TPC cluster");
886     fhQA[kCutChi2ITS][i]        ->SetXTitle("#chi^{2} per ITS cluster");
887     fhQA[kCutChi2TRD][i]        ->SetXTitle("#chi^{2} per TRD tracklet");
888     fhQA[kCutdEdxClusterTPC][i] ->SetXTitle("n dEdx TPC clusters");
889     fhQA[kCutCovElement11][i]   ->SetXTitle("cov 11 : #sigma_{y}^{2} (cm^{2})");
890     fhQA[kCutCovElement22][i]   ->SetXTitle("cov 22 : #sigma_{z}^{2} (cm^{2})");
891     fhQA[kCutCovElement33][i]   ->SetXTitle("cov 33 : #sigma_{sin(#phi)}^{2}");
892     fhQA[kCutCovElement44][i]   ->SetXTitle("cov 44 : #sigma_{tan(#theta_{dip})}^{2}");
893     fhQA[kCutCovElement55][i]   ->SetXTitle("cov 55 : #sigma_{1/p_{T}}^{2} ((c/GeV)^{2})");
894   }
895
896   for(Int_t i=0; i<kNHist; i++) fhQA[i][1]->SetLineColor(color);
897 }
898 //__________________________________________________________________________________
899 void AliCFTrackQualityCuts::FillHistograms(TObject* obj, Bool_t b)
900 {
901   //
902   // fill the QA histograms
903   //
904
905   if (!obj) return;
906   if (!obj->InheritsFrom("AliVParticle")) {
907     AliError("object must derived from AliVParticle !");
908     return;
909   }
910
911   Bool_t isESDTrack = strcmp(obj->ClassName(),"AliESDtrack") == 0 ? kTRUE : kFALSE ;
912   Bool_t isAODTrack = strcmp(obj->ClassName(),"AliAODTrack") == 0 ? kTRUE : kFALSE ;
913
914   AliESDtrack * esdTrack = 0x0 ;
915   AliAODTrack * aodTrack = 0x0 ; 
916   if (isESDTrack) esdTrack = dynamic_cast<AliESDtrack*>(obj);
917   if (isAODTrack) aodTrack = dynamic_cast<AliAODTrack*>(obj);
918
919   // b = 0: fill histograms before cuts
920   // b = 1: fill histograms after cuts
921
922 // // // remove following 5 lines when AliESDtrackCuts is updated
923    Int_t    nClustersTPC = 0;
924    Int_t    nClustersITS = 0 ;
925    Float_t  chi2PerClusterTPC =  0 ;
926    Float_t  chi2PerClusterITS = 0 ;
927    Double_t extCov[15]={0,0,0,0,0,0,0,0,0,0,0,0,0,0,0};
928
929   Int_t   nClustersTRD = 0 ;
930   Int_t   nTrackletsTRD = 0 ;
931   Float_t chi2PerTrackletTRD = 0 ;
932   Float_t fractionFoundClustersTPC = 0;
933
934   if (isESDTrack) {
935     nClustersTRD = esdTrack->GetTRDncls();
936     nTrackletsTRD = esdTrack->GetTRDntracklets();
937     if (nTrackletsTRD != 0) chi2PerTrackletTRD = esdTrack->GetTRDchi2() / Float_t(nTrackletsTRD);
938
939 // // // include following line when AliESDtrackCuts is updated
940 //     if (esdTrack->GetTPCNclsF() != 0) fractionFoundClustersTPC = fTrackCuts->GetCutVariable(2) / float(esdTrack->GetTPCNclsF());
941
942 // // // remove following 6 lines when AliESDtrackCuts is updated
943      nClustersTPC = esdTrack->GetTPCclusters(0x0);
944      nClustersITS = esdTrack->GetITSclusters(0x0);
945      if (nClustersTPC != 0) chi2PerClusterTPC = esdTrack->GetTPCchi2() / Float_t(nClustersTPC);
946      if (nClustersITS != 0) chi2PerClusterITS = esdTrack->GetITSchi2() / Float_t(nClustersITS);
947      esdTrack->GetExternalCovariance(extCov);
948      if (esdTrack->GetTPCNclsF() != 0) fractionFoundClustersTPC = float(nClustersTPC) / float(esdTrack->GetTPCNclsF());
949   }
950
951 // // // include following lines when AliESDtrackCuts is updated
952 //   fhQA[kCutClusterTPC][b]->Fill(fTrackCuts->GetCutVariable(2));
953 //   fhQA[kCutClusterITS][b]->Fill(fTrackCuts->GetCutVariable(3));
954 //   fhQA[kCutChi2TPC][b]->Fill(fTrackCuts->GetCutVariable(4));
955 //   fhQA[kCutChi2ITS][b]->Fill(fTrackCuts->GetCutVariable(5));
956 // // // remove following 4 lines when AliESDtrackCuts is updated
957    fhQA[kCutClusterTPC][b]->Fill((float)nClustersTPC);
958    fhQA[kCutChi2TPC][b]->Fill(chi2PerClusterTPC);
959    fhQA[kCutClusterITS][b]->Fill((float)nClustersITS);
960    fhQA[kCutChi2ITS][b]->Fill(chi2PerClusterITS);
961
962
963   fhQA[kCutClusterTRD][b]->Fill((float)nClustersTRD);
964   fhQA[kCutChi2TRD][b]->Fill(chi2PerTrackletTRD);
965 //   if (b==0 || (b==1 && fTrackCuts->GetCutVariable(2)>0)) fhQA[kCutMinFoundClusterTPC][b]->Fill((float)fractionFoundClustersTPC);
966   if (b==0 || (b==1 && nClustersTPC>0)) fhQA[kCutMinFoundClusterTPC][b]->Fill((float)fractionFoundClustersTPC);
967   fhQA[kCutTrackletTRD][b]->Fill((float)nTrackletsTRD);
968   if (isESDTrack) {
969    fhQA[kCutTrackletTRDpid][b]->Fill((float)esdTrack->GetTRDntrackletsPID());
970    fhQA[kCutdEdxClusterTPC][b]->Fill((float)esdTrack->GetTPCsignalN());
971   }
972 // // // include following lines when AliESDtrackCuts is updated
973 //   fhQA[kCutCovElement11][b]->Fill(fTrackCuts->GetCutVariable(6));
974 //   fhQA[kCutCovElement22][b]->Fill(fTrackCuts->GetCutVariable(7));
975 //   fhQA[kCutCovElement33][b]->Fill(fTrackCuts->GetCutVariable(8));
976 //   fhQA[kCutCovElement44][b]->Fill(fTrackCuts->GetCutVariable(9));
977 //   fhQA[kCutCovElement55][b]->Fill(fTrackCuts->GetCutVariable(10));
978 // // // remove following 5 lines when AliESDtrackCuts is updated
979    fhQA[kCutCovElement11][b]->Fill(extCov[0]);
980    fhQA[kCutCovElement22][b]->Fill(extCov[2]);
981    fhQA[kCutCovElement33][b]->Fill(extCov[5]);
982    fhQA[kCutCovElement44][b]->Fill(extCov[9]);
983    fhQA[kCutCovElement55][b]->Fill(extCov[14]);
984
985
986   // fill cut statistics and cut correlation histograms with information from the bitmap
987   if (b) return;
988
989   // Get the bitmap of the single cuts
990   if ( !obj ) return;
991   SelectionBitMap(obj);
992
993   // Number of single cuts in this class
994   UInt_t ncuts = fBitmap->GetNbits();
995   for(UInt_t bit=0; bit<ncuts;bit++) {
996     if (!fBitmap->TestBitNumber(bit)) {
997         fhCutStatistics->Fill(bit+1);
998         for (UInt_t bit2=bit; bit2<ncuts;bit2++) {
999           if (!fBitmap->TestBitNumber(bit2)) 
1000             fhCutCorrelation->Fill(bit+1,bit2+1);
1001         }
1002     }
1003   }
1004 }
1005 //__________________________________________________________________________________
1006 void AliCFTrackQualityCuts::SaveHistograms(const Char_t* dir) {
1007   //
1008   // saves the histograms in a directory (dir)
1009   //
1010   if(!fIsQAOn) return;
1011
1012   if (!dir)
1013     dir = GetName();
1014
1015   gDirectory->mkdir(dir);
1016   gDirectory->cd(dir);
1017
1018   gDirectory->mkdir("before_cuts");
1019   gDirectory->mkdir("after_cuts");
1020
1021   fhCutStatistics->Write();
1022   fhCutCorrelation->Write();
1023
1024   for (Int_t j=0; j<kNStepQA; j++) {
1025     if (j==0)
1026       gDirectory->cd("before_cuts");
1027     else
1028       gDirectory->cd("after_cuts");
1029
1030     for(Int_t i=0; i<kNHist; i++) fhQA[i][j]->Write();
1031
1032     gDirectory->cd("../");
1033   }
1034   gDirectory->cd("../");
1035 }
1036 //__________________________________________________________________________________
1037 void AliCFTrackQualityCuts::DrawHistograms(Bool_t drawLogScale)
1038 {
1039   //
1040   // draws some histograms
1041   //
1042   if(!fIsQAOn) return;
1043
1044   // pad margins
1045   Float_t right = 0.03;
1046   Float_t left = 0.175;
1047   Float_t top = 0.03;
1048   Float_t bottom = 0.175;
1049
1050   TCanvas* canvas1 = new TCanvas("Track_QA_Quality_1", "Track QA Quality 1", 800, 500);
1051   canvas1->Divide(2, 1);
1052
1053   canvas1->cd(1);
1054   fhCutStatistics->SetStats(kFALSE);
1055   fhCutStatistics->LabelsOption("v");
1056   gPad->SetLeftMargin(left);
1057   gPad->SetBottomMargin(0.25);
1058   gPad->SetRightMargin(right);
1059   gPad->SetTopMargin(0.1);
1060   fhCutStatistics->Draw();
1061
1062   canvas1->cd(2);
1063   fhCutCorrelation->SetStats(kFALSE);
1064   fhCutCorrelation->LabelsOption("v");
1065   gPad->SetLeftMargin(0.30);
1066   gPad->SetRightMargin(bottom);
1067   gPad->SetTopMargin(0.1);
1068   gPad->SetBottomMargin(0.25);
1069   fhCutCorrelation->Draw("COLZ");
1070
1071   canvas1->SaveAs(Form("%s.eps", canvas1->GetName()));
1072   canvas1->SaveAs(Form("%s.ps", canvas1->GetName()));
1073
1074   // -----
1075
1076   TCanvas* canvas2 = new TCanvas("Track_QA_Quality_2", "Track QA Quality 2", 1200, 800);
1077   canvas2->Divide(2, 2);
1078
1079   canvas2->cd(1);
1080   gPad->SetRightMargin(right);
1081   gPad->SetLeftMargin(left);
1082   gPad->SetTopMargin(top);
1083   gPad->SetBottomMargin(bottom);
1084   fhQA[kCutClusterTPC][0]->SetStats(kFALSE);
1085   fhQA[kCutClusterTPC][0]->Draw();
1086   fhQA[kCutClusterTPC][1]->Draw("same");
1087
1088   canvas2->cd(2);
1089   gPad->SetRightMargin(right);
1090   gPad->SetLeftMargin(left);
1091   gPad->SetTopMargin(top);
1092   gPad->SetBottomMargin(bottom);
1093   fhQA[kCutChi2TPC][0]->SetStats(kFALSE);
1094   fhQA[kCutChi2TPC][0]->Draw();
1095   fhQA[kCutChi2TPC][1]->Draw("same");
1096
1097   canvas2->cd(3);
1098   gPad->SetRightMargin(right);
1099   gPad->SetLeftMargin(left);
1100   gPad->SetTopMargin(top);
1101   gPad->SetBottomMargin(bottom);
1102   fhQA[kCutClusterITS][0]->SetStats(kFALSE);
1103   fhQA[kCutClusterITS][0]->Draw();
1104   fhQA[kCutClusterITS][1]->Draw("same");
1105
1106   canvas2->cd(4);
1107   gPad->SetRightMargin(right);
1108   gPad->SetLeftMargin(left);
1109   gPad->SetTopMargin(top);
1110   gPad->SetBottomMargin(bottom);
1111   fhQA[kCutChi2ITS][0]->SetStats(kFALSE);
1112   fhQA[kCutChi2ITS][0]->Draw();
1113   fhQA[kCutChi2ITS][1]->Draw("same");
1114
1115   canvas2->SaveAs(Form("%s.eps", canvas2->GetName()));
1116   canvas2->SaveAs(Form("%s.ps", canvas2->GetName()));
1117
1118   // -----
1119
1120   TCanvas* canvas3 = new TCanvas("Track_QA_Quality_3", "Track QA Quality 3", 1200, 800);
1121   canvas3->Divide(3, 2);
1122
1123   canvas3->cd(1);
1124   if(drawLogScale) gPad->SetLogy();
1125   gPad->SetRightMargin(right);
1126   gPad->SetLeftMargin(left);
1127   gPad->SetTopMargin(top);
1128   gPad->SetBottomMargin(bottom);
1129   fhQA[kCutCovElement11][0]->SetStats(kFALSE);
1130   fhQA[kCutCovElement11][0]->Draw();
1131   fhQA[kCutCovElement11][1]->Draw("same");
1132
1133   canvas3->cd(2);
1134   if(drawLogScale) gPad->SetLogy();
1135   gPad->SetRightMargin(right);
1136   gPad->SetLeftMargin(left);
1137   gPad->SetTopMargin(top);
1138   gPad->SetBottomMargin(bottom);
1139   fhQA[kCutCovElement22][0]->SetStats(kFALSE);
1140   fhQA[kCutCovElement22][0]->Draw();
1141   fhQA[kCutCovElement22][1]->Draw("same");
1142
1143   canvas3->cd(3);
1144   if(drawLogScale) gPad->SetLogy();
1145   gPad->SetRightMargin(right);
1146   gPad->SetLeftMargin(left);
1147   gPad->SetTopMargin(top);
1148   gPad->SetBottomMargin(bottom);
1149   fhQA[kCutCovElement33][0]->SetStats(kFALSE);
1150   fhQA[kCutCovElement33][0]->Draw();
1151   fhQA[kCutCovElement33][1]->Draw("same");
1152
1153   canvas3->cd(4);
1154   if(drawLogScale) gPad->SetLogy();
1155   gPad->SetRightMargin(right);
1156   gPad->SetLeftMargin(left);
1157   gPad->SetTopMargin(top);
1158   gPad->SetBottomMargin(bottom);
1159   fhQA[kCutCovElement44][0]->SetStats(kFALSE);
1160   fhQA[kCutCovElement44][0]->Draw();
1161   fhQA[kCutCovElement44][1]->Draw("same");
1162
1163   canvas3->cd(5);
1164   if(drawLogScale) gPad->SetLogy();
1165   gPad->SetRightMargin(right);
1166   gPad->SetLeftMargin(left);
1167   gPad->SetTopMargin(top);
1168   gPad->SetBottomMargin(bottom);
1169   fhQA[kCutCovElement55][0]->SetStats(kFALSE);
1170   fhQA[kCutCovElement55][0]->Draw();
1171   fhQA[kCutCovElement55][1]->Draw("same");
1172
1173   canvas3->SaveAs(Form("%s.eps", canvas3->GetName()));
1174   canvas3->SaveAs(Form("%s.ps", canvas3->GetName()));
1175 }
1176 //__________________________________________________________________________________
1177 void AliCFTrackQualityCuts::AddQAHistograms(TList *qaList) {
1178   //
1179   // saves the histograms in a TList
1180   //
1181   DefineHistograms();
1182
1183   qaList->Add(fhCutStatistics);
1184   qaList->Add(fhCutCorrelation);
1185
1186   for (Int_t j=0; j<kNStepQA; j++) {
1187     for(Int_t i=0; i<kNHist; i++)
1188         qaList->Add(fhQA[i][j]);
1189   }
1190 }