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