]> git.uio.no Git - u/mrichter/AliRoot.git/blob - CORRFW/AliCFTrackQualityCuts.cxx
Revert unwanted changes concerning PID in HLT
[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   AliESDtrack * esdTrack = 0x0 ;
442   AliAODTrack * aodTrack = 0x0 ; 
443   if (isESDTrack) esdTrack = dynamic_cast<AliESDtrack*>(obj);
444   if (isAODTrack) aodTrack = dynamic_cast<AliAODTrack*>(obj);
445
446   fTrackCuts->SetMinNClustersTPC(fMinNClusterTPC);
447   fTrackCuts->SetMinNClustersITS(fMinNClusterITS);
448   fTrackCuts->SetMaxChi2PerClusterTPC(fMaxChi2PerClusterTPC);
449   fTrackCuts->SetMaxChi2PerClusterITS(fMaxChi2PerClusterITS);
450   fTrackCuts->SetMaxCovDiagonalElements(fCovariance11Max,fCovariance22Max,fCovariance33Max,fCovariance44Max,fCovariance55Max);
451   if (isESDTrack) Bool_t ok = fTrackCuts->AcceptTrack(esdTrack);
452
453 // // // remove following 5 lines when AliESDtrackCuts is updated
454    Int_t    nClustersTPC = 0;
455    Int_t    nClustersITS = 0 ;
456    Float_t  chi2PerClusterTPC =  0 ;
457    Float_t  chi2PerClusterITS = 0 ;
458    Double_t extCov[15]={0,0,0,0,0,0,0,0,0,0,0,0,0,0,0};
459
460   // get cut quantities
461   Int_t   nClustersTRD = 0;
462   Int_t   nTrackletsTRD = 0;
463   Float_t chi2PerTrackletTRD = 0;
464   Float_t fractionFoundClustersTPC = 0;
465
466   if (isESDTrack) {
467     nClustersTRD = esdTrack->GetTRDncls();
468     nTrackletsTRD = esdTrack->GetTRDntracklets();
469     if (nTrackletsTRD != 0) chi2PerTrackletTRD = esdTrack->GetTRDchi2() / Float_t(nTrackletsTRD);
470 // // // include following line when AliESDtrackCuts is updated
471 //     if (esdTrack->GetTPCNclsF() != 0) fractionFoundClustersTPC = fTrackCuts->GetCutVariable(2) / float(esdTrack->GetTPCNclsF());
472
473 // // // remove following 6 lines when AliESDtrackCuts is updated
474      nClustersTPC = esdTrack->GetTPCclusters(0x0);
475      nClustersITS = esdTrack->GetITSclusters(0x0);
476      if (nClustersTPC != 0) chi2PerClusterTPC = esdTrack->GetTPCchi2() / Float_t(nClustersTPC);
477      if (nClustersITS != 0) chi2PerClusterITS = esdTrack->GetITSchi2() / Float_t(nClustersITS);
478      esdTrack->GetExternalCovariance(extCov);
479      if (esdTrack->GetTPCNclsF() != 0) fractionFoundClustersTPC = float(nClustersTPC) / float(esdTrack->GetTPCNclsF());
480   }
481
482   // fill the bitmap
483   Int_t iCutBit = 0;
484
485 // // // include following lines when AliESDtrackCuts is updated
486 //   fBitmap->SetBitNumber(iCutBit,fTrackCuts->GetCutDecision(2)); iCutBit++; // nClustersTPC
487 //   fBitmap->SetBitNumber(iCutBit,fTrackCuts->GetCutDecision(3)); iCutBit++; // nClustersITS
488 // // // remove following 6 lines when AliESDtrackCuts is updated
489    if (nClustersTPC >= fMinNClusterTPC)
490      fBitmap->SetBitNumber(iCutBit,kTRUE);
491    iCutBit++;
492    if (nClustersITS >= fMinNClusterITS)
493      fBitmap->SetBitNumber(iCutBit,kTRUE);
494    iCutBit++;
495
496   if (nClustersTRD >= fMinNClusterTRD)
497     fBitmap->SetBitNumber(iCutBit,kTRUE);
498   iCutBit++;
499 // // //   if ((fMinFoundClusterTPC <= 0) || (fTrackCuts->GetCutVariable(2) > 0 && (fractionFoundClustersTPC >= fMinFoundClusterTPC)))
500   if ((fMinFoundClusterTPC <= 0) || (nClustersTPC > 0 && (fractionFoundClustersTPC >= fMinFoundClusterTPC)))
501     fBitmap->SetBitNumber(iCutBit,kTRUE);
502   iCutBit++;
503   if (nTrackletsTRD >= fMinNTrackletTRD)
504     fBitmap->SetBitNumber(iCutBit,kTRUE);
505   iCutBit++;
506   if (esdTrack->GetTRDntrackletsPID() >= fMinNTrackletTRDpid)
507     fBitmap->SetBitNumber(iCutBit,kTRUE);
508   iCutBit++;
509
510 // // // include following lines when AliESDtrackCuts is updated
511 //   fBitmap->SetBitNumber(iCutBit,fTrackCuts->GetCutDecision(4)); iCutBit++; // chi2PerClusterTPC
512 //   fBitmap->SetBitNumber(iCutBit,fTrackCuts->GetCutDecision(5)); iCutBit++; // chi2PerClusterITS
513 // // // remove following 6 lines when AliESDtrackCuts is updated
514    if (chi2PerClusterTPC <= fMaxChi2PerClusterTPC)
515      fBitmap->SetBitNumber(iCutBit,kTRUE);
516    iCutBit++;
517    if (chi2PerClusterITS <= fMaxChi2PerClusterITS)
518      fBitmap->SetBitNumber(iCutBit,kTRUE);
519    iCutBit++;
520
521   if (chi2PerTrackletTRD <= fMaxChi2PerTrackletTRD)
522     fBitmap->SetBitNumber(iCutBit,kTRUE);
523   iCutBit++;
524   if (esdTrack->GetTPCsignalN() >= fMinNdEdxClusterTPC)
525     fBitmap->SetBitNumber(iCutBit,kTRUE);
526   iCutBit++;
527
528 // // // include following lines when AliESDtrackCuts is updated
529 //   fBitmap->SetBitNumber(iCutBit,fTrackCuts->GetCutDecision(6)); iCutBit++; // extCov[0]
530 //   fBitmap->SetBitNumber(iCutBit,fTrackCuts->GetCutDecision(7)); iCutBit++; // extCov[2]
531 //   fBitmap->SetBitNumber(iCutBit,fTrackCuts->GetCutDecision(8)); iCutBit++; // extCov[5]
532 //   fBitmap->SetBitNumber(iCutBit,fTrackCuts->GetCutDecision(9)); iCutBit++; // extCov[9]
533 //   fBitmap->SetBitNumber(iCutBit,fTrackCuts->GetCutDecision(10)); iCutBit++; // extCov[14]
534 // // // remove following lines when AliESDtrackCuts is updated
535    if (extCov[0]  <= fCovariance11Max)
536      fBitmap->SetBitNumber(iCutBit,kTRUE);
537    iCutBit++;
538    if (extCov[2]  <= fCovariance22Max)
539      fBitmap->SetBitNumber(iCutBit,kTRUE);
540    iCutBit++;
541    if (extCov[5]  <= fCovariance33Max)
542      fBitmap->SetBitNumber(iCutBit,kTRUE);
543    iCutBit++;
544    if (extCov[9]  <= fCovariance44Max)
545      fBitmap->SetBitNumber(iCutBit,kTRUE);
546    iCutBit++;
547    if (extCov[14] <= fCovariance55Max)
548      fBitmap->SetBitNumber(iCutBit,kTRUE);
549    iCutBit++;
550
551
552   if (isESDTrack) {
553     if ( (esdTrack->GetStatus() & fStatus) == fStatus ) fBitmap->SetBitNumber(iCutBit,kTRUE);
554   }
555   else {
556     if ( (aodTrack->GetStatus() & fStatus) == fStatus ) fBitmap->SetBitNumber(iCutBit,kTRUE);
557   }
558
559   return;
560 }
561 //__________________________________________________________________________________
562 Bool_t AliCFTrackQualityCuts::IsSelected(TObject* obj) {
563   //
564   // loops over decisions of single cuts and returns if the track is accepted
565   //
566   SelectionBitMap(obj);
567
568   if (fIsQAOn) FillHistograms(obj,0);
569   Bool_t isSelected = kTRUE;
570
571   for (UInt_t icut=0; icut<fBitmap->GetNbits();icut++) {
572     if(!fBitmap->TestBitNumber(icut)) {
573         isSelected = kFALSE;
574         break;
575     }
576   }
577   if (!isSelected) return kFALSE ;
578   if (fIsQAOn) FillHistograms(obj,1);
579   return kTRUE;
580 }
581 //__________________________________________________________________________________
582 void AliCFTrackQualityCuts::SetHistogramBins(Int_t index, Int_t nbins, Double_t *bins)
583 {
584   //
585   // variable bin size
586   //
587   if(!fIsQAOn) return;
588
589   if (index<0 || index>=kNHist) {
590     Error("SetHistogramBins","could not determine histogram from index %d",index);
591     return;
592   }
593
594   switch(index){
595   case kCutClusterTPC:
596     fhNBinsClusterTPC=nbins+1;
597     fhBinLimClusterTPC=new Double_t[nbins+1];
598     for(Int_t i=0;i<nbins+1;i++)fhBinLimClusterTPC[i]=bins[i];
599     break;
600
601   case kCutClusterITS:
602     fhNBinsClusterITS=nbins+1;
603     fhBinLimClusterITS=new Double_t[nbins+1];
604     for(Int_t i=0;i<nbins+1;i++)fhBinLimClusterITS[i]=bins[i];
605     break;
606
607   case kCutClusterTRD:
608     fhNBinsClusterTRD=nbins+1;
609     fhBinLimClusterTRD=new Double_t[nbins+1];
610     for(Int_t i=0;i<nbins+1;i++)fhBinLimClusterTRD[i]=bins[i];
611     break;
612
613   case kCutMinFoundClusterTPC:
614     fhNBinsFoundClusterTPC=nbins+1;
615     fhBinLimFoundClusterTPC=new Double_t[nbins+1];
616     for(Int_t i=0;i<nbins+1;i++)fhBinLimFoundClusterTPC[i]=bins[i];
617     break;
618
619   case kCutTrackletTRD:
620     fhNBinsTrackletTRD=nbins+1;
621     fhBinLimTrackletTRD=new Double_t[nbins+1];
622     for(Int_t i=0;i<nbins+1;i++)fhBinLimTrackletTRD[i]=bins[i];
623     break;
624
625   case kCutTrackletTRDpid:
626     fhNBinsTrackletTRDpid=nbins+1;
627     fhBinLimTrackletTRDpid=new Double_t[nbins+1];
628     for(Int_t i=0;i<nbins+1;i++)fhBinLimTrackletTRDpid[i]=bins[i];
629     break;
630
631   case kCutChi2TPC:
632     fhNBinsChi2TPC=nbins+1;
633     fhBinLimChi2TPC=new Double_t[nbins+1];
634     for(Int_t i=0;i<nbins+1;i++)fhBinLimChi2TPC[i]=bins[i];
635     break;
636
637   case kCutChi2ITS:
638     fhNBinsChi2ITS=nbins+1;
639     fhBinLimChi2ITS=new Double_t[nbins+1];
640     for(Int_t i=0;i<nbins+1;i++)fhBinLimChi2ITS[i]=bins[i];
641     break;
642
643   case kCutChi2TRD:
644     fhNBinsChi2TRD=nbins+1;
645     fhBinLimChi2TRD=new Double_t[nbins+1];
646     for(Int_t i=0;i<nbins+1;i++)fhBinLimChi2TRD[i]=bins[i];
647     break;
648
649   case kCutdEdxClusterTPC:
650     fhNBinsdEdxClusterTPC=nbins+1;
651     fhBinLimdEdxClusterTPC=new Double_t[nbins+1];
652     for(Int_t i=0;i<nbins+1;i++)fhBinLimdEdxClusterTPC[i]=bins[i];
653     break;
654
655   case kCutCovElement11:
656     fhNBinsCovariance11=nbins+1;
657     fhBinLimCovariance11=new Double_t[nbins+1];
658     for(Int_t i=0;i<nbins+1;i++)fhBinLimCovariance11[i]=bins[i];
659     break;
660
661   case kCutCovElement22:
662     fhNBinsCovariance22=nbins+1;
663     fhBinLimCovariance22=new Double_t[nbins+1];
664     for(Int_t i=0;i<nbins+1;i++)fhBinLimCovariance22[i]=bins[i];
665     break;
666
667   case kCutCovElement33:
668     fhNBinsCovariance33=nbins+1;
669     fhBinLimCovariance33=new Double_t[nbins+1];
670     for(Int_t i=0;i<nbins+1;i++)fhBinLimCovariance33[i]=bins[i];
671     break;
672
673   case kCutCovElement44:
674     fhNBinsCovariance44=nbins+1;
675     fhBinLimCovariance44=new Double_t[nbins+1];
676     for(Int_t i=0;i<nbins+1;i++)fhBinLimCovariance44[i]=bins[i];
677     break;
678
679   case kCutCovElement55:
680     fhNBinsCovariance55=nbins+1;
681     fhBinLimCovariance55=new Double_t[nbins+1];
682     for(Int_t i=0;i<nbins+1;i++)fhBinLimCovariance55[i]=bins[i];
683     break;
684  }
685 }
686 //__________________________________________________________________________________
687 void AliCFTrackQualityCuts::SetHistogramBins(Int_t index, Int_t nbins, Double_t xmin, Double_t xmax)
688 {
689   //
690   // fixed bin size
691   //
692   if (index<0 || index>=kNHist) {
693     Error("SetHistogramBins","could not determine histogram from index %d",index);
694     return;
695   }
696
697   switch(index){
698   case kCutClusterTPC:
699     fhNBinsClusterTPC=nbins+1;
700     fhBinLimClusterTPC=new Double_t[nbins+1];
701     for(Int_t i=0;i<nbins+1;i++)fhBinLimClusterTPC[i]=xmin+i*(xmax-xmin)/Double_t(nbins);
702     break;
703
704   case kCutClusterITS:
705     fhNBinsClusterITS=nbins+1;
706     fhBinLimClusterITS=new Double_t[nbins+1];
707     for(Int_t i=0;i<nbins+1;i++)fhBinLimClusterITS[i]=xmin+i*(xmax-xmin)/Double_t(nbins);
708     break;
709
710   case kCutClusterTRD:
711     fhNBinsClusterTRD=nbins+1;
712     fhBinLimClusterTRD=new Double_t[nbins+1];
713     for(Int_t i=0;i<nbins+1;i++)fhBinLimClusterTRD[i]=xmin+i*(xmax-xmin)/Double_t(nbins);
714     break;
715
716   case kCutMinFoundClusterTPC:
717     fhNBinsFoundClusterTPC=nbins+1;
718     fhBinLimFoundClusterTPC=new Double_t[nbins+1];
719     for(Int_t i=0;i<nbins+1;i++)fhBinLimFoundClusterTPC[i]=xmin+i*(xmax-xmin)/Double_t(nbins);
720     break;
721
722   case kCutTrackletTRD:
723     fhNBinsTrackletTRD=nbins+1;
724     fhBinLimTrackletTRD=new Double_t[nbins+1];
725     for(Int_t i=0;i<nbins+1;i++)fhBinLimTrackletTRD[i]=xmin+i*(xmax-xmin)/Double_t(nbins);
726     break;
727
728   case kCutTrackletTRDpid:
729     fhNBinsTrackletTRDpid=nbins+1;
730     fhBinLimTrackletTRDpid=new Double_t[nbins+1];
731     for(Int_t i=0;i<nbins+1;i++)fhBinLimTrackletTRDpid[i]=xmin+i*(xmax-xmin)/Double_t(nbins);
732     break;
733
734   case kCutChi2TPC:
735     fhNBinsChi2TPC=nbins+1;
736     fhBinLimChi2TPC=new Double_t[nbins+1];
737     for(Int_t i=0;i<nbins+1;i++)fhBinLimChi2TPC[i]=xmin+i*(xmax-xmin)/Double_t(nbins);
738     break;
739
740   case kCutChi2ITS:
741     fhNBinsChi2ITS=nbins+1;
742     fhBinLimChi2ITS=new Double_t[nbins+1];
743     for(Int_t i=0;i<nbins+1;i++)fhBinLimChi2ITS[i]=xmin+i*(xmax-xmin)/Double_t(nbins);
744     break;
745
746   case kCutChi2TRD:
747     fhNBinsChi2TRD=nbins+1;
748     fhBinLimChi2TRD=new Double_t[nbins+1];
749     for(Int_t i=0;i<nbins+1;i++)fhBinLimChi2TRD[i]=xmin+i*(xmax-xmin)/Double_t(nbins);
750     break;
751
752   case kCutdEdxClusterTPC:
753     fhNBinsdEdxClusterTPC=nbins+1;
754     fhBinLimdEdxClusterTPC=new Double_t[nbins+1];
755     for(Int_t i=0;i<nbins+1;i++)fhBinLimdEdxClusterTPC[i]=xmin+i*(xmax-xmin)/Double_t(nbins);
756     break;
757
758   case kCutCovElement11:
759     fhNBinsCovariance11=nbins+1;
760     fhBinLimCovariance11=new Double_t[nbins+1];
761     for(Int_t i=0;i<nbins+1;i++)fhBinLimCovariance11[i]=xmin+i*(xmax-xmin)/Double_t(nbins);
762     break;
763
764   case kCutCovElement22:
765     fhNBinsCovariance22=nbins+1;
766     fhBinLimCovariance22=new Double_t[nbins+1];
767     for(Int_t i=0;i<nbins+1;i++)fhBinLimCovariance22[i]=xmin+i*(xmax-xmin)/Double_t(nbins);
768     break;
769
770   case kCutCovElement33:
771     fhNBinsCovariance33=nbins+1;
772     fhBinLimCovariance33=new Double_t[nbins+1];
773     for(Int_t i=0;i<nbins+1;i++)fhBinLimCovariance33[i]=xmin+i*(xmax-xmin)/Double_t(nbins);
774     break;
775
776   case kCutCovElement44:
777     fhNBinsCovariance44=nbins+1;
778     fhBinLimCovariance44=new Double_t[nbins+1];
779     for(Int_t i=0;i<nbins+1;i++)fhBinLimCovariance44[i]=xmin+i*(xmax-xmin)/Double_t(nbins);
780     break;
781
782   case kCutCovElement55:
783     fhNBinsCovariance55=nbins+1;
784     fhBinLimCovariance55=new Double_t[nbins+1];
785     for(Int_t i=0;i<nbins+1;i++)fhBinLimCovariance55[i]=xmin+i*(xmax-xmin)/Double_t(nbins);
786     break;
787  }
788 }
789 //__________________________________________________________________________________
790  void AliCFTrackQualityCuts::DefineHistograms() {
791   //
792   // histograms for cut variables, cut statistics and cut correlations
793   //
794   Int_t color = 2;
795
796   // book cut statistics and cut correlation histograms
797   fhCutStatistics = new TH1F(Form("%s_cut_statistics",GetName()), Form("%s cut statistics",GetName()), kNCuts,0.5,kNCuts+0.5);
798   fhCutStatistics->SetLineWidth(2);
799   fhCutStatistics->GetXaxis()->SetBinLabel(1, "nClustersTPC");
800   fhCutStatistics->GetXaxis()->SetBinLabel(2, "nClustersITS");
801   fhCutStatistics->GetXaxis()->SetBinLabel(3, "nClustersTRD");
802   fhCutStatistics->GetXaxis()->SetBinLabel(4, "fractionClustersTPC");
803   fhCutStatistics->GetXaxis()->SetBinLabel(5, "ntrackletsTRD");
804   fhCutStatistics->GetXaxis()->SetBinLabel(6, "ntrackletsTRDpid");
805   fhCutStatistics->GetXaxis()->SetBinLabel(7, "chi2PerClusterTPC");
806   fhCutStatistics->GetXaxis()->SetBinLabel(8, "chi2PerClusterITS");
807   fhCutStatistics->GetXaxis()->SetBinLabel(9, "chi2PerTrackletTRD");
808   fhCutStatistics->GetXaxis()->SetBinLabel(10, "ndEdxClusterTPC");
809   fhCutStatistics->GetXaxis()->SetBinLabel(11, "covElement11");
810   fhCutStatistics->GetXaxis()->SetBinLabel(12, "covElement22");
811   fhCutStatistics->GetXaxis()->SetBinLabel(13, "covElement33");
812   fhCutStatistics->GetXaxis()->SetBinLabel(14, "covElement44");
813   fhCutStatistics->GetXaxis()->SetBinLabel(15, "covElement55");
814   fhCutStatistics->GetXaxis()->SetBinLabel(16, "status");
815
816   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);
817   fhCutCorrelation->SetLineWidth(2);
818   fhCutCorrelation->GetXaxis()->SetBinLabel(1, "nClustersTPC");
819   fhCutCorrelation->GetXaxis()->SetBinLabel(2, "nClustersITS");
820   fhCutCorrelation->GetXaxis()->SetBinLabel(3, "nClustersTRD");
821   fhCutCorrelation->GetXaxis()->SetBinLabel(4, "fractionClustersTPC");
822   fhCutCorrelation->GetXaxis()->SetBinLabel(5, "ntrackletsTRD");
823   fhCutCorrelation->GetXaxis()->SetBinLabel(6, "ntrackletsTRDpid");
824   fhCutCorrelation->GetXaxis()->SetBinLabel(7, "chi2PerClusterTPC");
825   fhCutCorrelation->GetXaxis()->SetBinLabel(8, "chi2PerClusterITS");
826   fhCutCorrelation->GetXaxis()->SetBinLabel(9, "chi2PerTrackletTRD");
827   fhCutCorrelation->GetXaxis()->SetBinLabel(10, "ndEdxClustersTPC");
828   fhCutCorrelation->GetXaxis()->SetBinLabel(11, "covElement11");
829   fhCutCorrelation->GetXaxis()->SetBinLabel(12, "covElement22");
830   fhCutCorrelation->GetXaxis()->SetBinLabel(13, "covElement33");
831   fhCutCorrelation->GetXaxis()->SetBinLabel(14, "covElement44");
832   fhCutCorrelation->GetXaxis()->SetBinLabel(15, "covElement55");
833   fhCutCorrelation->GetXaxis()->SetBinLabel(16, "status");
834
835   fhCutCorrelation->GetYaxis()->SetBinLabel(1, "nClustersTPC");
836   fhCutCorrelation->GetYaxis()->SetBinLabel(2, "nClustersITS");
837   fhCutCorrelation->GetYaxis()->SetBinLabel(3, "nClustersTRD");
838   fhCutCorrelation->GetYaxis()->SetBinLabel(4, "fractionClustersTPC");
839   fhCutCorrelation->GetYaxis()->SetBinLabel(5, "ntrackletsTRD");
840   fhCutCorrelation->GetYaxis()->SetBinLabel(6, "ntrackletsTRDpid");
841   fhCutCorrelation->GetYaxis()->SetBinLabel(7, "chi2PerClusterTPC");
842   fhCutCorrelation->GetYaxis()->SetBinLabel(8, "chi2PerClusterITS");
843   fhCutCorrelation->GetYaxis()->SetBinLabel(9, "chi2PerTrackletTRD");
844   fhCutCorrelation->GetYaxis()->SetBinLabel(10, "ndEdxClustersTPC");
845   fhCutCorrelation->GetYaxis()->SetBinLabel(11, "covElement11");
846   fhCutCorrelation->GetYaxis()->SetBinLabel(12, "covElement22");
847   fhCutCorrelation->GetYaxis()->SetBinLabel(13, "covElement33");
848   fhCutCorrelation->GetYaxis()->SetBinLabel(14, "covElement44");
849   fhCutCorrelation->GetYaxis()->SetBinLabel(15, "covElement55");
850   fhCutCorrelation->GetYaxis()->SetBinLabel(16, "status");
851
852
853   // book QA histograms
854   Char_t str[256];
855   for (Int_t i=0; i<kNStepQA; i++) {
856     if (i==0) sprintf(str," ");
857     else sprintf(str,"_cut");
858
859     fhQA[kCutClusterTPC][i]     = new TH1F(Form("%s_nClustersTPC%s",GetName(),str)     ,"",fhNBinsClusterTPC-1,fhBinLimClusterTPC);
860     fhQA[kCutClusterITS][i]     = new TH1F(Form("%s_nClustersITS%s",GetName(),str)     ,"",fhNBinsClusterITS-1,fhBinLimClusterITS);
861     fhQA[kCutClusterTRD][i]     = new TH1F(Form("%s_nClustersTRD%s",GetName(),str)     ,"",fhNBinsClusterTRD-1,fhBinLimClusterTRD);
862     fhQA[kCutMinFoundClusterTPC][i]     = new TH1F(Form("%s_fractionClustersTPC%s",GetName(),str)     ,"",fhNBinsFoundClusterTPC-1,fhBinLimFoundClusterTPC);
863     fhQA[kCutTrackletTRD][i]    = new TH1F(Form("%s_ntrackletsTRD%s",GetName(),str)     ,"",fhNBinsTrackletTRD-1,fhBinLimTrackletTRD);
864     fhQA[kCutTrackletTRDpid][i] = new TH1F(Form("%s_ntrackletsTRDpid%s",GetName(),str)     ,"",fhNBinsTrackletTRDpid-1,fhBinLimTrackletTRDpid);
865     fhQA[kCutChi2TPC][i]        = new TH1F(Form("%s_chi2PerClusterTPC%s",GetName(),str),"",fhNBinsChi2TPC-1,fhBinLimChi2TPC);
866     fhQA[kCutChi2ITS][i]        = new TH1F(Form("%s_chi2PerClusterITS%s",GetName(),str),"",fhNBinsChi2ITS-1,fhBinLimChi2ITS);
867     fhQA[kCutChi2TRD][i]        = new TH1F(Form("%s_chi2PerTrackletTRD%s",GetName(),str),"",fhNBinsChi2TRD-1,fhBinLimChi2TRD);
868     fhQA[kCutdEdxClusterTPC][i] = new TH1F(Form("%s_ndEdxClustersTPC%s",GetName(),str)     ,"",fhNBinsdEdxClusterTPC-1,fhBinLimdEdxClusterTPC);
869     fhQA[kCutCovElement11][i]   = new TH1F(Form("%s_covMatrixDiagonal11%s",GetName(),str),"",fhNBinsCovariance11-1,fhBinLimCovariance11);
870     fhQA[kCutCovElement22][i]   = new TH1F(Form("%s_covMatrixDiagonal22%s",GetName(),str),"",fhNBinsCovariance22-1,fhBinLimCovariance22);
871     fhQA[kCutCovElement33][i]   = new TH1F(Form("%s_covMatrixDiagonal33%s",GetName(),str),"",fhNBinsCovariance33-1,fhBinLimCovariance33);
872     fhQA[kCutCovElement44][i]   = new TH1F(Form("%s_covMatrixDiagonal44%s",GetName(),str),"",fhNBinsCovariance44-1,fhBinLimCovariance44);
873     fhQA[kCutCovElement55][i]   = new TH1F(Form("%s_covMatrixDiagonal55%s",GetName(),str),"",fhNBinsCovariance55-1,fhBinLimCovariance55);
874
875     fhQA[kCutClusterTPC][i]     ->SetXTitle("n TPC clusters");
876     fhQA[kCutClusterITS][i]     ->SetXTitle("n ITS clusters");
877     fhQA[kCutClusterTRD][i]     ->SetXTitle("n TRD clusters");
878     fhQA[kCutMinFoundClusterTPC][i]->SetXTitle("fraction TPC clusters");
879     fhQA[kCutTrackletTRD][i] ->SetXTitle("n tracklets TRD");
880     fhQA[kCutTrackletTRDpid][i]->SetXTitle("n tracklets TRD pid");
881     fhQA[kCutChi2TPC][i]        ->SetXTitle("#chi^{2} per TPC cluster");
882     fhQA[kCutChi2ITS][i]        ->SetXTitle("#chi^{2} per ITS cluster");
883     fhQA[kCutChi2TRD][i]        ->SetXTitle("#chi^{2} per TRD tracklet");
884     fhQA[kCutdEdxClusterTPC][i] ->SetXTitle("n dEdx TPC clusters");
885     fhQA[kCutCovElement11][i]   ->SetXTitle("cov 11 : #sigma_{y}^{2} (cm^{2})");
886     fhQA[kCutCovElement22][i]   ->SetXTitle("cov 22 : #sigma_{z}^{2} (cm^{2})");
887     fhQA[kCutCovElement33][i]   ->SetXTitle("cov 33 : #sigma_{sin(#phi)}^{2}");
888     fhQA[kCutCovElement44][i]   ->SetXTitle("cov 44 : #sigma_{tan(#theta_{dip})}^{2}");
889     fhQA[kCutCovElement55][i]   ->SetXTitle("cov 55 : #sigma_{1/p_{T}}^{2} ((c/GeV)^{2})");
890   }
891
892   for(Int_t i=0; i<kNHist; i++) fhQA[i][1]->SetLineColor(color);
893 }
894 //__________________________________________________________________________________
895 void AliCFTrackQualityCuts::FillHistograms(TObject* obj, Bool_t b)
896 {
897   //
898   // fill the QA histograms
899   //
900
901   if (!obj) return;
902   if (!obj->InheritsFrom("AliVParticle")) {
903     AliError("object must derived from AliVParticle !");
904     return;
905   }
906
907   Bool_t isESDTrack = strcmp(obj->ClassName(),"AliESDtrack") == 0 ? kTRUE : kFALSE ;
908   Bool_t isAODTrack = strcmp(obj->ClassName(),"AliAODTrack") == 0 ? kTRUE : kFALSE ;
909
910   AliESDtrack * esdTrack = 0x0 ;
911   AliAODTrack * aodTrack = 0x0 ; 
912   if (isESDTrack) esdTrack = dynamic_cast<AliESDtrack*>(obj);
913   if (isAODTrack) aodTrack = dynamic_cast<AliAODTrack*>(obj);
914
915   // b = 0: fill histograms before cuts
916   // b = 1: fill histograms after cuts
917
918 // // // remove following 5 lines when AliESDtrackCuts is updated
919    Int_t    nClustersTPC = 0;
920    Int_t    nClustersITS = 0 ;
921    Float_t  chi2PerClusterTPC =  0 ;
922    Float_t  chi2PerClusterITS = 0 ;
923    Double_t extCov[15]={0,0,0,0,0,0,0,0,0,0,0,0,0,0,0};
924
925   Int_t   nClustersTRD = 0 ;
926   Int_t   nTrackletsTRD = 0 ;
927   Float_t chi2PerTrackletTRD = 0 ;
928   Float_t fractionFoundClustersTPC = 0;
929
930   if (isESDTrack) {
931     nClustersTRD = esdTrack->GetTRDncls();
932     nTrackletsTRD = esdTrack->GetTRDntracklets();
933     if (nTrackletsTRD != 0) chi2PerTrackletTRD = esdTrack->GetTRDchi2() / Float_t(nTrackletsTRD);
934
935 // // // include following line when AliESDtrackCuts is updated
936 //     if (esdTrack->GetTPCNclsF() != 0) fractionFoundClustersTPC = fTrackCuts->GetCutVariable(2) / float(esdTrack->GetTPCNclsF());
937
938 // // // remove following 6 lines when AliESDtrackCuts is updated
939      nClustersTPC = esdTrack->GetTPCclusters(0x0);
940      nClustersITS = esdTrack->GetITSclusters(0x0);
941      if (nClustersTPC != 0) chi2PerClusterTPC = esdTrack->GetTPCchi2() / Float_t(nClustersTPC);
942      if (nClustersITS != 0) chi2PerClusterITS = esdTrack->GetITSchi2() / Float_t(nClustersITS);
943      esdTrack->GetExternalCovariance(extCov);
944      if (esdTrack->GetTPCNclsF() != 0) fractionFoundClustersTPC = float(nClustersTPC) / float(esdTrack->GetTPCNclsF());
945   }
946
947 // // // include following lines when AliESDtrackCuts is updated
948 //   fhQA[kCutClusterTPC][b]->Fill(fTrackCuts->GetCutVariable(2));
949 //   fhQA[kCutClusterITS][b]->Fill(fTrackCuts->GetCutVariable(3));
950 //   fhQA[kCutChi2TPC][b]->Fill(fTrackCuts->GetCutVariable(4));
951 //   fhQA[kCutChi2ITS][b]->Fill(fTrackCuts->GetCutVariable(5));
952 // // // remove following 4 lines when AliESDtrackCuts is updated
953    fhQA[kCutClusterTPC][b]->Fill((float)nClustersTPC);
954    fhQA[kCutChi2TPC][b]->Fill(chi2PerClusterTPC);
955    fhQA[kCutClusterITS][b]->Fill((float)nClustersITS);
956    fhQA[kCutChi2ITS][b]->Fill(chi2PerClusterITS);
957
958
959   fhQA[kCutClusterTRD][b]->Fill((float)nClustersTRD);
960   fhQA[kCutChi2TRD][b]->Fill(chi2PerTrackletTRD);
961 //   if (b==0 || (b==1 && fTrackCuts->GetCutVariable(2)>0)) fhQA[kCutMinFoundClusterTPC][b]->Fill((float)fractionFoundClustersTPC);
962   if (b==0 || (b==1 && nClustersTPC>0)) fhQA[kCutMinFoundClusterTPC][b]->Fill((float)fractionFoundClustersTPC);
963   fhQA[kCutTrackletTRD][b]->Fill((float)nTrackletsTRD);
964   fhQA[kCutTrackletTRDpid][b]->Fill((float)esdTrack->GetTRDntrackletsPID());
965   fhQA[kCutdEdxClusterTPC][b]->Fill((float)esdTrack->GetTPCsignalN());
966
967 // // // include following lines when AliESDtrackCuts is updated
968 //   fhQA[kCutCovElement11][b]->Fill(fTrackCuts->GetCutVariable(6));
969 //   fhQA[kCutCovElement22][b]->Fill(fTrackCuts->GetCutVariable(7));
970 //   fhQA[kCutCovElement33][b]->Fill(fTrackCuts->GetCutVariable(8));
971 //   fhQA[kCutCovElement44][b]->Fill(fTrackCuts->GetCutVariable(9));
972 //   fhQA[kCutCovElement55][b]->Fill(fTrackCuts->GetCutVariable(10));
973 // // // remove following 5 lines when AliESDtrackCuts is updated
974    fhQA[kCutCovElement11][b]->Fill(extCov[0]);
975    fhQA[kCutCovElement22][b]->Fill(extCov[2]);
976    fhQA[kCutCovElement33][b]->Fill(extCov[5]);
977    fhQA[kCutCovElement44][b]->Fill(extCov[9]);
978    fhQA[kCutCovElement55][b]->Fill(extCov[14]);
979
980
981   // fill cut statistics and cut correlation histograms with information from the bitmap
982   if (b) return;
983
984   // Get the bitmap of the single cuts
985   if ( !obj ) return;
986   SelectionBitMap(obj);
987
988   // Number of single cuts in this class
989   UInt_t ncuts = fBitmap->GetNbits();
990   for(UInt_t bit=0; bit<ncuts;bit++) {
991     if (!fBitmap->TestBitNumber(bit)) {
992         fhCutStatistics->Fill(bit+1);
993         for (UInt_t bit2=bit; bit2<ncuts;bit2++) {
994           if (!fBitmap->TestBitNumber(bit2)) 
995             fhCutCorrelation->Fill(bit+1,bit2+1);
996         }
997     }
998   }
999 }
1000 //__________________________________________________________________________________
1001 void AliCFTrackQualityCuts::SaveHistograms(const Char_t* dir) {
1002   //
1003   // saves the histograms in a directory (dir)
1004   //
1005   if(!fIsQAOn) return;
1006
1007   if (!dir)
1008     dir = GetName();
1009
1010   gDirectory->mkdir(dir);
1011   gDirectory->cd(dir);
1012
1013   gDirectory->mkdir("before_cuts");
1014   gDirectory->mkdir("after_cuts");
1015
1016   fhCutStatistics->Write();
1017   fhCutCorrelation->Write();
1018
1019   for (Int_t j=0; j<kNStepQA; j++) {
1020     if (j==0)
1021       gDirectory->cd("before_cuts");
1022     else
1023       gDirectory->cd("after_cuts");
1024
1025     for(Int_t i=0; i<kNHist; i++) fhQA[i][j]->Write();
1026
1027     gDirectory->cd("../");
1028   }
1029   gDirectory->cd("../");
1030 }
1031 //__________________________________________________________________________________
1032 void AliCFTrackQualityCuts::DrawHistograms(Bool_t drawLogScale)
1033 {
1034   //
1035   // draws some histograms
1036   //
1037   if(!fIsQAOn) return;
1038
1039   // pad margins
1040   Float_t right = 0.03;
1041   Float_t left = 0.175;
1042   Float_t top = 0.03;
1043   Float_t bottom = 0.175;
1044
1045   TCanvas* canvas1 = new TCanvas("Track_QA_Quality_1", "Track QA Quality 1", 800, 500);
1046   canvas1->Divide(2, 1);
1047
1048   canvas1->cd(1);
1049   fhCutStatistics->SetStats(kFALSE);
1050   fhCutStatistics->LabelsOption("v");
1051   gPad->SetLeftMargin(left);
1052   gPad->SetBottomMargin(0.25);
1053   gPad->SetRightMargin(right);
1054   gPad->SetTopMargin(0.1);
1055   fhCutStatistics->Draw();
1056
1057   canvas1->cd(2);
1058   fhCutCorrelation->SetStats(kFALSE);
1059   fhCutCorrelation->LabelsOption("v");
1060   gPad->SetLeftMargin(0.30);
1061   gPad->SetRightMargin(bottom);
1062   gPad->SetTopMargin(0.1);
1063   gPad->SetBottomMargin(0.25);
1064   fhCutCorrelation->Draw("COLZ");
1065
1066   canvas1->SaveAs(Form("%s.eps", canvas1->GetName()));
1067   canvas1->SaveAs(Form("%s.ps", canvas1->GetName()));
1068
1069   // -----
1070
1071   TCanvas* canvas2 = new TCanvas("Track_QA_Quality_2", "Track QA Quality 2", 1200, 800);
1072   canvas2->Divide(2, 2);
1073
1074   canvas2->cd(1);
1075   gPad->SetRightMargin(right);
1076   gPad->SetLeftMargin(left);
1077   gPad->SetTopMargin(top);
1078   gPad->SetBottomMargin(bottom);
1079   fhQA[kCutClusterTPC][0]->SetStats(kFALSE);
1080   fhQA[kCutClusterTPC][0]->Draw();
1081   fhQA[kCutClusterTPC][1]->Draw("same");
1082
1083   canvas2->cd(2);
1084   gPad->SetRightMargin(right);
1085   gPad->SetLeftMargin(left);
1086   gPad->SetTopMargin(top);
1087   gPad->SetBottomMargin(bottom);
1088   fhQA[kCutChi2TPC][0]->SetStats(kFALSE);
1089   fhQA[kCutChi2TPC][0]->Draw();
1090   fhQA[kCutChi2TPC][1]->Draw("same");
1091
1092   canvas2->cd(3);
1093   gPad->SetRightMargin(right);
1094   gPad->SetLeftMargin(left);
1095   gPad->SetTopMargin(top);
1096   gPad->SetBottomMargin(bottom);
1097   fhQA[kCutClusterITS][0]->SetStats(kFALSE);
1098   fhQA[kCutClusterITS][0]->Draw();
1099   fhQA[kCutClusterITS][1]->Draw("same");
1100
1101   canvas2->cd(4);
1102   gPad->SetRightMargin(right);
1103   gPad->SetLeftMargin(left);
1104   gPad->SetTopMargin(top);
1105   gPad->SetBottomMargin(bottom);
1106   fhQA[kCutChi2ITS][0]->SetStats(kFALSE);
1107   fhQA[kCutChi2ITS][0]->Draw();
1108   fhQA[kCutChi2ITS][1]->Draw("same");
1109
1110   canvas2->SaveAs(Form("%s.eps", canvas2->GetName()));
1111   canvas2->SaveAs(Form("%s.ps", canvas2->GetName()));
1112
1113   // -----
1114
1115   TCanvas* canvas3 = new TCanvas("Track_QA_Quality_3", "Track QA Quality 3", 1200, 800);
1116   canvas3->Divide(3, 2);
1117
1118   canvas3->cd(1);
1119   if(drawLogScale) gPad->SetLogy();
1120   gPad->SetRightMargin(right);
1121   gPad->SetLeftMargin(left);
1122   gPad->SetTopMargin(top);
1123   gPad->SetBottomMargin(bottom);
1124   fhQA[kCutCovElement11][0]->SetStats(kFALSE);
1125   fhQA[kCutCovElement11][0]->Draw();
1126   fhQA[kCutCovElement11][1]->Draw("same");
1127
1128   canvas3->cd(2);
1129   if(drawLogScale) gPad->SetLogy();
1130   gPad->SetRightMargin(right);
1131   gPad->SetLeftMargin(left);
1132   gPad->SetTopMargin(top);
1133   gPad->SetBottomMargin(bottom);
1134   fhQA[kCutCovElement22][0]->SetStats(kFALSE);
1135   fhQA[kCutCovElement22][0]->Draw();
1136   fhQA[kCutCovElement22][1]->Draw("same");
1137
1138   canvas3->cd(3);
1139   if(drawLogScale) gPad->SetLogy();
1140   gPad->SetRightMargin(right);
1141   gPad->SetLeftMargin(left);
1142   gPad->SetTopMargin(top);
1143   gPad->SetBottomMargin(bottom);
1144   fhQA[kCutCovElement33][0]->SetStats(kFALSE);
1145   fhQA[kCutCovElement33][0]->Draw();
1146   fhQA[kCutCovElement33][1]->Draw("same");
1147
1148   canvas3->cd(4);
1149   if(drawLogScale) gPad->SetLogy();
1150   gPad->SetRightMargin(right);
1151   gPad->SetLeftMargin(left);
1152   gPad->SetTopMargin(top);
1153   gPad->SetBottomMargin(bottom);
1154   fhQA[kCutCovElement44][0]->SetStats(kFALSE);
1155   fhQA[kCutCovElement44][0]->Draw();
1156   fhQA[kCutCovElement44][1]->Draw("same");
1157
1158   canvas3->cd(5);
1159   if(drawLogScale) gPad->SetLogy();
1160   gPad->SetRightMargin(right);
1161   gPad->SetLeftMargin(left);
1162   gPad->SetTopMargin(top);
1163   gPad->SetBottomMargin(bottom);
1164   fhQA[kCutCovElement55][0]->SetStats(kFALSE);
1165   fhQA[kCutCovElement55][0]->Draw();
1166   fhQA[kCutCovElement55][1]->Draw("same");
1167
1168   canvas3->SaveAs(Form("%s.eps", canvas3->GetName()));
1169   canvas3->SaveAs(Form("%s.ps", canvas3->GetName()));
1170 }
1171 //__________________________________________________________________________________
1172 void AliCFTrackQualityCuts::AddQAHistograms(TList *qaList) {
1173   //
1174   // saves the histograms in a TList
1175   //
1176   DefineHistograms();
1177
1178   qaList->Add(fhCutStatistics);
1179   qaList->Add(fhCutCorrelation);
1180
1181   for (Int_t j=0; j<kNStepQA; j++) {
1182     for(Int_t i=0; i<kNHist; i++)
1183         qaList->Add(fhQA[i][j]);
1184   }
1185 }