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 ITS
25 // - number of clusters in the TPC
26 // - number of clusters in the TRD
27 // - ratio of found / finable number of clusters in the TPC
28 // - number of tracklets in the TRD
29 // - number TRD tracklets used for pid
30 // - chi2 / cluster in the ITS
31 // - chi2 / cluster in the TPC
32 // - chi2 / tracklet in the TRD
33 // - number of clusters in the TPC used for dEdx calculation
34 // - successful TPC refit
35 // - successful ITS refit
36 // - covariance matrix diagonal elements
38 // The cut values for these cuts are set with the corresponding set functions.
39 // All cut classes provided by the correction framework are supposed to be
40 // added in the Analysis Framwork's class AliAnalysisFilter and applied by
41 // the filter via a loop.
43 // author: I. Kraus (Ingrid.Kraus@cern.ch)
45 // AliESDtrackCuts writte by Jan Fiete Grosse-Oetringhaus and
46 // AliRsnDaughterCut class written by A. Pulvirenti.
50 #include <TDirectory.h>
54 #include <AliESDtrack.h>
55 #include <AliESDtrackCuts.h>
57 #include "AliCFTrackQualityCuts.h"
58 #include "AliAODTrack.h"
60 ClassImp(AliCFTrackQualityCuts)
62 //__________________________________________________________________________________
63 AliCFTrackQualityCuts::AliCFTrackQualityCuts() :
68 fMinFoundClusterTPC(0),
70 fMinNTrackletTRDpid(0),
71 fMaxChi2PerClusterTPC(0),
72 fMaxChi2PerClusterITS(0),
73 fMaxChi2PerTrackletTRD(0),
74 fMinNdEdxClusterTPC(0),
88 fhNBinsFoundClusterTPC(0),
89 fhNBinsTrackletTRD(0),
90 fhNBinsTrackletTRDpid(0),
94 fhNBinsdEdxClusterTPC(0),
95 fhNBinsCovariance11(0),
96 fhNBinsCovariance22(0),
97 fhNBinsCovariance33(0),
98 fhNBinsCovariance44(0),
99 fhNBinsCovariance55(0),
100 fhBinLimClusterTPC(0x0),
101 fhBinLimClusterITS(0x0),
102 fhBinLimClusterTRD(0x0),
103 fhBinLimFoundClusterTPC(0x0),
104 fhBinLimTrackletTRD(0x0),
105 fhBinLimTrackletTRDpid(0x0),
106 fhBinLimChi2TPC(0x0),
107 fhBinLimChi2ITS(0x0),
108 fhBinLimChi2TRD(0x0),
109 fhBinLimdEdxClusterTPC(0x0),
110 fhBinLimCovariance11(0x0),
111 fhBinLimCovariance22(0x0),
112 fhBinLimCovariance33(0x0),
113 fhBinLimCovariance44(0x0),
114 fhBinLimCovariance55(0x0)
117 // Default constructor
121 //__________________________________________________________________________________
122 AliCFTrackQualityCuts::AliCFTrackQualityCuts(Char_t* name, Char_t* title) :
123 AliCFCutBase(name,title),
127 fMinFoundClusterTPC(0),
129 fMinNTrackletTRDpid(0),
130 fMaxChi2PerClusterTPC(0),
131 fMaxChi2PerClusterITS(0),
132 fMaxChi2PerTrackletTRD(0),
133 fMinNdEdxClusterTPC(0),
144 fhNBinsClusterTPC(0),
145 fhNBinsClusterITS(0),
146 fhNBinsClusterTRD(0),
147 fhNBinsFoundClusterTPC(0),
148 fhNBinsTrackletTRD(0),
149 fhNBinsTrackletTRDpid(0),
153 fhNBinsdEdxClusterTPC(0),
154 fhNBinsCovariance11(0),
155 fhNBinsCovariance22(0),
156 fhNBinsCovariance33(0),
157 fhNBinsCovariance44(0),
158 fhNBinsCovariance55(0),
159 fhBinLimClusterTPC(0x0),
160 fhBinLimClusterITS(0x0),
161 fhBinLimClusterTRD(0x0),
162 fhBinLimFoundClusterTPC(0x0),
163 fhBinLimTrackletTRD(0x0),
164 fhBinLimTrackletTRDpid(0x0),
165 fhBinLimChi2TPC(0x0),
166 fhBinLimChi2ITS(0x0),
167 fhBinLimChi2TRD(0x0),
168 fhBinLimdEdxClusterTPC(0x0),
169 fhBinLimCovariance11(0x0),
170 fhBinLimCovariance22(0x0),
171 fhBinLimCovariance33(0x0),
172 fhBinLimCovariance44(0x0),
173 fhBinLimCovariance55(0x0)
180 //__________________________________________________________________________________
181 AliCFTrackQualityCuts::AliCFTrackQualityCuts(const AliCFTrackQualityCuts& c) :
183 fMinNClusterTPC(c.fMinNClusterTPC),
184 fMinNClusterITS(c.fMinNClusterITS),
185 fMinNClusterTRD(c.fMinNClusterTRD),
186 fMinFoundClusterTPC(c.fMinFoundClusterTPC),
187 fMinNTrackletTRD(c.fMinNTrackletTRD),
188 fMinNTrackletTRDpid(c.fMinNTrackletTRDpid),
189 fMaxChi2PerClusterTPC(c.fMaxChi2PerClusterTPC),
190 fMaxChi2PerClusterITS(c.fMaxChi2PerClusterITS),
191 fMaxChi2PerTrackletTRD(c.fMaxChi2PerTrackletTRD),
192 fMinNdEdxClusterTPC(c.fMinNdEdxClusterTPC),
193 fCovariance11Max(c.fCovariance11Max),
194 fCovariance22Max(c.fCovariance22Max),
195 fCovariance33Max(c.fCovariance33Max),
196 fCovariance44Max(c.fCovariance44Max),
197 fCovariance55Max(c.fCovariance55Max),
199 fhCutStatistics(c.fhCutStatistics),
200 fhCutCorrelation(c.fhCutCorrelation),
202 fTrackCuts(c.fTrackCuts),
203 fhNBinsClusterTPC(c.fhNBinsClusterTPC),
204 fhNBinsClusterITS(c.fhNBinsClusterITS),
205 fhNBinsClusterTRD(c.fhNBinsClusterTRD),
206 fhNBinsFoundClusterTPC(c.fhNBinsFoundClusterTPC),
207 fhNBinsTrackletTRD(c.fhNBinsTrackletTRD),
208 fhNBinsTrackletTRDpid(c.fhNBinsTrackletTRDpid),
209 fhNBinsChi2TPC(c.fhNBinsChi2TPC),
210 fhNBinsChi2ITS(c.fhNBinsChi2ITS),
211 fhNBinsChi2TRD(c.fhNBinsChi2TRD),
212 fhNBinsdEdxClusterTPC(c.fhNBinsdEdxClusterTPC),
213 fhNBinsCovariance11(c.fhNBinsCovariance11),
214 fhNBinsCovariance22(c.fhNBinsCovariance22),
215 fhNBinsCovariance33(c.fhNBinsCovariance33),
216 fhNBinsCovariance44(c.fhNBinsCovariance44),
217 fhNBinsCovariance55(c.fhNBinsCovariance55),
218 fhBinLimClusterTPC(c.fhBinLimClusterTPC),
219 fhBinLimClusterITS(c.fhBinLimClusterITS),
220 fhBinLimClusterTRD(c.fhBinLimClusterTRD),
221 fhBinLimFoundClusterTPC(c.fhBinLimFoundClusterTPC),
222 fhBinLimTrackletTRD(c.fhBinLimTrackletTRD),
223 fhBinLimTrackletTRDpid(c.fhBinLimTrackletTRDpid),
224 fhBinLimChi2TPC(c.fhBinLimChi2TPC),
225 fhBinLimChi2ITS(c.fhBinLimChi2ITS),
226 fhBinLimChi2TRD(c.fhBinLimChi2TRD),
227 fhBinLimdEdxClusterTPC(c.fhBinLimdEdxClusterTPC),
228 fhBinLimCovariance11(c.fhBinLimCovariance11),
229 fhBinLimCovariance22(c.fhBinLimCovariance22),
230 fhBinLimCovariance33(c.fhBinLimCovariance33),
231 fhBinLimCovariance44(c.fhBinLimCovariance44),
232 fhBinLimCovariance55(c.fhBinLimCovariance55)
237 ((AliCFTrackQualityCuts &) c).Copy(*this);
239 //__________________________________________________________________________________
240 AliCFTrackQualityCuts& AliCFTrackQualityCuts::operator=(const AliCFTrackQualityCuts& c)
243 // Assignment operator
246 AliCFCutBase::operator=(c) ;
247 fMinNClusterTPC = c.fMinNClusterTPC ;
248 fMinNClusterITS = c.fMinNClusterITS ;
249 fMinNClusterTRD = c.fMinNClusterTRD ;
250 fMinFoundClusterTPC = c.fMinFoundClusterTPC ;
251 fMinNTrackletTRD = c.fMinNTrackletTRD ;
252 fMinNTrackletTRDpid = c.fMinNTrackletTRDpid ;
253 fMaxChi2PerClusterTPC = c.fMaxChi2PerClusterTPC ;
254 fMaxChi2PerClusterITS = c.fMaxChi2PerClusterITS ;
255 fMaxChi2PerTrackletTRD = c.fMaxChi2PerTrackletTRD ;
256 fMinNdEdxClusterTPC = c.fMinNdEdxClusterTPC;
257 fCovariance11Max = c.fCovariance11Max ;
258 fCovariance22Max = c.fCovariance22Max ;
259 fCovariance33Max = c.fCovariance33Max ;
260 fCovariance44Max = c.fCovariance44Max ;
261 fCovariance55Max = c.fCovariance55Max ;
262 fStatus = c.fStatus ;
263 fhCutStatistics = c.fhCutStatistics ;
264 fhCutCorrelation = c.fhCutCorrelation ;
265 fBitmap = c.fBitmap ;
266 fTrackCuts = c.fTrackCuts ;
267 fhNBinsClusterTPC = c.fhNBinsClusterTPC ;
268 fhNBinsClusterITS = c.fhNBinsClusterITS ;
269 fhNBinsClusterTRD = c.fhNBinsClusterTRD ;
270 fhNBinsFoundClusterTPC = c.fhNBinsFoundClusterTPC ;
271 fhNBinsTrackletTRD = c.fhNBinsTrackletTRD ;
272 fhNBinsTrackletTRDpid = c.fhNBinsTrackletTRDpid ;
273 fhNBinsChi2TPC = c.fhNBinsChi2TPC ;
274 fhNBinsChi2ITS = c.fhNBinsChi2ITS ;
275 fhNBinsChi2TRD = c.fhNBinsChi2TRD ;
276 fhNBinsdEdxClusterTPC = c.fhNBinsdEdxClusterTPC ;
277 fhNBinsCovariance11 = c.fhNBinsCovariance11 ;
278 fhNBinsCovariance22 = c.fhNBinsCovariance22 ;
279 fhNBinsCovariance33 = c.fhNBinsCovariance33 ;
280 fhNBinsCovariance44 = c.fhNBinsCovariance44 ;
281 fhNBinsCovariance55 = c.fhNBinsCovariance55 ;
282 fhBinLimClusterTPC = c.fhBinLimClusterTPC ;
283 fhBinLimClusterITS = c.fhBinLimClusterITS ;
284 fhBinLimClusterTRD = c.fhBinLimClusterTRD ;
285 fhBinLimFoundClusterTPC = c.fhBinLimFoundClusterTPC ;
286 fhBinLimTrackletTRD = c.fhBinLimTrackletTRD ;
287 fhBinLimTrackletTRDpid = c.fhBinLimTrackletTRDpid ;
288 fhBinLimChi2TPC = c.fhBinLimChi2TPC ;
289 fhBinLimChi2ITS = c.fhBinLimChi2ITS ;
290 fhBinLimChi2TRD = c.fhBinLimChi2TRD ;
291 fhBinLimdEdxClusterTPC = c.fhBinLimdEdxClusterTPC ;
292 fhBinLimCovariance11 = c.fhBinLimCovariance11 ;
293 fhBinLimCovariance22 = c.fhBinLimCovariance22 ;
294 fhBinLimCovariance33 = c.fhBinLimCovariance33 ;
295 fhBinLimCovariance44 = c.fhBinLimCovariance44 ;
296 fhBinLimCovariance55 = c.fhBinLimCovariance55 ;
298 for (Int_t i=0; i<c.kNHist; i++){
299 for (Int_t j=0; j<c.kNStepQA; j++){
300 if(c.fhQA[i][j]) fhQA[i][j] = (TH1F*)c.fhQA[i][j]->Clone();
303 ((AliCFTrackQualityCuts &) c).Copy(*this);
307 //__________________________________________________________________________________
308 AliCFTrackQualityCuts::~AliCFTrackQualityCuts()
313 if (fhCutStatistics) delete fhCutStatistics;
314 if (fhCutCorrelation) delete fhCutCorrelation;
316 for (Int_t i=0; i<kNHist; i++){
317 for (Int_t j=0; j<kNStepQA; j++){
318 if(fhQA[i][j]) delete fhQA[i][j];
321 if(fBitmap) delete fBitmap;
322 if(fTrackCuts) delete fTrackCuts;
323 if(fhBinLimClusterTPC) delete fhBinLimClusterTPC;
324 if(fhBinLimClusterITS) delete fhBinLimClusterITS;
325 if(fhBinLimClusterTRD) delete fhBinLimClusterTRD;
326 if(fhBinLimFoundClusterTPC) delete fhBinLimFoundClusterTPC;
327 if(fhBinLimTrackletTRD) delete fhBinLimTrackletTRD;
328 if(fhBinLimTrackletTRDpid) delete fhBinLimTrackletTRDpid;
329 if(fhBinLimChi2TPC) delete fhBinLimChi2TPC;
330 if(fhBinLimChi2ITS) delete fhBinLimChi2ITS;
331 if(fhBinLimChi2TRD) delete fhBinLimChi2TRD;
332 if(fhBinLimdEdxClusterTPC) delete fhBinLimdEdxClusterTPC;
333 if(fhBinLimCovariance11) delete fhBinLimCovariance11;
334 if(fhBinLimCovariance22) delete fhBinLimCovariance22;
335 if(fhBinLimCovariance33) delete fhBinLimCovariance33;
336 if(fhBinLimCovariance44) delete fhBinLimCovariance44;
337 if(fhBinLimCovariance55) delete fhBinLimCovariance55;
339 //__________________________________________________________________________________
340 void AliCFTrackQualityCuts::Initialise()
343 // sets everything to zero
348 fMinFoundClusterTPC = 0;
349 fMinNTrackletTRD = 0;
350 fMinNTrackletTRDpid = 0;
351 fMaxChi2PerClusterTPC = 0;
352 fMaxChi2PerClusterITS = 0;
353 fMaxChi2PerTrackletTRD = 0;
354 fMinNdEdxClusterTPC = 0;
355 fCovariance11Max = 0;
356 fCovariance22Max = 0;
357 fCovariance33Max = 0;
358 fCovariance44Max = 0;
359 fCovariance55Max = 0;
365 SetMinFoundClusterTPC();
366 SetMinNTrackletTRD();
367 SetMinNTrackletTRDpid();
368 SetMaxChi2PerClusterTPC();
369 SetMaxChi2PerClusterITS();
370 SetMaxChi2PerTrackletTRD();
371 SetMinNdEdxClusterTPC();
372 SetMaxCovDiagonalElements();
375 for (Int_t i=0; i<kNHist; i++){
376 for (Int_t j=0; j<kNStepQA; j++)
380 fhCutCorrelation = 0;
381 fBitmap=new TBits(0);
382 fTrackCuts=new AliESDtrackCuts("aliESDtrackCuts","aliESDtrackCuts");
384 //set default bining for QA histograms
385 SetHistogramBins(kCutClusterTPC,165,-0.5,164.5);
386 SetHistogramBins(kCutClusterITS,8,-0.5,7.5);
387 SetHistogramBins(kCutClusterTRD,120,-0.5,119.5);
388 SetHistogramBins(kCutMinFoundClusterTPC,110,-0.05,1.05);
389 SetHistogramBins(kCutTrackletTRD,7,-0.5,6.5);
390 SetHistogramBins(kCutTrackletTRDpid,7,-0.5,6.5);
391 SetHistogramBins(kCutChi2TPC,500,0.,10.);
392 SetHistogramBins(kCutChi2ITS,500,0.,10.);
393 SetHistogramBins(kCutChi2TRD,500,0.,10.);
394 SetHistogramBins(kCutdEdxClusterTPC,165,-0.5,164.5);
395 SetHistogramBins(kCutCovElement11,200,0.,20.);
396 SetHistogramBins(kCutCovElement22,200,0.,20.);
397 SetHistogramBins(kCutCovElement33,100,0.,1.);
398 SetHistogramBins(kCutCovElement44,100,0.,5.);
399 SetHistogramBins(kCutCovElement55,100,0.,5.);
401 //__________________________________________________________________________________
402 void AliCFTrackQualityCuts::Copy(TObject &c) const
407 AliCFTrackQualityCuts& target = (AliCFTrackQualityCuts &) c;
411 if (fhCutStatistics) target.fhCutStatistics = (TH1F*) fhCutStatistics->Clone();
412 if (fhCutCorrelation) target.fhCutCorrelation = (TH2F*) fhCutCorrelation->Clone();
414 for (Int_t i=0; i<kNHist; i++){
415 for (Int_t j=0; j<kNStepQA; j++){
416 if(fhQA[i][j]) target.fhQA[i][j] = (TH1F*)fhQA[i][j]->Clone();
421 //__________________________________________________________________________________
422 void AliCFTrackQualityCuts::SelectionBitMap(TObject* obj)
425 // test if the track passes the single cuts
426 // and store the information in a bitmap
429 // bitmap stores the decision of each single cut
430 for(Int_t i=0; i<kNCuts; i++)fBitmap->SetBitNumber(i,kFALSE);
433 if (!obj->InheritsFrom("AliVParticle")) {
434 AliError("object must derived from AliVParticle !");
438 Bool_t isESDTrack = strcmp(obj->ClassName(),"AliESDtrack") == 0 ? kTRUE : kFALSE ;
439 Bool_t isAODTrack = strcmp(obj->ClassName(),"AliAODTrack") == 0 ? kTRUE : kFALSE ;
441 AliESDtrack * esdTrack = 0x0 ;
442 AliAODTrack * aodTrack = 0x0 ;
443 if (isESDTrack) esdTrack = dynamic_cast<AliESDtrack*>(obj);
444 if (isAODTrack) aodTrack = dynamic_cast<AliAODTrack*>(obj);
446 fTrackCuts->SetMinNClustersTPC(fMinNClusterTPC);
447 fTrackCuts->SetMinNClustersITS(fMinNClusterITS);
448 fTrackCuts->SetMaxChi2PerClusterTPC(fMaxChi2PerClusterTPC);
449 fTrackCuts->SetMaxChi2PerClusterITS(fMaxChi2PerClusterITS);
450 fTrackCuts->SetMaxCovDiagonalElements(fCovariance11Max,fCovariance22Max,fCovariance33Max,fCovariance44Max,fCovariance55Max);
451 if (isESDTrack) Bool_t ok = fTrackCuts->AcceptTrack(esdTrack);
453 // // // remove following 5 lines when AliESDtrackCuts is updated
454 Int_t nClustersTPC = 0;
455 Int_t nClustersITS = 0 ;
456 Float_t chi2PerClusterTPC = 0 ;
457 Float_t chi2PerClusterITS = 0 ;
458 Double_t extCov[15]={0,0,0,0,0,0,0,0,0,0,0,0,0,0,0};
460 // get cut quantities
461 Int_t nClustersTRD = 0;
462 Int_t nTrackletsTRD = 0;
463 Float_t chi2PerTrackletTRD = 0;
464 Float_t fractionFoundClustersTPC = 0;
467 nClustersTRD = esdTrack->GetTRDncls();
468 nTrackletsTRD = esdTrack->GetTRDntracklets();
469 if (nTrackletsTRD != 0) chi2PerTrackletTRD = esdTrack->GetTRDchi2() / Float_t(nTrackletsTRD);
470 // // // include following line when AliESDtrackCuts is updated
471 // if (esdTrack->GetTPCNclsF() != 0) fractionFoundClustersTPC = fTrackCuts->GetCutVariable(2) / float(esdTrack->GetTPCNclsF());
473 // // // remove following 6 lines when AliESDtrackCuts is updated
474 nClustersTPC = esdTrack->GetTPCclusters(0x0);
475 nClustersITS = esdTrack->GetITSclusters(0x0);
476 if (nClustersTPC != 0) chi2PerClusterTPC = esdTrack->GetTPCchi2() / Float_t(nClustersTPC);
477 if (nClustersITS != 0) chi2PerClusterITS = esdTrack->GetITSchi2() / Float_t(nClustersITS);
478 esdTrack->GetExternalCovariance(extCov);
479 if (esdTrack->GetTPCNclsF() != 0) fractionFoundClustersTPC = float(nClustersTPC) / float(esdTrack->GetTPCNclsF());
485 // // // include following lines when AliESDtrackCuts is updated
486 // fBitmap->SetBitNumber(iCutBit,fTrackCuts->GetCutDecision(2)); iCutBit++; // nClustersTPC
487 // fBitmap->SetBitNumber(iCutBit,fTrackCuts->GetCutDecision(3)); iCutBit++; // nClustersITS
488 // // // remove following 6 lines when AliESDtrackCuts is updated
489 if (nClustersTPC >= fMinNClusterTPC)
490 fBitmap->SetBitNumber(iCutBit,kTRUE);
492 if (nClustersITS >= fMinNClusterITS)
493 fBitmap->SetBitNumber(iCutBit,kTRUE);
496 if (nClustersTRD >= fMinNClusterTRD)
497 fBitmap->SetBitNumber(iCutBit,kTRUE);
499 // // // if ((fMinFoundClusterTPC <= 0) || (fTrackCuts->GetCutVariable(2) > 0 && (fractionFoundClustersTPC >= fMinFoundClusterTPC)))
500 if ((fMinFoundClusterTPC <= 0) || (nClustersTPC > 0 && (fractionFoundClustersTPC >= fMinFoundClusterTPC)))
501 fBitmap->SetBitNumber(iCutBit,kTRUE);
503 if (nTrackletsTRD >= fMinNTrackletTRD)
504 fBitmap->SetBitNumber(iCutBit,kTRUE);
506 if (esdTrack->GetTRDntrackletsPID() >= fMinNTrackletTRDpid)
507 fBitmap->SetBitNumber(iCutBit,kTRUE);
510 // // // include following lines when AliESDtrackCuts is updated
511 // fBitmap->SetBitNumber(iCutBit,fTrackCuts->GetCutDecision(4)); iCutBit++; // chi2PerClusterTPC
512 // fBitmap->SetBitNumber(iCutBit,fTrackCuts->GetCutDecision(5)); iCutBit++; // chi2PerClusterITS
513 // // // remove following 6 lines when AliESDtrackCuts is updated
514 if (chi2PerClusterTPC <= fMaxChi2PerClusterTPC)
515 fBitmap->SetBitNumber(iCutBit,kTRUE);
517 if (chi2PerClusterITS <= fMaxChi2PerClusterITS)
518 fBitmap->SetBitNumber(iCutBit,kTRUE);
521 if (chi2PerTrackletTRD <= fMaxChi2PerTrackletTRD)
522 fBitmap->SetBitNumber(iCutBit,kTRUE);
524 if (esdTrack->GetTPCsignalN() >= fMinNdEdxClusterTPC)
525 fBitmap->SetBitNumber(iCutBit,kTRUE);
528 // // // include following lines when AliESDtrackCuts is updated
529 // fBitmap->SetBitNumber(iCutBit,fTrackCuts->GetCutDecision(6)); iCutBit++; // extCov[0]
530 // fBitmap->SetBitNumber(iCutBit,fTrackCuts->GetCutDecision(7)); iCutBit++; // extCov[2]
531 // fBitmap->SetBitNumber(iCutBit,fTrackCuts->GetCutDecision(8)); iCutBit++; // extCov[5]
532 // fBitmap->SetBitNumber(iCutBit,fTrackCuts->GetCutDecision(9)); iCutBit++; // extCov[9]
533 // fBitmap->SetBitNumber(iCutBit,fTrackCuts->GetCutDecision(10)); iCutBit++; // extCov[14]
534 // // // remove following lines when AliESDtrackCuts is updated
535 if (extCov[0] <= fCovariance11Max)
536 fBitmap->SetBitNumber(iCutBit,kTRUE);
538 if (extCov[2] <= fCovariance22Max)
539 fBitmap->SetBitNumber(iCutBit,kTRUE);
541 if (extCov[5] <= fCovariance33Max)
542 fBitmap->SetBitNumber(iCutBit,kTRUE);
544 if (extCov[9] <= fCovariance44Max)
545 fBitmap->SetBitNumber(iCutBit,kTRUE);
547 if (extCov[14] <= fCovariance55Max)
548 fBitmap->SetBitNumber(iCutBit,kTRUE);
553 if ( (esdTrack->GetStatus() & fStatus) == fStatus ) fBitmap->SetBitNumber(iCutBit,kTRUE);
556 if ( (aodTrack->GetStatus() & fStatus) == fStatus ) fBitmap->SetBitNumber(iCutBit,kTRUE);
561 //__________________________________________________________________________________
562 Bool_t AliCFTrackQualityCuts::IsSelected(TObject* obj) {
564 // loops over decisions of single cuts and returns if the track is accepted
566 SelectionBitMap(obj);
568 if (fIsQAOn) FillHistograms(obj,0);
569 Bool_t isSelected = kTRUE;
571 for (UInt_t icut=0; icut<fBitmap->GetNbits();icut++) {
572 if(!fBitmap->TestBitNumber(icut)) {
577 if (!isSelected) return kFALSE ;
578 if (fIsQAOn) FillHistograms(obj,1);
581 //__________________________________________________________________________________
582 void AliCFTrackQualityCuts::SetHistogramBins(Int_t index, Int_t nbins, Double_t *bins)
589 if (index<0 || index>=kNHist) {
590 Error("SetHistogramBins","could not determine histogram from index %d",index);
596 fhNBinsClusterTPC=nbins+1;
597 fhBinLimClusterTPC=new Double_t[nbins+1];
598 for(Int_t i=0;i<nbins+1;i++)fhBinLimClusterTPC[i]=bins[i];
602 fhNBinsClusterITS=nbins+1;
603 fhBinLimClusterITS=new Double_t[nbins+1];
604 for(Int_t i=0;i<nbins+1;i++)fhBinLimClusterITS[i]=bins[i];
608 fhNBinsClusterTRD=nbins+1;
609 fhBinLimClusterTRD=new Double_t[nbins+1];
610 for(Int_t i=0;i<nbins+1;i++)fhBinLimClusterTRD[i]=bins[i];
613 case kCutMinFoundClusterTPC:
614 fhNBinsFoundClusterTPC=nbins+1;
615 fhBinLimFoundClusterTPC=new Double_t[nbins+1];
616 for(Int_t i=0;i<nbins+1;i++)fhBinLimFoundClusterTPC[i]=bins[i];
619 case kCutTrackletTRD:
620 fhNBinsTrackletTRD=nbins+1;
621 fhBinLimTrackletTRD=new Double_t[nbins+1];
622 for(Int_t i=0;i<nbins+1;i++)fhBinLimTrackletTRD[i]=bins[i];
625 case kCutTrackletTRDpid:
626 fhNBinsTrackletTRDpid=nbins+1;
627 fhBinLimTrackletTRDpid=new Double_t[nbins+1];
628 for(Int_t i=0;i<nbins+1;i++)fhBinLimTrackletTRDpid[i]=bins[i];
632 fhNBinsChi2TPC=nbins+1;
633 fhBinLimChi2TPC=new Double_t[nbins+1];
634 for(Int_t i=0;i<nbins+1;i++)fhBinLimChi2TPC[i]=bins[i];
638 fhNBinsChi2ITS=nbins+1;
639 fhBinLimChi2ITS=new Double_t[nbins+1];
640 for(Int_t i=0;i<nbins+1;i++)fhBinLimChi2ITS[i]=bins[i];
644 fhNBinsChi2TRD=nbins+1;
645 fhBinLimChi2TRD=new Double_t[nbins+1];
646 for(Int_t i=0;i<nbins+1;i++)fhBinLimChi2TRD[i]=bins[i];
649 case kCutdEdxClusterTPC:
650 fhNBinsdEdxClusterTPC=nbins+1;
651 fhBinLimdEdxClusterTPC=new Double_t[nbins+1];
652 for(Int_t i=0;i<nbins+1;i++)fhBinLimdEdxClusterTPC[i]=bins[i];
655 case kCutCovElement11:
656 fhNBinsCovariance11=nbins+1;
657 fhBinLimCovariance11=new Double_t[nbins+1];
658 for(Int_t i=0;i<nbins+1;i++)fhBinLimCovariance11[i]=bins[i];
661 case kCutCovElement22:
662 fhNBinsCovariance22=nbins+1;
663 fhBinLimCovariance22=new Double_t[nbins+1];
664 for(Int_t i=0;i<nbins+1;i++)fhBinLimCovariance22[i]=bins[i];
667 case kCutCovElement33:
668 fhNBinsCovariance33=nbins+1;
669 fhBinLimCovariance33=new Double_t[nbins+1];
670 for(Int_t i=0;i<nbins+1;i++)fhBinLimCovariance33[i]=bins[i];
673 case kCutCovElement44:
674 fhNBinsCovariance44=nbins+1;
675 fhBinLimCovariance44=new Double_t[nbins+1];
676 for(Int_t i=0;i<nbins+1;i++)fhBinLimCovariance44[i]=bins[i];
679 case kCutCovElement55:
680 fhNBinsCovariance55=nbins+1;
681 fhBinLimCovariance55=new Double_t[nbins+1];
682 for(Int_t i=0;i<nbins+1;i++)fhBinLimCovariance55[i]=bins[i];
686 //__________________________________________________________________________________
687 void AliCFTrackQualityCuts::SetHistogramBins(Int_t index, Int_t nbins, Double_t xmin, Double_t xmax)
692 if (index<0 || index>=kNHist) {
693 Error("SetHistogramBins","could not determine histogram from index %d",index);
699 fhNBinsClusterTPC=nbins+1;
700 fhBinLimClusterTPC=new Double_t[nbins+1];
701 for(Int_t i=0;i<nbins+1;i++)fhBinLimClusterTPC[i]=xmin+i*(xmax-xmin)/Double_t(nbins);
705 fhNBinsClusterITS=nbins+1;
706 fhBinLimClusterITS=new Double_t[nbins+1];
707 for(Int_t i=0;i<nbins+1;i++)fhBinLimClusterITS[i]=xmin+i*(xmax-xmin)/Double_t(nbins);
711 fhNBinsClusterTRD=nbins+1;
712 fhBinLimClusterTRD=new Double_t[nbins+1];
713 for(Int_t i=0;i<nbins+1;i++)fhBinLimClusterTRD[i]=xmin+i*(xmax-xmin)/Double_t(nbins);
716 case kCutMinFoundClusterTPC:
717 fhNBinsFoundClusterTPC=nbins+1;
718 fhBinLimFoundClusterTPC=new Double_t[nbins+1];
719 for(Int_t i=0;i<nbins+1;i++)fhBinLimFoundClusterTPC[i]=xmin+i*(xmax-xmin)/Double_t(nbins);
722 case kCutTrackletTRD:
723 fhNBinsTrackletTRD=nbins+1;
724 fhBinLimTrackletTRD=new Double_t[nbins+1];
725 for(Int_t i=0;i<nbins+1;i++)fhBinLimTrackletTRD[i]=xmin+i*(xmax-xmin)/Double_t(nbins);
728 case kCutTrackletTRDpid:
729 fhNBinsTrackletTRDpid=nbins+1;
730 fhBinLimTrackletTRDpid=new Double_t[nbins+1];
731 for(Int_t i=0;i<nbins+1;i++)fhBinLimTrackletTRDpid[i]=xmin+i*(xmax-xmin)/Double_t(nbins);
735 fhNBinsChi2TPC=nbins+1;
736 fhBinLimChi2TPC=new Double_t[nbins+1];
737 for(Int_t i=0;i<nbins+1;i++)fhBinLimChi2TPC[i]=xmin+i*(xmax-xmin)/Double_t(nbins);
741 fhNBinsChi2ITS=nbins+1;
742 fhBinLimChi2ITS=new Double_t[nbins+1];
743 for(Int_t i=0;i<nbins+1;i++)fhBinLimChi2ITS[i]=xmin+i*(xmax-xmin)/Double_t(nbins);
747 fhNBinsChi2TRD=nbins+1;
748 fhBinLimChi2TRD=new Double_t[nbins+1];
749 for(Int_t i=0;i<nbins+1;i++)fhBinLimChi2TRD[i]=xmin+i*(xmax-xmin)/Double_t(nbins);
752 case kCutdEdxClusterTPC:
753 fhNBinsdEdxClusterTPC=nbins+1;
754 fhBinLimdEdxClusterTPC=new Double_t[nbins+1];
755 for(Int_t i=0;i<nbins+1;i++)fhBinLimdEdxClusterTPC[i]=xmin+i*(xmax-xmin)/Double_t(nbins);
758 case kCutCovElement11:
759 fhNBinsCovariance11=nbins+1;
760 fhBinLimCovariance11=new Double_t[nbins+1];
761 for(Int_t i=0;i<nbins+1;i++)fhBinLimCovariance11[i]=xmin+i*(xmax-xmin)/Double_t(nbins);
764 case kCutCovElement22:
765 fhNBinsCovariance22=nbins+1;
766 fhBinLimCovariance22=new Double_t[nbins+1];
767 for(Int_t i=0;i<nbins+1;i++)fhBinLimCovariance22[i]=xmin+i*(xmax-xmin)/Double_t(nbins);
770 case kCutCovElement33:
771 fhNBinsCovariance33=nbins+1;
772 fhBinLimCovariance33=new Double_t[nbins+1];
773 for(Int_t i=0;i<nbins+1;i++)fhBinLimCovariance33[i]=xmin+i*(xmax-xmin)/Double_t(nbins);
776 case kCutCovElement44:
777 fhNBinsCovariance44=nbins+1;
778 fhBinLimCovariance44=new Double_t[nbins+1];
779 for(Int_t i=0;i<nbins+1;i++)fhBinLimCovariance44[i]=xmin+i*(xmax-xmin)/Double_t(nbins);
782 case kCutCovElement55:
783 fhNBinsCovariance55=nbins+1;
784 fhBinLimCovariance55=new Double_t[nbins+1];
785 for(Int_t i=0;i<nbins+1;i++)fhBinLimCovariance55[i]=xmin+i*(xmax-xmin)/Double_t(nbins);
789 //__________________________________________________________________________________
790 void AliCFTrackQualityCuts::DefineHistograms() {
792 // histograms for cut variables, cut statistics and cut correlations
796 // book cut statistics and cut correlation histograms
797 fhCutStatistics = new TH1F(Form("%s_cut_statistics",GetName()), Form("%s cut statistics",GetName()), kNCuts,0.5,kNCuts+0.5);
798 fhCutStatistics->SetLineWidth(2);
799 fhCutStatistics->GetXaxis()->SetBinLabel(1, "nClustersTPC");
800 fhCutStatistics->GetXaxis()->SetBinLabel(2, "nClustersITS");
801 fhCutStatistics->GetXaxis()->SetBinLabel(3, "nClustersTRD");
802 fhCutStatistics->GetXaxis()->SetBinLabel(4, "fractionClustersTPC");
803 fhCutStatistics->GetXaxis()->SetBinLabel(5, "ntrackletsTRD");
804 fhCutStatistics->GetXaxis()->SetBinLabel(6, "ntrackletsTRDpid");
805 fhCutStatistics->GetXaxis()->SetBinLabel(7, "chi2PerClusterTPC");
806 fhCutStatistics->GetXaxis()->SetBinLabel(8, "chi2PerClusterITS");
807 fhCutStatistics->GetXaxis()->SetBinLabel(9, "chi2PerTrackletTRD");
808 fhCutStatistics->GetXaxis()->SetBinLabel(10, "ndEdxClusterTPC");
809 fhCutStatistics->GetXaxis()->SetBinLabel(11, "covElement11");
810 fhCutStatistics->GetXaxis()->SetBinLabel(12, "covElement22");
811 fhCutStatistics->GetXaxis()->SetBinLabel(13, "covElement33");
812 fhCutStatistics->GetXaxis()->SetBinLabel(14, "covElement44");
813 fhCutStatistics->GetXaxis()->SetBinLabel(15, "covElement55");
814 fhCutStatistics->GetXaxis()->SetBinLabel(16, "status");
816 fhCutCorrelation = new TH2F(Form("%s_cut_correlation",GetName()), Form("%s cut correlation",GetName()), kNCuts,0.5,kNCuts+0.5,kNCuts,0.5,kNCuts+0.5);
817 fhCutCorrelation->SetLineWidth(2);
818 fhCutCorrelation->GetXaxis()->SetBinLabel(1, "nClustersTPC");
819 fhCutCorrelation->GetXaxis()->SetBinLabel(2, "nClustersITS");
820 fhCutCorrelation->GetXaxis()->SetBinLabel(3, "nClustersTRD");
821 fhCutCorrelation->GetXaxis()->SetBinLabel(4, "fractionClustersTPC");
822 fhCutCorrelation->GetXaxis()->SetBinLabel(5, "ntrackletsTRD");
823 fhCutCorrelation->GetXaxis()->SetBinLabel(6, "ntrackletsTRDpid");
824 fhCutCorrelation->GetXaxis()->SetBinLabel(7, "chi2PerClusterTPC");
825 fhCutCorrelation->GetXaxis()->SetBinLabel(8, "chi2PerClusterITS");
826 fhCutCorrelation->GetXaxis()->SetBinLabel(9, "chi2PerTrackletTRD");
827 fhCutCorrelation->GetXaxis()->SetBinLabel(10, "ndEdxClustersTPC");
828 fhCutCorrelation->GetXaxis()->SetBinLabel(11, "covElement11");
829 fhCutCorrelation->GetXaxis()->SetBinLabel(12, "covElement22");
830 fhCutCorrelation->GetXaxis()->SetBinLabel(13, "covElement33");
831 fhCutCorrelation->GetXaxis()->SetBinLabel(14, "covElement44");
832 fhCutCorrelation->GetXaxis()->SetBinLabel(15, "covElement55");
833 fhCutCorrelation->GetXaxis()->SetBinLabel(16, "status");
835 fhCutCorrelation->GetYaxis()->SetBinLabel(1, "nClustersTPC");
836 fhCutCorrelation->GetYaxis()->SetBinLabel(2, "nClustersITS");
837 fhCutCorrelation->GetYaxis()->SetBinLabel(3, "nClustersTRD");
838 fhCutCorrelation->GetYaxis()->SetBinLabel(4, "fractionClustersTPC");
839 fhCutCorrelation->GetYaxis()->SetBinLabel(5, "ntrackletsTRD");
840 fhCutCorrelation->GetYaxis()->SetBinLabel(6, "ntrackletsTRDpid");
841 fhCutCorrelation->GetYaxis()->SetBinLabel(7, "chi2PerClusterTPC");
842 fhCutCorrelation->GetYaxis()->SetBinLabel(8, "chi2PerClusterITS");
843 fhCutCorrelation->GetYaxis()->SetBinLabel(9, "chi2PerTrackletTRD");
844 fhCutCorrelation->GetYaxis()->SetBinLabel(10, "ndEdxClustersTPC");
845 fhCutCorrelation->GetYaxis()->SetBinLabel(11, "covElement11");
846 fhCutCorrelation->GetYaxis()->SetBinLabel(12, "covElement22");
847 fhCutCorrelation->GetYaxis()->SetBinLabel(13, "covElement33");
848 fhCutCorrelation->GetYaxis()->SetBinLabel(14, "covElement44");
849 fhCutCorrelation->GetYaxis()->SetBinLabel(15, "covElement55");
850 fhCutCorrelation->GetYaxis()->SetBinLabel(16, "status");
853 // book QA histograms
855 for (Int_t i=0; i<kNStepQA; i++) {
856 if (i==0) sprintf(str," ");
857 else sprintf(str,"_cut");
859 fhQA[kCutClusterTPC][i] = new TH1F(Form("%s_nClustersTPC%s",GetName(),str) ,"",fhNBinsClusterTPC-1,fhBinLimClusterTPC);
860 fhQA[kCutClusterITS][i] = new TH1F(Form("%s_nClustersITS%s",GetName(),str) ,"",fhNBinsClusterITS-1,fhBinLimClusterITS);
861 fhQA[kCutClusterTRD][i] = new TH1F(Form("%s_nClustersTRD%s",GetName(),str) ,"",fhNBinsClusterTRD-1,fhBinLimClusterTRD);
862 fhQA[kCutMinFoundClusterTPC][i] = new TH1F(Form("%s_fractionClustersTPC%s",GetName(),str) ,"",fhNBinsFoundClusterTPC-1,fhBinLimFoundClusterTPC);
863 fhQA[kCutTrackletTRD][i] = new TH1F(Form("%s_ntrackletsTRD%s",GetName(),str) ,"",fhNBinsTrackletTRD-1,fhBinLimTrackletTRD);
864 fhQA[kCutTrackletTRDpid][i] = new TH1F(Form("%s_ntrackletsTRDpid%s",GetName(),str) ,"",fhNBinsTrackletTRDpid-1,fhBinLimTrackletTRDpid);
865 fhQA[kCutChi2TPC][i] = new TH1F(Form("%s_chi2PerClusterTPC%s",GetName(),str),"",fhNBinsChi2TPC-1,fhBinLimChi2TPC);
866 fhQA[kCutChi2ITS][i] = new TH1F(Form("%s_chi2PerClusterITS%s",GetName(),str),"",fhNBinsChi2ITS-1,fhBinLimChi2ITS);
867 fhQA[kCutChi2TRD][i] = new TH1F(Form("%s_chi2PerTrackletTRD%s",GetName(),str),"",fhNBinsChi2TRD-1,fhBinLimChi2TRD);
868 fhQA[kCutdEdxClusterTPC][i] = new TH1F(Form("%s_ndEdxClustersTPC%s",GetName(),str) ,"",fhNBinsdEdxClusterTPC-1,fhBinLimdEdxClusterTPC);
869 fhQA[kCutCovElement11][i] = new TH1F(Form("%s_covMatrixDiagonal11%s",GetName(),str),"",fhNBinsCovariance11-1,fhBinLimCovariance11);
870 fhQA[kCutCovElement22][i] = new TH1F(Form("%s_covMatrixDiagonal22%s",GetName(),str),"",fhNBinsCovariance22-1,fhBinLimCovariance22);
871 fhQA[kCutCovElement33][i] = new TH1F(Form("%s_covMatrixDiagonal33%s",GetName(),str),"",fhNBinsCovariance33-1,fhBinLimCovariance33);
872 fhQA[kCutCovElement44][i] = new TH1F(Form("%s_covMatrixDiagonal44%s",GetName(),str),"",fhNBinsCovariance44-1,fhBinLimCovariance44);
873 fhQA[kCutCovElement55][i] = new TH1F(Form("%s_covMatrixDiagonal55%s",GetName(),str),"",fhNBinsCovariance55-1,fhBinLimCovariance55);
875 fhQA[kCutClusterTPC][i] ->SetXTitle("n TPC clusters");
876 fhQA[kCutClusterITS][i] ->SetXTitle("n ITS clusters");
877 fhQA[kCutClusterTRD][i] ->SetXTitle("n TRD clusters");
878 fhQA[kCutMinFoundClusterTPC][i]->SetXTitle("fraction TPC clusters");
879 fhQA[kCutTrackletTRD][i] ->SetXTitle("n tracklets TRD");
880 fhQA[kCutTrackletTRDpid][i]->SetXTitle("n tracklets TRD pid");
881 fhQA[kCutChi2TPC][i] ->SetXTitle("#chi^{2} per TPC cluster");
882 fhQA[kCutChi2ITS][i] ->SetXTitle("#chi^{2} per ITS cluster");
883 fhQA[kCutChi2TRD][i] ->SetXTitle("#chi^{2} per TRD tracklet");
884 fhQA[kCutdEdxClusterTPC][i] ->SetXTitle("n dEdx TPC clusters");
885 fhQA[kCutCovElement11][i] ->SetXTitle("cov 11 : #sigma_{y}^{2} (cm^{2})");
886 fhQA[kCutCovElement22][i] ->SetXTitle("cov 22 : #sigma_{z}^{2} (cm^{2})");
887 fhQA[kCutCovElement33][i] ->SetXTitle("cov 33 : #sigma_{sin(#phi)}^{2}");
888 fhQA[kCutCovElement44][i] ->SetXTitle("cov 44 : #sigma_{tan(#theta_{dip})}^{2}");
889 fhQA[kCutCovElement55][i] ->SetXTitle("cov 55 : #sigma_{1/p_{T}}^{2} ((c/GeV)^{2})");
892 for(Int_t i=0; i<kNHist; i++) fhQA[i][1]->SetLineColor(color);
894 //__________________________________________________________________________________
895 void AliCFTrackQualityCuts::FillHistograms(TObject* obj, Bool_t b)
898 // fill the QA histograms
902 if (!obj->InheritsFrom("AliVParticle")) {
903 AliError("object must derived from AliVParticle !");
907 Bool_t isESDTrack = strcmp(obj->ClassName(),"AliESDtrack") == 0 ? kTRUE : kFALSE ;
908 Bool_t isAODTrack = strcmp(obj->ClassName(),"AliAODTrack") == 0 ? kTRUE : kFALSE ;
910 AliESDtrack * esdTrack = 0x0 ;
911 AliAODTrack * aodTrack = 0x0 ;
912 if (isESDTrack) esdTrack = dynamic_cast<AliESDtrack*>(obj);
913 if (isAODTrack) aodTrack = dynamic_cast<AliAODTrack*>(obj);
915 // b = 0: fill histograms before cuts
916 // b = 1: fill histograms after cuts
918 // // // remove following 5 lines when AliESDtrackCuts is updated
919 Int_t nClustersTPC = 0;
920 Int_t nClustersITS = 0 ;
921 Float_t chi2PerClusterTPC = 0 ;
922 Float_t chi2PerClusterITS = 0 ;
923 Double_t extCov[15]={0,0,0,0,0,0,0,0,0,0,0,0,0,0,0};
925 Int_t nClustersTRD = 0 ;
926 Int_t nTrackletsTRD = 0 ;
927 Float_t chi2PerTrackletTRD = 0 ;
928 Float_t fractionFoundClustersTPC = 0;
931 nClustersTRD = esdTrack->GetTRDncls();
932 nTrackletsTRD = esdTrack->GetTRDntracklets();
933 if (nTrackletsTRD != 0) chi2PerTrackletTRD = esdTrack->GetTRDchi2() / Float_t(nTrackletsTRD);
934 if (esdTrack->GetTPCNclsF() != 0) fractionFoundClustersTPC = fTrackCuts->GetCutVariable(2) / float(esdTrack->GetTPCNclsF());
936 // // // include following line when AliESDtrackCuts is updated
937 // if (esdTrack->GetTPCNclsF() != 0) fractionFoundClustersTPC = fTrackCuts->GetCutVariable(2) / float(esdTrack->GetTPCNclsF());
939 // // // remove following 6 lines when AliESDtrackCuts is updated
940 nClustersTPC = esdTrack->GetTPCclusters(0x0);
941 nClustersITS = esdTrack->GetITSclusters(0x0);
942 if (nClustersTPC != 0) chi2PerClusterTPC = esdTrack->GetTPCchi2() / Float_t(nClustersTPC);
943 if (nClustersITS != 0) chi2PerClusterITS = esdTrack->GetITSchi2() / Float_t(nClustersITS);
944 esdTrack->GetExternalCovariance(extCov);
945 if (esdTrack->GetTPCNclsF() != 0) fractionFoundClustersTPC = float(nClustersTPC) / float(esdTrack->GetTPCNclsF());
948 // // // include following lines when AliESDtrackCuts is updated
949 // fhQA[kCutClusterTPC][b]->Fill(fTrackCuts->GetCutVariable(2));
950 // fhQA[kCutClusterITS][b]->Fill(fTrackCuts->GetCutVariable(3));
951 // fhQA[kCutChi2TPC][b]->Fill(fTrackCuts->GetCutVariable(4));
952 // fhQA[kCutChi2ITS][b]->Fill(fTrackCuts->GetCutVariable(5));
953 // // // remove following 4 lines when AliESDtrackCuts is updated
954 fhQA[kCutClusterTPC][b]->Fill((float)nClustersTPC);
955 fhQA[kCutChi2TPC][b]->Fill(chi2PerClusterTPC);
956 fhQA[kCutClusterITS][b]->Fill((float)nClustersITS);
957 fhQA[kCutChi2ITS][b]->Fill(chi2PerClusterITS);
960 fhQA[kCutClusterTRD][b]->Fill((float)nClustersTRD);
961 fhQA[kCutChi2TRD][b]->Fill(chi2PerTrackletTRD);
962 // if (b==0 || (b==1 && fTrackCuts->GetCutVariable(2)>0)) fhQA[kCutMinFoundClusterTPC][b]->Fill((float)fractionFoundClustersTPC);
963 if (b==0 || (b==1 && nClustersTPC>0)) fhQA[kCutMinFoundClusterTPC][b]->Fill((float)fractionFoundClustersTPC);
964 fhQA[kCutTrackletTRD][b]->Fill((float)nTrackletsTRD);
965 fhQA[kCutTrackletTRDpid][b]->Fill((float)esdTrack->GetTRDntrackletsPID());
966 fhQA[kCutdEdxClusterTPC][b]->Fill((float)esdTrack->GetTPCsignalN());
968 // // // include following lines when AliESDtrackCuts is updated
969 // fhQA[kCutCovElement11][b]->Fill(fTrackCuts->GetCutVariable(6));
970 // fhQA[kCutCovElement22][b]->Fill(fTrackCuts->GetCutVariable(7));
971 // fhQA[kCutCovElement33][b]->Fill(fTrackCuts->GetCutVariable(8));
972 // fhQA[kCutCovElement44][b]->Fill(fTrackCuts->GetCutVariable(9));
973 // fhQA[kCutCovElement55][b]->Fill(fTrackCuts->GetCutVariable(10));
974 // // // remove following 5 lines when AliESDtrackCuts is updated
975 fhQA[kCutCovElement11][b]->Fill(extCov[0]);
976 fhQA[kCutCovElement22][b]->Fill(extCov[2]);
977 fhQA[kCutCovElement33][b]->Fill(extCov[5]);
978 fhQA[kCutCovElement44][b]->Fill(extCov[9]);
979 fhQA[kCutCovElement55][b]->Fill(extCov[14]);
982 // fill cut statistics and cut correlation histograms with information from the bitmap
985 // Get the bitmap of the single cuts
987 SelectionBitMap(obj);
989 // Number of single cuts in this class
990 UInt_t ncuts = fBitmap->GetNbits();
991 for(UInt_t bit=0; bit<ncuts;bit++) {
992 if (!fBitmap->TestBitNumber(bit)) {
993 fhCutStatistics->Fill(bit+1);
994 for (UInt_t bit2=bit; bit2<ncuts;bit2++) {
995 if (!fBitmap->TestBitNumber(bit2))
996 fhCutCorrelation->Fill(bit+1,bit2+1);
1001 //__________________________________________________________________________________
1002 void AliCFTrackQualityCuts::SaveHistograms(const Char_t* dir) {
1004 // saves the histograms in a directory (dir)
1006 if(!fIsQAOn) return;
1011 gDirectory->mkdir(dir);
1012 gDirectory->cd(dir);
1014 gDirectory->mkdir("before_cuts");
1015 gDirectory->mkdir("after_cuts");
1017 fhCutStatistics->Write();
1018 fhCutCorrelation->Write();
1020 for (Int_t j=0; j<kNStepQA; j++) {
1022 gDirectory->cd("before_cuts");
1024 gDirectory->cd("after_cuts");
1026 for(Int_t i=0; i<kNHist; i++) fhQA[i][j]->Write();
1028 gDirectory->cd("../");
1030 gDirectory->cd("../");
1032 //__________________________________________________________________________________
1033 void AliCFTrackQualityCuts::DrawHistograms(Bool_t drawLogScale)
1036 // draws some histograms
1038 if(!fIsQAOn) return;
1041 Float_t right = 0.03;
1042 Float_t left = 0.175;
1044 Float_t bottom = 0.175;
1046 TCanvas* canvas1 = new TCanvas("Track_QA_Quality_1", "Track QA Quality 1", 800, 500);
1047 canvas1->Divide(2, 1);
1050 fhCutStatistics->SetStats(kFALSE);
1051 fhCutStatistics->LabelsOption("v");
1052 gPad->SetLeftMargin(left);
1053 gPad->SetBottomMargin(0.25);
1054 gPad->SetRightMargin(right);
1055 gPad->SetTopMargin(0.1);
1056 fhCutStatistics->Draw();
1059 fhCutCorrelation->SetStats(kFALSE);
1060 fhCutCorrelation->LabelsOption("v");
1061 gPad->SetLeftMargin(0.30);
1062 gPad->SetRightMargin(bottom);
1063 gPad->SetTopMargin(0.1);
1064 gPad->SetBottomMargin(0.25);
1065 fhCutCorrelation->Draw("COLZ");
1067 canvas1->SaveAs(Form("%s.eps", canvas1->GetName()));
1068 canvas1->SaveAs(Form("%s.ps", canvas1->GetName()));
1072 TCanvas* canvas2 = new TCanvas("Track_QA_Quality_2", "Track QA Quality 2", 1200, 800);
1073 canvas2->Divide(2, 2);
1076 gPad->SetRightMargin(right);
1077 gPad->SetLeftMargin(left);
1078 gPad->SetTopMargin(top);
1079 gPad->SetBottomMargin(bottom);
1080 fhQA[kCutClusterTPC][0]->SetStats(kFALSE);
1081 fhQA[kCutClusterTPC][0]->Draw();
1082 fhQA[kCutClusterTPC][1]->Draw("same");
1085 gPad->SetRightMargin(right);
1086 gPad->SetLeftMargin(left);
1087 gPad->SetTopMargin(top);
1088 gPad->SetBottomMargin(bottom);
1089 fhQA[kCutChi2TPC][0]->SetStats(kFALSE);
1090 fhQA[kCutChi2TPC][0]->Draw();
1091 fhQA[kCutChi2TPC][1]->Draw("same");
1094 gPad->SetRightMargin(right);
1095 gPad->SetLeftMargin(left);
1096 gPad->SetTopMargin(top);
1097 gPad->SetBottomMargin(bottom);
1098 fhQA[kCutClusterITS][0]->SetStats(kFALSE);
1099 fhQA[kCutClusterITS][0]->Draw();
1100 fhQA[kCutClusterITS][1]->Draw("same");
1103 gPad->SetRightMargin(right);
1104 gPad->SetLeftMargin(left);
1105 gPad->SetTopMargin(top);
1106 gPad->SetBottomMargin(bottom);
1107 fhQA[kCutChi2ITS][0]->SetStats(kFALSE);
1108 fhQA[kCutChi2ITS][0]->Draw();
1109 fhQA[kCutChi2ITS][1]->Draw("same");
1111 canvas2->SaveAs(Form("%s.eps", canvas2->GetName()));
1112 canvas2->SaveAs(Form("%s.ps", canvas2->GetName()));
1116 TCanvas* canvas3 = new TCanvas("Track_QA_Quality_3", "Track QA Quality 3", 1200, 800);
1117 canvas3->Divide(3, 2);
1120 if(drawLogScale) gPad->SetLogy();
1121 gPad->SetRightMargin(right);
1122 gPad->SetLeftMargin(left);
1123 gPad->SetTopMargin(top);
1124 gPad->SetBottomMargin(bottom);
1125 fhQA[kCutCovElement11][0]->SetStats(kFALSE);
1126 fhQA[kCutCovElement11][0]->Draw();
1127 fhQA[kCutCovElement11][1]->Draw("same");
1130 if(drawLogScale) gPad->SetLogy();
1131 gPad->SetRightMargin(right);
1132 gPad->SetLeftMargin(left);
1133 gPad->SetTopMargin(top);
1134 gPad->SetBottomMargin(bottom);
1135 fhQA[kCutCovElement22][0]->SetStats(kFALSE);
1136 fhQA[kCutCovElement22][0]->Draw();
1137 fhQA[kCutCovElement22][1]->Draw("same");
1140 if(drawLogScale) gPad->SetLogy();
1141 gPad->SetRightMargin(right);
1142 gPad->SetLeftMargin(left);
1143 gPad->SetTopMargin(top);
1144 gPad->SetBottomMargin(bottom);
1145 fhQA[kCutCovElement33][0]->SetStats(kFALSE);
1146 fhQA[kCutCovElement33][0]->Draw();
1147 fhQA[kCutCovElement33][1]->Draw("same");
1150 if(drawLogScale) gPad->SetLogy();
1151 gPad->SetRightMargin(right);
1152 gPad->SetLeftMargin(left);
1153 gPad->SetTopMargin(top);
1154 gPad->SetBottomMargin(bottom);
1155 fhQA[kCutCovElement44][0]->SetStats(kFALSE);
1156 fhQA[kCutCovElement44][0]->Draw();
1157 fhQA[kCutCovElement44][1]->Draw("same");
1160 if(drawLogScale) gPad->SetLogy();
1161 gPad->SetRightMargin(right);
1162 gPad->SetLeftMargin(left);
1163 gPad->SetTopMargin(top);
1164 gPad->SetBottomMargin(bottom);
1165 fhQA[kCutCovElement55][0]->SetStats(kFALSE);
1166 fhQA[kCutCovElement55][0]->Draw();
1167 fhQA[kCutCovElement55][1]->Draw("same");
1169 canvas3->SaveAs(Form("%s.eps", canvas3->GetName()));
1170 canvas3->SaveAs(Form("%s.ps", canvas3->GetName()));
1172 //__________________________________________________________________________________
1173 void AliCFTrackQualityCuts::AddQAHistograms(TList *qaList) {
1175 // saves the histograms in a TList
1179 qaList->Add(fhCutStatistics);
1180 qaList->Add(fhCutCorrelation);
1182 for (Int_t j=0; j<kNStepQA; j++) {
1183 for(Int_t i=0; i<kNHist; i++)
1184 qaList->Add(fhQA[i][j]);