1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
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
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
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.
37 // author: I. Kraus (Ingrid.Kraus@cern.ch)
39 // AliESDtrackCuts writte by Jan Fiete Grosse-Oetringhaus and
40 // AliRsnDaughterCut class written by A. Pulvirenti.
44 #include <TDirectory.h>
48 #include <AliESDtrack.h>
50 #include "AliCFTrackQualityCuts.h"
51 #include "AliAODTrack.h"
53 ClassImp(AliCFTrackQualityCuts)
55 //__________________________________________________________________________________
56 AliCFTrackQualityCuts::AliCFTrackQualityCuts() :
60 fMaxChi2PerClusterTPC(0),
61 fMaxChi2PerClusterITS(0),
75 fhNBinsCovariance11(0),
76 fhNBinsCovariance22(0),
77 fhNBinsCovariance33(0),
78 fhNBinsCovariance44(0),
79 fhNBinsCovariance55(0),
80 fhBinLimClusterTPC(0x0),
81 fhBinLimClusterITS(0x0),
84 fhBinLimCovariance11(0x0),
85 fhBinLimCovariance22(0x0),
86 fhBinLimCovariance33(0x0),
87 fhBinLimCovariance44(0x0),
88 fhBinLimCovariance55(0x0)
91 // Default constructor
95 //__________________________________________________________________________________
96 AliCFTrackQualityCuts::AliCFTrackQualityCuts(Char_t* name, Char_t* title) :
97 AliCFCutBase(name,title),
100 fMaxChi2PerClusterTPC(0),
101 fMaxChi2PerClusterITS(0),
111 fhNBinsClusterTPC(0),
112 fhNBinsClusterITS(0),
115 fhNBinsCovariance11(0),
116 fhNBinsCovariance22(0),
117 fhNBinsCovariance33(0),
118 fhNBinsCovariance44(0),
119 fhNBinsCovariance55(0),
120 fhBinLimClusterTPC(0x0),
121 fhBinLimClusterITS(0x0),
122 fhBinLimChi2TPC(0x0),
123 fhBinLimChi2ITS(0x0),
124 fhBinLimCovariance11(0x0),
125 fhBinLimCovariance22(0x0),
126 fhBinLimCovariance33(0x0),
127 fhBinLimCovariance44(0x0),
128 fhBinLimCovariance55(0x0)
135 //__________________________________________________________________________________
136 AliCFTrackQualityCuts::AliCFTrackQualityCuts(const AliCFTrackQualityCuts& c) :
138 fMinNClusterTPC(c.fMinNClusterTPC),
139 fMinNClusterITS(c.fMinNClusterITS),
140 fMaxChi2PerClusterTPC(c.fMaxChi2PerClusterTPC),
141 fMaxChi2PerClusterITS(c.fMaxChi2PerClusterITS),
142 fCovariance11Max(c.fCovariance11Max),
143 fCovariance22Max(c.fCovariance22Max),
144 fCovariance33Max(c.fCovariance33Max),
145 fCovariance44Max(c.fCovariance44Max),
146 fCovariance55Max(c.fCovariance55Max),
148 fhCutStatistics(c.fhCutStatistics),
149 fhCutCorrelation(c.fhCutCorrelation),
151 fhNBinsClusterTPC(c.fhNBinsClusterTPC),
152 fhNBinsClusterITS(c.fhNBinsClusterITS),
153 fhNBinsChi2TPC(c.fhNBinsChi2TPC),
154 fhNBinsChi2ITS(c.fhNBinsChi2ITS),
155 fhNBinsCovariance11(c.fhNBinsCovariance11),
156 fhNBinsCovariance22(c.fhNBinsCovariance22),
157 fhNBinsCovariance33(c.fhNBinsCovariance33),
158 fhNBinsCovariance44(c.fhNBinsCovariance44),
159 fhNBinsCovariance55(c.fhNBinsCovariance55),
160 fhBinLimClusterTPC(c.fhBinLimClusterTPC),
161 fhBinLimClusterITS(c.fhBinLimClusterITS),
162 fhBinLimChi2TPC(c.fhBinLimChi2TPC),
163 fhBinLimChi2ITS(c.fhBinLimChi2ITS),
164 fhBinLimCovariance11(c.fhBinLimCovariance11),
165 fhBinLimCovariance22(c.fhBinLimCovariance22),
166 fhBinLimCovariance33(c.fhBinLimCovariance33),
167 fhBinLimCovariance44(c.fhBinLimCovariance44),
168 fhBinLimCovariance55(c.fhBinLimCovariance55)
173 ((AliCFTrackQualityCuts &) c).Copy(*this);
175 //__________________________________________________________________________________
176 AliCFTrackQualityCuts& AliCFTrackQualityCuts::operator=(const AliCFTrackQualityCuts& c)
179 // Assignment operator
182 AliCFCutBase::operator=(c) ;
183 fMinNClusterTPC = c.fMinNClusterTPC ;
184 fMinNClusterITS = c.fMinNClusterITS ;
185 fMaxChi2PerClusterTPC = c.fMaxChi2PerClusterTPC ;
186 fMaxChi2PerClusterITS = c.fMaxChi2PerClusterITS ;
187 fCovariance11Max = c.fCovariance11Max ;
188 fCovariance22Max = c.fCovariance22Max ;
189 fCovariance33Max = c.fCovariance33Max ;
190 fCovariance44Max = c.fCovariance44Max ;
191 fCovariance55Max = c.fCovariance55Max ;
192 fStatus = c.fStatus ;
193 fhCutStatistics = c.fhCutStatistics ;
194 fhCutCorrelation = c.fhCutCorrelation ;
195 fBitmap = c.fBitmap ;
196 fhNBinsClusterTPC = c.fhNBinsClusterTPC ;
197 fhNBinsClusterITS = c.fhNBinsClusterITS ;
198 fhNBinsChi2TPC = c.fhNBinsChi2TPC ;
199 fhNBinsChi2ITS = c.fhNBinsChi2ITS ;
200 fhNBinsCovariance11 = c.fhNBinsCovariance11 ;
201 fhNBinsCovariance22 = c.fhNBinsCovariance22 ;
202 fhNBinsCovariance33 = c.fhNBinsCovariance33 ;
203 fhNBinsCovariance44 = c.fhNBinsCovariance44 ;
204 fhNBinsCovariance55 = c.fhNBinsCovariance55 ;
205 fhBinLimClusterTPC = c.fhBinLimClusterTPC ;
206 fhBinLimClusterITS = c.fhBinLimClusterITS ;
207 fhBinLimChi2TPC = c.fhBinLimChi2TPC ;
208 fhBinLimChi2ITS = c.fhBinLimChi2ITS ;
209 fhBinLimCovariance11 = c.fhBinLimCovariance11 ;
210 fhBinLimCovariance22 = c.fhBinLimCovariance22 ;
211 fhBinLimCovariance33 = c.fhBinLimCovariance33 ;
212 fhBinLimCovariance44 = c.fhBinLimCovariance44 ;
213 fhBinLimCovariance55 = c.fhBinLimCovariance55 ;
215 for (Int_t i=0; i<c.kNHist; i++){
216 for (Int_t j=0; j<c.kNStepQA; j++){
217 if(c.fhQA[i][j]) fhQA[i][j] = (TH1F*)c.fhQA[i][j]->Clone();
220 ((AliCFTrackQualityCuts &) c).Copy(*this);
224 //__________________________________________________________________________________
225 AliCFTrackQualityCuts::~AliCFTrackQualityCuts()
230 if (fhCutStatistics) delete fhCutStatistics;
231 if (fhCutCorrelation) delete fhCutCorrelation;
233 for (Int_t i=0; i<kNHist; i++){
234 for (Int_t j=0; j<kNStepQA; j++){
235 if(fhQA[i][j]) delete fhQA[i][j];
238 if(fBitmap) delete fBitmap;
239 if(fhBinLimClusterTPC) delete fhBinLimClusterTPC;
240 if(fhBinLimClusterITS) delete fhBinLimClusterITS;
241 if(fhBinLimChi2TPC) delete fhBinLimChi2TPC;
242 if(fhBinLimChi2ITS) delete fhBinLimChi2ITS;
243 if(fhBinLimCovariance11) delete fhBinLimCovariance11;
244 if(fhBinLimCovariance22) delete fhBinLimCovariance22;
245 if(fhBinLimCovariance33) delete fhBinLimCovariance33;
246 if(fhBinLimCovariance44) delete fhBinLimCovariance44;
247 if(fhBinLimCovariance55) delete fhBinLimCovariance55;
249 //__________________________________________________________________________________
250 void AliCFTrackQualityCuts::Initialise()
253 // sets everything to zero
257 fMaxChi2PerClusterTPC = 0;
258 fMaxChi2PerClusterITS = 0;
259 fCovariance11Max = 0;
260 fCovariance22Max = 0;
261 fCovariance33Max = 0;
262 fCovariance44Max = 0;
263 fCovariance55Max = 0;
268 SetMaxChi2PerClusterTPC();
269 SetMaxChi2PerClusterITS();
270 SetMaxCovDiagonalElements();
273 for (Int_t i=0; i<kNHist; i++){
274 for (Int_t j=0; j<kNStepQA; j++)
278 fhCutCorrelation = 0;
279 fBitmap=new TBits(0);
281 //set default bining for QA histograms
282 SetHistogramBins(kCutClusterTPC,165,-0.5,164.5);
283 SetHistogramBins(kCutClusterITS,8,-0.5,7.5);
284 SetHistogramBins(kCutChi2TPC,500,0.,10.);
285 SetHistogramBins(kCutChi2ITS,500,0.,10.);
286 SetHistogramBins(kCutCovElement11,200,0.,20.);
287 SetHistogramBins(kCutCovElement22,200,0.,20.);
288 SetHistogramBins(kCutCovElement33,100,0.,1.);
289 SetHistogramBins(kCutCovElement44,100,0.,5.);
290 SetHistogramBins(kCutCovElement55,100,0.,5.);
292 //__________________________________________________________________________________
293 void AliCFTrackQualityCuts::Copy(TObject &c) const
298 AliCFTrackQualityCuts& target = (AliCFTrackQualityCuts &) c;
302 if (fhCutStatistics) target.fhCutStatistics = (TH1F*) fhCutStatistics->Clone();
303 if (fhCutCorrelation) target.fhCutCorrelation = (TH2F*) fhCutCorrelation->Clone();
305 for (Int_t i=0; i<kNHist; i++){
306 for (Int_t j=0; j<kNStepQA; j++){
307 if(fhQA[i][j]) target.fhQA[i][j] = (TH1F*)fhQA[i][j]->Clone();
312 //__________________________________________________________________________________
313 void AliCFTrackQualityCuts::SelectionBitMap(TObject* obj)
316 // test if the track passes the single cuts
317 // and store the information in a bitmap
320 // bitmap stores the decision of each single cut
321 for(Int_t i=0; i<kNCuts; i++)fBitmap->SetBitNumber(i,kFALSE);
324 if (!obj->InheritsFrom("AliVParticle")) {
325 AliError("object must derived from AliVParticle !");
329 Bool_t isESDTrack = strcmp(obj->ClassName(),"AliESDtrack") == 0 ? kTRUE : kFALSE ;
330 Bool_t isAODTrack = strcmp(obj->ClassName(),"AliAODTrack") == 0 ? kTRUE : kFALSE ;
332 AliESDtrack * esdTrack = 0x0 ;
333 AliAODTrack * aodTrack = 0x0 ;
334 if (isESDTrack) esdTrack = dynamic_cast<AliESDtrack*>(obj);
335 if (isAODTrack) aodTrack = dynamic_cast<AliAODTrack*>(obj);
337 // get cut quantities
339 Int_t nClustersTPC = 0;
340 Int_t nClustersITS = 0 ;
341 Float_t chi2PerClusterTPC = 0 ;
342 Float_t chi2PerClusterITS = 0 ;
343 Double_t extCov[15]={0,0,0,0,0,0,0,0,0,0,0,0,0,0,0};
346 nClustersTPC = esdTrack->GetTPCclusters(fIdxInt);
347 nClustersITS = esdTrack->GetITSclusters(fIdxInt);
348 if (nClustersTPC != 0) chi2PerClusterTPC = esdTrack->GetTPCchi2() / Float_t(nClustersTPC);
349 if (nClustersITS != 0) chi2PerClusterITS = esdTrack->GetITSchi2() / Float_t(nClustersITS);
350 esdTrack->GetExternalCovariance(extCov);
356 if (nClustersTPC >= fMinNClusterTPC)
357 fBitmap->SetBitNumber(iCutBit,kTRUE);
359 if (nClustersITS >= fMinNClusterITS)
360 fBitmap->SetBitNumber(iCutBit,kTRUE);
362 if (chi2PerClusterTPC <= fMaxChi2PerClusterTPC)
363 fBitmap->SetBitNumber(iCutBit,kTRUE);
365 if (chi2PerClusterITS <= fMaxChi2PerClusterITS)
366 fBitmap->SetBitNumber(iCutBit,kTRUE);
368 if (extCov[0] <= fCovariance11Max)
369 fBitmap->SetBitNumber(iCutBit,kTRUE);
371 if (extCov[2] <= fCovariance22Max)
372 fBitmap->SetBitNumber(iCutBit,kTRUE);
374 if (extCov[5] <= fCovariance33Max)
375 fBitmap->SetBitNumber(iCutBit,kTRUE);
377 if (extCov[9] <= fCovariance44Max)
378 fBitmap->SetBitNumber(iCutBit,kTRUE);
380 if (extCov[14] <= fCovariance55Max)
381 fBitmap->SetBitNumber(iCutBit,kTRUE);
385 if ( (esdTrack->GetStatus() & fStatus) == fStatus ) fBitmap->SetBitNumber(iCutBit,kTRUE);
388 if ( (aodTrack->GetStatus() & fStatus) == fStatus ) fBitmap->SetBitNumber(iCutBit,kTRUE);
393 //__________________________________________________________________________________
394 Bool_t AliCFTrackQualityCuts::IsSelected(TObject* obj) {
396 // loops over decisions of single cuts and returns if the track is accepted
398 SelectionBitMap(obj);
400 if (fIsQAOn) FillHistograms(obj,0);
401 Bool_t isSelected = kTRUE;
403 for (UInt_t icut=0; icut<fBitmap->GetNbits();icut++) {
404 if(!fBitmap->TestBitNumber(icut)) {
409 if (!isSelected) return kFALSE ;
410 if (fIsQAOn) FillHistograms(obj,1);
413 //__________________________________________________________________________________
414 void AliCFTrackQualityCuts::SetHistogramBins(Int_t index, Int_t nbins, Double_t *bins)
421 if (index<0 || index>=kNHist) {
422 Error("SetHistogramBins","could not determine histogram from index %d",index);
428 fhNBinsClusterTPC=nbins+1;
429 fhBinLimClusterTPC=new Double_t[nbins+1];
430 for(Int_t i=0;i<nbins+1;i++)fhBinLimClusterTPC[i]=bins[i];
434 fhNBinsClusterITS=nbins+1;
435 fhBinLimClusterITS=new Double_t[nbins+1];
436 for(Int_t i=0;i<nbins+1;i++)fhBinLimClusterITS[i]=bins[i];
440 fhNBinsChi2TPC=nbins+1;
441 fhBinLimChi2TPC=new Double_t[nbins+1];
442 for(Int_t i=0;i<nbins+1;i++)fhBinLimChi2TPC[i]=bins[i];
446 fhNBinsChi2ITS=nbins+1;
447 fhBinLimChi2ITS=new Double_t[nbins+1];
448 for(Int_t i=0;i<nbins+1;i++)fhBinLimChi2ITS[i]=bins[i];
451 case kCutCovElement11:
452 fhNBinsCovariance11=nbins+1;
453 fhBinLimCovariance11=new Double_t[nbins+1];
454 for(Int_t i=0;i<nbins+1;i++)fhBinLimCovariance11[i]=bins[i];
457 case kCutCovElement22:
458 fhNBinsCovariance22=nbins+1;
459 fhBinLimCovariance22=new Double_t[nbins+1];
460 for(Int_t i=0;i<nbins+1;i++)fhBinLimCovariance22[i]=bins[i];
463 case kCutCovElement33:
464 fhNBinsCovariance33=nbins+1;
465 fhBinLimCovariance33=new Double_t[nbins+1];
466 for(Int_t i=0;i<nbins+1;i++)fhBinLimCovariance33[i]=bins[i];
469 case kCutCovElement44:
470 fhNBinsCovariance44=nbins+1;
471 fhBinLimCovariance44=new Double_t[nbins+1];
472 for(Int_t i=0;i<nbins+1;i++)fhBinLimCovariance44[i]=bins[i];
475 case kCutCovElement55:
476 fhNBinsCovariance55=nbins+1;
477 fhBinLimCovariance55=new Double_t[nbins+1];
478 for(Int_t i=0;i<nbins+1;i++)fhBinLimCovariance55[i]=bins[i];
482 //__________________________________________________________________________________
483 void AliCFTrackQualityCuts::SetHistogramBins(Int_t index, Int_t nbins, Double_t xmin, Double_t xmax)
488 if (index<0 || index>=kNHist) {
489 Error("SetHistogramBins","could not determine histogram from index %d",index);
495 fhNBinsClusterTPC=nbins+1;
496 fhBinLimClusterTPC=new Double_t[nbins+1];
497 for(Int_t i=0;i<nbins+1;i++)fhBinLimClusterTPC[i]=xmin+i*(xmax-xmin)/Double_t(nbins);
501 fhNBinsClusterITS=nbins+1;
502 fhBinLimClusterITS=new Double_t[nbins+1];
503 for(Int_t i=0;i<nbins+1;i++)fhBinLimClusterITS[i]=xmin+i*(xmax-xmin)/Double_t(nbins);
507 fhNBinsChi2TPC=nbins+1;
508 fhBinLimChi2TPC=new Double_t[nbins+1];
509 for(Int_t i=0;i<nbins+1;i++)fhBinLimChi2TPC[i]=xmin+i*(xmax-xmin)/Double_t(nbins);
513 fhNBinsChi2ITS=nbins+1;
514 fhBinLimChi2ITS=new Double_t[nbins+1];
515 for(Int_t i=0;i<nbins+1;i++)fhBinLimChi2ITS[i]=xmin+i*(xmax-xmin)/Double_t(nbins);
518 case kCutCovElement11:
519 fhNBinsCovariance11=nbins+1;
520 fhBinLimCovariance11=new Double_t[nbins+1];
521 for(Int_t i=0;i<nbins+1;i++)fhBinLimCovariance11[i]=xmin+i*(xmax-xmin)/Double_t(nbins);
524 case kCutCovElement22:
525 fhNBinsCovariance22=nbins+1;
526 fhBinLimCovariance22=new Double_t[nbins+1];
527 for(Int_t i=0;i<nbins+1;i++)fhBinLimCovariance22[i]=xmin+i*(xmax-xmin)/Double_t(nbins);
530 case kCutCovElement33:
531 fhNBinsCovariance33=nbins+1;
532 fhBinLimCovariance33=new Double_t[nbins+1];
533 for(Int_t i=0;i<nbins+1;i++)fhBinLimCovariance33[i]=xmin+i*(xmax-xmin)/Double_t(nbins);
536 case kCutCovElement44:
537 fhNBinsCovariance44=nbins+1;
538 fhBinLimCovariance44=new Double_t[nbins+1];
539 for(Int_t i=0;i<nbins+1;i++)fhBinLimCovariance44[i]=xmin+i*(xmax-xmin)/Double_t(nbins);
542 case kCutCovElement55:
543 fhNBinsCovariance55=nbins+1;
544 fhBinLimCovariance55=new Double_t[nbins+1];
545 for(Int_t i=0;i<nbins+1;i++)fhBinLimCovariance55[i]=xmin+i*(xmax-xmin)/Double_t(nbins);
549 //__________________________________________________________________________________
550 void AliCFTrackQualityCuts::DefineHistograms() {
552 // histograms for cut variables, cut statistics and cut correlations
556 // book cut statistics and cut correlation histograms
557 fhCutStatistics = new TH1F(Form("%s_cut_statistics",GetName()), Form("%s cut statistics",GetName()), kNCuts,0.5,kNCuts+0.5);
558 fhCutStatistics->SetLineWidth(2);
559 fhCutStatistics->GetXaxis()->SetBinLabel(1,"nClustersTPC");
560 fhCutStatistics->GetXaxis()->SetBinLabel(2,"nClustersITS");
561 fhCutStatistics->GetXaxis()->SetBinLabel(3,"chi2PerClusterTPC");
562 fhCutStatistics->GetXaxis()->SetBinLabel(4,"chi2PerClusterITS");
563 fhCutStatistics->GetXaxis()->SetBinLabel(5,"covElement11");
564 fhCutStatistics->GetXaxis()->SetBinLabel(6,"covElement22");
565 fhCutStatistics->GetXaxis()->SetBinLabel(7,"covElement33");
566 fhCutStatistics->GetXaxis()->SetBinLabel(8,"covElement44");
567 fhCutStatistics->GetXaxis()->SetBinLabel(9,"covElement55");
569 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);
570 fhCutCorrelation->SetLineWidth(2);
571 fhCutCorrelation->GetXaxis()->SetBinLabel(1,"nClustersTPC");
572 fhCutCorrelation->GetXaxis()->SetBinLabel(2,"nClustersITS");
573 fhCutCorrelation->GetXaxis()->SetBinLabel(3,"chi2PerClusterTPC");
574 fhCutCorrelation->GetXaxis()->SetBinLabel(4,"chi2PerClusterITS");
575 fhCutCorrelation->GetXaxis()->SetBinLabel(5,"covElement11");
576 fhCutCorrelation->GetXaxis()->SetBinLabel(6,"covElement22");
577 fhCutCorrelation->GetXaxis()->SetBinLabel(7,"covElement33");
578 fhCutCorrelation->GetXaxis()->SetBinLabel(8,"covElement44");
579 fhCutCorrelation->GetXaxis()->SetBinLabel(9,"covElement55");
581 fhCutCorrelation->GetYaxis()->SetBinLabel(1,"nClustersTPC");
582 fhCutCorrelation->GetYaxis()->SetBinLabel(2,"nClustersITS");
583 fhCutCorrelation->GetYaxis()->SetBinLabel(3,"chi2PerClusterTPC");
584 fhCutCorrelation->GetYaxis()->SetBinLabel(4,"chi2PerClusterITS");
585 fhCutCorrelation->GetYaxis()->SetBinLabel(5,"covElement11");
586 fhCutCorrelation->GetYaxis()->SetBinLabel(6,"covElement22");
587 fhCutCorrelation->GetYaxis()->SetBinLabel(7,"covElement33");
588 fhCutCorrelation->GetYaxis()->SetBinLabel(8,"covElement44");
589 fhCutCorrelation->GetYaxis()->SetBinLabel(9,"covElement55");
592 // book QA histograms
594 for (Int_t i=0; i<kNStepQA; i++) {
595 if (i==0) sprintf(str," ");
596 else sprintf(str,"_cut");
598 fhQA[kCutClusterTPC][i] = new TH1F(Form("%s_nClustersTPC%s",GetName(),str) ,"",fhNBinsClusterTPC-1,fhBinLimClusterTPC);
599 fhQA[kCutClusterITS][i] = new TH1F(Form("%s_nClustersITS%s",GetName(),str) ,"",fhNBinsClusterITS-1,fhBinLimClusterITS);
600 fhQA[kCutChi2TPC][i] = new TH1F(Form("%s_chi2PerClusterTPC%s",GetName(),str),"",fhNBinsChi2TPC-1,fhBinLimChi2TPC);
601 fhQA[kCutChi2ITS][i] = new TH1F(Form("%s_chi2PerClusterITS%s",GetName(),str),"",fhNBinsChi2ITS-1,fhBinLimChi2ITS);
602 fhQA[kCutCovElement11][i] = new TH1F(Form("%s_covMatrixDiagonal11%s",GetName(),str),"",fhNBinsCovariance11-1,fhBinLimCovariance11);
603 fhQA[kCutCovElement22][i] = new TH1F(Form("%s_covMatrixDiagonal22%s",GetName(),str),"",fhNBinsCovariance22-1,fhBinLimCovariance22);
604 fhQA[kCutCovElement33][i] = new TH1F(Form("%s_covMatrixDiagonal33%s",GetName(),str),"",fhNBinsCovariance33-1,fhBinLimCovariance33);
605 fhQA[kCutCovElement44][i] = new TH1F(Form("%s_covMatrixDiagonal44%s",GetName(),str),"",fhNBinsCovariance44-1,fhBinLimCovariance44);
606 fhQA[kCutCovElement55][i] = new TH1F(Form("%s_covMatrixDiagonal55%s",GetName(),str),"",fhNBinsCovariance55-1,fhBinLimCovariance55);
608 fhQA[kCutClusterTPC][i] ->SetXTitle("n TPC clusters");
609 fhQA[kCutClusterITS][i] ->SetXTitle("n ITS clusters");
610 fhQA[kCutChi2TPC][i] ->SetXTitle("#chi^{2} per TPC cluster");
611 fhQA[kCutChi2ITS][i] ->SetXTitle("#chi^{2} per ITS cluster");
612 fhQA[kCutCovElement11][i] ->SetXTitle("cov 11 : #sigma_{y}^{2} (cm^{2})");
613 fhQA[kCutCovElement22][i] ->SetXTitle("cov 22 : #sigma_{z}^{2} (cm^{2})");
614 fhQA[kCutCovElement33][i] ->SetXTitle("cov 33 : #sigma_{sin(#phi)}^{2}");
615 fhQA[kCutCovElement44][i] ->SetXTitle("cov 44 : #sigma_{tan(#theta_{dip})}^{2}");
616 fhQA[kCutCovElement55][i] ->SetXTitle("cov 55 : #sigma_{1/p_{T}}^{2} ((c/GeV)^{2})");
619 for(Int_t i=0; i<kNHist; i++) fhQA[i][1]->SetLineColor(color);
621 //__________________________________________________________________________________
622 void AliCFTrackQualityCuts::FillHistograms(TObject* obj, Bool_t b)
625 // fill the QA histograms
629 if (!obj->InheritsFrom("AliVParticle")) {
630 AliError("object must derived from AliVParticle !");
634 Bool_t isESDTrack = strcmp(obj->ClassName(),"AliESDtrack") == 0 ? kTRUE : kFALSE ;
635 Bool_t isAODTrack = strcmp(obj->ClassName(),"AliAODTrack") == 0 ? kTRUE : kFALSE ;
637 AliESDtrack * esdTrack = 0x0 ;
638 AliAODTrack * aodTrack = 0x0 ;
639 if (isESDTrack) esdTrack = dynamic_cast<AliESDtrack*>(obj);
640 if (isAODTrack) aodTrack = dynamic_cast<AliAODTrack*>(obj);
642 // b = 0: fill histograms before cuts
643 // b = 1: fill histograms after cuts
646 Int_t nClustersTPC = 0;
647 Int_t nClustersITS = 0 ;
648 Float_t chi2PerClusterTPC = 0 ;
649 Float_t chi2PerClusterITS = 0 ;
650 Double_t extCov[15]={0,0,0,0,0,0,0,0,0,0,0,0,0,0,0};
653 nClustersTPC = esdTrack->GetTPCclusters(fIdxInt);
654 nClustersITS = esdTrack->GetITSclusters(fIdxInt);
655 if (nClustersTPC != 0) chi2PerClusterTPC = esdTrack->GetTPCchi2() / Float_t(nClustersTPC);
656 if (nClustersITS!=0) chi2PerClusterITS = esdTrack->GetITSchi2() / Float_t(nClustersITS);
657 esdTrack->GetExternalCovariance(extCov);
659 fhQA[kCutClusterTPC][b]->Fill((float)nClustersTPC);
660 fhQA[kCutChi2TPC][b]->Fill(chi2PerClusterTPC);
661 fhQA[kCutClusterITS][b]->Fill((float)nClustersITS);
662 fhQA[kCutChi2ITS][b]->Fill(chi2PerClusterITS);
663 fhQA[kCutCovElement11][b]->Fill(extCov[0]);
664 fhQA[kCutCovElement22][b]->Fill(extCov[2]);
665 fhQA[kCutCovElement33][b]->Fill(extCov[5]);
666 fhQA[kCutCovElement44][b]->Fill(extCov[9]);
667 fhQA[kCutCovElement55][b]->Fill(extCov[14]);
669 // fill cut statistics and cut correlation histograms with information from the bitmap
672 // Get the bitmap of the single cuts
674 SelectionBitMap(obj);
676 // Number of single cuts in this class
677 UInt_t ncuts = fBitmap->GetNbits();
678 for(UInt_t bit=0; bit<ncuts;bit++) {
679 if (!fBitmap->TestBitNumber(bit)) {
680 fhCutStatistics->Fill(bit+1);
681 for (UInt_t bit2=bit; bit2<ncuts;bit2++) {
682 if (!fBitmap->TestBitNumber(bit2))
683 fhCutCorrelation->Fill(bit+1,bit2+1);
688 //__________________________________________________________________________________
689 void AliCFTrackQualityCuts::SaveHistograms(const Char_t* dir) {
691 // saves the histograms in a directory (dir)
698 gDirectory->mkdir(dir);
701 gDirectory->mkdir("before_cuts");
702 gDirectory->mkdir("after_cuts");
704 fhCutStatistics->Write();
705 fhCutCorrelation->Write();
707 for (Int_t j=0; j<kNStepQA; j++) {
709 gDirectory->cd("before_cuts");
711 gDirectory->cd("after_cuts");
713 for(Int_t i=0; i<kNHist; i++) fhQA[i][j]->Write();
715 gDirectory->cd("../");
717 gDirectory->cd("../");
719 //__________________________________________________________________________________
720 void AliCFTrackQualityCuts::DrawHistograms(Bool_t drawLogScale)
723 // draws some histograms
728 Float_t right = 0.03;
729 Float_t left = 0.175;
731 Float_t bottom = 0.175;
733 TCanvas* canvas1 = new TCanvas("Track_QA_Quality_1", "Track QA Quality 1", 800, 500);
734 canvas1->Divide(2, 1);
737 fhCutStatistics->SetStats(kFALSE);
738 fhCutStatistics->LabelsOption("v");
739 gPad->SetLeftMargin(left);
740 gPad->SetBottomMargin(0.25);
741 gPad->SetRightMargin(right);
742 gPad->SetTopMargin(0.1);
743 fhCutStatistics->Draw();
746 fhCutCorrelation->SetStats(kFALSE);
747 fhCutCorrelation->LabelsOption("v");
748 gPad->SetLeftMargin(0.30);
749 gPad->SetRightMargin(bottom);
750 gPad->SetTopMargin(0.1);
751 gPad->SetBottomMargin(0.25);
752 fhCutCorrelation->Draw("COLZ");
754 canvas1->SaveAs(Form("%s.eps", canvas1->GetName()));
755 canvas1->SaveAs(Form("%s.ps", canvas1->GetName()));
759 TCanvas* canvas2 = new TCanvas("Track_QA_Quality_2", "Track QA Quality 2", 1200, 800);
760 canvas2->Divide(2, 2);
763 gPad->SetRightMargin(right);
764 gPad->SetLeftMargin(left);
765 gPad->SetTopMargin(top);
766 gPad->SetBottomMargin(bottom);
767 fhQA[kCutClusterTPC][0]->SetStats(kFALSE);
768 fhQA[kCutClusterTPC][0]->Draw();
769 fhQA[kCutClusterTPC][1]->Draw("same");
772 gPad->SetRightMargin(right);
773 gPad->SetLeftMargin(left);
774 gPad->SetTopMargin(top);
775 gPad->SetBottomMargin(bottom);
776 fhQA[kCutChi2TPC][0]->SetStats(kFALSE);
777 fhQA[kCutChi2TPC][0]->Draw();
778 fhQA[kCutChi2TPC][1]->Draw("same");
781 gPad->SetRightMargin(right);
782 gPad->SetLeftMargin(left);
783 gPad->SetTopMargin(top);
784 gPad->SetBottomMargin(bottom);
785 fhQA[kCutClusterITS][0]->SetStats(kFALSE);
786 fhQA[kCutClusterITS][0]->Draw();
787 fhQA[kCutClusterITS][1]->Draw("same");
790 gPad->SetRightMargin(right);
791 gPad->SetLeftMargin(left);
792 gPad->SetTopMargin(top);
793 gPad->SetBottomMargin(bottom);
794 fhQA[kCutChi2ITS][0]->SetStats(kFALSE);
795 fhQA[kCutChi2ITS][0]->Draw();
796 fhQA[kCutChi2ITS][1]->Draw("same");
798 canvas2->SaveAs(Form("%s.eps", canvas2->GetName()));
799 canvas2->SaveAs(Form("%s.ps", canvas2->GetName()));
803 TCanvas* canvas3 = new TCanvas("Track_QA_Quality_3", "Track QA Quality 3", 1200, 800);
804 canvas3->Divide(3, 2);
807 if(drawLogScale) gPad->SetLogy();
808 gPad->SetRightMargin(right);
809 gPad->SetLeftMargin(left);
810 gPad->SetTopMargin(top);
811 gPad->SetBottomMargin(bottom);
812 fhQA[kCutCovElement11][0]->SetStats(kFALSE);
813 fhQA[kCutCovElement11][0]->Draw();
814 fhQA[kCutCovElement11][1]->Draw("same");
817 if(drawLogScale) gPad->SetLogy();
818 gPad->SetRightMargin(right);
819 gPad->SetLeftMargin(left);
820 gPad->SetTopMargin(top);
821 gPad->SetBottomMargin(bottom);
822 fhQA[kCutCovElement22][0]->SetStats(kFALSE);
823 fhQA[kCutCovElement22][0]->Draw();
824 fhQA[kCutCovElement22][1]->Draw("same");
827 if(drawLogScale) gPad->SetLogy();
828 gPad->SetRightMargin(right);
829 gPad->SetLeftMargin(left);
830 gPad->SetTopMargin(top);
831 gPad->SetBottomMargin(bottom);
832 fhQA[kCutCovElement33][0]->SetStats(kFALSE);
833 fhQA[kCutCovElement33][0]->Draw();
834 fhQA[kCutCovElement33][1]->Draw("same");
837 if(drawLogScale) gPad->SetLogy();
838 gPad->SetRightMargin(right);
839 gPad->SetLeftMargin(left);
840 gPad->SetTopMargin(top);
841 gPad->SetBottomMargin(bottom);
842 fhQA[kCutCovElement44][0]->SetStats(kFALSE);
843 fhQA[kCutCovElement44][0]->Draw();
844 fhQA[kCutCovElement44][1]->Draw("same");
847 if(drawLogScale) gPad->SetLogy();
848 gPad->SetRightMargin(right);
849 gPad->SetLeftMargin(left);
850 gPad->SetTopMargin(top);
851 gPad->SetBottomMargin(bottom);
852 fhQA[kCutCovElement55][0]->SetStats(kFALSE);
853 fhQA[kCutCovElement55][0]->Draw();
854 fhQA[kCutCovElement55][1]->Draw("same");
856 canvas3->SaveAs(Form("%s.eps", canvas3->GetName()));
857 canvas3->SaveAs(Form("%s.ps", canvas3->GetName()));
859 //__________________________________________________________________________________
860 void AliCFTrackQualityCuts::AddQAHistograms(TList *qaList) {
862 // saves the histograms in a TList
866 qaList->Add(fhCutStatistics);
867 qaList->Add(fhCutCorrelation);
869 for (Int_t j=0; j<kNStepQA; j++) {
870 for(Int_t i=0; i<kNHist; i++)
871 qaList->Add(fhQA[i][j]);