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 // - number of clusters in the TRD
27 // - chi2 / cluster in the TPC
28 // - chi2 / cluster in the ITS
29 // - chi2 / cluster in the TRD
30 // - successful TPC refit
31 // - successful ITS refit
32 // - covariance matrix diagonal elements
34 // The cut values for these cuts are set with the corresponding set functions.
35 // All cut classes provided by the correction framework are supposed to be
36 // added in the Analysis Framwork's class AliAnalysisFilter and applied by
37 // the filter via a loop.
39 // author: I. Kraus (Ingrid.Kraus@cern.ch)
41 // AliESDtrackCuts writte by Jan Fiete Grosse-Oetringhaus and
42 // AliRsnDaughterCut class written by A. Pulvirenti.
46 #include <TDirectory.h>
50 #include <AliESDtrack.h>
52 #include "AliCFTrackQualityCuts.h"
53 #include "AliAODTrack.h"
55 ClassImp(AliCFTrackQualityCuts)
57 //__________________________________________________________________________________
58 AliCFTrackQualityCuts::AliCFTrackQualityCuts() :
63 fMaxChi2PerClusterTPC(0),
64 fMaxChi2PerClusterITS(0),
65 fMaxChi2PerClusterTRD(0),
81 fhNBinsCovariance11(0),
82 fhNBinsCovariance22(0),
83 fhNBinsCovariance33(0),
84 fhNBinsCovariance44(0),
85 fhNBinsCovariance55(0),
86 fhBinLimClusterTPC(0x0),
87 fhBinLimClusterITS(0x0),
88 fhBinLimClusterTRD(0x0),
92 fhBinLimCovariance11(0x0),
93 fhBinLimCovariance22(0x0),
94 fhBinLimCovariance33(0x0),
95 fhBinLimCovariance44(0x0),
96 fhBinLimCovariance55(0x0)
99 // Default constructor
103 //__________________________________________________________________________________
104 AliCFTrackQualityCuts::AliCFTrackQualityCuts(Char_t* name, Char_t* title) :
105 AliCFCutBase(name,title),
109 fMaxChi2PerClusterTPC(0),
110 fMaxChi2PerClusterITS(0),
111 fMaxChi2PerClusterTRD(0),
121 fhNBinsClusterTPC(0),
122 fhNBinsClusterITS(0),
123 fhNBinsClusterTRD(0),
127 fhNBinsCovariance11(0),
128 fhNBinsCovariance22(0),
129 fhNBinsCovariance33(0),
130 fhNBinsCovariance44(0),
131 fhNBinsCovariance55(0),
132 fhBinLimClusterTPC(0x0),
133 fhBinLimClusterITS(0x0),
134 fhBinLimClusterTRD(0x0),
135 fhBinLimChi2TPC(0x0),
136 fhBinLimChi2ITS(0x0),
137 fhBinLimChi2TRD(0x0),
138 fhBinLimCovariance11(0x0),
139 fhBinLimCovariance22(0x0),
140 fhBinLimCovariance33(0x0),
141 fhBinLimCovariance44(0x0),
142 fhBinLimCovariance55(0x0)
149 //__________________________________________________________________________________
150 AliCFTrackQualityCuts::AliCFTrackQualityCuts(const AliCFTrackQualityCuts& c) :
152 fMinNClusterTPC(c.fMinNClusterTPC),
153 fMinNClusterITS(c.fMinNClusterITS),
154 fMinNClusterTRD(c.fMinNClusterTRD),
155 fMaxChi2PerClusterTPC(c.fMaxChi2PerClusterTPC),
156 fMaxChi2PerClusterITS(c.fMaxChi2PerClusterITS),
157 fMaxChi2PerClusterTRD(c.fMaxChi2PerClusterTRD),
158 fCovariance11Max(c.fCovariance11Max),
159 fCovariance22Max(c.fCovariance22Max),
160 fCovariance33Max(c.fCovariance33Max),
161 fCovariance44Max(c.fCovariance44Max),
162 fCovariance55Max(c.fCovariance55Max),
164 fhCutStatistics(c.fhCutStatistics),
165 fhCutCorrelation(c.fhCutCorrelation),
167 fhNBinsClusterTPC(c.fhNBinsClusterTPC),
168 fhNBinsClusterITS(c.fhNBinsClusterITS),
169 fhNBinsClusterTRD(c.fhNBinsClusterTRD),
170 fhNBinsChi2TPC(c.fhNBinsChi2TPC),
171 fhNBinsChi2ITS(c.fhNBinsChi2ITS),
172 fhNBinsChi2TRD(c.fhNBinsChi2TRD),
173 fhNBinsCovariance11(c.fhNBinsCovariance11),
174 fhNBinsCovariance22(c.fhNBinsCovariance22),
175 fhNBinsCovariance33(c.fhNBinsCovariance33),
176 fhNBinsCovariance44(c.fhNBinsCovariance44),
177 fhNBinsCovariance55(c.fhNBinsCovariance55),
178 fhBinLimClusterTPC(c.fhBinLimClusterTPC),
179 fhBinLimClusterITS(c.fhBinLimClusterITS),
180 fhBinLimClusterTRD(c.fhBinLimClusterTRD),
181 fhBinLimChi2TPC(c.fhBinLimChi2TPC),
182 fhBinLimChi2ITS(c.fhBinLimChi2ITS),
183 fhBinLimChi2TRD(c.fhBinLimChi2TRD),
184 fhBinLimCovariance11(c.fhBinLimCovariance11),
185 fhBinLimCovariance22(c.fhBinLimCovariance22),
186 fhBinLimCovariance33(c.fhBinLimCovariance33),
187 fhBinLimCovariance44(c.fhBinLimCovariance44),
188 fhBinLimCovariance55(c.fhBinLimCovariance55)
193 ((AliCFTrackQualityCuts &) c).Copy(*this);
195 //__________________________________________________________________________________
196 AliCFTrackQualityCuts& AliCFTrackQualityCuts::operator=(const AliCFTrackQualityCuts& c)
199 // Assignment operator
202 AliCFCutBase::operator=(c) ;
203 fMinNClusterTPC = c.fMinNClusterTPC ;
204 fMinNClusterITS = c.fMinNClusterITS ;
205 fMinNClusterTRD = c.fMinNClusterTRD ;
206 fMaxChi2PerClusterTPC = c.fMaxChi2PerClusterTPC ;
207 fMaxChi2PerClusterITS = c.fMaxChi2PerClusterITS ;
208 fMaxChi2PerClusterTRD = c.fMaxChi2PerClusterTRD ;
209 fCovariance11Max = c.fCovariance11Max ;
210 fCovariance22Max = c.fCovariance22Max ;
211 fCovariance33Max = c.fCovariance33Max ;
212 fCovariance44Max = c.fCovariance44Max ;
213 fCovariance55Max = c.fCovariance55Max ;
214 fStatus = c.fStatus ;
215 fhCutStatistics = c.fhCutStatistics ;
216 fhCutCorrelation = c.fhCutCorrelation ;
217 fBitmap = c.fBitmap ;
218 fhNBinsClusterTPC = c.fhNBinsClusterTPC ;
219 fhNBinsClusterITS = c.fhNBinsClusterITS ;
220 fhNBinsClusterTRD = c.fhNBinsClusterTRD ;
221 fhNBinsChi2TPC = c.fhNBinsChi2TPC ;
222 fhNBinsChi2ITS = c.fhNBinsChi2ITS ;
223 fhNBinsChi2TRD = c.fhNBinsChi2TRD ;
224 fhNBinsCovariance11 = c.fhNBinsCovariance11 ;
225 fhNBinsCovariance22 = c.fhNBinsCovariance22 ;
226 fhNBinsCovariance33 = c.fhNBinsCovariance33 ;
227 fhNBinsCovariance44 = c.fhNBinsCovariance44 ;
228 fhNBinsCovariance55 = c.fhNBinsCovariance55 ;
229 fhBinLimClusterTPC = c.fhBinLimClusterTPC ;
230 fhBinLimClusterITS = c.fhBinLimClusterITS ;
231 fhBinLimClusterTRD = c.fhBinLimClusterTRD ;
232 fhBinLimChi2TPC = c.fhBinLimChi2TPC ;
233 fhBinLimChi2ITS = c.fhBinLimChi2ITS ;
234 fhBinLimChi2TRD = c.fhBinLimChi2TRD ;
235 fhBinLimCovariance11 = c.fhBinLimCovariance11 ;
236 fhBinLimCovariance22 = c.fhBinLimCovariance22 ;
237 fhBinLimCovariance33 = c.fhBinLimCovariance33 ;
238 fhBinLimCovariance44 = c.fhBinLimCovariance44 ;
239 fhBinLimCovariance55 = c.fhBinLimCovariance55 ;
241 for (Int_t i=0; i<c.kNHist; i++){
242 for (Int_t j=0; j<c.kNStepQA; j++){
243 if(c.fhQA[i][j]) fhQA[i][j] = (TH1F*)c.fhQA[i][j]->Clone();
246 ((AliCFTrackQualityCuts &) c).Copy(*this);
250 //__________________________________________________________________________________
251 AliCFTrackQualityCuts::~AliCFTrackQualityCuts()
256 if (fhCutStatistics) delete fhCutStatistics;
257 if (fhCutCorrelation) delete fhCutCorrelation;
259 for (Int_t i=0; i<kNHist; i++){
260 for (Int_t j=0; j<kNStepQA; j++){
261 if(fhQA[i][j]) delete fhQA[i][j];
264 if(fBitmap) delete fBitmap;
265 if(fhBinLimClusterTPC) delete fhBinLimClusterTPC;
266 if(fhBinLimClusterITS) delete fhBinLimClusterITS;
267 if(fhBinLimClusterTRD) delete fhBinLimClusterTRD;
268 if(fhBinLimChi2TPC) delete fhBinLimChi2TPC;
269 if(fhBinLimChi2ITS) delete fhBinLimChi2ITS;
270 if(fhBinLimChi2TRD) delete fhBinLimChi2TRD;
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;
277 //__________________________________________________________________________________
278 void AliCFTrackQualityCuts::Initialise()
281 // sets everything to zero
286 fMaxChi2PerClusterTPC = 0;
287 fMaxChi2PerClusterITS = 0;
288 fMaxChi2PerClusterTRD = 0;
289 fCovariance11Max = 0;
290 fCovariance22Max = 0;
291 fCovariance33Max = 0;
292 fCovariance44Max = 0;
293 fCovariance55Max = 0;
299 SetMaxChi2PerClusterTPC();
300 SetMaxChi2PerClusterITS();
301 SetMaxChi2PerClusterTRD();
302 SetMaxCovDiagonalElements();
305 for (Int_t i=0; i<kNHist; i++){
306 for (Int_t j=0; j<kNStepQA; j++)
310 fhCutCorrelation = 0;
311 fBitmap=new TBits(0);
313 //set default bining for QA histograms
314 SetHistogramBins(kCutClusterTPC,165,-0.5,164.5);
315 SetHistogramBins(kCutClusterITS,8,-0.5,7.5);
316 SetHistogramBins(kCutClusterTRD,100,-0.5,99.5);
317 SetHistogramBins(kCutChi2TPC,500,0.,10.);
318 SetHistogramBins(kCutChi2ITS,500,0.,10.);
319 SetHistogramBins(kCutChi2TRD,500,0.,10.);
320 SetHistogramBins(kCutCovElement11,200,0.,20.);
321 SetHistogramBins(kCutCovElement22,200,0.,20.);
322 SetHistogramBins(kCutCovElement33,100,0.,1.);
323 SetHistogramBins(kCutCovElement44,100,0.,5.);
324 SetHistogramBins(kCutCovElement55,100,0.,5.);
326 //__________________________________________________________________________________
327 void AliCFTrackQualityCuts::Copy(TObject &c) const
332 AliCFTrackQualityCuts& target = (AliCFTrackQualityCuts &) c;
336 if (fhCutStatistics) target.fhCutStatistics = (TH1F*) fhCutStatistics->Clone();
337 if (fhCutCorrelation) target.fhCutCorrelation = (TH2F*) fhCutCorrelation->Clone();
339 for (Int_t i=0; i<kNHist; i++){
340 for (Int_t j=0; j<kNStepQA; j++){
341 if(fhQA[i][j]) target.fhQA[i][j] = (TH1F*)fhQA[i][j]->Clone();
346 //__________________________________________________________________________________
347 void AliCFTrackQualityCuts::SelectionBitMap(TObject* obj)
350 // test if the track passes the single cuts
351 // and store the information in a bitmap
354 // bitmap stores the decision of each single cut
355 for(Int_t i=0; i<kNCuts; i++)fBitmap->SetBitNumber(i,kFALSE);
358 if (!obj->InheritsFrom("AliVParticle")) {
359 AliError("object must derived from AliVParticle !");
363 Bool_t isESDTrack = strcmp(obj->ClassName(),"AliESDtrack") == 0 ? kTRUE : kFALSE ;
364 Bool_t isAODTrack = strcmp(obj->ClassName(),"AliAODTrack") == 0 ? kTRUE : kFALSE ;
366 AliESDtrack * esdTrack = 0x0 ;
367 AliAODTrack * aodTrack = 0x0 ;
368 if (isESDTrack) esdTrack = dynamic_cast<AliESDtrack*>(obj);
369 if (isAODTrack) aodTrack = dynamic_cast<AliAODTrack*>(obj);
371 // get cut quantities
373 Int_t nClustersTPC = 0;
374 Int_t nClustersITS = 0 ;
375 Int_t nClustersTRD = 0 ;
376 Float_t chi2PerClusterTPC = 0 ;
377 Float_t chi2PerClusterITS = 0 ;
378 Float_t chi2PerClusterTRD = 0 ;
379 Double_t extCov[15]={0,0,0,0,0,0,0,0,0,0,0,0,0,0,0};
382 nClustersTPC = esdTrack->GetTPCclusters(fIdxInt);
383 nClustersITS = esdTrack->GetITSclusters(fIdxInt);
384 nClustersTRD = esdTrack->GetTRDclusters(fIdxInt);
385 if (nClustersTPC != 0) chi2PerClusterTPC = esdTrack->GetTPCchi2() / Float_t(nClustersTPC);
386 if (nClustersITS != 0) chi2PerClusterITS = esdTrack->GetITSchi2() / Float_t(nClustersITS);
387 if (nClustersTRD != 0) chi2PerClusterTRD = esdTrack->GetTRDchi2() / Float_t(nClustersTRD);
388 esdTrack->GetExternalCovariance(extCov);
394 if (nClustersTPC >= fMinNClusterTPC)
395 fBitmap->SetBitNumber(iCutBit,kTRUE);
397 if (nClustersITS >= fMinNClusterITS)
398 fBitmap->SetBitNumber(iCutBit,kTRUE);
400 if (nClustersTRD >= fMinNClusterTRD)
401 fBitmap->SetBitNumber(iCutBit,kTRUE);
403 if (chi2PerClusterTPC <= fMaxChi2PerClusterTPC)
404 fBitmap->SetBitNumber(iCutBit,kTRUE);
406 if (chi2PerClusterITS <= fMaxChi2PerClusterITS)
407 fBitmap->SetBitNumber(iCutBit,kTRUE);
409 if (chi2PerClusterTRD <= fMaxChi2PerClusterTRD)
410 fBitmap->SetBitNumber(iCutBit,kTRUE);
412 if (extCov[0] <= fCovariance11Max)
413 fBitmap->SetBitNumber(iCutBit,kTRUE);
415 if (extCov[2] <= fCovariance22Max)
416 fBitmap->SetBitNumber(iCutBit,kTRUE);
418 if (extCov[5] <= fCovariance33Max)
419 fBitmap->SetBitNumber(iCutBit,kTRUE);
421 if (extCov[9] <= fCovariance44Max)
422 fBitmap->SetBitNumber(iCutBit,kTRUE);
424 if (extCov[14] <= fCovariance55Max)
425 fBitmap->SetBitNumber(iCutBit,kTRUE);
429 if ( (esdTrack->GetStatus() & fStatus) == fStatus ) fBitmap->SetBitNumber(iCutBit,kTRUE);
432 if ( (aodTrack->GetStatus() & fStatus) == fStatus ) fBitmap->SetBitNumber(iCutBit,kTRUE);
437 //__________________________________________________________________________________
438 Bool_t AliCFTrackQualityCuts::IsSelected(TObject* obj) {
440 // loops over decisions of single cuts and returns if the track is accepted
442 SelectionBitMap(obj);
444 if (fIsQAOn) FillHistograms(obj,0);
445 Bool_t isSelected = kTRUE;
447 for (UInt_t icut=0; icut<fBitmap->GetNbits();icut++) {
448 if(!fBitmap->TestBitNumber(icut)) {
453 if (!isSelected) return kFALSE ;
454 if (fIsQAOn) FillHistograms(obj,1);
457 //__________________________________________________________________________________
458 void AliCFTrackQualityCuts::SetHistogramBins(Int_t index, Int_t nbins, Double_t *bins)
465 if (index<0 || index>=kNHist) {
466 Error("SetHistogramBins","could not determine histogram from index %d",index);
472 fhNBinsClusterTPC=nbins+1;
473 fhBinLimClusterTPC=new Double_t[nbins+1];
474 for(Int_t i=0;i<nbins+1;i++)fhBinLimClusterTPC[i]=bins[i];
478 fhNBinsClusterITS=nbins+1;
479 fhBinLimClusterITS=new Double_t[nbins+1];
480 for(Int_t i=0;i<nbins+1;i++)fhBinLimClusterITS[i]=bins[i];
484 fhNBinsClusterTRD=nbins+1;
485 fhBinLimClusterTRD=new Double_t[nbins+1];
486 for(Int_t i=0;i<nbins+1;i++)fhBinLimClusterTRD[i]=bins[i];
490 fhNBinsChi2TPC=nbins+1;
491 fhBinLimChi2TPC=new Double_t[nbins+1];
492 for(Int_t i=0;i<nbins+1;i++)fhBinLimChi2TPC[i]=bins[i];
496 fhNBinsChi2ITS=nbins+1;
497 fhBinLimChi2ITS=new Double_t[nbins+1];
498 for(Int_t i=0;i<nbins+1;i++)fhBinLimChi2ITS[i]=bins[i];
502 fhNBinsChi2TRD=nbins+1;
503 fhBinLimChi2TRD=new Double_t[nbins+1];
504 for(Int_t i=0;i<nbins+1;i++)fhBinLimChi2TRD[i]=bins[i];
507 case kCutCovElement11:
508 fhNBinsCovariance11=nbins+1;
509 fhBinLimCovariance11=new Double_t[nbins+1];
510 for(Int_t i=0;i<nbins+1;i++)fhBinLimCovariance11[i]=bins[i];
513 case kCutCovElement22:
514 fhNBinsCovariance22=nbins+1;
515 fhBinLimCovariance22=new Double_t[nbins+1];
516 for(Int_t i=0;i<nbins+1;i++)fhBinLimCovariance22[i]=bins[i];
519 case kCutCovElement33:
520 fhNBinsCovariance33=nbins+1;
521 fhBinLimCovariance33=new Double_t[nbins+1];
522 for(Int_t i=0;i<nbins+1;i++)fhBinLimCovariance33[i]=bins[i];
525 case kCutCovElement44:
526 fhNBinsCovariance44=nbins+1;
527 fhBinLimCovariance44=new Double_t[nbins+1];
528 for(Int_t i=0;i<nbins+1;i++)fhBinLimCovariance44[i]=bins[i];
531 case kCutCovElement55:
532 fhNBinsCovariance55=nbins+1;
533 fhBinLimCovariance55=new Double_t[nbins+1];
534 for(Int_t i=0;i<nbins+1;i++)fhBinLimCovariance55[i]=bins[i];
538 //__________________________________________________________________________________
539 void AliCFTrackQualityCuts::SetHistogramBins(Int_t index, Int_t nbins, Double_t xmin, Double_t xmax)
544 if (index<0 || index>=kNHist) {
545 Error("SetHistogramBins","could not determine histogram from index %d",index);
551 fhNBinsClusterTPC=nbins+1;
552 fhBinLimClusterTPC=new Double_t[nbins+1];
553 for(Int_t i=0;i<nbins+1;i++)fhBinLimClusterTPC[i]=xmin+i*(xmax-xmin)/Double_t(nbins);
557 fhNBinsClusterITS=nbins+1;
558 fhBinLimClusterITS=new Double_t[nbins+1];
559 for(Int_t i=0;i<nbins+1;i++)fhBinLimClusterITS[i]=xmin+i*(xmax-xmin)/Double_t(nbins);
563 fhNBinsClusterTRD=nbins+1;
564 fhBinLimClusterTRD=new Double_t[nbins+1];
565 for(Int_t i=0;i<nbins+1;i++)fhBinLimClusterTRD[i]=xmin+i*(xmax-xmin)/Double_t(nbins);
569 fhNBinsChi2TPC=nbins+1;
570 fhBinLimChi2TPC=new Double_t[nbins+1];
571 for(Int_t i=0;i<nbins+1;i++)fhBinLimChi2TPC[i]=xmin+i*(xmax-xmin)/Double_t(nbins);
575 fhNBinsChi2ITS=nbins+1;
576 fhBinLimChi2ITS=new Double_t[nbins+1];
577 for(Int_t i=0;i<nbins+1;i++)fhBinLimChi2ITS[i]=xmin+i*(xmax-xmin)/Double_t(nbins);
581 fhNBinsChi2TRD=nbins+1;
582 fhBinLimChi2TRD=new Double_t[nbins+1];
583 for(Int_t i=0;i<nbins+1;i++)fhBinLimChi2TRD[i]=xmin+i*(xmax-xmin)/Double_t(nbins);
586 case kCutCovElement11:
587 fhNBinsCovariance11=nbins+1;
588 fhBinLimCovariance11=new Double_t[nbins+1];
589 for(Int_t i=0;i<nbins+1;i++)fhBinLimCovariance11[i]=xmin+i*(xmax-xmin)/Double_t(nbins);
592 case kCutCovElement22:
593 fhNBinsCovariance22=nbins+1;
594 fhBinLimCovariance22=new Double_t[nbins+1];
595 for(Int_t i=0;i<nbins+1;i++)fhBinLimCovariance22[i]=xmin+i*(xmax-xmin)/Double_t(nbins);
598 case kCutCovElement33:
599 fhNBinsCovariance33=nbins+1;
600 fhBinLimCovariance33=new Double_t[nbins+1];
601 for(Int_t i=0;i<nbins+1;i++)fhBinLimCovariance33[i]=xmin+i*(xmax-xmin)/Double_t(nbins);
604 case kCutCovElement44:
605 fhNBinsCovariance44=nbins+1;
606 fhBinLimCovariance44=new Double_t[nbins+1];
607 for(Int_t i=0;i<nbins+1;i++)fhBinLimCovariance44[i]=xmin+i*(xmax-xmin)/Double_t(nbins);
610 case kCutCovElement55:
611 fhNBinsCovariance55=nbins+1;
612 fhBinLimCovariance55=new Double_t[nbins+1];
613 for(Int_t i=0;i<nbins+1;i++)fhBinLimCovariance55[i]=xmin+i*(xmax-xmin)/Double_t(nbins);
617 //__________________________________________________________________________________
618 void AliCFTrackQualityCuts::DefineHistograms() {
620 // histograms for cut variables, cut statistics and cut correlations
624 // book cut statistics and cut correlation histograms
625 fhCutStatistics = new TH1F(Form("%s_cut_statistics",GetName()), Form("%s cut statistics",GetName()), kNCuts,0.5,kNCuts+0.5);
626 fhCutStatistics->SetLineWidth(2);
627 fhCutStatistics->GetXaxis()->SetBinLabel(1, "nClustersTPC");
628 fhCutStatistics->GetXaxis()->SetBinLabel(2, "nClustersITS");
629 fhCutStatistics->GetXaxis()->SetBinLabel(3, "nClustersTRD");
630 fhCutStatistics->GetXaxis()->SetBinLabel(4, "chi2PerClusterTPC");
631 fhCutStatistics->GetXaxis()->SetBinLabel(5, "chi2PerClusterITS");
632 fhCutStatistics->GetXaxis()->SetBinLabel(6, "chi2PerClusterTRD");
633 fhCutStatistics->GetXaxis()->SetBinLabel(7, "covElement11");
634 fhCutStatistics->GetXaxis()->SetBinLabel(8, "covElement22");
635 fhCutStatistics->GetXaxis()->SetBinLabel(9, "covElement33");
636 fhCutStatistics->GetXaxis()->SetBinLabel(10,"covElement44");
637 fhCutStatistics->GetXaxis()->SetBinLabel(11,"covElement55");
638 fhCutStatistics->GetXaxis()->SetBinLabel(12,"status");
640 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);
641 fhCutCorrelation->SetLineWidth(2);
642 fhCutCorrelation->GetXaxis()->SetBinLabel(1, "nClustersTPC");
643 fhCutCorrelation->GetXaxis()->SetBinLabel(2, "nClustersITS");
644 fhCutCorrelation->GetXaxis()->SetBinLabel(3, "nClustersTRD");
645 fhCutCorrelation->GetXaxis()->SetBinLabel(4, "chi2PerClusterTPC");
646 fhCutCorrelation->GetXaxis()->SetBinLabel(5, "chi2PerClusterITS");
647 fhCutCorrelation->GetXaxis()->SetBinLabel(6, "chi2PerClusterTRD");
648 fhCutCorrelation->GetXaxis()->SetBinLabel(7, "covElement11");
649 fhCutCorrelation->GetXaxis()->SetBinLabel(8, "covElement22");
650 fhCutCorrelation->GetXaxis()->SetBinLabel(9, "covElement33");
651 fhCutCorrelation->GetXaxis()->SetBinLabel(10,"covElement44");
652 fhCutCorrelation->GetXaxis()->SetBinLabel(11,"covElement55");
653 fhCutStatistics->GetXaxis()->SetBinLabel(12,"status");
655 fhCutCorrelation->GetYaxis()->SetBinLabel(1, "nClustersTPC");
656 fhCutCorrelation->GetYaxis()->SetBinLabel(2, "nClustersITS");
657 fhCutCorrelation->GetYaxis()->SetBinLabel(3, "nClustersTRD");
658 fhCutCorrelation->GetYaxis()->SetBinLabel(4, "chi2PerClusterTPC");
659 fhCutCorrelation->GetYaxis()->SetBinLabel(5, "chi2PerClusterITS");
660 fhCutCorrelation->GetYaxis()->SetBinLabel(6, "chi2PerClusterTRD");
661 fhCutCorrelation->GetYaxis()->SetBinLabel(7, "covElement11");
662 fhCutCorrelation->GetYaxis()->SetBinLabel(8, "covElement22");
663 fhCutCorrelation->GetYaxis()->SetBinLabel(9, "covElement33");
664 fhCutCorrelation->GetYaxis()->SetBinLabel(10,"covElement44");
665 fhCutCorrelation->GetYaxis()->SetBinLabel(11,"covElement55");
666 fhCutCorrelation->GetXaxis()->SetBinLabel(12,"status");
669 // book QA histograms
671 for (Int_t i=0; i<kNStepQA; i++) {
672 if (i==0) sprintf(str," ");
673 else sprintf(str,"_cut");
675 fhQA[kCutClusterTPC][i] = new TH1F(Form("%s_nClustersTPC%s",GetName(),str) ,"",fhNBinsClusterTPC-1,fhBinLimClusterTPC);
676 fhQA[kCutClusterITS][i] = new TH1F(Form("%s_nClustersITS%s",GetName(),str) ,"",fhNBinsClusterITS-1,fhBinLimClusterITS);
677 fhQA[kCutClusterTRD][i] = new TH1F(Form("%s_nClustersTRD%s",GetName(),str) ,"",fhNBinsClusterTRD-1,fhBinLimClusterTRD);
678 fhQA[kCutChi2TPC][i] = new TH1F(Form("%s_chi2PerClusterTPC%s",GetName(),str),"",fhNBinsChi2TPC-1,fhBinLimChi2TPC);
679 fhQA[kCutChi2ITS][i] = new TH1F(Form("%s_chi2PerClusterITS%s",GetName(),str),"",fhNBinsChi2ITS-1,fhBinLimChi2ITS);
680 fhQA[kCutChi2TRD][i] = new TH1F(Form("%s_chi2PerClusterTRD%s",GetName(),str),"",fhNBinsChi2TRD-1,fhBinLimChi2TRD);
681 fhQA[kCutCovElement11][i] = new TH1F(Form("%s_covMatrixDiagonal11%s",GetName(),str),"",fhNBinsCovariance11-1,fhBinLimCovariance11);
682 fhQA[kCutCovElement22][i] = new TH1F(Form("%s_covMatrixDiagonal22%s",GetName(),str),"",fhNBinsCovariance22-1,fhBinLimCovariance22);
683 fhQA[kCutCovElement33][i] = new TH1F(Form("%s_covMatrixDiagonal33%s",GetName(),str),"",fhNBinsCovariance33-1,fhBinLimCovariance33);
684 fhQA[kCutCovElement44][i] = new TH1F(Form("%s_covMatrixDiagonal44%s",GetName(),str),"",fhNBinsCovariance44-1,fhBinLimCovariance44);
685 fhQA[kCutCovElement55][i] = new TH1F(Form("%s_covMatrixDiagonal55%s",GetName(),str),"",fhNBinsCovariance55-1,fhBinLimCovariance55);
687 fhQA[kCutClusterTPC][i] ->SetXTitle("n TPC clusters");
688 fhQA[kCutClusterITS][i] ->SetXTitle("n ITS clusters");
689 fhQA[kCutClusterTRD][i] ->SetXTitle("n TRD clusters");
690 fhQA[kCutChi2TPC][i] ->SetXTitle("#chi^{2} per TPC cluster");
691 fhQA[kCutChi2ITS][i] ->SetXTitle("#chi^{2} per ITS cluster");
692 fhQA[kCutChi2TRD][i] ->SetXTitle("#chi^{2} per TRD cluster");
693 fhQA[kCutCovElement11][i] ->SetXTitle("cov 11 : #sigma_{y}^{2} (cm^{2})");
694 fhQA[kCutCovElement22][i] ->SetXTitle("cov 22 : #sigma_{z}^{2} (cm^{2})");
695 fhQA[kCutCovElement33][i] ->SetXTitle("cov 33 : #sigma_{sin(#phi)}^{2}");
696 fhQA[kCutCovElement44][i] ->SetXTitle("cov 44 : #sigma_{tan(#theta_{dip})}^{2}");
697 fhQA[kCutCovElement55][i] ->SetXTitle("cov 55 : #sigma_{1/p_{T}}^{2} ((c/GeV)^{2})");
700 for(Int_t i=0; i<kNHist; i++) fhQA[i][1]->SetLineColor(color);
702 //__________________________________________________________________________________
703 void AliCFTrackQualityCuts::FillHistograms(TObject* obj, Bool_t b)
706 // fill the QA histograms
710 if (!obj->InheritsFrom("AliVParticle")) {
711 AliError("object must derived from AliVParticle !");
715 Bool_t isESDTrack = strcmp(obj->ClassName(),"AliESDtrack") == 0 ? kTRUE : kFALSE ;
716 Bool_t isAODTrack = strcmp(obj->ClassName(),"AliAODTrack") == 0 ? kTRUE : kFALSE ;
718 AliESDtrack * esdTrack = 0x0 ;
719 AliAODTrack * aodTrack = 0x0 ;
720 if (isESDTrack) esdTrack = dynamic_cast<AliESDtrack*>(obj);
721 if (isAODTrack) aodTrack = dynamic_cast<AliAODTrack*>(obj);
723 // b = 0: fill histograms before cuts
724 // b = 1: fill histograms after cuts
727 Int_t nClustersTPC = 0;
728 Int_t nClustersITS = 0 ;
729 Int_t nClustersTRD = 0 ;
730 Float_t chi2PerClusterTPC = 0 ;
731 Float_t chi2PerClusterITS = 0 ;
732 Float_t chi2PerClusterTRD = 0 ;
733 Double_t extCov[15]={0,0,0,0,0,0,0,0,0,0,0,0,0,0,0};
736 nClustersTPC = esdTrack->GetTPCclusters(fIdxInt);
737 nClustersITS = esdTrack->GetITSclusters(fIdxInt);
738 nClustersTRD = esdTrack->GetTRDclusters(fIdxInt);
739 if (nClustersTPC != 0) chi2PerClusterTPC = esdTrack->GetTPCchi2() / Float_t(nClustersTPC);
740 if (nClustersITS != 0) chi2PerClusterITS = esdTrack->GetITSchi2() / Float_t(nClustersITS);
741 if (nClustersTRD != 0) chi2PerClusterTRD = esdTrack->GetTRDchi2() / Float_t(nClustersTRD);
742 esdTrack->GetExternalCovariance(extCov);
744 fhQA[kCutClusterTPC][b]->Fill((float)nClustersTPC);
745 fhQA[kCutChi2TPC][b]->Fill(chi2PerClusterTPC);
746 fhQA[kCutClusterITS][b]->Fill((float)nClustersITS);
747 fhQA[kCutChi2ITS][b]->Fill(chi2PerClusterITS);
748 fhQA[kCutClusterTRD][b]->Fill((float)nClustersTRD);
749 fhQA[kCutChi2TRD][b]->Fill(chi2PerClusterTRD);
750 fhQA[kCutCovElement11][b]->Fill(extCov[0]);
751 fhQA[kCutCovElement22][b]->Fill(extCov[2]);
752 fhQA[kCutCovElement33][b]->Fill(extCov[5]);
753 fhQA[kCutCovElement44][b]->Fill(extCov[9]);
754 fhQA[kCutCovElement55][b]->Fill(extCov[14]);
756 // fill cut statistics and cut correlation histograms with information from the bitmap
759 // Get the bitmap of the single cuts
761 SelectionBitMap(obj);
763 // Number of single cuts in this class
764 UInt_t ncuts = fBitmap->GetNbits();
765 for(UInt_t bit=0; bit<ncuts;bit++) {
766 if (!fBitmap->TestBitNumber(bit)) {
767 fhCutStatistics->Fill(bit+1);
768 for (UInt_t bit2=bit; bit2<ncuts;bit2++) {
769 if (!fBitmap->TestBitNumber(bit2))
770 fhCutCorrelation->Fill(bit+1,bit2+1);
775 //__________________________________________________________________________________
776 void AliCFTrackQualityCuts::SaveHistograms(const Char_t* dir) {
778 // saves the histograms in a directory (dir)
785 gDirectory->mkdir(dir);
788 gDirectory->mkdir("before_cuts");
789 gDirectory->mkdir("after_cuts");
791 fhCutStatistics->Write();
792 fhCutCorrelation->Write();
794 for (Int_t j=0; j<kNStepQA; j++) {
796 gDirectory->cd("before_cuts");
798 gDirectory->cd("after_cuts");
800 for(Int_t i=0; i<kNHist; i++) fhQA[i][j]->Write();
802 gDirectory->cd("../");
804 gDirectory->cd("../");
806 //__________________________________________________________________________________
807 void AliCFTrackQualityCuts::DrawHistograms(Bool_t drawLogScale)
810 // draws some histograms
815 Float_t right = 0.03;
816 Float_t left = 0.175;
818 Float_t bottom = 0.175;
820 TCanvas* canvas1 = new TCanvas("Track_QA_Quality_1", "Track QA Quality 1", 800, 500);
821 canvas1->Divide(2, 1);
824 fhCutStatistics->SetStats(kFALSE);
825 fhCutStatistics->LabelsOption("v");
826 gPad->SetLeftMargin(left);
827 gPad->SetBottomMargin(0.25);
828 gPad->SetRightMargin(right);
829 gPad->SetTopMargin(0.1);
830 fhCutStatistics->Draw();
833 fhCutCorrelation->SetStats(kFALSE);
834 fhCutCorrelation->LabelsOption("v");
835 gPad->SetLeftMargin(0.30);
836 gPad->SetRightMargin(bottom);
837 gPad->SetTopMargin(0.1);
838 gPad->SetBottomMargin(0.25);
839 fhCutCorrelation->Draw("COLZ");
841 canvas1->SaveAs(Form("%s.eps", canvas1->GetName()));
842 canvas1->SaveAs(Form("%s.ps", canvas1->GetName()));
846 TCanvas* canvas2 = new TCanvas("Track_QA_Quality_2", "Track QA Quality 2", 1200, 800);
847 canvas2->Divide(2, 2);
850 gPad->SetRightMargin(right);
851 gPad->SetLeftMargin(left);
852 gPad->SetTopMargin(top);
853 gPad->SetBottomMargin(bottom);
854 fhQA[kCutClusterTPC][0]->SetStats(kFALSE);
855 fhQA[kCutClusterTPC][0]->Draw();
856 fhQA[kCutClusterTPC][1]->Draw("same");
859 gPad->SetRightMargin(right);
860 gPad->SetLeftMargin(left);
861 gPad->SetTopMargin(top);
862 gPad->SetBottomMargin(bottom);
863 fhQA[kCutChi2TPC][0]->SetStats(kFALSE);
864 fhQA[kCutChi2TPC][0]->Draw();
865 fhQA[kCutChi2TPC][1]->Draw("same");
868 gPad->SetRightMargin(right);
869 gPad->SetLeftMargin(left);
870 gPad->SetTopMargin(top);
871 gPad->SetBottomMargin(bottom);
872 fhQA[kCutClusterITS][0]->SetStats(kFALSE);
873 fhQA[kCutClusterITS][0]->Draw();
874 fhQA[kCutClusterITS][1]->Draw("same");
877 gPad->SetRightMargin(right);
878 gPad->SetLeftMargin(left);
879 gPad->SetTopMargin(top);
880 gPad->SetBottomMargin(bottom);
881 fhQA[kCutChi2ITS][0]->SetStats(kFALSE);
882 fhQA[kCutChi2ITS][0]->Draw();
883 fhQA[kCutChi2ITS][1]->Draw("same");
885 canvas2->SaveAs(Form("%s.eps", canvas2->GetName()));
886 canvas2->SaveAs(Form("%s.ps", canvas2->GetName()));
890 TCanvas* canvas3 = new TCanvas("Track_QA_Quality_3", "Track QA Quality 3", 1200, 800);
891 canvas3->Divide(3, 2);
894 if(drawLogScale) gPad->SetLogy();
895 gPad->SetRightMargin(right);
896 gPad->SetLeftMargin(left);
897 gPad->SetTopMargin(top);
898 gPad->SetBottomMargin(bottom);
899 fhQA[kCutCovElement11][0]->SetStats(kFALSE);
900 fhQA[kCutCovElement11][0]->Draw();
901 fhQA[kCutCovElement11][1]->Draw("same");
904 if(drawLogScale) gPad->SetLogy();
905 gPad->SetRightMargin(right);
906 gPad->SetLeftMargin(left);
907 gPad->SetTopMargin(top);
908 gPad->SetBottomMargin(bottom);
909 fhQA[kCutCovElement22][0]->SetStats(kFALSE);
910 fhQA[kCutCovElement22][0]->Draw();
911 fhQA[kCutCovElement22][1]->Draw("same");
914 if(drawLogScale) gPad->SetLogy();
915 gPad->SetRightMargin(right);
916 gPad->SetLeftMargin(left);
917 gPad->SetTopMargin(top);
918 gPad->SetBottomMargin(bottom);
919 fhQA[kCutCovElement33][0]->SetStats(kFALSE);
920 fhQA[kCutCovElement33][0]->Draw();
921 fhQA[kCutCovElement33][1]->Draw("same");
924 if(drawLogScale) gPad->SetLogy();
925 gPad->SetRightMargin(right);
926 gPad->SetLeftMargin(left);
927 gPad->SetTopMargin(top);
928 gPad->SetBottomMargin(bottom);
929 fhQA[kCutCovElement44][0]->SetStats(kFALSE);
930 fhQA[kCutCovElement44][0]->Draw();
931 fhQA[kCutCovElement44][1]->Draw("same");
934 if(drawLogScale) gPad->SetLogy();
935 gPad->SetRightMargin(right);
936 gPad->SetLeftMargin(left);
937 gPad->SetTopMargin(top);
938 gPad->SetBottomMargin(bottom);
939 fhQA[kCutCovElement55][0]->SetStats(kFALSE);
940 fhQA[kCutCovElement55][0]->Draw();
941 fhQA[kCutCovElement55][1]->Draw("same");
943 canvas3->SaveAs(Form("%s.eps", canvas3->GetName()));
944 canvas3->SaveAs(Form("%s.ps", canvas3->GetName()));
946 //__________________________________________________________________________________
947 void AliCFTrackQualityCuts::AddQAHistograms(TList *qaList) {
949 // saves the histograms in a TList
953 qaList->Add(fhCutStatistics);
954 qaList->Add(fhCutCorrelation);
956 for (Int_t j=0; j<kNStepQA; j++) {
957 for(Int_t i=0; i<kNHist; i++)
958 qaList->Add(fhQA[i][j]);