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 AliCFTrackIsPrimaryCut is designed to select reconstructed tracks
17 // with a small impact parameter and tracks which are (not) daughters of kink
18 // decays and to provide corresponding QA histograms.
19 // This class inherits from the Analysis' Framework abstract base class
20 // AliAnalysisCuts and is a part of the Correction Framework.
21 // This class acts on single, reconstructed tracks, it is applicable on
23 // It mainly consists of a IsSelected function that returns a boolean.
24 // This function checks whether the considered track passes a set of cuts:
25 // - distance to main vertex in units of sigma (resolution)
26 // - require that the dca calculation doesn't fail
27 // - accept or not accept daughter tracks of kink decays
29 // The cut values for these cuts are set with the corresponding set functions.
30 // All cut classes provided by the correction framework are supposed to be
31 // added in the Analysis Framwork's class AliAnalysisFilter and applied by
32 // the filter via a loop.
34 // author: I. Kraus (Ingrid.Kraus@cern.ch)
36 // AliESDtrackCuts writte by Jan Fiete Grosse-Oetringhaus and
37 // AliRsnDaughterCut class written by A. Pulvirenti.
40 #include <TDirectory.h>
43 #include <AliESDtrack.h>
45 #include "AliCFTrackIsPrimaryCuts.h"
47 ClassImp(AliCFTrackIsPrimaryCuts)
49 //__________________________________________________________________________________
50 AliCFTrackIsPrimaryCuts::AliCFTrackIsPrimaryCuts() :
52 fNSigmaToVertexMax(0),
53 fRequireSigmaToVertex(0),
54 fAcceptKinkDaughters(0),
59 fhNBinsRequireSigma(0),
66 fhBinLimRequireSigma(0x0),
67 fhBinLimAcceptKink(0x0),
70 fhBinLimDcaXYnorm(0x0),
74 // Default constructor
79 //__________________________________________________________________________________
80 AliCFTrackIsPrimaryCuts::AliCFTrackIsPrimaryCuts(Char_t* name, Char_t* title) :
81 AliCFCutBase(name,title),
82 fNSigmaToVertexMax(0),
83 fRequireSigmaToVertex(0),
84 fAcceptKinkDaughters(0),
89 fhNBinsRequireSigma(0),
96 fhBinLimRequireSigma(0x0),
97 fhBinLimAcceptKink(0x0),
100 fhBinLimDcaXYnorm(0x0),
101 fhBinLimDcaZnorm(0x0)
106 fBitmap=new TBits(0);
109 //__________________________________________________________________________________
110 AliCFTrackIsPrimaryCuts::AliCFTrackIsPrimaryCuts(const AliCFTrackIsPrimaryCuts& c) :
112 fNSigmaToVertexMax(c.fNSigmaToVertexMax),
113 fRequireSigmaToVertex(c.fRequireSigmaToVertex),
114 fAcceptKinkDaughters(c.fAcceptKinkDaughters),
115 fhCutStatistics(c.fhCutStatistics),
116 fhCutCorrelation(c.fhCutCorrelation),
118 fhNBinsNSigma(c.fhNBinsNSigma),
119 fhNBinsRequireSigma(c.fhNBinsRequireSigma),
120 fhNBinsAcceptKink(c.fhNBinsAcceptKink),
121 fhNBinsDcaXY(c.fhNBinsDcaXY),
122 fhNBinsDcaZ(c.fhNBinsDcaZ),
123 fhNBinsDcaXYnorm(c.fhNBinsDcaXYnorm),
124 fhNBinsDcaZnorm(c.fhNBinsDcaZnorm),
125 fhBinLimNSigma(c.fhBinLimNSigma),
126 fhBinLimRequireSigma(c.fhBinLimRequireSigma),
127 fhBinLimAcceptKink(c.fhBinLimAcceptKink),
128 fhBinLimDcaXY(c.fhBinLimDcaXY),
129 fhBinLimDcaZ(c.fhBinLimDcaZ),
130 fhBinLimDcaXYnorm(c.fhBinLimDcaXYnorm),
131 fhBinLimDcaZnorm(c.fhBinLimDcaZnorm)
136 ((AliCFTrackIsPrimaryCuts &) c).Copy(*this);
138 //__________________________________________________________________________________
139 AliCFTrackIsPrimaryCuts& AliCFTrackIsPrimaryCuts::operator=(const AliCFTrackIsPrimaryCuts& c)
142 // Assignment operator
145 AliCFCutBase::operator=(c) ;
146 fNSigmaToVertexMax = c.fNSigmaToVertexMax ;
147 fRequireSigmaToVertex = c.fRequireSigmaToVertex ;
148 fAcceptKinkDaughters = c.fAcceptKinkDaughters ;
149 fhCutStatistics = c.fhCutStatistics ;
150 fhCutCorrelation = c.fhCutCorrelation ;
152 fhNBinsNSigma = c.fhNBinsNSigma;
153 fhNBinsRequireSigma = c.fhNBinsRequireSigma;
154 fhNBinsAcceptKink = c.fhNBinsAcceptKink;
155 fhNBinsDcaXY = c.fhNBinsDcaXY;
156 fhNBinsDcaZ = c.fhNBinsDcaZ;
157 fhNBinsDcaXYnorm = c.fhNBinsDcaXYnorm;
158 fhNBinsDcaZnorm = c.fhNBinsDcaZnorm;
159 fhBinLimNSigma = c.fhBinLimNSigma;
160 fhBinLimRequireSigma = c.fhBinLimRequireSigma;
161 fhBinLimAcceptKink = c.fhBinLimAcceptKink;
162 fhBinLimDcaXY = c.fhBinLimDcaXY;
163 fhBinLimDcaZ = c.fhBinLimDcaZ;
164 fhBinLimDcaXYnorm = c.fhBinLimDcaXYnorm;
165 fhBinLimDcaZnorm = c.fhBinLimDcaZnorm;
167 for (Int_t i=0; i<c.kNHist; i++){
168 for (Int_t j=0; j<c.kNStepQA; j++){
169 if(c.fhQA[i][j]) fhQA[i][j] = (TH1F*)c.fhQA[i][j]->Clone();
170 if(c.fhDcaXYvsDcaZ[j]) fhDcaXYvsDcaZ[j] = (TH2F*)c.fhDcaXYvsDcaZ[j]->Clone();
171 if(c.fhDcaXYvsDcaZnorm[j]) fhDcaXYvsDcaZnorm[j] = (TH2F*)c.fhDcaXYvsDcaZnorm[j]->Clone();
175 ((AliCFTrackIsPrimaryCuts &) c).Copy(*this);
179 //__________________________________________________________________________________
180 AliCFTrackIsPrimaryCuts::~AliCFTrackIsPrimaryCuts()
185 if (fhCutStatistics) delete fhCutStatistics;
186 if (fhCutCorrelation) delete fhCutCorrelation;
188 for (Int_t j=0; j<kNStepQA; j++){
189 for (Int_t i=0; i<kNHist; i++){
190 if(fhQA[i][j]) delete fhQA[i][j];
192 if(fhDcaXYvsDcaZ[j]) delete fhDcaXYvsDcaZ[j];
193 if(fhDcaXYvsDcaZnorm[j]) delete fhDcaXYvsDcaZnorm[j];
196 if (fBitmap) delete fBitmap;
198 if(fhBinLimNSigma) delete fhBinLimNSigma;
199 if(fhBinLimRequireSigma) delete fhBinLimRequireSigma;
200 if(fhBinLimAcceptKink) delete fhBinLimAcceptKink;
201 if(fhBinLimDcaXY) delete fhBinLimDcaXY;
202 if(fhBinLimDcaZ) delete fhBinLimDcaZ;
203 if(fhBinLimDcaXYnorm) delete fhBinLimDcaXYnorm;
204 if(fhBinLimDcaZnorm) delete fhBinLimDcaZnorm;
206 //__________________________________________________________________________________
207 void AliCFTrackIsPrimaryCuts::Initialise()
210 // sets everything to zero
212 fNSigmaToVertexMax = 0;
213 fRequireSigmaToVertex = 0;
214 fAcceptKinkDaughters = 0;
216 SetMaxNSigmaToVertex();
217 SetRequireSigmaToVertex();
219 for (Int_t j=0; j<kNStepQA; j++) {
220 for (Int_t i=0; i<kNHist; i++){
223 fhDcaXYvsDcaZ[j] = 0x0;
224 fhDcaXYvsDcaZnorm[j] = 0x0;
227 fhCutCorrelation = 0;
229 //set default bining for QA histograms
230 SetHistogramBins(kCutNSigmaToVertex,500,0,50);
231 SetHistogramBins(kCutRequireSigmaToVertex,2,-0.5,1.5);
232 SetHistogramBins(kCutAcceptKinkDaughters,2,-0.5,1.5);
233 SetHistogramBins(kDcaXY,500,-10,10);
234 SetHistogramBins(kDcaZ,500,-10,10);
235 SetHistogramBins(kDcaXYnorm,500,-10,10);
236 SetHistogramBins(kDcaZnorm,500,-10,10);
238 //__________________________________________________________________________________
239 void AliCFTrackIsPrimaryCuts::Copy(TObject &c) const
244 AliCFTrackIsPrimaryCuts& target = (AliCFTrackIsPrimaryCuts &) c;
248 if (fhCutStatistics) target.fhCutStatistics = (TH1F*) fhCutStatistics->Clone();
249 if (fhCutCorrelation) target.fhCutCorrelation = (TH2F*) fhCutCorrelation->Clone();
251 for (Int_t j=0; j<kNStepQA; j++){
252 for (Int_t i=0; i<kNHist; i++){
253 if(fhQA[i][j]) target.fhQA[i][j] = (TH1F*)fhQA[i][j]->Clone();
256 if(fhDcaXYvsDcaZ[j]) target.fhDcaXYvsDcaZ[j] = (TH2F*)fhDcaXYvsDcaZ[j]->Clone();
257 if(fhDcaXYvsDcaZnorm[j]) target.fhDcaXYvsDcaZnorm[j] = (TH2F*)fhDcaXYvsDcaZnorm[j]->Clone();
262 //____________________________________________________________________
263 Float_t AliCFTrackIsPrimaryCuts::GetSigmaToVertex(AliESDtrack* esdTrack) const
266 // Calculates the number of sigma to the vertex.
271 esdTrack->GetImpactParameters(b,bCov);
272 if (bCov[0]<=0 || bCov[2]<=0) {
273 AliDebug(1, "Estimated b resolution lower or equal zero!");
274 bCov[0]=0; bCov[2]=0;
276 bRes[0] = TMath::Sqrt(bCov[0]);
277 bRes[1] = TMath::Sqrt(bCov[2]);
279 // -----------------------------------
280 // How to get to a n-sigma cut?
282 // The accumulated statistics from 0 to d is
284 // -> Erf(d/Sqrt(2)) for a 1-dim gauss (d = n_sigma)
285 // -> 1 - Exp(-d**2) for a 2-dim gauss (d*d = dx*dx + dy*dy != n_sigma)
287 // It means that for a 2-dim gauss: n_sigma(d) = Sqrt(2)*ErfInv(1 - Exp((-x**2)/2)
288 // Can this be expressed in a different way?
290 if (bRes[0] == 0 || bRes[1] ==0)
293 Float_t d = TMath::Sqrt(TMath::Power(b[0]/bRes[0],2) + TMath::Power(b[1]/bRes[1],2));
295 // stupid rounding problem screws up everything:
296 // if d is too big, TMath::Exp(...) gets 0, and TMath::ErfInverse(1) that should be infinite, gets 0 :(
297 if (TMath::Exp(-d * d / 2) < 1e-10)
300 d = TMath::ErfInverse(1 - TMath::Exp(-d * d / 2)) * TMath::Sqrt(2);
303 //__________________________________________________________________________________
304 void AliCFTrackIsPrimaryCuts::GetBitMap(TObject* obj, TBits *bitmap) {
306 // retrieve the pointer to the bitmap
308 bitmap = SelectionBitMap(obj);
310 //__________________________________________________________________________________
311 TBits* AliCFTrackIsPrimaryCuts::SelectionBitMap(TObject* obj)
314 // test if the track passes the single cuts
315 // and store the information in a bitmap
318 // bitmap stores the decision of each single cut
319 for(Int_t i=0; i<kNCuts; i++)fBitmap->SetBitNumber(i,kFALSE);
321 // cast TObject into ESDtrack
322 AliESDtrack* esdTrack = dynamic_cast<AliESDtrack *>(obj);
323 if ( !esdTrack ) return fBitmap ;
326 for(Int_t i=0; i<kNCuts; i++)fBitmap->SetBitNumber(i,kTRUE);
328 // get the track to vertex parameter
329 Float_t nSigmaToVertex = GetSigmaToVertex(esdTrack);
333 if (nSigmaToVertex > fNSigmaToVertexMax)
334 fBitmap->SetBitNumber(iCutBit,kFALSE);
336 if (nSigmaToVertex<0 && fRequireSigmaToVertex)
337 fBitmap->SetBitNumber(iCutBit,kFALSE);
339 if (!fAcceptKinkDaughters && esdTrack->GetKinkIndex(0)>0)
340 fBitmap->SetBitNumber(iCutBit,kFALSE);
344 //__________________________________________________________________________________
345 Bool_t AliCFTrackIsPrimaryCuts::IsSelected(TObject* obj) {
347 // loops over decisions of single cuts and returns if the track is accepted
349 TBits* bitmap = SelectionBitMap(obj);
351 Bool_t isSelected = kTRUE;
353 for (UInt_t icut=0; icut<bitmap->GetNbits();icut++)
354 if(!bitmap->TestBitNumber(icut)) isSelected = kFALSE;
358 //__________________________________________________________________________________
359 void AliCFTrackIsPrimaryCuts::Init() {
361 // initialises all histograms and the TList which holds the histograms
366 //__________________________________________________________________________________
367 void AliCFTrackIsPrimaryCuts::SetHistogramBins(Int_t index, Int_t nbins, Double_t *bins)
375 case kCutNSigmaToVertex:
377 fhBinLimNSigma=new Double_t[nbins+1];
378 for(Int_t i=0;i<nbins+1;i++)fhBinLimNSigma[i]=bins[i];
381 case kCutRequireSigmaToVertex:
382 fhNBinsRequireSigma=nbins;
383 fhBinLimRequireSigma=new Double_t[nbins+1];
384 for(Int_t i=0;i<nbins+1;i++)fhBinLimRequireSigma[i]=bins[i];
387 case kCutAcceptKinkDaughters:
388 fhNBinsAcceptKink=nbins;
389 fhBinLimAcceptKink=new Double_t[nbins+1];
390 for(Int_t i=0;i<nbins+1;i++)fhBinLimAcceptKink[i]=bins[i];
395 fhBinLimDcaXY=new Double_t[nbins+1];
396 for(Int_t i=0;i<nbins+1;i++)fhBinLimDcaXY[i]=bins[i];
401 fhBinLimDcaZ=new Double_t[nbins+1];
402 for(Int_t i=0;i<nbins+1;i++)fhBinLimDcaZ[i]=bins[i];
406 fhNBinsDcaXYnorm=nbins;
407 fhBinLimDcaXYnorm=new Double_t[nbins+1];
408 for(Int_t i=0;i<nbins+1;i++)fhBinLimDcaXYnorm[i]=bins[i];
412 fhNBinsDcaZnorm=nbins;
413 fhBinLimDcaZnorm=new Double_t[nbins+1];
414 for(Int_t i=0;i<nbins+1;i++)fhBinLimDcaZnorm[i]=bins[i];
418 //__________________________________________________________________________________
419 void AliCFTrackIsPrimaryCuts::SetHistogramBins(Int_t index, Int_t nbins, Double_t xmin, Double_t xmax)
427 case kCutNSigmaToVertex:
429 fhBinLimNSigma=new Double_t[nbins+1];
430 for(Int_t i=0;i<nbins+1;i++)fhBinLimNSigma[i]=xmin+i*(xmax-xmin)/Double_t(nbins);
433 case kCutRequireSigmaToVertex:
434 fhNBinsRequireSigma=nbins;
435 fhBinLimRequireSigma=new Double_t[nbins+1];
436 for(Int_t i=0;i<nbins+1;i++)fhBinLimRequireSigma[i]=xmin+i*(xmax-xmin)/Double_t(nbins);
439 case kCutAcceptKinkDaughters:
440 fhNBinsAcceptKink=nbins;
441 fhBinLimAcceptKink=new Double_t[nbins+1];
442 for(Int_t i=0;i<nbins+1;i++)fhBinLimAcceptKink[i]=xmin+i*(xmax-xmin)/Double_t(nbins);
447 fhBinLimDcaXY=new Double_t[nbins+1];
448 for(Int_t i=0;i<nbins+1;i++)fhBinLimDcaXY[i]=xmin+i*(xmax-xmin)/Double_t(nbins);
453 fhBinLimDcaZ=new Double_t[nbins+1];
454 for(Int_t i=0;i<nbins+1;i++)fhBinLimDcaZ[i]=xmin+i*(xmax-xmin)/Double_t(nbins);
458 fhNBinsDcaXYnorm=nbins;
459 fhBinLimDcaXYnorm=new Double_t[nbins+1];
460 for(Int_t i=0;i<nbins+1;i++)fhBinLimDcaXYnorm[i]=xmin+i*(xmax-xmin)/Double_t(nbins);
464 fhNBinsDcaZnorm=nbins;
465 fhBinLimDcaZnorm=new Double_t[nbins+1];
466 for(Int_t i=0;i<nbins+1;i++)fhBinLimDcaZnorm[i]=xmin+i*(xmax-xmin)/Double_t(nbins);
470 //__________________________________________________________________________________
471 void AliCFTrackIsPrimaryCuts::DefineHistograms() {
473 // histograms for cut variables, cut statistics and cut correlations
478 // book cut statistics and cut correlation histograms
479 fhCutStatistics = new TH1F(Form("%s_cut_statistics",GetName()), Form("%s cut statistics",GetName()), kNCuts,0.5,kNCuts+0.5);
480 fhCutStatistics->SetLineWidth(2);
481 fhCutStatistics->GetXaxis()->SetBinLabel(1,"n dca");
482 fhCutStatistics->GetXaxis()->SetBinLabel(2,"require dca");
483 fhCutStatistics->GetXaxis()->SetBinLabel(3,"kink daughter");
485 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);
486 fhCutCorrelation->SetLineWidth(2);
487 fhCutCorrelation->GetXaxis()->SetBinLabel(1,"n dca");
488 fhCutCorrelation->GetXaxis()->SetBinLabel(2,"require dca");
489 fhCutCorrelation->GetXaxis()->SetBinLabel(3,"kink daughter");
491 fhCutCorrelation->GetYaxis()->SetBinLabel(1,"n dca");
492 fhCutCorrelation->GetYaxis()->SetBinLabel(2,"require dca");
493 fhCutCorrelation->GetYaxis()->SetBinLabel(3,"kink daughter");
495 // book QA histograms
497 for (Int_t i=0; i<kNStepQA; i++) {
498 if (i==0) sprintf(str," ");
499 else sprintf(str,"_cut");
501 fhDcaXYvsDcaZ[i] = new TH2F(Form("%s_dcaXYvsDcaZ%s",GetName(),str),"",200,-10,10,200,-10,10);
502 fhDcaXYvsDcaZnorm[i] = new TH2F(Form("%s_dcaXYvsDcaZnorm%s",GetName(),str),"",200,-10,10,200,-10,10);
504 fhQA[kCutNSigmaToVertex][i] = new TH1F(Form("%s_nSigmaToVertex%s",GetName(),str),"",fhNBinsNSigma,fhBinLimNSigma);
505 fhQA[kCutRequireSigmaToVertex][i] = new TH1F(Form("%s_requireSigmaToVertex%s",GetName(),str),"",fhNBinsRequireSigma,fhBinLimRequireSigma);
506 fhQA[kCutAcceptKinkDaughters][i] = new TH1F(Form("%s_acceptKinkDaughters%s",GetName(),str),"",fhNBinsAcceptKink,fhBinLimAcceptKink);
507 fhQA[kDcaXY][i] = new TH1F(Form("%s_dcaXY%s",GetName(),str),"",fhNBinsDcaXY,fhBinLimDcaXY);
508 fhQA[kDcaZ][i] = new TH1F(Form("%s_dcaZ%s",GetName(),str),"",fhNBinsDcaZ,fhBinLimDcaZ);
509 fhQA[kDcaXYnorm][i] = new TH1F(Form("%s_dcaXYnorm%s",GetName(),str),"",fhNBinsDcaXYnorm,fhBinLimDcaXYnorm);
510 fhQA[kDcaZnorm][i] = new TH1F(Form("%s_dcaZnorm%s",GetName(),str),"",fhNBinsDcaZnorm,fhBinLimDcaZnorm);
512 fhDcaXYvsDcaZ[i]->SetXTitle("impact par. d_{z}");
513 fhDcaXYvsDcaZ[i]->SetYTitle("impact par. d_{xy}");
514 fhDcaXYvsDcaZnorm[i]->SetXTitle("norm. impact par. d_{z} / #sigma_{z}");
515 fhDcaXYvsDcaZnorm[i]->SetYTitle("norm. impact par. d_{xy} / #sigma_{xy}");
517 fhQA[kCutNSigmaToVertex][i]->SetXTitle("n #sigma to vertex");
518 fhQA[kCutRequireSigmaToVertex][i]->SetXTitle("require #sigma to vertex");
519 fhQA[kCutAcceptKinkDaughters][i]->SetXTitle("accept kink daughters");
520 fhQA[kDcaXY][i]->SetXTitle("impact par. d_{xy}");
521 fhQA[kDcaZ][i]->SetXTitle("impact par. d_{z}");
522 fhQA[kDcaXYnorm][i]->SetXTitle("norm. impact par. d_{xy} / #sigma_{xy}");
523 fhQA[kDcaZnorm][i]->SetXTitle("norm. impact par. d_{z} / #sigma_{z}");
526 for(Int_t i=0; i<kNHist; i++) fhQA[i][1]->SetLineColor(color);
529 //__________________________________________________________________________________
530 void AliCFTrackIsPrimaryCuts::FillHistograms(TObject* obj, Bool_t f)
533 // fill the QA histograms
537 // cast TObject into ESDtrack
538 AliESDtrack* esdTrack = dynamic_cast<AliESDtrack *>(obj);
539 if ( !esdTrack ) return;
541 // index = 0: fill histograms before cuts
542 // index = 1: fill histograms after cuts
544 index = ((f) ? 1 : 0);
549 esdTrack->GetImpactParameters(b,bCov);
550 if (bCov[0]<=0 || bCov[2]<=0) {
551 AliDebug(1, "Estimated b resolution lower or equal zero!");
552 bCov[0]=0; bCov[2]=0;
554 bRes[0] = TMath::Sqrt(bCov[0]);
555 bRes[1] = TMath::Sqrt(bCov[2]);
557 fhQA[kDcaZ][index]->Fill(b[1]);
558 fhQA[kDcaXY][index]->Fill(b[0]);
559 fhDcaXYvsDcaZ[index]->Fill(b[1],b[0]);
561 if (bRes[0]!=0 && bRes[1]!=0) {
562 fhQA[kDcaZnorm][index]->Fill(b[1]/bRes[1]);
563 fhQA[kDcaXYnorm][index]->Fill(b[0]/bRes[0]);
564 fhDcaXYvsDcaZnorm[index]->Fill(b[1]/bRes[1], b[0]/bRes[0]);
567 // getting the track to vertex parameters
568 Float_t nSigmaToVertex = GetSigmaToVertex(esdTrack);
570 fhQA[kCutNSigmaToVertex][index]->Fill(nSigmaToVertex);
571 if (nSigmaToVertex<0 && fRequireSigmaToVertex) fhQA[kCutRequireSigmaToVertex][index]->Fill(0.);
572 if (!(nSigmaToVertex<0 && fRequireSigmaToVertex)) fhQA[kCutRequireSigmaToVertex][index]->Fill(1.);
575 if (!fAcceptKinkDaughters && esdTrack->GetKinkIndex(0)>0) fhQA[kCutAcceptKinkDaughters][index]->Fill(0.);
576 if (!(!fAcceptKinkDaughters && esdTrack->GetKinkIndex(0)>0)) fhQA[kCutAcceptKinkDaughters][index]->Fill(0.);
578 // fill cut statistics and cut correlation histograms with information from the bitmap
581 // Get the bitmap of the single cuts
583 TBits* bitmap = SelectionBitMap(obj);
585 // Number of single cuts in this class
586 UInt_t ncuts = bitmap->GetNbits();
587 for(UInt_t bit=0; bit<ncuts;bit++) {
588 if (!bitmap->TestBitNumber(bit)) {
589 fhCutStatistics->Fill(bit+1);
590 for (UInt_t bit2=bit; bit2<ncuts;bit2++) {
591 if (!bitmap->TestBitNumber(bit2))
592 fhCutCorrelation->Fill(bit+1,bit2+1);
597 //__________________________________________________________________________________
598 void AliCFTrackIsPrimaryCuts::SaveHistograms(const Char_t* dir) {
600 // saves the histograms in a directory (dir)
607 gDirectory->mkdir(dir);
610 gDirectory->mkdir("before_cuts");
611 gDirectory->mkdir("after_cuts");
613 fhCutStatistics->Write();
614 fhCutCorrelation->Write();
616 for (Int_t j=0; j<kNStepQA; j++) {
618 gDirectory->cd("before_cuts");
620 gDirectory->cd("after_cuts");
622 fhDcaXYvsDcaZ[j] ->Write();
623 fhDcaXYvsDcaZnorm[j]->Write();
625 for(Int_t i=0; i<kNHist; i++) fhQA[i][j]->Write();
627 gDirectory->cd("../");
629 gDirectory->cd("../");
631 //__________________________________________________________________________________
632 void AliCFTrackIsPrimaryCuts::DrawHistograms()
635 // draws some histograms
640 Float_t right = 0.03;
641 Float_t left = 0.175;
643 Float_t bottom = 0.175;
645 TCanvas* canvas1 = new TCanvas("Track_QA_Primary_1", "Track QA Primary 1", 800, 500);
646 canvas1->Divide(2, 1);
649 fhCutStatistics->SetStats(kFALSE);
650 fhCutStatistics->LabelsOption("v");
651 gPad->SetLeftMargin(left);
652 gPad->SetBottomMargin(0.25);
653 gPad->SetRightMargin(right);
654 gPad->SetTopMargin(0.1);
655 fhCutStatistics->Draw();
658 fhCutCorrelation->SetStats(kFALSE);
659 fhCutCorrelation->LabelsOption("v");
660 gPad->SetLeftMargin(0.30);
661 gPad->SetRightMargin(bottom);
662 gPad->SetTopMargin(0.1);
663 gPad->SetBottomMargin(0.25);
664 fhCutCorrelation->Draw("COLZ");
666 canvas1->SaveAs(Form("%s.eps", canvas1->GetName()));
667 canvas1->SaveAs(Form("%s.ps", canvas1->GetName()));
671 TCanvas* canvas2 = new TCanvas("Track_QA_Primary_2", "Track QA Primary 2", 800, 800);
672 canvas2->Divide(2, 2);
675 gPad->SetRightMargin(right);
676 gPad->SetLeftMargin(left);
677 gPad->SetTopMargin(top);
678 gPad->SetBottomMargin(bottom);
679 fhQA[kDcaXY][0]->SetStats(kFALSE);
680 fhQA[kDcaXY][0]->Draw();
681 fhQA[kDcaXY][1]->Draw("same");
684 gPad->SetRightMargin(right);
685 gPad->SetLeftMargin(left);
686 gPad->SetTopMargin(top);
687 gPad->SetBottomMargin(bottom);
688 fhQA[kDcaZ][0]->SetStats(kFALSE);
689 fhQA[kDcaZ][0]->Draw();
690 fhQA[kDcaZ][1]->Draw("same");
693 // fhDXYvsDZ[0]->SetStats(kFALSE);
695 gPad->SetLeftMargin(bottom);
696 gPad->SetTopMargin(0.1);
697 gPad->SetBottomMargin(bottom);
698 gPad->SetRightMargin(0.2);
699 fhDcaXYvsDcaZ[0]->Draw("COLZ");
702 // fhDXYvsDZ[1]->SetStats(kFALSE);
704 gPad->SetLeftMargin(bottom);
705 gPad->SetTopMargin(0.1);
706 gPad->SetBottomMargin(bottom);
707 gPad->SetRightMargin(0.2);
708 fhDcaXYvsDcaZ[1]->Draw("COLZ");
710 canvas2->SaveAs(Form("%s.eps", canvas2->GetName()));
711 canvas2->SaveAs(Form("%s.ps", canvas2->GetName()));
715 TCanvas* canvas3 = new TCanvas("Track_QA_Primary_3", "Track QA Primary 3", 800, 800);
716 canvas3->Divide(2, 2);
719 gPad->SetRightMargin(right);
720 gPad->SetLeftMargin(left);
721 gPad->SetTopMargin(top);
722 gPad->SetBottomMargin(bottom);
723 fhQA[kDcaXYnorm][0]->SetStats(kFALSE);
724 fhQA[kDcaXYnorm][0]->Draw();
725 fhQA[kDcaXYnorm][1]->Draw("same");
728 gPad->SetRightMargin(right);
729 gPad->SetLeftMargin(left);
730 gPad->SetTopMargin(top);
731 gPad->SetBottomMargin(bottom);
732 fhQA[kDcaZnorm][0]->SetStats(kFALSE);
733 fhQA[kDcaZnorm][0]->Draw();
734 fhQA[kDcaZnorm][1]->Draw("same");
737 // fhDXYvsDZ[0]->SetStats(kFALSE);
739 gPad->SetLeftMargin(bottom);
740 gPad->SetTopMargin(0.1);
741 gPad->SetBottomMargin(bottom);
742 gPad->SetRightMargin(0.2);
743 fhDcaXYvsDcaZnorm[0]->Draw("COLZ");
746 // fhDXYvsDZ[1]->SetStats(kFALSE);
748 gPad->SetLeftMargin(bottom);
749 gPad->SetTopMargin(0.1);
750 gPad->SetBottomMargin(bottom);
751 gPad->SetRightMargin(0.2);
752 fhDcaXYvsDcaZnorm[1]->Draw("COLZ");
754 canvas3->SaveAs(Form("%s.eps", canvas3->GetName()));
755 canvas3->SaveAs(Form("%s.ps", canvas3->GetName()));
759 TCanvas* canvas4 = new TCanvas("Track_QA_Primary_4", "Track QA Primary 4", 1200, 500);
760 canvas4->Divide(3, 1);
763 gPad->SetRightMargin(right);
764 gPad->SetLeftMargin(left);
765 gPad->SetTopMargin(top);
766 gPad->SetBottomMargin(bottom);
767 fhQA[kCutNSigmaToVertex][0]->SetStats(kFALSE);
768 fhQA[kCutNSigmaToVertex][0]->Draw();
769 fhQA[kCutNSigmaToVertex][1]->Draw("same");
772 gPad->SetRightMargin(right);
773 gPad->SetLeftMargin(left);
774 gPad->SetTopMargin(top);
775 gPad->SetBottomMargin(bottom);
776 fhQA[kCutRequireSigmaToVertex][0]->SetStats(kFALSE);
777 fhQA[kCutRequireSigmaToVertex][0]->Draw();
778 fhQA[kCutRequireSigmaToVertex][1]->Draw("same");
781 gPad->SetRightMargin(right);
782 gPad->SetLeftMargin(left);
783 gPad->SetTopMargin(top);
784 gPad->SetBottomMargin(bottom);
785 fhQA[kCutAcceptKinkDaughters][0]->SetStats(kFALSE);
786 fhQA[kCutAcceptKinkDaughters][0]->Draw();
787 fhQA[kCutAcceptKinkDaughters][1]->Draw("same");
789 canvas4->SaveAs(Form("%s.eps", canvas4->GetName()));
790 canvas4->SaveAs(Form("%s.ps", canvas4->GetName()));
792 //__________________________________________________________________________________
793 void AliCFTrackIsPrimaryCuts::AddQAHistograms(TList *qaList) const {
795 // saves the histograms in a TList
799 qaList->Add(fhCutStatistics);
800 qaList->Add(fhCutCorrelation);
802 for (Int_t j=0; j<kNStepQA; j++) {
803 qaList->Add(fhDcaXYvsDcaZ[j]);
804 qaList->Add(fhDcaXYvsDcaZnorm[j]);
805 for(Int_t i=0; i<kNHist; i++)
806 qaList->Add(fhQA[i][j]);