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