]> git.uio.no Git - u/mrichter/AliRoot.git/blame - CORRFW/AliCFTrackQualityCuts.cxx
QA added
[u/mrichter/AliRoot.git] / CORRFW / AliCFTrackQualityCuts.cxx
CommitLineData
563113d0 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 TPC
25// - number of clusters in the ITS
26// - chi2 / cluster in the TPC
27// - chi2 / cluster in the ITS
28// - successful TPC refit
29// - successful ITS refit
30// - covariance matrix diagonal elements
31//
32// The cut values for these cuts are set with the corresponding set functions.
33// All cut classes provided by the correction framework are supposed to be
34// added in the Analysis Framwork's class AliAnalysisFilter and applied by
35// the filter via a loop.
36//
37// author: I. Kraus (Ingrid.Kraus@cern.ch)
38// idea taken form
39// AliESDtrackCuts writte by Jan Fiete Grosse-Oetringhaus and
40// AliRsnDaughterCut class written by A. Pulvirenti.
41
42
43#include <TCanvas.h>
44#include <TDirectory.h>
563113d0 45#include <TH2.h>
6cb5e116 46#include <TBits.h>
563113d0 47
48#include <AliESDtrack.h>
49#include <AliLog.h>
50#include "AliCFTrackQualityCuts.h"
51
52ClassImp(AliCFTrackQualityCuts)
53
54//__________________________________________________________________________________
55AliCFTrackQualityCuts::AliCFTrackQualityCuts() :
56 AliCFCutBase(),
57 fMinNClusterTPC(0),
58 fMinNClusterITS(0),
59 fMaxChi2PerClusterTPC(0),
60 fMaxChi2PerClusterITS(0),
61 fRequireTPCRefit(0),
62 fRequireITSRefit(0),
63 fCovariance11Max(0),
64 fCovariance22Max(0),
65 fCovariance33Max(0),
66 fCovariance44Max(0),
67 fCovariance55Max(0),
68 fhCutStatistics(0),
69 fhCutCorrelation(0),
70 fBitmap(0x0),
71 fhNBinsClusterTPC(0),
72 fhNBinsClusterITS(0),
73 fhNBinsChi2TPC(0),
74 fhNBinsChi2ITS(0),
75 fhNBinsRefitTPC(0),
76 fhNBinsRefitITS(0),
77 fhNBinsCovariance11(0),
78 fhNBinsCovariance22(0),
79 fhNBinsCovariance33(0),
80 fhNBinsCovariance44(0),
81 fhNBinsCovariance55(0),
82 fhBinLimClusterTPC(0x0),
83 fhBinLimClusterITS(0x0),
84 fhBinLimChi2TPC(0x0),
85 fhBinLimChi2ITS(0x0),
86 fhBinLimRefitTPC(0x0),
87 fhBinLimRefitITS(0x0),
88 fhBinLimCovariance11(0x0),
89 fhBinLimCovariance22(0x0),
90 fhBinLimCovariance33(0x0),
91 fhBinLimCovariance44(0x0),
92 fhBinLimCovariance55(0x0)
6cb5e116 93{
563113d0 94 //
95 // Default constructor
96 //
563113d0 97 Initialise();
98}
99//__________________________________________________________________________________
100AliCFTrackQualityCuts::AliCFTrackQualityCuts(Char_t* name, Char_t* title) :
101 AliCFCutBase(name,title),
102 fMinNClusterTPC(0),
103 fMinNClusterITS(0),
104 fMaxChi2PerClusterTPC(0),
105 fMaxChi2PerClusterITS(0),
106 fRequireTPCRefit(0),
107 fRequireITSRefit(0),
108 fCovariance11Max(0),
109 fCovariance22Max(0),
110 fCovariance33Max(0),
111 fCovariance44Max(0),
112 fCovariance55Max(0),
113 fhCutStatistics(0),
114 fhCutCorrelation(0),
115 fBitmap(0x0),
116 fhNBinsClusterTPC(0),
117 fhNBinsClusterITS(0),
118 fhNBinsChi2TPC(0),
119 fhNBinsChi2ITS(0),
120 fhNBinsRefitTPC(0),
121 fhNBinsRefitITS(0),
122 fhNBinsCovariance11(0),
123 fhNBinsCovariance22(0),
124 fhNBinsCovariance33(0),
125 fhNBinsCovariance44(0),
126 fhNBinsCovariance55(0),
127 fhBinLimClusterTPC(0x0),
128 fhBinLimClusterITS(0x0),
129 fhBinLimChi2TPC(0x0),
130 fhBinLimChi2ITS(0x0),
131 fhBinLimRefitTPC(0x0),
132 fhBinLimRefitITS(0x0),
133 fhBinLimCovariance11(0x0),
134 fhBinLimCovariance22(0x0),
135 fhBinLimCovariance33(0x0),
136 fhBinLimCovariance44(0x0),
137 fhBinLimCovariance55(0x0)
563113d0 138{
139 //
140 // Constructor
141 //
563113d0 142 Initialise();
143}
144//__________________________________________________________________________________
145AliCFTrackQualityCuts::AliCFTrackQualityCuts(const AliCFTrackQualityCuts& c) :
146 AliCFCutBase(c),
147 fMinNClusterTPC(c.fMinNClusterTPC),
148 fMinNClusterITS(c.fMinNClusterITS),
149 fMaxChi2PerClusterTPC(c.fMaxChi2PerClusterTPC),
150 fMaxChi2PerClusterITS(c.fMaxChi2PerClusterITS),
151 fRequireTPCRefit(c.fRequireTPCRefit),
152 fRequireITSRefit(c.fRequireITSRefit),
153 fCovariance11Max(c.fCovariance11Max),
154 fCovariance22Max(c.fCovariance22Max),
155 fCovariance33Max(c.fCovariance33Max),
156 fCovariance44Max(c.fCovariance44Max),
157 fCovariance55Max(c.fCovariance55Max),
158 fhCutStatistics(c.fhCutStatistics),
159 fhCutCorrelation(c.fhCutCorrelation),
160 fBitmap(c.fBitmap),
161 fhNBinsClusterTPC(c.fhNBinsClusterTPC),
162 fhNBinsClusterITS(c.fhNBinsClusterITS),
163 fhNBinsChi2TPC(c.fhNBinsChi2TPC),
164 fhNBinsChi2ITS(c.fhNBinsChi2ITS),
165 fhNBinsRefitTPC(c.fhNBinsRefitTPC),
166 fhNBinsRefitITS(c.fhNBinsRefitITS),
167 fhNBinsCovariance11(c.fhNBinsCovariance11),
168 fhNBinsCovariance22(c.fhNBinsCovariance22),
169 fhNBinsCovariance33(c.fhNBinsCovariance33),
170 fhNBinsCovariance44(c.fhNBinsCovariance44),
171 fhNBinsCovariance55(c.fhNBinsCovariance55),
172 fhBinLimClusterTPC(c.fhBinLimClusterTPC),
173 fhBinLimClusterITS(c.fhBinLimClusterITS),
174 fhBinLimChi2TPC(c.fhBinLimChi2TPC),
175 fhBinLimChi2ITS(c.fhBinLimChi2ITS),
176 fhBinLimRefitTPC(c.fhBinLimRefitTPC),
177 fhBinLimRefitITS(c.fhBinLimRefitITS),
178 fhBinLimCovariance11(c.fhBinLimCovariance11),
179 fhBinLimCovariance22(c.fhBinLimCovariance22),
180 fhBinLimCovariance33(c.fhBinLimCovariance33),
181 fhBinLimCovariance44(c.fhBinLimCovariance44),
182 fhBinLimCovariance55(c.fhBinLimCovariance55)
183{
184 //
185 // copy constructor
186 //
187 ((AliCFTrackQualityCuts &) c).Copy(*this);
188}
189//__________________________________________________________________________________
190AliCFTrackQualityCuts& AliCFTrackQualityCuts::operator=(const AliCFTrackQualityCuts& c)
191{
192 //
193 // Assignment operator
194 //
195 if (this != &c) {
196 AliCFCutBase::operator=(c) ;
197 fMinNClusterTPC = c.fMinNClusterTPC ;
198 fMinNClusterITS = c.fMinNClusterITS ;
199 fMaxChi2PerClusterTPC = c.fMaxChi2PerClusterTPC ;
200 fMaxChi2PerClusterITS = c.fMaxChi2PerClusterITS ;
201 fRequireTPCRefit = c.fRequireTPCRefit ;
202 fRequireITSRefit = c.fRequireITSRefit ;
203 fCovariance11Max = c.fCovariance11Max ;
204 fCovariance22Max = c.fCovariance22Max ;
205 fCovariance33Max = c.fCovariance33Max ;
206 fCovariance44Max = c.fCovariance44Max ;
207 fCovariance55Max = c.fCovariance55Max ;
208 fhCutStatistics = c.fhCutStatistics ;
209 fhCutCorrelation = c.fhCutCorrelation ;
210 fBitmap = c.fBitmap ;
563113d0 211 fhNBinsClusterTPC = c.fhNBinsClusterTPC ;
212 fhNBinsClusterITS = c.fhNBinsClusterITS ;
213 fhNBinsChi2TPC = c.fhNBinsChi2TPC ;
214 fhNBinsChi2ITS = c.fhNBinsChi2ITS ;
215 fhNBinsRefitTPC = c.fhNBinsRefitTPC ;
216 fhNBinsRefitITS = c.fhNBinsRefitITS ;
217 fhNBinsCovariance11 = c.fhNBinsCovariance11 ;
218 fhNBinsCovariance22 = c.fhNBinsCovariance22 ;
219 fhNBinsCovariance33 = c.fhNBinsCovariance33 ;
220 fhNBinsCovariance44 = c.fhNBinsCovariance44 ;
221 fhNBinsCovariance55 = c.fhNBinsCovariance55 ;
222 fhBinLimClusterTPC = c.fhBinLimClusterTPC ;
223 fhBinLimClusterITS = c.fhBinLimClusterITS ;
224 fhBinLimChi2TPC = c.fhBinLimChi2TPC ;
225 fhBinLimChi2ITS = c.fhBinLimChi2ITS ;
226 fhBinLimRefitTPC = c.fhBinLimRefitTPC ;
227 fhBinLimRefitITS = c.fhBinLimRefitITS ;
228 fhBinLimCovariance11 = c.fhBinLimCovariance11 ;
229 fhBinLimCovariance22 = c.fhBinLimCovariance22 ;
230 fhBinLimCovariance33 = c.fhBinLimCovariance33 ;
231 fhBinLimCovariance44 = c.fhBinLimCovariance44 ;
232 fhBinLimCovariance55 = c.fhBinLimCovariance55 ;
6cb5e116 233
563113d0 234 for (Int_t i=0; i<c.kNHist; i++){
235 for (Int_t j=0; j<c.kNStepQA; j++){
236 if(c.fhQA[i][j]) fhQA[i][j] = (TH1F*)c.fhQA[i][j]->Clone();
237 }
238 }
563113d0 239 ((AliCFTrackQualityCuts &) c).Copy(*this);
6cb5e116 240 }
563113d0 241 return *this;
242}
243//__________________________________________________________________________________
244AliCFTrackQualityCuts::~AliCFTrackQualityCuts()
245{
246 //
247 // destructor
248 //
249 if (fhCutStatistics) delete fhCutStatistics;
250 if (fhCutCorrelation) delete fhCutCorrelation;
251
252 for (Int_t i=0; i<kNHist; i++){
253 for (Int_t j=0; j<kNStepQA; j++){
254 if(fhQA[i][j]) delete fhQA[i][j];
255 }
256 }
6cb5e116 257 if(fBitmap) delete fBitmap;
563113d0 258 if(fhBinLimClusterTPC) delete fhBinLimClusterTPC;
259 if(fhBinLimClusterITS) delete fhBinLimClusterITS;
260 if(fhBinLimChi2TPC) delete fhBinLimChi2TPC;
261 if(fhBinLimChi2ITS) delete fhBinLimChi2ITS;
262 if(fhBinLimRefitTPC) delete fhBinLimRefitTPC;
263 if(fhBinLimRefitITS) delete fhBinLimRefitITS;
264 if(fhBinLimCovariance11) delete fhBinLimCovariance11;
265 if(fhBinLimCovariance22) delete fhBinLimCovariance22;
266 if(fhBinLimCovariance33) delete fhBinLimCovariance33;
267 if(fhBinLimCovariance44) delete fhBinLimCovariance44;
268 if(fhBinLimCovariance55) delete fhBinLimCovariance55;
563113d0 269}
270//__________________________________________________________________________________
271void AliCFTrackQualityCuts::Initialise()
272{
273 //
274 // sets everything to zero
275 //
276 fMinNClusterTPC = 0;
277 fMinNClusterITS = 0;
278 fMaxChi2PerClusterTPC = 0;
279 fMaxChi2PerClusterITS = 0;
280 fRequireTPCRefit = 0;
281 fRequireITSRefit = 0;
282 fCovariance11Max = 0;
283 fCovariance22Max = 0;
284 fCovariance33Max = 0;
285 fCovariance44Max = 0;
286 fCovariance55Max = 0;
287
288 SetMinNClusterTPC();
289 SetMinNClusterITS();
290 SetMaxChi2PerClusterTPC();
291 SetMaxChi2PerClusterITS();
292 SetRequireTPCRefit();
293 SetRequireITSRefit();
294 SetMaxCovDiagonalElements();
295
296 for (Int_t i=0; i<kNHist; i++){
6cb5e116 297 for (Int_t j=0; j<kNStepQA; j++)
563113d0 298 fhQA[i][j] = 0x0;
563113d0 299 }
300 fhCutStatistics = 0;
301 fhCutCorrelation = 0;
6cb5e116 302 fBitmap=new TBits(0);
563113d0 303
6cb5e116 304 //set default bining for QA histograms
563113d0 305 SetHistogramBins(kCutClusterTPC,165,-0.5,164.5);
306 SetHistogramBins(kCutClusterITS,8,-0.5,7.5);
6cb5e116 307 SetHistogramBins(kCutChi2TPC,500,0.,10.);
308 SetHistogramBins(kCutChi2ITS,500,0.,10.);
309 SetHistogramBins(kCutRefitTPC,5,-0.75,1.75);
310 SetHistogramBins(kCutRefitITS,5,-0.75,1.75);
311 SetHistogramBins(kCutCovElement11,200,0.,20.);
312 SetHistogramBins(kCutCovElement22,200,0.,20.);
313 SetHistogramBins(kCutCovElement33,100,0.,1.);
314 SetHistogramBins(kCutCovElement44,100,0.,5.);
315 SetHistogramBins(kCutCovElement55,100,0.,5.);
563113d0 316}
317//__________________________________________________________________________________
318void AliCFTrackQualityCuts::Copy(TObject &c) const
319{
320 //
321 // Copy function
322 //
323 AliCFTrackQualityCuts& target = (AliCFTrackQualityCuts &) c;
324
325 target.Initialise();
326
327 if (fhCutStatistics) target.fhCutStatistics = (TH1F*) fhCutStatistics->Clone();
328 if (fhCutCorrelation) target.fhCutCorrelation = (TH2F*) fhCutCorrelation->Clone();
329
330 for (Int_t i=0; i<kNHist; i++){
331 for (Int_t j=0; j<kNStepQA; j++){
332 if(fhQA[i][j]) target.fhQA[i][j] = (TH1F*)fhQA[i][j]->Clone();
333 }
334 }
563113d0 335 TNamed::Copy(c);
336}
337//__________________________________________________________________________________
338void AliCFTrackQualityCuts::GetBitMap(TObject* obj, TBits *bitmap) {
339 //
340 // retrieve the pointer to the bitmap
341 //
6cb5e116 342 TBits *bm = SelectionBitMap(obj);
343 *bitmap = *bm;
563113d0 344}
345//__________________________________________________________________________________
346TBits* AliCFTrackQualityCuts::SelectionBitMap(TObject* obj)
347{
348 //
349 // test if the track passes the single cuts
350 // and store the information in a bitmap
351 //
352
353 // bitmap stores the decision of each single cut
354 for(Int_t i=0; i<kNCuts; i++)fBitmap->SetBitNumber(i,kFALSE);
355
356 // cast TObject into ESDtrack
357 AliESDtrack* esdTrack = dynamic_cast<AliESDtrack *>(obj);
358 if ( !esdTrack ) return fBitmap ;
359
360 for(Int_t i=0; i<kNCuts; i++)fBitmap->SetBitNumber(i,kTRUE);
361
362 // get cut quantities
363 Int_t fIdxInt[200];
364 Int_t nClustersTPC = esdTrack->GetTPCclusters(fIdxInt);
365 Int_t nClustersITS = esdTrack->GetITSclusters(fIdxInt);
366
367 Float_t chi2PerClusterTPC = -1;
368 Float_t chi2PerClusterITS = -1;
369
370 if (nClustersTPC != 0)
371 chi2PerClusterTPC = esdTrack->GetTPCchi2() / Float_t(nClustersTPC);
372 if (nClustersITS!=0)
373 chi2PerClusterITS = esdTrack->GetITSchi2() / Float_t(nClustersITS);
374
375 Double_t extCov[15];
376 esdTrack->GetExternalCovariance(extCov);
377
563113d0 378 // fill the bitmap
379 Int_t iCutBit = 0;
380
381 if (nClustersTPC < fMinNClusterTPC)
382 fBitmap->SetBitNumber(iCutBit,kFALSE);
383 iCutBit++;
384 if (nClustersITS < fMinNClusterITS)
385 fBitmap->SetBitNumber(iCutBit,kFALSE);
386 iCutBit++;
387 if (chi2PerClusterTPC > fMaxChi2PerClusterTPC)
388 fBitmap->SetBitNumber(iCutBit,kFALSE);
389 iCutBit++;
390 if (chi2PerClusterITS > fMaxChi2PerClusterITS)
391 fBitmap->SetBitNumber(iCutBit,kFALSE);
392 iCutBit++;
393 if (fRequireTPCRefit && (esdTrack->GetStatus()&AliESDtrack::kTPCrefit)==0)
394 fBitmap->SetBitNumber(iCutBit,kFALSE);
395 iCutBit++;
396 if (fRequireITSRefit && (esdTrack->GetStatus()&AliESDtrack::kITSrefit)==0)
397 fBitmap->SetBitNumber(iCutBit,kFALSE);
398 iCutBit++;
399 if (extCov[0] > fCovariance11Max)
400 fBitmap->SetBitNumber(iCutBit,kFALSE);
401 iCutBit++;
402 if (extCov[2] > fCovariance22Max)
403 fBitmap->SetBitNumber(iCutBit,kFALSE);
404 iCutBit++;
405 if (extCov[5] > fCovariance33Max)
406 fBitmap->SetBitNumber(iCutBit,kFALSE);
407 iCutBit++;
408 if (extCov[9] > fCovariance44Max)
409 fBitmap->SetBitNumber(iCutBit,kFALSE);
410 iCutBit++;
411 if (extCov[14] > fCovariance55Max)
412 fBitmap->SetBitNumber(iCutBit,kFALSE);
413
414 return fBitmap;
415}
416//__________________________________________________________________________________
417Bool_t AliCFTrackQualityCuts::IsSelected(TObject* obj) {
418 //
419 // loops over decisions of single cuts and returns if the track is accepted
420 //
421 TBits* bitmap = SelectionBitMap(obj);
422
423 Bool_t isSelected = kTRUE;
424
425 for (UInt_t icut=0; icut<bitmap->GetNbits();icut++)
426 if(!bitmap->TestBitNumber(icut)) isSelected = kFALSE;
427
428 return isSelected;
429}
430//__________________________________________________________________________________
431void AliCFTrackQualityCuts::Init() {
432 //
433 // initialises all histograms and the TList which holds the histograms
434 //
435 if(fIsQAOn)
436 DefineHistograms();
437}
438//__________________________________________________________________________________
439void AliCFTrackQualityCuts::SetHistogramBins(Int_t index, Int_t nbins, Double_t *bins)
440{
441 //
442 // variable bin size
443 //
444 if(!fIsQAOn) return;
445
446 switch(index){
447 case kCutClusterTPC:
db6722a5 448 fhNBinsClusterTPC=nbins+1;
563113d0 449 fhBinLimClusterTPC=new Double_t[nbins+1];
450 for(Int_t i=0;i<nbins+1;i++)fhBinLimClusterTPC[i]=bins[i];
451 break;
6cb5e116 452
563113d0 453 case kCutClusterITS:
db6722a5 454 fhNBinsClusterITS=nbins+1;
563113d0 455 fhBinLimClusterITS=new Double_t[nbins+1];
456 for(Int_t i=0;i<nbins+1;i++)fhBinLimClusterITS[i]=bins[i];
457 break;
6cb5e116 458
563113d0 459 case kCutChi2TPC:
db6722a5 460 fhNBinsChi2TPC=nbins+1;
563113d0 461 fhBinLimChi2TPC=new Double_t[nbins+1];
462 for(Int_t i=0;i<nbins+1;i++)fhBinLimChi2TPC[i]=bins[i];
463 break;
6cb5e116 464
563113d0 465 case kCutChi2ITS:
db6722a5 466 fhNBinsChi2ITS=nbins+1;
563113d0 467 fhBinLimChi2ITS=new Double_t[nbins+1];
468 for(Int_t i=0;i<nbins+1;i++)fhBinLimChi2ITS[i]=bins[i];
469 break;
6cb5e116 470
563113d0 471 case kCutRefitTPC:
db6722a5 472 fhNBinsRefitTPC=nbins+1;
563113d0 473 fhBinLimRefitTPC=new Double_t[nbins+1];
474 for(Int_t i=0;i<nbins+1;i++)fhBinLimRefitTPC[i]=bins[i];
475 break;
6cb5e116 476
563113d0 477 case kCutRefitITS:
db6722a5 478 fhNBinsRefitITS=nbins+1;
563113d0 479 fhBinLimRefitITS=new Double_t[nbins+1];
480 for(Int_t i=0;i<nbins+1;i++)fhBinLimRefitITS[i]=bins[i];
481 break;
6cb5e116 482
563113d0 483 case kCutCovElement11:
db6722a5 484 fhNBinsCovariance11=nbins+1;
563113d0 485 fhBinLimCovariance11=new Double_t[nbins+1];
486 for(Int_t i=0;i<nbins+1;i++)fhBinLimCovariance11[i]=bins[i];
487 break;
6cb5e116 488
563113d0 489 case kCutCovElement22:
db6722a5 490 fhNBinsCovariance22=nbins+1;
563113d0 491 fhBinLimCovariance22=new Double_t[nbins+1];
492 for(Int_t i=0;i<nbins+1;i++)fhBinLimCovariance22[i]=bins[i];
493 break;
6cb5e116 494
563113d0 495 case kCutCovElement33:
db6722a5 496 fhNBinsCovariance33=nbins+1;
563113d0 497 fhBinLimCovariance33=new Double_t[nbins+1];
498 for(Int_t i=0;i<nbins+1;i++)fhBinLimCovariance33[i]=bins[i];
499 break;
6cb5e116 500
563113d0 501 case kCutCovElement44:
db6722a5 502 fhNBinsCovariance44=nbins+1;
563113d0 503 fhBinLimCovariance44=new Double_t[nbins+1];
504 for(Int_t i=0;i<nbins+1;i++)fhBinLimCovariance44[i]=bins[i];
505 break;
6cb5e116 506
563113d0 507 case kCutCovElement55:
db6722a5 508 fhNBinsCovariance55=nbins+1;
563113d0 509 fhBinLimCovariance55=new Double_t[nbins+1];
510 for(Int_t i=0;i<nbins+1;i++)fhBinLimCovariance55[i]=bins[i];
511 break;
512 }
513}
514//__________________________________________________________________________________
515void AliCFTrackQualityCuts::SetHistogramBins(Int_t index, Int_t nbins, Double_t xmin, Double_t xmax)
516{
517 //
518 // fixed bin size
519 //
563113d0 520 switch(index){
521 case kCutClusterTPC:
db6722a5 522 fhNBinsClusterTPC=nbins+1;
563113d0 523 fhBinLimClusterTPC=new Double_t[nbins+1];
524 for(Int_t i=0;i<nbins+1;i++)fhBinLimClusterTPC[i]=xmin+i*(xmax-xmin)/Double_t(nbins);
525 break;
6cb5e116 526
563113d0 527 case kCutClusterITS:
db6722a5 528 fhNBinsClusterITS=nbins+1;
563113d0 529 fhBinLimClusterITS=new Double_t[nbins+1];
530 for(Int_t i=0;i<nbins+1;i++)fhBinLimClusterITS[i]=xmin+i*(xmax-xmin)/Double_t(nbins);
531 break;
6cb5e116 532
563113d0 533 case kCutChi2TPC:
db6722a5 534 fhNBinsChi2TPC=nbins+1;
563113d0 535 fhBinLimChi2TPC=new Double_t[nbins+1];
536 for(Int_t i=0;i<nbins+1;i++)fhBinLimChi2TPC[i]=xmin+i*(xmax-xmin)/Double_t(nbins);
537 break;
6cb5e116 538
563113d0 539 case kCutChi2ITS:
db6722a5 540 fhNBinsChi2ITS=nbins+1;
563113d0 541 fhBinLimChi2ITS=new Double_t[nbins+1];
542 for(Int_t i=0;i<nbins+1;i++)fhBinLimChi2ITS[i]=xmin+i*(xmax-xmin)/Double_t(nbins);
543 break;
6cb5e116 544
563113d0 545 case kCutRefitTPC:
db6722a5 546 fhNBinsRefitTPC=nbins+1;
563113d0 547 fhBinLimRefitTPC=new Double_t[nbins+1];
548 for(Int_t i=0;i<nbins+1;i++)fhBinLimRefitTPC[i]=xmin+i*(xmax-xmin)/Double_t(nbins);
549 break;
6cb5e116 550
563113d0 551 case kCutRefitITS:
db6722a5 552 fhNBinsRefitITS=nbins+1;
563113d0 553 fhBinLimRefitITS=new Double_t[nbins+1];
554 for(Int_t i=0;i<nbins+1;i++)fhBinLimRefitITS[i]=xmin+i*(xmax-xmin)/Double_t(nbins);
555 break;
6cb5e116 556
563113d0 557 case kCutCovElement11:
db6722a5 558 fhNBinsCovariance11=nbins+1;
563113d0 559 fhBinLimCovariance11=new Double_t[nbins+1];
560 for(Int_t i=0;i<nbins+1;i++)fhBinLimCovariance11[i]=xmin+i*(xmax-xmin)/Double_t(nbins);
561 break;
6cb5e116 562
563113d0 563 case kCutCovElement22:
db6722a5 564 fhNBinsCovariance22=nbins+1;
563113d0 565 fhBinLimCovariance22=new Double_t[nbins+1];
566 for(Int_t i=0;i<nbins+1;i++)fhBinLimCovariance22[i]=xmin+i*(xmax-xmin)/Double_t(nbins);
567 break;
6cb5e116 568
563113d0 569 case kCutCovElement33:
db6722a5 570 fhNBinsCovariance33=nbins+1;
563113d0 571 fhBinLimCovariance33=new Double_t[nbins+1];
572 for(Int_t i=0;i<nbins+1;i++)fhBinLimCovariance33[i]=xmin+i*(xmax-xmin)/Double_t(nbins);
573 break;
6cb5e116 574
563113d0 575 case kCutCovElement44:
db6722a5 576 fhNBinsCovariance44=nbins+1;
563113d0 577 fhBinLimCovariance44=new Double_t[nbins+1];
578 for(Int_t i=0;i<nbins+1;i++)fhBinLimCovariance44[i]=xmin+i*(xmax-xmin)/Double_t(nbins);
579 break;
6cb5e116 580
563113d0 581 case kCutCovElement55:
db6722a5 582 fhNBinsCovariance55=nbins+1;
563113d0 583 fhBinLimCovariance55=new Double_t[nbins+1];
584 for(Int_t i=0;i<nbins+1;i++)fhBinLimCovariance55[i]=xmin+i*(xmax-xmin)/Double_t(nbins);
585 break;
586 }
587}
588//__________________________________________________________________________________
589 void AliCFTrackQualityCuts::DefineHistograms() {
590 //
591 // histograms for cut variables, cut statistics and cut correlations
592 //
563113d0 593 Int_t color = 2;
594
595 // book cut statistics and cut correlation histograms
596 fhCutStatistics = new TH1F(Form("%s_cut_statistics",GetName()), Form("%s cut statistics",GetName()), kNCuts,0.5,kNCuts+0.5);
597 fhCutStatistics->SetLineWidth(2);
598 fhCutStatistics->GetXaxis()->SetBinLabel(1,"nClustersTPC");
599 fhCutStatistics->GetXaxis()->SetBinLabel(2,"nClustersITS");
600 fhCutStatistics->GetXaxis()->SetBinLabel(3,"chi2PerClusterTPC");
601 fhCutStatistics->GetXaxis()->SetBinLabel(4,"chi2PerClusterITS");
602 fhCutStatistics->GetXaxis()->SetBinLabel(5,"refitTPC");
603 fhCutStatistics->GetXaxis()->SetBinLabel(6,"refitITS");
604 fhCutStatistics->GetXaxis()->SetBinLabel(7,"covElement11");
605 fhCutStatistics->GetXaxis()->SetBinLabel(8,"covElement22");
606 fhCutStatistics->GetXaxis()->SetBinLabel(9,"covElement33");
607 fhCutStatistics->GetXaxis()->SetBinLabel(10,"covElement44");
608 fhCutStatistics->GetXaxis()->SetBinLabel(11,"covElement55");
609
610 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);
611 fhCutCorrelation->SetLineWidth(2);
612 fhCutCorrelation->GetXaxis()->SetBinLabel(1,"nClustersTPC");
613 fhCutCorrelation->GetXaxis()->SetBinLabel(2,"nClustersITS");
614 fhCutCorrelation->GetXaxis()->SetBinLabel(3,"chi2PerClusterTPC");
615 fhCutCorrelation->GetXaxis()->SetBinLabel(4,"chi2PerClusterITS");
616 fhCutCorrelation->GetXaxis()->SetBinLabel(5,"refitTPC");
617 fhCutCorrelation->GetXaxis()->SetBinLabel(6,"refitITS");
618 fhCutCorrelation->GetXaxis()->SetBinLabel(7,"covElement11");
619 fhCutCorrelation->GetXaxis()->SetBinLabel(8,"covElement22");
620 fhCutCorrelation->GetXaxis()->SetBinLabel(9,"covElement33");
621 fhCutCorrelation->GetXaxis()->SetBinLabel(10,"covElement44");
622 fhCutCorrelation->GetXaxis()->SetBinLabel(11,"covElement55");
623
624 fhCutCorrelation->GetYaxis()->SetBinLabel(1,"nClustersTPC");
625 fhCutCorrelation->GetYaxis()->SetBinLabel(2,"nClustersITS");
626 fhCutCorrelation->GetYaxis()->SetBinLabel(3,"chi2PerClusterTPC");
627 fhCutCorrelation->GetYaxis()->SetBinLabel(4,"chi2PerClusterITS");
628 fhCutCorrelation->GetYaxis()->SetBinLabel(5,"refitTPC");
629 fhCutCorrelation->GetYaxis()->SetBinLabel(6,"refitITS");
630 fhCutCorrelation->GetYaxis()->SetBinLabel(7,"covElement11");
631 fhCutCorrelation->GetYaxis()->SetBinLabel(8,"covElement22");
632 fhCutCorrelation->GetYaxis()->SetBinLabel(9,"covElement33");
633 fhCutCorrelation->GetYaxis()->SetBinLabel(10,"covElement44");
634 fhCutCorrelation->GetYaxis()->SetBinLabel(11,"covElement55");
635
636
637 // book QA histograms
638 Char_t str[256];
639 for (Int_t i=0; i<kNStepQA; i++) {
640 if (i==0) sprintf(str," ");
641 else sprintf(str,"_cut");
6cb5e116 642
db6722a5 643 fhQA[kCutClusterTPC][i] = new TH1F(Form("%s_nClustersTPC%s",GetName(),str) ,"",fhNBinsClusterTPC-1,fhBinLimClusterTPC);
644 fhQA[kCutClusterITS][i] = new TH1F(Form("%s_nClustersITS%s",GetName(),str) ,"",fhNBinsClusterITS-1,fhBinLimClusterITS);
645 fhQA[kCutChi2TPC][i] = new TH1F(Form("%s_chi2PerClusterTPC%s",GetName(),str),"",fhNBinsChi2TPC-1,fhBinLimChi2TPC);
646 fhQA[kCutChi2ITS][i] = new TH1F(Form("%s_chi2PerClusterITS%s",GetName(),str),"",fhNBinsChi2ITS-1,fhBinLimChi2ITS);
647 fhQA[kCutRefitTPC][i] = new TH1F(Form("%s_refitTPC%s",GetName(),str) ,"",fhNBinsRefitTPC-1,fhBinLimRefitTPC);
648 fhQA[kCutRefitITS][i] = new TH1F(Form("%s_refitITS%s",GetName(),str) ,"",fhNBinsRefitITS-1,fhBinLimRefitITS);
649 fhQA[kCutCovElement11][i] = new TH1F(Form("%s_covMatrixDiagonal11%s",GetName(),str),"",fhNBinsCovariance11-1,fhBinLimCovariance11);
650 fhQA[kCutCovElement22][i] = new TH1F(Form("%s_covMatrixDiagonal22%s",GetName(),str),"",fhNBinsCovariance22-1,fhBinLimCovariance22);
651 fhQA[kCutCovElement33][i] = new TH1F(Form("%s_covMatrixDiagonal33%s",GetName(),str),"",fhNBinsCovariance33-1,fhBinLimCovariance33);
652 fhQA[kCutCovElement44][i] = new TH1F(Form("%s_covMatrixDiagonal44%s",GetName(),str),"",fhNBinsCovariance44-1,fhBinLimCovariance44);
653 fhQA[kCutCovElement55][i] = new TH1F(Form("%s_covMatrixDiagonal55%s",GetName(),str),"",fhNBinsCovariance55-1,fhBinLimCovariance55);
563113d0 654
655 fhQA[kCutClusterTPC][i] ->SetXTitle("n TPC clusters");
656 fhQA[kCutClusterITS][i] ->SetXTitle("n ITS clusters");
657 fhQA[kCutChi2TPC][i] ->SetXTitle("#chi^{2} per TPC cluster");
658 fhQA[kCutChi2ITS][i] ->SetXTitle("#chi^{2} per ITS cluster");
659 fhQA[kCutRefitTPC][i] ->SetXTitle("refit TPC");
660 fhQA[kCutRefitITS][i] ->SetXTitle("refit ITS");
661 fhQA[kCutCovElement11][i] ->SetXTitle("cov 11 : #sigma_{y}^{2} (cm^{2})");
662 fhQA[kCutCovElement22][i] ->SetXTitle("cov 22 : #sigma_{z}^{2} (cm^{2})");
663 fhQA[kCutCovElement33][i] ->SetXTitle("cov 33 : #sigma_{sin(#phi)}^{2}");
664 fhQA[kCutCovElement44][i] ->SetXTitle("cov 44 : #sigma_{tan(#theta_{dip})}^{2}");
6cb5e116 665 fhQA[kCutCovElement55][i] ->SetXTitle("cov 55 : #sigma_{1/p_{T}}^{2} ((c/GeV)^{2})");
563113d0 666 }
667
668 for(Int_t i=0; i<kNHist; i++) fhQA[i][1]->SetLineColor(color);
563113d0 669}
670//__________________________________________________________________________________
671void AliCFTrackQualityCuts::FillHistograms(TObject* obj, Bool_t b)
672{
673 //
674 // fill the QA histograms
675 //
676 if(!fIsQAOn) return;
677
678 // cast TObject into ESDtrack
679 AliESDtrack* esdTrack = dynamic_cast<AliESDtrack *>(obj);
680 if ( !esdTrack ) return;
681
682 // index = 0: fill histograms before cuts
683 // index = 1: fill histograms after cuts
684 Int_t index = -1;
685 index = ((b) ? 1 : 0);
686
687 Int_t fIdxInt[200];
688 Int_t nClustersTPC = esdTrack->GetTPCclusters(fIdxInt);
689 fhQA[kCutClusterTPC][index]->Fill((float)nClustersTPC);
690 Float_t chi2PerClusterTPC = -1.;
691 if (nClustersTPC!=0)
692 chi2PerClusterTPC = esdTrack->GetTPCchi2()/((float)nClustersTPC);
693 fhQA[kCutChi2TPC][index]->Fill(chi2PerClusterTPC);
6cb5e116 694
563113d0 695 Int_t nClustersITS = esdTrack->GetITSclusters(fIdxInt);
696 fhQA[kCutClusterITS][index]->Fill((float)nClustersITS);
697 Float_t chi2PerClusterITS = -1.;
698 if (nClustersITS!=0)
699 chi2PerClusterITS = esdTrack->GetITSchi2()/((float)nClustersITS);
700 fhQA[kCutChi2ITS][index]->Fill(chi2PerClusterITS);
6cb5e116 701
563113d0 702 if ((esdTrack->GetStatus()&AliESDtrack::kTPCrefit)==0)
703 fhQA[kCutRefitTPC][index]->Fill(0.);
704 if (!((esdTrack->GetStatus()&AliESDtrack::kTPCrefit)==0))
705 fhQA[kCutRefitTPC][index]->Fill(1.);
6cb5e116 706
563113d0 707 if ((esdTrack->GetStatus()&AliESDtrack::kITSrefit)==0)
708 fhQA[kCutRefitITS][index]->Fill(0.);
709 if (!((esdTrack->GetStatus()&AliESDtrack::kITSrefit)==0))
710 fhQA[kCutRefitITS][index]->Fill(1.);
6cb5e116 711
563113d0 712 Double_t extCov[15];
713 esdTrack->GetExternalCovariance(extCov);
6cb5e116 714
563113d0 715 fhQA[kCutCovElement11][index]->Fill(extCov[0]);
716 fhQA[kCutCovElement22][index]->Fill(extCov[2]);
717 fhQA[kCutCovElement33][index]->Fill(extCov[5]);
718 fhQA[kCutCovElement44][index]->Fill(extCov[9]);
719 fhQA[kCutCovElement55][index]->Fill(extCov[14]);
6cb5e116 720
563113d0 721 // fill cut statistics and cut correlation histograms with information from the bitmap
722 if (b) return;
723
724 // Get the bitmap of the single cuts
725 if ( !obj ) return;
726 TBits* bitmap = SelectionBitMap(obj);
727
728 // Number of single cuts in this class
729 UInt_t ncuts = bitmap->GetNbits();
730 for(UInt_t bit=0; bit<ncuts;bit++) {
731 if (!bitmap->TestBitNumber(bit)) {
732 fhCutStatistics->Fill(bit+1);
733 for (UInt_t bit2=bit; bit2<ncuts;bit2++) {
734 if (!bitmap->TestBitNumber(bit2))
735 fhCutCorrelation->Fill(bit+1,bit2+1);
736 }
737 }
738 }
739}
740//__________________________________________________________________________________
741void AliCFTrackQualityCuts::SaveHistograms(const Char_t* dir) {
742 //
743 // saves the histograms in a directory (dir)
744 //
745 if(!fIsQAOn) return;
746
747 if (!dir)
748 dir = GetName();
749
750 gDirectory->mkdir(dir);
751 gDirectory->cd(dir);
752
753 gDirectory->mkdir("before_cuts");
754 gDirectory->mkdir("after_cuts");
755
756 fhCutStatistics->Write();
757 fhCutCorrelation->Write();
758
759 for (Int_t j=0; j<kNStepQA; j++) {
760 if (j==0)
761 gDirectory->cd("before_cuts");
762 else
763 gDirectory->cd("after_cuts");
764
765 for(Int_t i=0; i<kNHist; i++) fhQA[i][j]->Write();
766
767 gDirectory->cd("../");
768 }
769 gDirectory->cd("../");
770}
771//__________________________________________________________________________________
772void AliCFTrackQualityCuts::DrawHistograms(Bool_t drawLogScale)
773{
774 //
775 // draws some histograms
776 //
777 if(!fIsQAOn) return;
778
779 // pad margins
780 Float_t right = 0.03;
781 Float_t left = 0.175;
782 Float_t top = 0.03;
783 Float_t bottom = 0.175;
784
785 TCanvas* canvas1 = new TCanvas("Track_QA_Quality_1", "Track QA Quality 1", 800, 500);
786 canvas1->Divide(2, 1);
787
788 canvas1->cd(1);
789 fhCutStatistics->SetStats(kFALSE);
790 fhCutStatistics->LabelsOption("v");
791 gPad->SetLeftMargin(left);
792 gPad->SetBottomMargin(0.25);
793 gPad->SetRightMargin(right);
794 gPad->SetTopMargin(0.1);
795 fhCutStatistics->Draw();
796
797 canvas1->cd(2);
798 fhCutCorrelation->SetStats(kFALSE);
799 fhCutCorrelation->LabelsOption("v");
800 gPad->SetLeftMargin(0.30);
801 gPad->SetRightMargin(bottom);
802 gPad->SetTopMargin(0.1);
803 gPad->SetBottomMargin(0.25);
804 fhCutCorrelation->Draw("COLZ");
805
806 canvas1->SaveAs(Form("%s.eps", canvas1->GetName()));
807 canvas1->SaveAs(Form("%s.ps", canvas1->GetName()));
808
809 // -----
810
811 TCanvas* canvas2 = new TCanvas("Track_QA_Quality_2", "Track QA Quality 2", 1200, 800);
812 canvas2->Divide(3, 2);
813
814 canvas2->cd(1);
815 gPad->SetRightMargin(right);
816 gPad->SetLeftMargin(left);
817 gPad->SetTopMargin(top);
818 gPad->SetBottomMargin(bottom);
819 fhQA[kCutClusterTPC][0]->SetStats(kFALSE);
820 fhQA[kCutClusterTPC][0]->Draw();
821 fhQA[kCutClusterTPC][1]->Draw("same");
822
823 canvas2->cd(2);
824 gPad->SetRightMargin(right);
825 gPad->SetLeftMargin(left);
826 gPad->SetTopMargin(top);
827 gPad->SetBottomMargin(bottom);
828 fhQA[kCutChi2TPC][0]->SetStats(kFALSE);
829 fhQA[kCutChi2TPC][0]->Draw();
830 fhQA[kCutChi2TPC][1]->Draw("same");
831
832 canvas2->cd(3);
833 gPad->SetRightMargin(right);
834 gPad->SetLeftMargin(left);
835 gPad->SetTopMargin(top);
836 gPad->SetBottomMargin(bottom);
6cb5e116 837 fhQA[kCutRefitTPC][0]->SetStats(kFALSE);
563113d0 838 fhQA[kCutRefitTPC][0]->Draw();
839 fhQA[kCutRefitTPC][1]->Draw("same");
840
841 canvas2->cd(4);
842 gPad->SetRightMargin(right);
843 gPad->SetLeftMargin(left);
844 gPad->SetTopMargin(top);
845 gPad->SetBottomMargin(bottom);
846 fhQA[kCutClusterITS][0]->SetStats(kFALSE);
847 fhQA[kCutClusterITS][0]->Draw();
848 fhQA[kCutClusterITS][1]->Draw("same");
849
850 canvas2->cd(5);
851 gPad->SetRightMargin(right);
852 gPad->SetLeftMargin(left);
853 gPad->SetTopMargin(top);
854 gPad->SetBottomMargin(bottom);
855 fhQA[kCutChi2ITS][0]->SetStats(kFALSE);
856 fhQA[kCutChi2ITS][0]->Draw();
857 fhQA[kCutChi2ITS][1]->Draw("same");
858
859 canvas2->cd(6);
860 gPad->SetRightMargin(right);
861 gPad->SetLeftMargin(left);
862 gPad->SetTopMargin(top);
863 gPad->SetBottomMargin(bottom);
6cb5e116 864 fhQA[kCutRefitITS][0]->SetStats(kFALSE);
563113d0 865 fhQA[kCutRefitITS][0]->Draw();
866 fhQA[kCutRefitITS][1]->Draw("same");
867
868 canvas2->SaveAs(Form("%s.eps", canvas2->GetName()));
869 canvas2->SaveAs(Form("%s.ps", canvas2->GetName()));
870
871 // -----
872
873 TCanvas* canvas3 = new TCanvas("Track_QA_Quality_3", "Track QA Quality 3", 1200, 800);
874 canvas3->Divide(3, 2);
875
876 canvas3->cd(1);
877 if(drawLogScale) gPad->SetLogy();
878 gPad->SetRightMargin(right);
879 gPad->SetLeftMargin(left);
880 gPad->SetTopMargin(top);
881 gPad->SetBottomMargin(bottom);
882 fhQA[kCutCovElement11][0]->SetStats(kFALSE);
883 fhQA[kCutCovElement11][0]->Draw();
884 fhQA[kCutCovElement11][1]->Draw("same");
885
886 canvas3->cd(2);
887 if(drawLogScale) gPad->SetLogy();
888 gPad->SetRightMargin(right);
889 gPad->SetLeftMargin(left);
890 gPad->SetTopMargin(top);
891 gPad->SetBottomMargin(bottom);
892 fhQA[kCutCovElement22][0]->SetStats(kFALSE);
893 fhQA[kCutCovElement22][0]->Draw();
894 fhQA[kCutCovElement22][1]->Draw("same");
895
896 canvas3->cd(3);
897 if(drawLogScale) gPad->SetLogy();
898 gPad->SetRightMargin(right);
899 gPad->SetLeftMargin(left);
900 gPad->SetTopMargin(top);
901 gPad->SetBottomMargin(bottom);
902 fhQA[kCutCovElement33][0]->SetStats(kFALSE);
903 fhQA[kCutCovElement33][0]->Draw();
904 fhQA[kCutCovElement33][1]->Draw("same");
905
906 canvas3->cd(4);
907 if(drawLogScale) gPad->SetLogy();
908 gPad->SetRightMargin(right);
909 gPad->SetLeftMargin(left);
910 gPad->SetTopMargin(top);
911 gPad->SetBottomMargin(bottom);
912 fhQA[kCutCovElement44][0]->SetStats(kFALSE);
913 fhQA[kCutCovElement44][0]->Draw();
914 fhQA[kCutCovElement44][1]->Draw("same");
915
916 canvas3->cd(5);
917 if(drawLogScale) gPad->SetLogy();
918 gPad->SetRightMargin(right);
919 gPad->SetLeftMargin(left);
920 gPad->SetTopMargin(top);
921 gPad->SetBottomMargin(bottom);
922 fhQA[kCutCovElement55][0]->SetStats(kFALSE);
923 fhQA[kCutCovElement55][0]->Draw();
924 fhQA[kCutCovElement55][1]->Draw("same");
925
926 canvas3->SaveAs(Form("%s.eps", canvas3->GetName()));
927 canvas3->SaveAs(Form("%s.ps", canvas3->GetName()));
928}
929//__________________________________________________________________________________
930void AliCFTrackQualityCuts::AddQAHistograms(TList *qaList) const {
931 //
932 // saves the histograms in a TList
933 //
934 if(!fIsQAOn) return;
935
936 qaList->Add(fhCutStatistics);
937 qaList->Add(fhCutCorrelation);
938
939 for (Int_t j=0; j<kNStepQA; j++) {
940 for(Int_t i=0; i<kNHist; i++)
941 qaList->Add(fhQA[i][j]);
942 }
943}