]> git.uio.no Git - u/mrichter/AliRoot.git/blob - CORRFW/AliCFTrackIsPrimaryCuts.cxx
added option to calculate dca relative to TPC vertex
[u/mrichter/AliRoot.git] / CORRFW / AliCFTrackIsPrimaryCuts.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
7  * Permission to use, copy, modify and distribute this software and its   *
8  * documentation strictly for non-commercial purposes is hereby granted   *
9  * without fee, provided that the above copyright notice appears in all   *
10  * copies and that both the copyright notice and this permission notice   *
11  * appear in the supporting documentation. The authors make no claims     *
12  * about the suitability of this software for any purpose. It is          *
13  * provided "as is" without express or implied warranty.                  *
14  **************************************************************************/
15
16 // The class 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
22 // ESD and AOD data.
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 // - min. and max. distance to main vertex in transverse plane (xy)
26 // - min. and max. longitudinal distance to main vertex (z)
27 // - min. and max. distance to main vertex as ellpise in xy - z plane
28 // - all above cuts on absolute values or in units of sigma (resolution)
29 // - min. and max. distance to main vertex in units of sigma (resolution)
30 // - max. transverse (xy) and longitudinal (z) impact parameter resolution
31 // - require that the dca calculation doesn't fail
32 // - accept or not accept daughter tracks of kink decays
33 //
34 // By default, the distance to 'vertex calculated from tracks' is used.
35 // Optionally the TPC (TPC only tracks based) vertex can be used.
36 //
37 // The cut values for these cuts are set with the corresponding set functions.
38 // All cut classes provided by the correction framework are supposed to be
39 // added in the Analysis Framwork's class AliAnalysisFilter and applied by
40 // the filter via a loop.
41 //
42 // author: I. Kraus (Ingrid.Kraus@cern.ch)
43 // idea taken form
44 // AliESDtrackCuts writte by Jan Fiete Grosse-Oetringhaus and
45 // AliRsnDaughterCut class written by A. Pulvirenti.
46
47 #include <TCanvas.h>
48 #include <TDirectory.h>
49 #include <TH2.h>
50 #include <TBits.h>
51
52 #include <AliESDtrack.h>
53 #include <AliLog.h>
54 #include <AliExternalTrackParam.h>
55 #include "AliCFTrackIsPrimaryCuts.h"
56
57 ClassImp(AliCFTrackIsPrimaryCuts)
58
59 //__________________________________________________________________________________
60 AliCFTrackIsPrimaryCuts::AliCFTrackIsPrimaryCuts() :
61   AliCFCutBase(),
62   fUseTPCvertex(0),
63   fMinDCAToVertexXY(0),
64   fMinDCAToVertexZ(0),
65   fMaxDCAToVertexXY(0),
66   fMaxDCAToVertexZ(0),
67   fDCAToVertex2D(0),
68   fAbsDCAToVertex(0),
69   fNSigmaToVertexMin(0),
70   fNSigmaToVertexMax(0),
71   fSigmaDCAxy(0),
72   fSigmaDCAz(0),
73   fRequireSigmaToVertex(0),
74   fAODType(AliAODTrack::kUndef),
75   fAcceptKinkDaughters(0),
76   fhCutStatistics(0),
77   fhCutCorrelation(0),
78   fBitmap(0x0),
79   fhNBinsNSigma(0),
80   fhNBinsRequireSigma(0),
81   fhNBinsAcceptKink(0),
82   fhNBinsDcaXY(0),
83   fhNBinsDcaZ(0),
84   fhNBinsDcaXYnorm(0),
85   fhNBinsDcaZnorm(0),
86   fhNBinsSigmaDcaXY(0),
87   fhNBinsSigmaDcaZ(0),
88   fhBinLimNSigma(0x0),
89   fhBinLimRequireSigma(0x0),
90   fhBinLimAcceptKink(0x0),
91   fhBinLimDcaXY(0x0),
92   fhBinLimDcaZ(0x0),
93   fhBinLimDcaXYnorm(0x0),
94   fhBinLimDcaZnorm(0x0),
95   fhBinLimSigmaDcaXY(0x0),
96   fhBinLimSigmaDcaZ(0x0)
97 {
98   //
99   // Default constructor
100   //
101   Initialise();
102 }
103 //__________________________________________________________________________________
104 AliCFTrackIsPrimaryCuts::AliCFTrackIsPrimaryCuts(Char_t* name, Char_t* title) :
105   AliCFCutBase(name,title),
106   fUseTPCvertex(0),
107   fMinDCAToVertexXY(0),
108   fMinDCAToVertexZ(0),
109   fMaxDCAToVertexXY(0),
110   fMaxDCAToVertexZ(0),
111   fDCAToVertex2D(0),
112   fAbsDCAToVertex(0),
113   fNSigmaToVertexMin(0),
114   fNSigmaToVertexMax(0),
115   fSigmaDCAxy(0),
116   fSigmaDCAz(0),
117   fRequireSigmaToVertex(0),
118   fAODType(AliAODTrack::kUndef),
119   fAcceptKinkDaughters(0),
120   fhCutStatistics(0),
121   fhCutCorrelation(0),
122   fBitmap(0x0),
123   fhNBinsNSigma(0),
124   fhNBinsRequireSigma(0),
125   fhNBinsAcceptKink(0),
126   fhNBinsDcaXY(0),
127   fhNBinsDcaZ(0),
128   fhNBinsDcaXYnorm(0),
129   fhNBinsDcaZnorm(0),
130   fhNBinsSigmaDcaXY(0),
131   fhNBinsSigmaDcaZ(0),
132   fhBinLimNSigma(0x0),
133   fhBinLimRequireSigma(0x0),
134   fhBinLimAcceptKink(0x0),
135   fhBinLimDcaXY(0x0),
136   fhBinLimDcaZ(0x0),
137   fhBinLimDcaXYnorm(0x0),
138   fhBinLimDcaZnorm(0x0),
139   fhBinLimSigmaDcaXY(0x0),
140   fhBinLimSigmaDcaZ(0x0)
141 {
142   //
143   // Constructor
144   //
145   Initialise();
146 }
147 //__________________________________________________________________________________
148 AliCFTrackIsPrimaryCuts::AliCFTrackIsPrimaryCuts(const AliCFTrackIsPrimaryCuts& c) :
149   AliCFCutBase(c),
150   fUseTPCvertex(c.fUseTPCvertex),
151   fMinDCAToVertexXY(c.fMinDCAToVertexXY),
152   fMinDCAToVertexZ(c.fMinDCAToVertexZ),
153   fMaxDCAToVertexXY(c.fMaxDCAToVertexXY),
154   fMaxDCAToVertexZ(c.fMaxDCAToVertexZ),
155   fDCAToVertex2D(c.fDCAToVertex2D),
156   fAbsDCAToVertex(c.fAbsDCAToVertex),
157   fNSigmaToVertexMin(c.fNSigmaToVertexMin),
158   fNSigmaToVertexMax(c.fNSigmaToVertexMax),
159   fSigmaDCAxy(c.fSigmaDCAxy),
160   fSigmaDCAz(c.fSigmaDCAz),
161   fRequireSigmaToVertex(c.fRequireSigmaToVertex),
162   fAODType(c.fAODType),
163   fAcceptKinkDaughters(c.fAcceptKinkDaughters),
164   fhCutStatistics(c.fhCutStatistics),
165   fhCutCorrelation(c.fhCutCorrelation),
166   fBitmap(c.fBitmap),
167   fhNBinsNSigma(c.fhNBinsNSigma),
168   fhNBinsRequireSigma(c.fhNBinsRequireSigma),
169   fhNBinsAcceptKink(c.fhNBinsAcceptKink),
170   fhNBinsDcaXY(c.fhNBinsDcaXY),
171   fhNBinsDcaZ(c.fhNBinsDcaZ),
172   fhNBinsDcaXYnorm(c.fhNBinsDcaXYnorm),
173   fhNBinsDcaZnorm(c.fhNBinsDcaZnorm),
174   fhNBinsSigmaDcaXY(c.fhNBinsSigmaDcaXY),
175   fhNBinsSigmaDcaZ(c.fhNBinsSigmaDcaZ),
176   fhBinLimNSigma(c.fhBinLimNSigma),
177   fhBinLimRequireSigma(c.fhBinLimRequireSigma),
178   fhBinLimAcceptKink(c.fhBinLimAcceptKink),
179   fhBinLimDcaXY(c.fhBinLimDcaXY),
180   fhBinLimDcaZ(c.fhBinLimDcaZ),
181   fhBinLimDcaXYnorm(c.fhBinLimDcaXYnorm),
182   fhBinLimDcaZnorm(c.fhBinLimDcaZnorm),
183   fhBinLimSigmaDcaXY(c.fhBinLimSigmaDcaXY),
184   fhBinLimSigmaDcaZ(c.fhBinLimSigmaDcaZ)
185 {
186   //
187   // copy constructor
188   //
189   ((AliCFTrackIsPrimaryCuts &) c).Copy(*this);
190 }
191 //__________________________________________________________________________________
192 AliCFTrackIsPrimaryCuts& AliCFTrackIsPrimaryCuts::operator=(const AliCFTrackIsPrimaryCuts& c)
193 {
194   //
195   // Assignment operator
196   //
197   if (this != &c) {
198     AliCFCutBase::operator=(c) ;
199     fUseTPCvertex = c.fUseTPCvertex;
200     fMinDCAToVertexXY = c.fMinDCAToVertexXY;
201     fMinDCAToVertexZ = c.fMinDCAToVertexZ;
202     fMaxDCAToVertexXY = c.fMaxDCAToVertexXY;
203     fMaxDCAToVertexZ = c.fMaxDCAToVertexZ;
204     fDCAToVertex2D = c.fDCAToVertex2D;
205     fAbsDCAToVertex = c.fAbsDCAToVertex;
206     fNSigmaToVertexMin = c.fNSigmaToVertexMin ;
207     fNSigmaToVertexMax = c.fNSigmaToVertexMax ;
208     fSigmaDCAxy = c.fSigmaDCAxy ;
209     fSigmaDCAz = c.fSigmaDCAz ;
210     fRequireSigmaToVertex = c.fRequireSigmaToVertex ;
211     fAODType = c.fAODType ;
212     fAcceptKinkDaughters = c.fAcceptKinkDaughters ;
213     fhCutStatistics = c.fhCutStatistics ;
214     fhCutCorrelation = c.fhCutCorrelation ;
215     fBitmap =  c.fBitmap;
216     fhNBinsNSigma = c.fhNBinsNSigma;
217     fhNBinsRequireSigma = c.fhNBinsRequireSigma;
218     fhNBinsAcceptKink = c.fhNBinsAcceptKink;
219     fhNBinsDcaXY = c.fhNBinsDcaXY;
220     fhNBinsDcaZ = c.fhNBinsDcaZ;
221     fhNBinsDcaXYnorm = c.fhNBinsDcaXYnorm;
222     fhNBinsDcaZnorm = c.fhNBinsDcaZnorm;
223     fhNBinsSigmaDcaXY = c.fhNBinsSigmaDcaXY;
224     fhNBinsSigmaDcaZ = c.fhNBinsSigmaDcaZ;
225     fhBinLimNSigma = c.fhBinLimNSigma;
226     fhBinLimRequireSigma = c.fhBinLimRequireSigma;
227     fhBinLimAcceptKink = c.fhBinLimAcceptKink;
228     fhBinLimDcaXY = c.fhBinLimDcaXY;
229     fhBinLimDcaZ = c.fhBinLimDcaZ;
230     fhBinLimDcaXYnorm = c.fhBinLimDcaXYnorm;
231     fhBinLimDcaZnorm = c.fhBinLimDcaZnorm;
232     fhBinLimSigmaDcaXY = c.fhBinLimSigmaDcaXY;
233     fhBinLimSigmaDcaZ = c.fhBinLimSigmaDcaZ;
234
235     for (Int_t j=0; j<6; j++) fDCA[j] = c.fDCA[j];
236     for (Int_t j=0; j<c.kNStepQA; j++){
237       if(c.fhDcaXYvsDcaZ[j]) fhDcaXYvsDcaZ[j] = (TH2F*)c.fhDcaXYvsDcaZ[j]->Clone();
238       for (Int_t i=0; i<c.kNHist; i++){
239         if(c.fhQA[i][j]) fhQA[i][j] = (TH1F*)c.fhQA[i][j]->Clone();
240       }
241     }
242     ((AliCFTrackIsPrimaryCuts &) c).Copy(*this);
243  }
244   return *this;
245 }
246 //__________________________________________________________________________________
247 AliCFTrackIsPrimaryCuts::~AliCFTrackIsPrimaryCuts()
248 {
249   //
250   // destructor
251   //
252   if (fhCutStatistics)                  delete fhCutStatistics;
253   if (fhCutCorrelation)                 delete fhCutCorrelation;
254
255   for (Int_t j=0; j<kNStepQA; j++){
256     if(fhDcaXYvsDcaZ[j])        delete fhDcaXYvsDcaZ[j];
257     for (Int_t i=0; i<kNHist; i++)
258       if(fhQA[i][j])            delete fhQA[i][j];
259   }
260   if(fBitmap)                   delete fBitmap;
261   if(fhBinLimNSigma)            delete fhBinLimNSigma;
262   if(fhBinLimRequireSigma)      delete fhBinLimRequireSigma;
263   if(fhBinLimAcceptKink)        delete fhBinLimAcceptKink;
264   if(fhBinLimDcaXY)             delete fhBinLimDcaXY;
265   if(fhBinLimDcaZ)              delete fhBinLimDcaZ;
266   if(fhBinLimDcaXYnorm)         delete fhBinLimDcaXYnorm;
267   if(fhBinLimDcaZnorm)          delete fhBinLimDcaZnorm;
268   if(fhBinLimSigmaDcaXY)        delete fhBinLimSigmaDcaXY;
269   if(fhBinLimSigmaDcaZ)         delete fhBinLimSigmaDcaZ;
270 }
271 //__________________________________________________________________________________
272 void AliCFTrackIsPrimaryCuts::Initialise()
273 {
274   //
275   // sets everything to zero
276   //
277   fUseTPCvertex = 0;
278   fMinDCAToVertexXY = 0;
279   fMinDCAToVertexZ = 0;
280   fMaxDCAToVertexXY = 0;
281   fMaxDCAToVertexZ = 0;
282   fDCAToVertex2D = 0;
283   fAbsDCAToVertex = 0;
284   fNSigmaToVertexMin = 0;
285   fNSigmaToVertexMax = 0;
286   fSigmaDCAxy = 0;
287   fSigmaDCAz = 0;
288   fRequireSigmaToVertex = 0;
289   fAcceptKinkDaughters = 0;
290   fAODType = AliAODTrack::kUndef;
291
292   SetMinDCAToVertexXY();
293   SetMinDCAToVertexZ();
294   SetMaxDCAToVertexXY();
295   SetMaxDCAToVertexZ();
296   SetDCAToVertex2D();
297   SetAbsDCAToVertex();
298   SetMinNSigmaToVertex();
299   SetMaxNSigmaToVertex();
300   SetMaxSigmaDCAxy();
301   SetMaxSigmaDCAz();
302   SetRequireSigmaToVertex();
303   SetAcceptKinkDaughters();
304   SetAODType();
305
306   for (Int_t j=0; j<6; j++) fDCA[j] = 0.;
307   for (Int_t j=0; j<kNStepQA; j++)  {
308     fhDcaXYvsDcaZ[j] = 0x0;
309     for (Int_t i=0; i<kNHist; i++)
310       fhQA[i][j] = 0x0;
311   }
312   fhCutStatistics = 0;
313   fhCutCorrelation = 0;
314   fBitmap=new TBits(0);
315
316   //set default bining for QA histograms
317   SetHistogramBins(kCutNSigmaToVertex,100,0.,10.);
318   SetHistogramBins(kCutRequireSigmaToVertex,5,-0.75,1.75);
319   SetHistogramBins(kCutAcceptKinkDaughters,5,-0.75,1.75);
320   SetHistogramBins(kDcaXY,500,-10.,10.);
321   SetHistogramBins(kDcaZ,500,-10.,10.);
322   SetHistogramBins(kDcaXYnorm,500,-10.,10.);
323   SetHistogramBins(kDcaZnorm,500,-10.,10.);
324   SetHistogramBins(kSigmaDcaXY,500,-0.1,0.9);
325   SetHistogramBins(kSigmaDcaZ,500,-0.1,0.9);
326 }
327 //__________________________________________________________________________________
328 void AliCFTrackIsPrimaryCuts::Copy(TObject &c) const
329 {
330   //
331   // Copy function
332   //
333   AliCFTrackIsPrimaryCuts& target = (AliCFTrackIsPrimaryCuts &) c;
334
335   target.Initialise();
336
337   if (fhCutStatistics)  target.fhCutStatistics = (TH1F*) fhCutStatistics->Clone();
338   if (fhCutCorrelation) target.fhCutCorrelation = (TH2F*) fhCutCorrelation->Clone();
339
340   for (Int_t j=0; j<6; j++) target.fDCA[j] = fDCA[j];
341   for (Int_t j=0; j<kNStepQA; j++){
342     if(fhDcaXYvsDcaZ[j]) target.fhDcaXYvsDcaZ[j] = (TH2F*)fhDcaXYvsDcaZ[j]->Clone();
343     for (Int_t i=0; i<kNHist; i++)
344       if(fhQA[i][j]) target.fhQA[i][j] = (TH1F*)fhQA[i][j]->Clone();
345   }
346   TNamed::Copy(c);
347 }
348 //__________________________________________________________________________________
349 void AliCFTrackIsPrimaryCuts::GetDCA(AliESDtrack* esdTrack)
350 {
351   if (!esdTrack) return;
352
353   Float_t b[2] = {0.,0.};
354   Float_t bCov[3] = {0.,0.,0.};
355   if(!fUseTPCvertex) esdTrack->GetImpactParameters(b,bCov);
356   if( fUseTPCvertex) esdTrack->GetImpactParametersTPC(b,bCov);
357
358   if (bCov[0]<=0 || bCov[2]<=0) {
359       bCov[0]=0; bCov[2]=0;
360   }
361   fDCA[0] = b[0]; // impact parameter xy
362   fDCA[1] = b[1]; // impact parameter z
363   fDCA[2] = TMath::Sqrt(bCov[0]); // resolution xy
364   fDCA[3] = TMath::Sqrt(bCov[2]); // resolution z
365
366   if (!fAbsDCAToVertex) {
367         if (fDCA[2] > 0) fDCA[0] = fDCA[0]/fDCA[2]; // normalised impact parameter xy
368         if (fDCA[3] > 0) fDCA[1] = fDCA[1]/fDCA[3]; // normalised impact parameter z
369   }
370
371   // get n_sigma
372   if(!fUseTPCvertex)
373         fDCA[5] = AliESDtrackCuts::GetSigmaToVertex(esdTrack);
374
375   if(fUseTPCvertex) {
376         fDCA[5] = -1;
377         if (fDCA[2]==0 || fDCA[3]==0)
378           return;
379         fDCA[5] = 1000.;
380         Float_t d = TMath::Sqrt(TMath::Power(b[0]/fDCA[2],2) + TMath::Power(b[1]/fDCA[3],2));
381         if (TMath::Exp(-d * d / 2) < 1e-15)
382           return;
383         fDCA[5] = TMath::ErfInverse(1 - TMath::Exp(-d * d / 2)) * TMath::Sqrt(2);
384   }
385   return;
386 }
387 //__________________________________________________________________________________
388 void AliCFTrackIsPrimaryCuts::SelectionBitMap(TObject* obj)
389 {
390   //
391   // test if the track passes the single cuts
392   // and store the information in a bitmap
393   //
394
395   // bitmap stores the decision of each single cut
396   for(Int_t i=0; i<kNCuts; i++)fBitmap->SetBitNumber(i,kFALSE);
397
398   // check TObject and cast into ESDtrack
399   if (!obj) return;
400   if (!obj->InheritsFrom("AliVParticle")) {
401     AliError("object must derived from AliVParticle !");
402     return;
403   }
404
405   Bool_t isESDTrack = strcmp(obj->ClassName(),"AliESDtrack") == 0 ? kTRUE : kFALSE ;
406   Bool_t isAODTrack = strcmp(obj->ClassName(),"AliAODTrack") == 0 ? kTRUE : kFALSE ;
407
408   AliESDtrack * esdTrack = 0x0 ;
409   AliAODTrack * aodTrack = 0x0 ; 
410   if (isESDTrack) esdTrack = dynamic_cast<AliESDtrack*>(obj);
411   if (isAODTrack) aodTrack = dynamic_cast<AliAODTrack*>(obj);
412
413   // get the track to vertex parameter for ESD track
414   if (isESDTrack) GetDCA(esdTrack);
415
416   Float_t bxy = 0, bz = 0;
417   if (isESDTrack) {
418         bxy = TMath::Abs(fDCA[0]);
419         bz  = TMath::Abs(fDCA[1]);
420  }
421   Float_t b2Dmin = 0, b2Dmax = 0;
422   if (isESDTrack) {
423     if (fMinDCAToVertexXY>0 && fMinDCAToVertexZ>0)
424         b2Dmin = fDCA[0]*fDCA[0]/fMinDCAToVertexXY/fMinDCAToVertexXY + fDCA[1]*fDCA[1]/fMinDCAToVertexZ/fMinDCAToVertexZ;
425     if (fMaxDCAToVertexXY>0 && fMaxDCAToVertexZ>0) 
426         b2Dmax = fDCA[0]*fDCA[0]/fMaxDCAToVertexXY/fMaxDCAToVertexXY + fDCA[1]*fDCA[1]/fMaxDCAToVertexZ/fMaxDCAToVertexZ;
427   }
428
429   // fill the bitmap
430   Int_t iCutBit = 0;
431
432   if (isESDTrack) {
433     if (fDCAToVertex2D || (!fDCAToVertex2D && bxy >= fMinDCAToVertexXY && bxy <= fMaxDCAToVertexXY)) {
434       fBitmap->SetBitNumber(iCutBit,kTRUE);
435     }
436   }
437   else fBitmap->SetBitNumber(iCutBit,kTRUE);
438
439   iCutBit++;
440
441   if (isESDTrack) {
442     if (fDCAToVertex2D || (!fDCAToVertex2D && bz  >= fMinDCAToVertexZ && bz  <= fMaxDCAToVertexZ)) {
443       fBitmap->SetBitNumber(iCutBit,kTRUE);
444     }
445   }
446   else fBitmap->SetBitNumber(iCutBit,kTRUE);
447
448   iCutBit++;
449
450   if (isESDTrack) {
451     if (!fDCAToVertex2D || (fDCAToVertex2D && b2Dmin > 1  && b2Dmax < 1)) {
452       fBitmap->SetBitNumber(iCutBit,kTRUE);
453     }
454   }
455   else fBitmap->SetBitNumber(iCutBit,kTRUE);
456
457   iCutBit++;
458
459   if (isESDTrack) {
460     if (fDCA[5] >= fNSigmaToVertexMin && fDCA[5] <= fNSigmaToVertexMax) {
461       fBitmap->SetBitNumber(iCutBit,kTRUE);
462     }
463   }
464   else fBitmap->SetBitNumber(iCutBit,kTRUE);
465
466   iCutBit++;
467
468   if (isESDTrack) {
469     if (fDCA[2] < fSigmaDCAxy) {
470       fBitmap->SetBitNumber(iCutBit,kTRUE);
471     }
472   }
473   else fBitmap->SetBitNumber(iCutBit,kTRUE);
474
475   iCutBit++;
476
477   if (isESDTrack) {
478     if (fDCA[3] < fSigmaDCAz) {
479       fBitmap->SetBitNumber(iCutBit,kTRUE);
480     }
481   }
482   else fBitmap->SetBitNumber(iCutBit,kTRUE);
483
484   iCutBit++;
485
486   if (isESDTrack) {
487     if (!fRequireSigmaToVertex || (fDCA[5]>=0 && fRequireSigmaToVertex)) {
488       fBitmap->SetBitNumber(iCutBit,kTRUE);
489     }
490   }
491   else fBitmap->SetBitNumber(iCutBit,kTRUE);
492
493   iCutBit++;
494
495   if (esdTrack) {
496     if (fAcceptKinkDaughters || (!fAcceptKinkDaughters && esdTrack->GetKinkIndex(0)<=0)) {
497       fBitmap->SetBitNumber(iCutBit,kTRUE);
498     }
499   }
500   else fBitmap->SetBitNumber(iCutBit,kTRUE);
501   
502   iCutBit++;
503   
504   if (isAODTrack) {
505     if (fAODType==AliAODTrack::kUndef || fAODType == aodTrack->GetType()) {
506       fBitmap->SetBitNumber(iCutBit,kTRUE);
507     }
508   }
509   else fBitmap->SetBitNumber(iCutBit,kTRUE);
510   
511   return;
512 }
513 //__________________________________________________________________________________
514 Bool_t AliCFTrackIsPrimaryCuts::IsSelected(TObject* obj) {
515   //
516   // loops over decisions of single cuts and returns if the track is accepted
517   //
518   SelectionBitMap(obj);
519
520   if (fIsQAOn) FillHistograms(obj,0);
521   Bool_t isSelected = kTRUE;
522
523   for (UInt_t icut=0; icut<fBitmap->GetNbits();icut++) {
524     if(!fBitmap->TestBitNumber(icut)) {
525         isSelected = kFALSE;
526         break;
527     }
528   }
529   if (!isSelected) return kFALSE ;
530   if (fIsQAOn) FillHistograms(obj,1);
531   return kTRUE;
532 }
533 //__________________________________________________________________________________
534 void AliCFTrackIsPrimaryCuts::SetHistogramBins(Int_t index, Int_t nbins, Double_t *bins)
535 {
536   //
537   // variable bin size
538   //
539   if(!fIsQAOn) return;
540
541   switch(index){
542   case kCutNSigmaToVertex:
543     fhNBinsNSigma=nbins+1;
544     fhBinLimNSigma=new Double_t[nbins+1];
545     for(Int_t i=0;i<nbins+1;i++)fhBinLimNSigma[i]=bins[i];
546     break;
547
548   case kCutRequireSigmaToVertex:
549     fhNBinsRequireSigma=nbins+1;
550     fhBinLimRequireSigma=new Double_t[nbins+1];
551     for(Int_t i=0;i<nbins+1;i++)fhBinLimRequireSigma[i]=bins[i];
552     break;
553
554   case kCutAcceptKinkDaughters:
555     fhNBinsAcceptKink=nbins+1;
556     fhBinLimAcceptKink=new Double_t[nbins+1];
557     for(Int_t i=0;i<nbins+1;i++)fhBinLimAcceptKink[i]=bins[i];
558     break;
559
560   case kDcaXY:
561     fhNBinsDcaXY=nbins+1;
562     fhBinLimDcaXY=new Double_t[nbins+1];
563     for(Int_t i=0;i<nbins+1;i++)fhBinLimDcaXY[i]=bins[i];
564     break;
565
566   case kDcaZ:
567     fhNBinsDcaZ=nbins+1;
568     fhBinLimDcaZ=new Double_t[nbins+1];
569     for(Int_t i=0;i<nbins+1;i++)fhBinLimDcaZ[i]=bins[i];
570     break;
571
572   case kDcaXYnorm:
573     fhNBinsDcaXYnorm=nbins+1;
574     fhBinLimDcaXYnorm=new Double_t[nbins+1];
575     for(Int_t i=0;i<nbins+1;i++)fhBinLimDcaXYnorm[i]=bins[i];
576     break;
577
578   case kDcaZnorm:
579     fhNBinsDcaZnorm=nbins+1;
580     fhBinLimDcaZnorm=new Double_t[nbins+1];
581     for(Int_t i=0;i<nbins+1;i++)fhBinLimDcaZnorm[i]=bins[i];
582     break;
583
584   case kSigmaDcaXY:
585     fhNBinsSigmaDcaXY=nbins+1;
586     fhBinLimSigmaDcaXY=new Double_t[nbins+1];
587     for(Int_t i=0;i<nbins+1;i++)fhBinLimSigmaDcaXY[i]=bins[i];
588     break;
589
590   case kSigmaDcaZ:
591     fhNBinsSigmaDcaZ=nbins+1;
592     fhBinLimSigmaDcaZ=new Double_t[nbins+1];
593     for(Int_t i=0;i<nbins+1;i++)fhBinLimSigmaDcaZ[i]=bins[i];
594     break;
595   }
596 }
597 //__________________________________________________________________________________
598 void AliCFTrackIsPrimaryCuts::SetHistogramBins(Int_t index, Int_t nbins, Double_t xmin, Double_t xmax)
599 {
600   //
601   // fixed bin size
602   //
603   switch(index){
604   case kCutNSigmaToVertex:
605     fhNBinsNSigma=nbins+1;
606     fhBinLimNSigma=new Double_t[nbins+1];
607     for(Int_t i=0;i<nbins+1;i++)fhBinLimNSigma[i]=xmin+i*(xmax-xmin)/Double_t(nbins);
608     break;
609
610   case kCutRequireSigmaToVertex:
611     fhNBinsRequireSigma=nbins+1;
612     fhBinLimRequireSigma=new Double_t[nbins+1];
613     for(Int_t i=0;i<nbins+1;i++)fhBinLimRequireSigma[i]=xmin+i*(xmax-xmin)/Double_t(nbins);
614     break;
615
616   case kCutAcceptKinkDaughters:
617     fhNBinsAcceptKink=nbins+1;
618     fhBinLimAcceptKink=new Double_t[nbins+1];
619     for(Int_t i=0;i<nbins+1;i++)fhBinLimAcceptKink[i]=xmin+i*(xmax-xmin)/Double_t(nbins);
620     break;
621
622   case kDcaXY:
623     fhNBinsDcaXY=nbins+1;
624     fhBinLimDcaXY=new Double_t[nbins+1];
625     for(Int_t i=0;i<nbins+1;i++)fhBinLimDcaXY[i]=xmin+i*(xmax-xmin)/Double_t(nbins);
626     break;
627
628   case kDcaZ:
629     fhNBinsDcaZ=nbins+1;
630     fhBinLimDcaZ=new Double_t[nbins+1];
631     for(Int_t i=0;i<nbins+1;i++)fhBinLimDcaZ[i]=xmin+i*(xmax-xmin)/Double_t(nbins);
632     break;
633
634   case kDcaXYnorm:
635     fhNBinsDcaXYnorm=nbins+1;
636     fhBinLimDcaXYnorm=new Double_t[nbins+1];
637     for(Int_t i=0;i<nbins+1;i++)fhBinLimDcaXYnorm[i]=xmin+i*(xmax-xmin)/Double_t(nbins);
638     break;
639
640   case kDcaZnorm:
641     fhNBinsDcaZnorm=nbins+1;
642     fhBinLimDcaZnorm=new Double_t[nbins+1];
643     for(Int_t i=0;i<nbins+1;i++)fhBinLimDcaZnorm[i]=xmin+i*(xmax-xmin)/Double_t(nbins);
644     break;
645
646   case kSigmaDcaXY:
647     fhNBinsSigmaDcaXY=nbins+1;
648     fhBinLimSigmaDcaXY=new Double_t[nbins+1];
649     for(Int_t i=0;i<nbins+1;i++)fhBinLimSigmaDcaXY[i]=xmin+i*(xmax-xmin)/Double_t(nbins);
650     break;
651
652   case kSigmaDcaZ:
653     fhNBinsSigmaDcaZ=nbins+1;
654     fhBinLimSigmaDcaZ=new Double_t[nbins+1];
655     for(Int_t i=0;i<nbins+1;i++)fhBinLimSigmaDcaZ[i]=xmin+i*(xmax-xmin)/Double_t(nbins);
656     break;
657   }
658 }
659 //__________________________________________________________________________________
660  void AliCFTrackIsPrimaryCuts::DefineHistograms() {
661   //
662   // histograms for cut variables, cut statistics and cut correlations
663   //
664   Int_t color = 2;
665
666   // book cut statistics and cut correlation histograms
667   fhCutStatistics = new TH1F(Form("%s_cut_statistics",GetName()), Form("%s cut statistics",GetName()), kNCuts,0.5,kNCuts+0.5);
668   fhCutStatistics->SetLineWidth(2);
669   fhCutStatistics->GetXaxis()->SetBinLabel(1,"dca xy");
670   fhCutStatistics->GetXaxis()->SetBinLabel(2,"dca z");
671   fhCutStatistics->GetXaxis()->SetBinLabel(3,"dca ellipse");
672   fhCutStatistics->GetXaxis()->SetBinLabel(4,"n dca");
673   fhCutStatistics->GetXaxis()->SetBinLabel(5,"sigma dca xy");
674   fhCutStatistics->GetXaxis()->SetBinLabel(6,"sigma dca z");
675   fhCutStatistics->GetXaxis()->SetBinLabel(7,"require dca");
676   fhCutStatistics->GetXaxis()->SetBinLabel(8,"kink daughter");
677   fhCutStatistics->GetXaxis()->SetBinLabel(9,"AOD type");
678
679   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);
680   fhCutCorrelation->SetLineWidth(2);
681   fhCutCorrelation->GetXaxis()->SetBinLabel(1,"dca xy");
682   fhCutCorrelation->GetXaxis()->SetBinLabel(2,"dca z");
683   fhCutCorrelation->GetXaxis()->SetBinLabel(3,"dca ellipse");
684   fhCutCorrelation->GetXaxis()->SetBinLabel(4,"n dca");
685   fhCutCorrelation->GetXaxis()->SetBinLabel(5,"sigma dca xy");
686   fhCutCorrelation->GetXaxis()->SetBinLabel(6,"sigma dca z");
687   fhCutCorrelation->GetXaxis()->SetBinLabel(7,"require dca");
688   fhCutCorrelation->GetXaxis()->SetBinLabel(8,"kink daughter");
689   fhCutCorrelation->GetXaxis()->SetBinLabel(9,"AOD type");
690
691   fhCutCorrelation->GetYaxis()->SetBinLabel(1,"dca xy");
692   fhCutCorrelation->GetYaxis()->SetBinLabel(2,"dca z");
693   fhCutCorrelation->GetYaxis()->SetBinLabel(3,"dca ellipse");
694   fhCutCorrelation->GetYaxis()->SetBinLabel(4,"n dca");
695   fhCutCorrelation->GetYaxis()->SetBinLabel(5,"sigma dca xy");
696   fhCutCorrelation->GetYaxis()->SetBinLabel(6,"sigma dca z");
697   fhCutCorrelation->GetYaxis()->SetBinLabel(7,"require dca");
698   fhCutCorrelation->GetYaxis()->SetBinLabel(8,"kink daughter");
699   fhCutCorrelation->GetYaxis()->SetBinLabel(9,"AOD type");
700
701   // book QA histograms
702   Char_t str[256];
703   for (Int_t i=0; i<kNStepQA; i++) {
704     if (i==0) sprintf(str," ");
705     else sprintf(str,"_cut");
706
707     fhDcaXYvsDcaZ[i]            = new  TH2F(Form("%s_dcaXYvsDcaZ%s",GetName(),str),"",200,-10,10,200,-10,10);
708     fhQA[kCutNSigmaToVertex][i] = new TH1F(Form("%s_nSigmaToVertex%s",GetName(),str),"",fhNBinsNSigma-1,fhBinLimNSigma);
709     fhQA[kCutRequireSigmaToVertex][i] = new TH1F(Form("%s_requireSigmaToVertex%s",GetName(),str),"",fhNBinsRequireSigma-1,fhBinLimRequireSigma);
710     fhQA[kCutAcceptKinkDaughters][i] = new TH1F(Form("%s_acceptKinkDaughters%s",GetName(),str),"",fhNBinsAcceptKink-1,fhBinLimAcceptKink);
711     fhQA[kDcaXY][i]             = new TH1F(Form("%s_dcaXY%s",GetName(),str),"",fhNBinsDcaXY-1,fhBinLimDcaXY);
712     fhQA[kDcaZ][i]              = new TH1F(Form("%s_dcaZ%s",GetName(),str),"",fhNBinsDcaZ-1,fhBinLimDcaZ);
713     fhQA[kDcaXYnorm][i]         = new TH1F(Form("%s_dcaXYnorm%s",GetName(),str),"",fhNBinsDcaXYnorm-1,fhBinLimDcaXYnorm);
714     fhQA[kDcaZnorm][i]          = new TH1F(Form("%s_dcaZnorm%s",GetName(),str),"",fhNBinsDcaZnorm-1,fhBinLimDcaZnorm);
715     fhQA[kSigmaDcaXY][i]        = new TH1F(Form("%s_sigmaDcaXY%s",GetName(),str),"",fhNBinsSigmaDcaXY-1,fhBinLimSigmaDcaXY);
716     fhQA[kSigmaDcaZ][i]         = new TH1F(Form("%s_sigmaDcaZ%s",GetName(),str),"",fhNBinsSigmaDcaZ-1,fhBinLimSigmaDcaZ);
717
718     fhDcaXYvsDcaZ[i]->SetXTitle("impact par. d_{z}");
719     fhDcaXYvsDcaZ[i]->SetYTitle("impact par. d_{xy}");
720
721     fhQA[kCutNSigmaToVertex][i]->SetXTitle("n #sigma to vertex");
722     fhQA[kCutRequireSigmaToVertex][i]->SetXTitle("require #sigma to vertex");
723     fhQA[kCutAcceptKinkDaughters][i]->SetXTitle("accept kink daughters");
724     fhQA[kDcaXY][i]->SetXTitle("impact par. d_{xy}");
725     fhQA[kDcaZ][i]->SetXTitle("impact par. d_{z}");
726     fhQA[kDcaXYnorm][i]->SetXTitle("norm. impact par. d_{xy} / #sigma_{xy}");
727     fhQA[kDcaZnorm][i]->SetXTitle("norm. impact par. d_{z} / #sigma_{z}");
728     fhQA[kSigmaDcaXY][i]->SetXTitle("impact par. resolution #sigma_{xy}");
729     fhQA[kSigmaDcaZ][i]->SetXTitle("impact par. resolution #sigma_{z}");
730   }
731   for(Int_t i=0; i<kNHist; i++) fhQA[i][1]->SetLineColor(color);
732 }
733 //__________________________________________________________________________________
734 void AliCFTrackIsPrimaryCuts::FillHistograms(TObject* obj, Bool_t f)
735 {
736   //
737   // fill the QA histograms
738   //
739
740   if (!obj) return;
741
742   Bool_t isESDTrack = strcmp(obj->ClassName(),"AliESDtrack") == 0 ? kTRUE : kFALSE ;
743   Bool_t isAODTrack = strcmp(obj->ClassName(),"AliAODTrack") == 0 ? kTRUE : kFALSE ;
744
745   AliESDtrack * esdTrack = 0x0 ;
746   AliAODTrack * aodTrack = 0x0 ; 
747   if (isESDTrack) esdTrack = dynamic_cast<AliESDtrack*>(obj);
748   if (isAODTrack) aodTrack = dynamic_cast<AliAODTrack*>(obj);
749
750   // f = 0: fill histograms before cuts
751   // f = 1: fill histograms after cuts
752
753   // get the track to vertex parameter for ESD track
754   if (esdTrack) {
755
756     fhQA[kDcaZ][f]->Fill(fDCA[1]);
757     fhQA[kDcaXY][f]->Fill(fDCA[0]);
758     fhDcaXYvsDcaZ[f]->Fill(fDCA[1],fDCA[0]);
759     fhQA[kSigmaDcaXY][f]->Fill(fDCA[2]);
760     fhQA[kSigmaDcaZ][f]->Fill(fDCA[3]);
761 // // // // // // //  delete histograms
762       fhQA[kDcaZnorm][f]->Fill(fDCA[1]);
763       fhQA[kDcaXYnorm][f]->Fill(fDCA[0]);
764
765     fhQA[kCutNSigmaToVertex][f]->Fill(fDCA[5]);
766     if (fDCA[5]<0 && fRequireSigmaToVertex) fhQA[kCutRequireSigmaToVertex][f]->Fill(0.);
767     if (!(fDCA[5]<0 && fRequireSigmaToVertex)) fhQA[kCutRequireSigmaToVertex][f]->Fill(1.);
768
769     if (!fAcceptKinkDaughters && esdTrack->GetKinkIndex(0)>0) fhQA[kCutAcceptKinkDaughters][f]->Fill(0.);
770     if (!(!fAcceptKinkDaughters && esdTrack->GetKinkIndex(0)>0)) fhQA[kCutAcceptKinkDaughters][f]->Fill(0.);
771   }
772
773   // fill cut statistics and cut correlation histograms with information from the bitmap
774   if (f) return;
775
776   // Get the bitmap of the single cuts
777   if ( !obj ) return;
778   SelectionBitMap(obj);
779
780   // Number of single cuts in this class
781   UInt_t ncuts = fBitmap->GetNbits();
782   for(UInt_t bit=0; bit<ncuts;bit++) {
783     if (!fBitmap->TestBitNumber(bit)) {
784         fhCutStatistics->Fill(bit+1);
785         for (UInt_t bit2=bit; bit2<ncuts;bit2++) {
786           if (!fBitmap->TestBitNumber(bit2)) 
787             fhCutCorrelation->Fill(bit+1,bit2+1);
788         }
789     }
790   }
791 }
792 //__________________________________________________________________________________
793 void AliCFTrackIsPrimaryCuts::SaveHistograms(const Char_t* dir) {
794   //
795   // saves the histograms in a directory (dir)
796   //
797   if(!fIsQAOn) return;
798
799   if (!dir)
800     dir = GetName();
801
802   gDirectory->mkdir(dir);
803   gDirectory->cd(dir);
804
805   gDirectory->mkdir("before_cuts");
806   gDirectory->mkdir("after_cuts");
807
808   fhCutStatistics->Write();
809   fhCutCorrelation->Write();
810
811   for (Int_t j=0; j<kNStepQA; j++) {
812     if (j==0)
813       gDirectory->cd("before_cuts");
814     else
815       gDirectory->cd("after_cuts");
816
817     fhDcaXYvsDcaZ[j]    ->Write();
818
819     for(Int_t i=0; i<kNHist; i++) fhQA[i][j]->Write();
820
821     gDirectory->cd("../");
822   }
823   gDirectory->cd("../");
824 }
825 //__________________________________________________________________________________
826 void AliCFTrackIsPrimaryCuts::DrawHistograms()
827 {
828   //
829   // draws some histograms
830   //
831   if(!fIsQAOn) return;
832
833   // pad margins
834   Float_t right = 0.03;
835   Float_t left = 0.175;
836   Float_t top = 0.03;
837   Float_t bottom = 0.175;
838
839   TCanvas* canvas1 = new TCanvas("Track_QA_Primary_1", "Track QA Primary 1", 800, 500);
840   canvas1->Divide(2, 1);
841
842   canvas1->cd(1);
843   fhCutStatistics->SetStats(kFALSE);
844   fhCutStatistics->LabelsOption("v");
845   gPad->SetLeftMargin(left);
846   gPad->SetBottomMargin(0.25);
847   gPad->SetRightMargin(right);
848   gPad->SetTopMargin(0.1);
849   fhCutStatistics->Draw();
850
851   canvas1->cd(2);
852   fhCutCorrelation->SetStats(kFALSE);
853   fhCutCorrelation->LabelsOption("v");
854   gPad->SetLeftMargin(0.30);
855   gPad->SetRightMargin(bottom);
856   gPad->SetTopMargin(0.1);
857   gPad->SetBottomMargin(0.25);
858   fhCutCorrelation->Draw("COLZ");
859
860   canvas1->SaveAs(Form("%s.eps", canvas1->GetName()));
861   canvas1->SaveAs(Form("%s.ps", canvas1->GetName()));
862
863   // -----
864
865   TCanvas* canvas2 = new TCanvas("Track_QA_Primary_2", "Track QA Primary 2", 800, 800);
866   canvas2->Divide(2, 2);
867
868   canvas2->cd(1);
869   gPad->SetRightMargin(right);
870   gPad->SetLeftMargin(left);
871   gPad->SetTopMargin(top);
872   gPad->SetBottomMargin(bottom);
873   fhQA[kDcaXY][0]->SetStats(kFALSE);
874   fhQA[kDcaXY][0]->Draw();
875   fhQA[kDcaXY][1]->Draw("same");
876
877   canvas2->cd(2);
878   gPad->SetRightMargin(right);
879   gPad->SetLeftMargin(left);
880   gPad->SetTopMargin(top);
881   gPad->SetBottomMargin(bottom);
882   fhQA[kDcaZ][0]->SetStats(kFALSE);
883   fhQA[kDcaZ][0]->Draw();
884   fhQA[kDcaZ][1]->Draw("same");
885
886   canvas2->cd(3);
887 //   fhDXYvsDZ[0]->SetStats(kFALSE);
888   gPad->SetLogz();
889   gPad->SetLeftMargin(bottom);
890   gPad->SetTopMargin(0.1);
891   gPad->SetBottomMargin(bottom);
892   gPad->SetRightMargin(0.2);
893   fhDcaXYvsDcaZ[0]->Draw("COLZ");
894
895   canvas2->cd(4);
896 //   fhDXYvsDZ[1]->SetStats(kFALSE);
897   gPad->SetLogz();
898   gPad->SetLeftMargin(bottom);
899   gPad->SetTopMargin(0.1);
900   gPad->SetBottomMargin(bottom);
901   gPad->SetRightMargin(0.2);
902   fhDcaXYvsDcaZ[1]->Draw("COLZ");
903
904   canvas2->SaveAs(Form("%s.eps", canvas2->GetName()));
905   canvas2->SaveAs(Form("%s.ps", canvas2->GetName()));
906
907   // -----
908
909   TCanvas* canvas3 = new TCanvas("Track_QA_Primary_3", "Track QA Primary 3", 800, 400);
910   canvas3->Divide(2, 1);
911
912   canvas3->cd(1);
913   gPad->SetRightMargin(right);
914   gPad->SetLeftMargin(left);
915   gPad->SetTopMargin(top);
916   gPad->SetBottomMargin(bottom);
917   fhQA[kDcaXYnorm][0]->SetStats(kFALSE);
918   fhQA[kDcaXYnorm][0]->Draw();
919   fhQA[kDcaXYnorm][1]->Draw("same");
920
921   canvas3->cd(2);
922   gPad->SetRightMargin(right);
923   gPad->SetLeftMargin(left);
924   gPad->SetTopMargin(top);
925   gPad->SetBottomMargin(bottom);
926   fhQA[kDcaZnorm][0]->SetStats(kFALSE);
927   fhQA[kDcaZnorm][0]->Draw();
928   fhQA[kDcaZnorm][1]->Draw("same");
929
930   canvas3->SaveAs(Form("%s.eps", canvas3->GetName()));
931   canvas3->SaveAs(Form("%s.ps", canvas3->GetName()));
932
933   // -----
934
935   TCanvas* canvas4 = new TCanvas("Track_QA_Primary_4", "Track QA Primary 4", 1200, 800);
936   canvas4->Divide(3, 2);
937
938   canvas4->cd(1);
939   gPad->SetRightMargin(right);
940   gPad->SetLeftMargin(left);
941   gPad->SetTopMargin(top);
942   gPad->SetBottomMargin(bottom);
943   fhQA[kSigmaDcaXY][0]->SetStats(kFALSE);
944   fhQA[kSigmaDcaXY][0]->Draw();
945   fhQA[kSigmaDcaXY][1]->Draw("same");
946
947   canvas4->cd(2);
948   gPad->SetRightMargin(right);
949   gPad->SetLeftMargin(left);
950   gPad->SetTopMargin(top);
951   gPad->SetBottomMargin(bottom);
952   fhQA[kSigmaDcaZ][0]->SetStats(kFALSE);
953   fhQA[kSigmaDcaZ][0]->Draw();
954   fhQA[kSigmaDcaZ][1]->Draw("same");
955
956   canvas4->cd(4);
957   gPad->SetRightMargin(right);
958   gPad->SetLeftMargin(left);
959   gPad->SetTopMargin(top);
960   gPad->SetBottomMargin(bottom);
961   fhQA[kCutNSigmaToVertex][0]->SetStats(kFALSE);
962   fhQA[kCutNSigmaToVertex][0]->Draw();
963   fhQA[kCutNSigmaToVertex][1]->Draw("same");
964
965   canvas4->cd(5);
966   gPad->SetRightMargin(right);
967   gPad->SetLeftMargin(left);
968   gPad->SetTopMargin(top);
969   gPad->SetBottomMargin(bottom);
970   fhQA[kCutRequireSigmaToVertex][0]->SetStats(kFALSE);
971   fhQA[kCutRequireSigmaToVertex][0]->Draw();
972   fhQA[kCutRequireSigmaToVertex][1]->Draw("same");
973
974   canvas4->cd(6);
975   gPad->SetRightMargin(right);
976   gPad->SetLeftMargin(left);
977   gPad->SetTopMargin(top);
978   gPad->SetBottomMargin(bottom);
979   fhQA[kCutAcceptKinkDaughters][0]->SetStats(kFALSE);
980   fhQA[kCutAcceptKinkDaughters][0]->Draw();
981   fhQA[kCutAcceptKinkDaughters][1]->Draw("same");
982
983   canvas4->SaveAs(Form("%s.eps", canvas4->GetName()));
984   canvas4->SaveAs(Form("%s.ps", canvas4->GetName()));
985 }
986 //__________________________________________________________________________________
987 void AliCFTrackIsPrimaryCuts::AddQAHistograms(TList *qaList) {
988   //
989   // saves the histograms in a TList
990   //
991   DefineHistograms();
992
993   qaList->Add(fhCutStatistics);
994   qaList->Add(fhCutCorrelation);
995
996   for (Int_t j=0; j<kNStepQA; j++) {
997     qaList->Add(fhDcaXYvsDcaZ[j]);
998     for(Int_t i=0; i<kNHist; i++)
999         qaList->Add(fhQA[i][j]);
1000   }
1001 }