]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ANALYSIS/AliESDtrackCuts.cxx
Moved from PWG0
[u/mrichter/AliRoot.git] / ANALYSIS / AliESDtrackCuts.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 /* $Id: AliESDtrackCuts.cxx 24534 2008-03-16 22:22:11Z fca $ */
17
18 #include "AliESDtrackCuts.h"
19
20 #include <AliESDtrack.h>
21 #include <AliESD.h>
22 #include <AliESDEvent.h>
23 #include <AliLog.h>
24
25 #include <TTree.h>
26 #include <TCanvas.h>
27 #include <TDirectory.h>
28
29 //____________________________________________________________________
30 ClassImp(AliESDtrackCuts)
31
32 // Cut names
33 const Char_t* AliESDtrackCuts::fgkCutNames[kNCuts] = {
34  "require TPC refit",
35  "require ITS refit",
36  "n clusters TPC",
37  "n clusters ITS",
38  "#Chi^{2}/clusters TPC",
39  "#Chi^{2}/clusters ITS",
40  "cov 11",
41  "cov 22",
42  "cov 33",
43  "cov 44",
44  "cov 55",
45  "trk-to-vtx",
46  "trk-to-vtx failed",
47  "kink daughters",
48  "p",
49  "p_{T}",
50  "p_{x}",
51  "p_{y}",
52  "p_{z}",
53  "y",
54  "eta"
55 };
56
57 //____________________________________________________________________
58 AliESDtrackCuts::AliESDtrackCuts(const Char_t* name, const Char_t* title) : AliAnalysisCuts(name,title),
59   fCutMinNClusterTPC(0),
60   fCutMinNClusterITS(0),
61   fCutMaxChi2PerClusterTPC(0),
62   fCutMaxChi2PerClusterITS(0),
63   fCutMaxC11(0),
64   fCutMaxC22(0),
65   fCutMaxC33(0),
66   fCutMaxC44(0),
67   fCutMaxC55(0),
68   fCutAcceptKinkDaughters(0),
69   fCutRequireTPCRefit(0),
70   fCutRequireITSRefit(0),
71   fCutNsigmaToVertex(0),
72   fCutSigmaToVertexRequired(0),
73   fPMin(0),
74   fPMax(0),
75   fPtMin(0),
76   fPtMax(0),
77   fPxMin(0),
78   fPxMax(0),
79   fPyMin(0),
80   fPyMax(0),
81   fPzMin(0),
82   fPzMax(0),
83   fEtaMin(0),
84   fEtaMax(0),
85   fRapMin(0),
86   fRapMax(0),
87   fHistogramsOn(0),
88   ffDTheoretical(0),
89   fhCutStatistics(0),         
90   fhCutCorrelation(0)
91 {
92   //
93   // constructor
94   //
95
96   Init();
97
98   //##############################################################################
99   // setting default cuts
100   SetMinNClustersTPC();
101   SetMinNClustersITS();
102   SetMaxChi2PerClusterTPC();
103   SetMaxChi2PerClusterITS();                                
104   SetMaxCovDiagonalElements();                                      
105   SetRequireTPCRefit();
106   SetRequireITSRefit();
107   SetAcceptKingDaughters();
108   SetMinNsigmaToVertex();
109   SetRequireSigmaToVertex();
110   SetPRange();
111   SetPtRange();
112   SetPxRange();
113   SetPyRange();
114   SetPzRange();
115   SetEtaRange();
116   SetRapRange();
117
118   SetHistogramsOn();
119 }
120
121 //_____________________________________________________________________________
122 AliESDtrackCuts::AliESDtrackCuts(const AliESDtrackCuts &c) : AliAnalysisCuts(c),
123   fCutMinNClusterTPC(0),
124   fCutMinNClusterITS(0),
125   fCutMaxChi2PerClusterTPC(0),
126   fCutMaxChi2PerClusterITS(0),
127   fCutMaxC11(0),
128   fCutMaxC22(0),
129   fCutMaxC33(0),
130   fCutMaxC44(0),
131   fCutMaxC55(0),
132   fCutAcceptKinkDaughters(0),
133   fCutRequireTPCRefit(0),
134   fCutRequireITSRefit(0),
135   fCutNsigmaToVertex(0),
136   fCutSigmaToVertexRequired(0),
137   fPMin(0),
138   fPMax(0),
139   fPtMin(0),
140   fPtMax(0),
141   fPxMin(0),
142   fPxMax(0),
143   fPyMin(0),
144   fPyMax(0),
145   fPzMin(0),
146   fPzMax(0),
147   fEtaMin(0),
148   fEtaMax(0),
149   fRapMin(0),
150   fRapMax(0),
151   fHistogramsOn(0),
152   ffDTheoretical(0),                                 
153   fhCutStatistics(0),         
154   fhCutCorrelation(0)
155 {
156   //
157   // copy constructor
158   //
159
160   ((AliESDtrackCuts &) c).Copy(*this);
161 }
162
163 AliESDtrackCuts::~AliESDtrackCuts()
164 {
165   //
166   // destructor
167   //
168
169   for (Int_t i=0; i<2; i++) {
170     
171     if (fhNClustersITS[i])
172       delete fhNClustersITS[i];            
173     if (fhNClustersTPC[i])
174       delete fhNClustersTPC[i];                
175     if (fhChi2PerClusterITS[i])
176       delete fhChi2PerClusterITS[i];       
177     if (fhChi2PerClusterTPC[i])
178       delete fhChi2PerClusterTPC[i];       
179     if (fhC11[i])
180       delete fhC11[i];                     
181     if (fhC22[i])
182       delete fhC22[i];                     
183     if (fhC33[i])
184       delete fhC33[i];                     
185     if (fhC44[i])
186       delete fhC44[i];                     
187     if (fhC55[i])
188     delete fhC55[i];                     
189     
190     if (fhDXY[i])
191       delete fhDXY[i];                     
192     if (fhDZ[i])
193       delete fhDZ[i];                      
194     if (fhDXYvsDZ[i])
195       delete fhDXYvsDZ[i];                 
196     
197     if (fhDXYNormalized[i])
198       delete fhDXYNormalized[i];           
199     if (fhDZNormalized[i])
200       delete fhDZNormalized[i];
201     if (fhDXYvsDZNormalized[i])
202       delete fhDXYvsDZNormalized[i];       
203     if (fhNSigmaToVertex[i])
204       delete fhNSigmaToVertex[i];
205     if (fhPt[i])
206       delete fhPt[i];
207     if (fhEta[i])
208       delete fhEta[i];
209   }
210
211   if (ffDTheoretical)
212     delete ffDTheoretical;
213
214   if (fhCutStatistics)
215     delete fhCutStatistics;             
216   if (fhCutCorrelation)
217     delete fhCutCorrelation;            
218 }
219
220 void AliESDtrackCuts::Init()
221 {
222   //
223   // sets everything to zero
224   //
225
226   fCutMinNClusterTPC = 0;
227   fCutMinNClusterITS = 0;
228
229   fCutMaxChi2PerClusterTPC = 0;
230   fCutMaxChi2PerClusterITS = 0;
231
232   fCutMaxC11 = 0;
233   fCutMaxC22 = 0;
234   fCutMaxC33 = 0;
235   fCutMaxC44 = 0;
236   fCutMaxC55 = 0;
237
238   fCutAcceptKinkDaughters = 0;
239   fCutRequireTPCRefit = 0;
240   fCutRequireITSRefit = 0;
241
242   fCutNsigmaToVertex = 0;
243   fCutSigmaToVertexRequired = 0;
244
245   fPMin = 0;
246   fPMax = 0;
247   fPtMin = 0;
248   fPtMax = 0;
249   fPxMin = 0;
250   fPxMax = 0;
251   fPyMin = 0;
252   fPyMax = 0;
253   fPzMin = 0;
254   fPzMax = 0;
255   fEtaMin = 0;
256   fEtaMax = 0;
257   fRapMin = 0;
258   fRapMax = 0;
259
260   fHistogramsOn = kFALSE;
261
262   for (Int_t i=0; i<2; ++i)
263   {
264     fhNClustersITS[i] = 0;
265     fhNClustersTPC[i] = 0;
266
267     fhChi2PerClusterITS[i] = 0;
268     fhChi2PerClusterTPC[i] = 0;
269
270     fhC11[i] = 0;
271     fhC22[i] = 0;
272     fhC33[i] = 0;
273     fhC44[i] = 0;
274     fhC55[i] = 0;
275
276     fhDXY[i] = 0;
277     fhDZ[i] = 0;
278     fhDXYvsDZ[i] = 0;
279
280     fhDXYNormalized[i] = 0;
281     fhDZNormalized[i] = 0;
282     fhDXYvsDZNormalized[i] = 0;
283     fhNSigmaToVertex[i] = 0;
284     
285     fhPt[i] = 0;
286     fhEta[i] = 0;
287   }
288   ffDTheoretical = 0;
289
290   fhCutStatistics = 0;
291   fhCutCorrelation = 0;
292 }
293
294 //_____________________________________________________________________________
295 AliESDtrackCuts &AliESDtrackCuts::operator=(const AliESDtrackCuts &c)
296 {
297   //
298   // Assignment operator
299   //
300
301   if (this != &c) ((AliESDtrackCuts &) c).Copy(*this);
302   return *this;
303 }
304
305 //_____________________________________________________________________________
306 void AliESDtrackCuts::Copy(TObject &c) const
307 {
308   //
309   // Copy function
310   //
311
312   AliESDtrackCuts& target = (AliESDtrackCuts &) c;
313
314   target.Init();
315
316   target.fCutMinNClusterTPC = fCutMinNClusterTPC;
317   target.fCutMinNClusterITS = fCutMinNClusterITS;
318
319   target.fCutMaxChi2PerClusterTPC = fCutMaxChi2PerClusterTPC;
320   target.fCutMaxChi2PerClusterITS = fCutMaxChi2PerClusterITS;
321
322   target.fCutMaxC11 = fCutMaxC11;
323   target.fCutMaxC22 = fCutMaxC22;
324   target.fCutMaxC33 = fCutMaxC33;
325   target.fCutMaxC44 = fCutMaxC44;
326   target.fCutMaxC55 = fCutMaxC55;
327
328   target.fCutAcceptKinkDaughters = fCutAcceptKinkDaughters;
329   target.fCutRequireTPCRefit = fCutRequireTPCRefit;
330   target.fCutRequireITSRefit = fCutRequireITSRefit;
331
332   target.fCutNsigmaToVertex = fCutNsigmaToVertex;
333   target.fCutSigmaToVertexRequired = fCutSigmaToVertexRequired;
334
335   target.fPMin = fPMin;
336   target.fPMax = fPMax;
337   target.fPtMin = fPtMin;
338   target.fPtMax = fPtMax;
339   target.fPxMin = fPxMin;
340   target.fPxMax = fPxMax;
341   target.fPyMin = fPyMin;
342   target.fPyMax = fPyMax;
343   target.fPzMin = fPzMin;
344   target.fPzMax = fPzMax;
345   target.fEtaMin = fEtaMin;
346   target.fEtaMax = fEtaMax;
347   target.fRapMin = fRapMin;
348   target.fRapMax = fRapMax;
349
350   target.fHistogramsOn = fHistogramsOn;
351
352   for (Int_t i=0; i<2; ++i)
353   {
354     if (fhNClustersITS[i]) target.fhNClustersITS[i] = (TH1F*) fhNClustersITS[i]->Clone();
355     if (fhNClustersTPC[i]) target.fhNClustersTPC[i] = (TH1F*) fhNClustersTPC[i]->Clone();
356
357     if (fhChi2PerClusterITS[i]) target.fhChi2PerClusterITS[i] = (TH1F*) fhChi2PerClusterITS[i]->Clone();
358     if (fhChi2PerClusterTPC[i]) target.fhChi2PerClusterTPC[i] = (TH1F*) fhChi2PerClusterTPC[i]->Clone();
359
360     if (fhC11[i]) target.fhC11[i] = (TH1F*) fhC11[i]->Clone();
361     if (fhC22[i]) target.fhC22[i] = (TH1F*) fhC22[i]->Clone();
362     if (fhC33[i]) target.fhC33[i] = (TH1F*) fhC33[i]->Clone();
363     if (fhC44[i]) target.fhC44[i] = (TH1F*) fhC44[i]->Clone();
364     if (fhC55[i]) target.fhC55[i] = (TH1F*) fhC55[i]->Clone();
365
366     if (fhDXY[i]) target.fhDXY[i] = (TH1F*) fhDXY[i]->Clone();
367     if (fhDZ[i]) target.fhDZ[i] = (TH1F*) fhDZ[i]->Clone();
368     if (fhDXYvsDZ[i]) target.fhDXYvsDZ[i] = (TH2F*) fhDXYvsDZ[i]->Clone();
369
370     if (fhDXYNormalized[i]) target.fhDXYNormalized[i] = (TH1F*) fhDXYNormalized[i]->Clone();
371     if (fhDZNormalized[i]) target.fhDZNormalized[i] = (TH1F*) fhDZNormalized[i]->Clone();
372     if (fhDXYvsDZNormalized[i]) target.fhDXYvsDZNormalized[i] = (TH2F*) fhDXYvsDZNormalized[i]->Clone();
373     if (fhNSigmaToVertex[i]) target.fhNSigmaToVertex[i] = (TH1F*) fhNSigmaToVertex[i]->Clone();
374     
375     if (fhPt[i]) target.fhPt[i] = (TH1F*) fhPt[i]->Clone();
376     if (fhEta[i]) target.fhEta[i] = (TH1F*) fhEta[i]->Clone();
377   }
378   if (ffDTheoretical) target.ffDTheoretical = (TF1*) ffDTheoretical->Clone();
379
380   if (fhCutStatistics) target.fhCutStatistics = (TH1F*) fhCutStatistics->Clone();
381   if (fhCutCorrelation) target.fhCutCorrelation = (TH2F*) fhCutCorrelation->Clone();
382
383   TNamed::Copy(c);
384 }
385
386 //_____________________________________________________________________________
387 Long64_t AliESDtrackCuts::Merge(TCollection* list) {
388   // Merge a list of AliESDtrackCuts objects with this (needed for PROOF)
389   // Returns the number of merged objects (including this)
390
391   if (!list)
392     return 0;
393   
394   if (list->IsEmpty())
395     return 1;
396
397   if (!fHistogramsOn)
398     return 0;
399
400   TIterator* iter = list->MakeIterator();
401   TObject* obj;
402
403
404   // collection of measured and generated histograms
405   Int_t count = 0;
406   while ((obj = iter->Next())) {
407
408     AliESDtrackCuts* entry = dynamic_cast<AliESDtrackCuts*>(obj);
409     if (entry == 0)
410       continue;
411
412     if (!entry->fHistogramsOn)
413       continue;
414     
415     for (Int_t i=0; i<2; i++) {
416       
417       fhNClustersITS[i]      ->Add(entry->fhNClustersITS[i]     );      
418       fhNClustersTPC[i]      ->Add(entry->fhNClustersTPC[i]     ); 
419                                                                     
420       fhChi2PerClusterITS[i] ->Add(entry->fhChi2PerClusterITS[i]); 
421       fhChi2PerClusterTPC[i] ->Add(entry->fhChi2PerClusterTPC[i]); 
422                                                                     
423       fhC11[i]               ->Add(entry->fhC11[i]              ); 
424       fhC22[i]               ->Add(entry->fhC22[i]              ); 
425       fhC33[i]               ->Add(entry->fhC33[i]              ); 
426       fhC44[i]               ->Add(entry->fhC44[i]              ); 
427       fhC55[i]               ->Add(entry->fhC55[i]              ); 
428                                                                     
429       fhDXY[i]               ->Add(entry->fhDXY[i]              ); 
430       fhDZ[i]                ->Add(entry->fhDZ[i]               ); 
431       fhDXYvsDZ[i]           ->Add(entry->fhDXYvsDZ[i]          ); 
432                                                                     
433       fhDXYNormalized[i]     ->Add(entry->fhDXYNormalized[i]    ); 
434       fhDZNormalized[i]      ->Add(entry->fhDZNormalized[i]     );
435       fhDXYvsDZNormalized[i] ->Add(entry->fhDXYvsDZNormalized[i]); 
436       fhNSigmaToVertex[i]    ->Add(entry->fhNSigmaToVertex[i]); 
437
438       fhPt[i]                ->Add(entry->fhPt[i]); 
439       fhEta[i]               ->Add(entry->fhEta[i]); 
440     }      
441
442     fhCutStatistics  ->Add(entry->fhCutStatistics);        
443     fhCutCorrelation ->Add(entry->fhCutCorrelation);      
444
445     count++;
446   }
447
448   return count+1;
449 }
450
451
452 //____________________________________________________________________
453 Float_t AliESDtrackCuts::GetSigmaToVertex(AliESDtrack* esdTrack)
454 {
455   // Calculates the number of sigma to the vertex.
456
457   Float_t b[2];
458   Float_t bRes[2];
459   Float_t bCov[3];
460   esdTrack->GetImpactParameters(b,bCov);
461   if (bCov[0]<=0 || bCov[2]<=0) {
462     AliDebug(1, "Estimated b resolution lower or equal zero!");
463     bCov[0]=0; bCov[2]=0;
464   }
465   bRes[0] = TMath::Sqrt(bCov[0]);
466   bRes[1] = TMath::Sqrt(bCov[2]);
467
468   // -----------------------------------
469   // How to get to a n-sigma cut?
470   //
471   // The accumulated statistics from 0 to d is
472   //
473   // ->  Erf(d/Sqrt(2)) for a 1-dim gauss (d = n_sigma)
474   // ->  1 - Exp(-d**2) for a 2-dim gauss (d*d = dx*dx + dy*dy != n_sigma)
475   //
476   // It means that for a 2-dim gauss: n_sigma(d) = Sqrt(2)*ErfInv(1 - Exp((-x**2)/2)
477   // Can this be expressed in a different way?
478
479   if (bRes[0] == 0 || bRes[1] ==0)
480     return -1;
481
482   Float_t d = TMath::Sqrt(TMath::Power(b[0]/bRes[0],2) + TMath::Power(b[1]/bRes[1],2));
483
484   // stupid rounding problem screws up everything:
485   // if d is too big, TMath::Exp(...) gets 0, and TMath::ErfInverse(1) that should be infinite, gets 0 :(
486   if (TMath::Exp(-d * d / 2) < 1e-10)
487     return 1000;
488
489   d = TMath::ErfInverse(1 - TMath::Exp(-d * d / 2)) * TMath::Sqrt(2);
490   return d;
491 }
492
493 void AliESDtrackCuts::EnableNeededBranches(TTree* tree)
494 {
495   // enables the branches needed by AcceptTrack, for a list see comment of AcceptTrack
496
497   tree->SetBranchStatus("fTracks.fFlags", 1);
498   tree->SetBranchStatus("fTracks.fITSncls", 1);
499   tree->SetBranchStatus("fTracks.fTPCncls", 1);
500   tree->SetBranchStatus("fTracks.fITSchi2", 1);
501   tree->SetBranchStatus("fTracks.fTPCchi2", 1);
502   tree->SetBranchStatus("fTracks.fC*", 1);
503   tree->SetBranchStatus("fTracks.fD", 1);
504   tree->SetBranchStatus("fTracks.fZ", 1);
505   tree->SetBranchStatus("fTracks.fCdd", 1);
506   tree->SetBranchStatus("fTracks.fCdz", 1);
507   tree->SetBranchStatus("fTracks.fCzz", 1);
508   tree->SetBranchStatus("fTracks.fP*", 1);
509   tree->SetBranchStatus("fTracks.fR*", 1);
510   tree->SetBranchStatus("fTracks.fKinkIndexes*", 1);
511 }
512
513 //____________________________________________________________________
514 Bool_t
515 AliESDtrackCuts::AcceptTrack(AliESDtrack* esdTrack) {
516   // 
517   // figure out if the tracks survives all the track cuts defined
518   //
519   // the different quality parameter and kinematic values are first
520   // retrieved from the track. then it is found out what cuts the
521   // track did not survive and finally the cuts are imposed.
522
523   // this function needs the following branches:
524   // fTracks.fFlags
525   // fTracks.fITSncls
526   // fTracks.fTPCncls
527   // fTracks.fITSchi2
528   // fTracks.fTPCchi2
529   // fTracks.fC   //GetExternalCovariance
530   // fTracks.fD   //GetImpactParameters
531   // fTracks.fZ   //GetImpactParameters
532   // fTracks.fCdd //GetImpactParameters
533   // fTracks.fCdz //GetImpactParameters
534   // fTracks.fCzz //GetImpactParameters
535   // fTracks.fP   //GetPxPyPz
536   // fTracks.fR   //GetMass
537   // fTracks.fP   //GetMass
538   // fTracks.fKinkIndexes
539
540   UInt_t status = esdTrack->GetStatus();
541
542   // dummy array
543   Int_t  fIdxInt[200];
544
545   // getting quality parameters from the ESD track
546   Int_t nClustersITS = esdTrack->GetITSclusters(fIdxInt);
547   Int_t nClustersTPC = esdTrack->GetTPCclusters(fIdxInt);
548   
549
550
551   Float_t chi2PerClusterITS = -1;
552   Float_t chi2PerClusterTPC = -1;
553   if (nClustersITS!=0)
554     chi2PerClusterITS = esdTrack->GetITSchi2()/Float_t(nClustersITS);
555   if (nClustersTPC!=0)
556     chi2PerClusterTPC = esdTrack->GetTPCchi2()/Float_t(nClustersTPC);
557
558   Double_t extCov[15];
559   esdTrack->GetExternalCovariance(extCov);
560
561   // getting the track to vertex parameters
562   Float_t nSigmaToVertex = GetSigmaToVertex(esdTrack);
563
564   // getting the kinematic variables of the track
565   // (assuming the mass is known)
566   Double_t p[3];
567   esdTrack->GetPxPyPz(p);
568   Float_t momentum = TMath::Sqrt(TMath::Power(p[0],2) + TMath::Power(p[1],2) + TMath::Power(p[2],2));
569   Float_t pt       = TMath::Sqrt(TMath::Power(p[0],2) + TMath::Power(p[1],2));
570   Float_t energy   = TMath::Sqrt(TMath::Power(esdTrack->GetMass(),2) + TMath::Power(momentum,2));
571
572
573   //y-eta related calculations
574   Float_t eta = -100.;
575   Float_t y   = -100.;
576   if((momentum != TMath::Abs(p[2]))&&(momentum != 0))
577     eta = 0.5*TMath::Log((momentum + p[2])/(momentum - p[2]));
578   if((energy != TMath::Abs(p[2]))&&(momentum != 0))
579     y = 0.5*TMath::Log((energy + p[2])/(energy - p[2]));
580
581   
582   //########################################################################
583   // cut the track?
584   
585   Bool_t cuts[kNCuts];
586   for (Int_t i=0; i<kNCuts; i++) cuts[i]=kFALSE;
587   
588   // track quality cuts
589   if (fCutRequireTPCRefit && (status&AliESDtrack::kTPCrefit)==0)
590     cuts[0]=kTRUE;
591   if (fCutRequireITSRefit && (status&AliESDtrack::kITSrefit)==0)
592     cuts[1]=kTRUE;
593   if (nClustersTPC<fCutMinNClusterTPC)
594     cuts[2]=kTRUE;
595   if (nClustersITS<fCutMinNClusterITS) 
596     cuts[3]=kTRUE;
597   if (chi2PerClusterTPC>fCutMaxChi2PerClusterTPC) 
598     cuts[4]=kTRUE; 
599   if (chi2PerClusterITS>fCutMaxChi2PerClusterITS) 
600     cuts[5]=kTRUE;
601   if (extCov[0]  > fCutMaxC11) 
602     cuts[6]=kTRUE;  
603   if (extCov[2]  > fCutMaxC22) 
604     cuts[7]=kTRUE;  
605   if (extCov[5]  > fCutMaxC33) 
606     cuts[8]=kTRUE;  
607   if (extCov[9]  > fCutMaxC44) 
608     cuts[9]=kTRUE;  
609   if (extCov[14]  > fCutMaxC55) 
610     cuts[10]=kTRUE;  
611   if (nSigmaToVertex > fCutNsigmaToVertex && fCutSigmaToVertexRequired)
612     cuts[11] = kTRUE;
613   // if n sigma could not be calculated
614   if (nSigmaToVertex<0 && fCutSigmaToVertexRequired)
615     cuts[12]=kTRUE;
616   if (!fCutAcceptKinkDaughters && esdTrack->GetKinkIndex(0)>0)
617     cuts[13]=kTRUE;
618   // track kinematics cut
619   if((momentum < fPMin) || (momentum > fPMax)) 
620     cuts[14]=kTRUE;
621   if((pt < fPtMin) || (pt > fPtMax)) 
622     cuts[15] = kTRUE;
623   if((p[0] < fPxMin) || (p[0] > fPxMax)) 
624     cuts[16] = kTRUE;
625   if((p[1] < fPyMin) || (p[1] > fPyMax)) 
626     cuts[17] = kTRUE;
627   if((p[2] < fPzMin) || (p[2] > fPzMax))
628     cuts[18] = kTRUE;
629   if((eta < fEtaMin) || (eta > fEtaMax)) 
630     cuts[19] = kTRUE;
631   if((y < fRapMin) || (y > fRapMax)) 
632     cuts[20] = kTRUE;
633
634   Bool_t cut=kFALSE;
635   for (Int_t i=0; i<kNCuts; i++) 
636     if (cuts[i]) cut = kTRUE;
637   
638   //########################################################################
639   // filling histograms
640   if (fHistogramsOn) {
641     fhCutStatistics->Fill(fhCutStatistics->GetBinCenter(fhCutStatistics->GetXaxis()->FindBin("n tracks")));
642     
643     if (cut)
644       fhCutStatistics->Fill(fhCutStatistics->GetBinCenter(fhCutStatistics->GetXaxis()->FindBin("n cut tracks")));
645     
646     for (Int_t i=0; i<kNCuts; i++) {
647       if (cuts[i])
648         fhCutStatistics->Fill(fhCutStatistics->GetBinCenter(fhCutStatistics->GetXaxis()->FindBin(fgkCutNames[i])));
649       
650       for (Int_t j=i; j<kNCuts; j++) {
651         if (cuts[i] && cuts[j]) {
652           Float_t x = fhCutCorrelation->GetXaxis()->GetBinCenter(fhCutCorrelation->GetXaxis()->FindBin(fgkCutNames[i]));
653           Float_t y = fhCutCorrelation->GetYaxis()->GetBinCenter(fhCutCorrelation->GetYaxis()->FindBin(fgkCutNames[j]));
654           fhCutCorrelation->Fill(x,y);
655         }
656       }
657     }
658     
659     fhNClustersITS[0]->Fill(nClustersITS);
660     fhNClustersTPC[0]->Fill(nClustersTPC);
661     fhChi2PerClusterITS[0]->Fill(chi2PerClusterITS);
662     fhChi2PerClusterTPC[0]->Fill(chi2PerClusterTPC);
663
664     fhC11[0]->Fill(extCov[0]);
665     fhC22[0]->Fill(extCov[2]);
666     fhC33[0]->Fill(extCov[5]);
667     fhC44[0]->Fill(extCov[9]);
668     fhC55[0]->Fill(extCov[14]);
669     
670     fhPt[0]->Fill(pt);
671     fhEta[0]->Fill(eta);
672
673     Float_t b[2];
674     Float_t bRes[2];
675     Float_t bCov[3];
676     esdTrack->GetImpactParameters(b,bCov);
677     if (bCov[0]<=0 || bCov[2]<=0) {
678       AliDebug(1, "Estimated b resolution lower or equal zero!");
679       bCov[0]=0; bCov[2]=0;
680     }
681     bRes[0] = TMath::Sqrt(bCov[0]);
682     bRes[1] = TMath::Sqrt(bCov[2]);
683
684     fhDZ[0]->Fill(b[1]);
685     fhDXY[0]->Fill(b[0]);
686     fhDXYvsDZ[0]->Fill(b[1],b[0]);
687
688     if (bRes[0]!=0 && bRes[1]!=0) {
689       fhDZNormalized[0]->Fill(b[1]/bRes[1]);
690       fhDXYNormalized[0]->Fill(b[0]/bRes[0]);
691       fhDXYvsDZNormalized[0]->Fill(b[1]/bRes[1], b[0]/bRes[0]);
692       fhNSigmaToVertex[0]->Fill(nSigmaToVertex);
693     }
694   }
695
696   //########################################################################
697   // cut the track!
698   if (cut) return kFALSE;
699
700   //########################################################################
701   // filling histograms after cut
702   if (fHistogramsOn) {
703     fhNClustersITS[1]->Fill(nClustersITS);
704     fhNClustersTPC[1]->Fill(nClustersTPC);
705     fhChi2PerClusterITS[1]->Fill(chi2PerClusterITS);
706     fhChi2PerClusterTPC[1]->Fill(chi2PerClusterTPC);
707
708     fhC11[1]->Fill(extCov[0]);
709     fhC22[1]->Fill(extCov[2]);
710     fhC33[1]->Fill(extCov[5]);
711     fhC44[1]->Fill(extCov[9]);
712     fhC55[1]->Fill(extCov[14]);
713
714     fhPt[1]->Fill(pt);
715     fhEta[1]->Fill(eta);
716     
717     Float_t b[2];
718     Float_t bRes[2];
719     Float_t bCov[3];
720     esdTrack->GetImpactParameters(b,bCov);
721     if (bCov[0]<=0 || bCov[2]<=0) {
722       AliDebug(1, "Estimated b resolution lower or equal zero!");
723       bCov[0]=0; bCov[2]=0;
724     }
725     bRes[0] = TMath::Sqrt(bCov[0]);
726     bRes[1] = TMath::Sqrt(bCov[2]);
727
728     fhDZ[1]->Fill(b[1]);
729     fhDXY[1]->Fill(b[0]);
730     fhDXYvsDZ[1]->Fill(b[1],b[0]);
731
732     if (bRes[0]!=0 && bRes[1]!=0)
733     {
734       fhDZNormalized[1]->Fill(b[1]/bRes[1]);
735       fhDXYNormalized[1]->Fill(b[0]/bRes[0]);
736       fhDXYvsDZNormalized[1]->Fill(b[1]/bRes[1], b[0]/bRes[0]);
737       fhNSigmaToVertex[1]->Fill(nSigmaToVertex);
738     }
739   }
740
741   return kTRUE;
742 }
743
744 //____________________________________________________________________
745 TObjArray* AliESDtrackCuts::GetAcceptedTracks(AliESD* esd)
746 {
747   //
748   // returns an array of all tracks that pass the cuts
749   //
750
751   TObjArray* acceptedTracks = new TObjArray();
752
753   // loop over esd tracks
754   for (Int_t iTrack = 0; iTrack < esd->GetNumberOfTracks(); iTrack++) {
755     AliESDtrack* track = esd->GetTrack(iTrack);
756
757     if (AcceptTrack(track))
758       acceptedTracks->Add(track);
759   }
760
761   return acceptedTracks;
762 }
763
764 //____________________________________________________________________
765 Int_t AliESDtrackCuts::CountAcceptedTracks(AliESD* esd)
766 {
767   //
768   // returns an the number of tracks that pass the cuts
769   //
770
771   Int_t count = 0;
772
773   // loop over esd tracks
774   for (Int_t iTrack = 0; iTrack < esd->GetNumberOfTracks(); iTrack++) {
775     AliESDtrack* track = esd->GetTrack(iTrack);
776
777     if (AcceptTrack(track))
778       count++;
779   }
780
781   return count;
782 }
783
784 //____________________________________________________________________
785 TObjArray* AliESDtrackCuts::GetAcceptedTracks(AliESDEvent* esd)
786 {
787   //
788   // returns an array of all tracks that pass the cuts
789   //
790
791   TObjArray* acceptedTracks = new TObjArray();
792
793   // loop over esd tracks
794   for (Int_t iTrack = 0; iTrack < esd->GetNumberOfTracks(); iTrack++) {
795     AliESDtrack* track = esd->GetTrack(iTrack);
796
797     if (AcceptTrack(track))
798       acceptedTracks->Add(track);
799   }
800
801   return acceptedTracks;
802 }
803
804 //____________________________________________________________________
805 Int_t AliESDtrackCuts::CountAcceptedTracks(AliESDEvent* esd)
806 {
807   //
808   // returns an the number of tracks that pass the cuts
809   //
810
811   Int_t count = 0;
812
813   // loop over esd tracks
814   for (Int_t iTrack = 0; iTrack < esd->GetNumberOfTracks(); iTrack++) {
815     AliESDtrack* track = esd->GetTrack(iTrack);
816
817     if (AcceptTrack(track))
818       count++;
819   }
820
821   return count;
822 }
823
824 //____________________________________________________________________
825  void AliESDtrackCuts::DefineHistograms(Int_t color) {
826    // 
827    // diagnostics histograms are defined
828    // 
829
830    fHistogramsOn=kTRUE;
831    
832    Bool_t oldStatus = TH1::AddDirectoryStatus();
833    TH1::AddDirectory(kFALSE);
834    
835    //###################################################################################
836    // defining histograms
837
838    fhCutStatistics = new TH1F("cut_statistics","cut statistics",kNCuts+4,-0.5,kNCuts+3.5);
839
840    fhCutStatistics->GetXaxis()->SetBinLabel(1,"n tracks");
841    fhCutStatistics->GetXaxis()->SetBinLabel(2,"n cut tracks");
842
843    fhCutCorrelation = new TH2F("cut_correlation","cut correlation",kNCuts,-0.5,kNCuts-0.5,kNCuts,-0.5,kNCuts-0.5);;
844   
845    for (Int_t i=0; i<kNCuts; i++) {
846      fhCutStatistics->GetXaxis()->SetBinLabel(i+4,fgkCutNames[i]);
847      fhCutCorrelation->GetXaxis()->SetBinLabel(i+1,fgkCutNames[i]);
848      fhCutCorrelation->GetYaxis()->SetBinLabel(i+1,fgkCutNames[i]);
849    } 
850
851   fhCutStatistics  ->SetLineColor(color);
852   fhCutCorrelation ->SetLineColor(color);
853   fhCutStatistics  ->SetLineWidth(2);
854   fhCutCorrelation ->SetLineWidth(2);
855
856   Char_t str[256];
857   for (Int_t i=0; i<2; i++) {
858     if (i==0) sprintf(str," ");
859     else sprintf(str,"_cut");
860
861     fhNClustersITS[i]        = new TH1F(Form("nClustersITS%s",str)     ,"",8,-0.5,7.5);
862     fhNClustersTPC[i]        = new TH1F(Form("nClustersTPC%s",str)     ,"",165,-0.5,164.5);
863     fhChi2PerClusterITS[i]   = new TH1F(Form("chi2PerClusterITS%s",str),"",500,0,10);
864     fhChi2PerClusterTPC[i]   = new TH1F(Form("chi2PerClusterTPC%s",str),"",500,0,10);
865
866     fhC11[i]                 = new  TH1F(Form("covMatrixDiagonal11%s",str),"",2000,0,20);
867     fhC22[i]                 = new  TH1F(Form("covMatrixDiagonal22%s",str),"",2000,0,20);
868     fhC33[i]                 = new  TH1F(Form("covMatrixDiagonal33%s",str),"",1000,0,1);
869     fhC44[i]                 = new  TH1F(Form("covMatrixDiagonal44%s",str),"",1000,0,5);
870     fhC55[i]                 = new  TH1F(Form("covMatrixDiagonal55%s",str),"",1000,0,5);
871
872     fhDXY[i]                 = new  TH1F(Form("dXY%s",str)    ,"",500,-10,10);
873     fhDZ[i]                  = new  TH1F(Form("dZ%s",str)     ,"",500,-10,10);
874     fhDXYvsDZ[i]             = new  TH2F(Form("dXYvsDZ%s",str),"",200,-10,10,200,-10,10);
875
876     fhDXYNormalized[i]       = new  TH1F(Form("dXYNormalized%s",str)    ,"",500,-10,10);
877     fhDZNormalized[i]        = new  TH1F(Form("dZNormalized%s",str)     ,"",500,-10,10);
878     fhDXYvsDZNormalized[i]   = new  TH2F(Form("dXYvsDZNormalized%s",str),"",200,-10,10,200,-10,10);
879
880     fhNSigmaToVertex[i]      = new  TH1F(Form("nSigmaToVertex%s",str),"",500,0,50);
881
882     fhPt[i]                  = new TH1F(Form("pt%s",str)     ,"p_{T} distribution;p_{T} (GeV/c)",500,0.0,100.0);
883     fhEta[i]                 = new TH1F(Form("eta%s",str)     ,"#eta distribution;#eta",40,-2.0,2.0);
884     
885     fhNClustersITS[i]->SetTitle("n ITS clusters");
886     fhNClustersTPC[i]->SetTitle("n TPC clusters");
887     fhChi2PerClusterITS[i]->SetTitle("#Chi^{2} per ITS cluster");
888     fhChi2PerClusterTPC[i]->SetTitle("#Chi^{2} per TPC cluster");
889
890     fhC11[i]->SetTitle("cov 11 : #sigma_{y}^{2} [cm^{2}]");
891     fhC22[i]->SetTitle("cov 22 : #sigma_{z}^{2} [cm^{2}]");
892     fhC33[i]->SetTitle("cov 33 : #sigma_{sin(#phi)}^{2}");
893     fhC44[i]->SetTitle("cov 44 : #sigma_{tan(#theta_{dip})}^{2}");
894     fhC55[i]->SetTitle("cov 55 : #sigma_{1/p_{T}}^{2} [(c/GeV)^2]");
895
896     fhDXY[i]->SetTitle("transverse impact parameter");
897     fhDZ[i]->SetTitle("longitudinal impact parameter");
898     fhDXYvsDZ[i]->SetTitle("longitudinal impact parameter");
899     fhDXYvsDZ[i]->SetYTitle("transverse impact parameter");
900
901     fhDXYNormalized[i]->SetTitle("normalized trans impact par");
902     fhDZNormalized[i]->SetTitle("normalized long impact par");
903     fhDXYvsDZNormalized[i]->SetTitle("normalized long impact par");
904     fhDXYvsDZNormalized[i]->SetYTitle("normalized trans impact par");
905     fhNSigmaToVertex[i]->SetTitle("n #sigma to vertex");
906
907     fhNClustersITS[i]->SetLineColor(color);   fhNClustersITS[i]->SetLineWidth(2);
908     fhNClustersTPC[i]->SetLineColor(color);   fhNClustersTPC[i]->SetLineWidth(2);
909     fhChi2PerClusterITS[i]->SetLineColor(color);   fhChi2PerClusterITS[i]->SetLineWidth(2);
910     fhChi2PerClusterTPC[i]->SetLineColor(color);   fhChi2PerClusterTPC[i]->SetLineWidth(2);
911
912     fhC11[i]->SetLineColor(color);   fhC11[i]->SetLineWidth(2);
913     fhC22[i]->SetLineColor(color);   fhC22[i]->SetLineWidth(2);
914     fhC33[i]->SetLineColor(color);   fhC33[i]->SetLineWidth(2);
915     fhC44[i]->SetLineColor(color);   fhC44[i]->SetLineWidth(2);
916     fhC55[i]->SetLineColor(color);   fhC55[i]->SetLineWidth(2);
917
918     fhDXY[i]->SetLineColor(color);   fhDXY[i]->SetLineWidth(2);
919     fhDZ[i]->SetLineColor(color);   fhDZ[i]->SetLineWidth(2);
920
921     fhDXYNormalized[i]->SetLineColor(color);   fhDXYNormalized[i]->SetLineWidth(2);
922     fhDZNormalized[i]->SetLineColor(color);    fhDZNormalized[i]->SetLineWidth(2);
923     fhNSigmaToVertex[i]->SetLineColor(color);  fhNSigmaToVertex[i]->SetLineWidth(2); 
924   }
925
926   // The number of sigmas to the vertex is per definition gaussian
927   ffDTheoretical = new TF1("nSigmaToVertexTheoretical","([0]/2.506628274)*exp(-(x**2)/2)",0,50);
928   ffDTheoretical->SetParameter(0,1);
929   
930   TH1::AddDirectory(oldStatus);
931 }
932
933 //____________________________________________________________________
934 Bool_t AliESDtrackCuts::LoadHistograms(const Char_t* dir)
935 {
936   //
937   // loads the histograms from a file
938   // if dir is empty a directory with the name of this object is taken (like in SaveHistogram)
939   //
940
941   if (!dir)
942     dir = GetName();
943
944   if (!gDirectory->cd(dir))
945     return kFALSE;
946
947   ffDTheoretical = dynamic_cast<TF1*> (gDirectory->Get("nSigmaToVertexTheory"));
948
949   fhCutStatistics = dynamic_cast<TH1F*> (gDirectory->Get("cut_statistics"));
950   fhCutCorrelation = dynamic_cast<TH2F*> (gDirectory->Get("cut_correlation"));
951
952   Char_t str[5];
953   for (Int_t i=0; i<2; i++) {
954     if (i==0)
955     {
956       gDirectory->cd("before_cuts");
957       str[0] = 0;
958     }
959     else
960     {
961       gDirectory->cd("after_cuts");
962       sprintf(str,"_cut");
963     }
964
965     fhNClustersITS[i]      = dynamic_cast<TH1F*> (gDirectory->Get(Form("nClustersITS%s",str)     ));
966     fhNClustersTPC[i]      = dynamic_cast<TH1F*> (gDirectory->Get(Form("nClustersTPC%s",str)     ));
967     fhChi2PerClusterITS[i] = dynamic_cast<TH1F*> (gDirectory->Get(Form("chi2PerClusterITS%s",str)));
968     fhChi2PerClusterTPC[i] = dynamic_cast<TH1F*> (gDirectory->Get(Form("chi2PerClusterTPC%s",str)));
969
970     fhC11[i] = dynamic_cast<TH1F*> (gDirectory->Get(Form("covMatrixDiagonal11%s",str)));
971     fhC22[i] = dynamic_cast<TH1F*> (gDirectory->Get(Form("covMatrixDiagonal22%s",str)));
972     fhC33[i] = dynamic_cast<TH1F*> (gDirectory->Get(Form("covMatrixDiagonal33%s",str)));
973     fhC44[i] = dynamic_cast<TH1F*> (gDirectory->Get(Form("covMatrixDiagonal44%s",str)));
974     fhC55[i] = dynamic_cast<TH1F*> (gDirectory->Get(Form("covMatrixDiagonal55%s",str)));
975
976     fhDXY[i] =     dynamic_cast<TH1F*> (gDirectory->Get(Form("dXY%s",str)    ));
977     fhDZ[i] =      dynamic_cast<TH1F*> (gDirectory->Get(Form("dZ%s",str)     ));
978     fhDXYvsDZ[i] = dynamic_cast<TH2F*> (gDirectory->Get(Form("dXYvsDZ%s",str)));
979
980     fhDXYNormalized[i] =     dynamic_cast<TH1F*> (gDirectory->Get(Form("dXYNormalized%s",str)    ));
981     fhDZNormalized[i] =      dynamic_cast<TH1F*> (gDirectory->Get(Form("dZNormalized%s",str)     ));
982     fhDXYvsDZNormalized[i] = dynamic_cast<TH2F*> (gDirectory->Get(Form("dXYvsDZNormalized%s",str)));
983     fhNSigmaToVertex[i] = dynamic_cast<TH1F*> (gDirectory->Get(Form("nSigmaToVertex%s",str)));
984
985     fhPt[i] = dynamic_cast<TH1F*> (gDirectory->Get(Form("pt%s",str)));
986     fhEta[i] = dynamic_cast<TH1F*> (gDirectory->Get(Form("eta%s",str)));
987
988     gDirectory->cd("../");
989   }
990
991   gDirectory->cd("..");
992
993   return kTRUE;
994 }
995
996 //____________________________________________________________________
997 void AliESDtrackCuts::SaveHistograms(const Char_t* dir) {
998   //
999   // saves the histograms in a directory (dir)
1000   //
1001
1002   if (!fHistogramsOn) {
1003     AliDebug(0, "Histograms not on - cannot save histograms!!!");
1004     return;
1005   }
1006
1007   if (!dir)
1008     dir = GetName();
1009
1010   gDirectory->mkdir(dir);
1011   gDirectory->cd(dir);
1012
1013   gDirectory->mkdir("before_cuts");
1014   gDirectory->mkdir("after_cuts");
1015  
1016   // a factor of 2 is needed since n sigma is positive
1017   ffDTheoretical->SetParameter(0,2*fhNSigmaToVertex[0]->Integral("width"));
1018   ffDTheoretical->Write("nSigmaToVertexTheory");
1019
1020   fhCutStatistics->Write();
1021   fhCutCorrelation->Write();
1022
1023   for (Int_t i=0; i<2; i++) {
1024     if (i==0)
1025       gDirectory->cd("before_cuts");
1026     else
1027       gDirectory->cd("after_cuts");
1028
1029     fhNClustersITS[i]        ->Write();
1030     fhNClustersTPC[i]        ->Write();
1031     fhChi2PerClusterITS[i]   ->Write();
1032     fhChi2PerClusterTPC[i]   ->Write();
1033
1034     fhC11[i]                 ->Write();
1035     fhC22[i]                 ->Write();
1036     fhC33[i]                 ->Write();
1037     fhC44[i]                 ->Write();
1038     fhC55[i]                 ->Write();
1039
1040     fhDXY[i]                 ->Write();
1041     fhDZ[i]                  ->Write();
1042     fhDXYvsDZ[i]             ->Write();
1043
1044     fhDXYNormalized[i]       ->Write();
1045     fhDZNormalized[i]        ->Write();
1046     fhDXYvsDZNormalized[i]   ->Write();
1047     fhNSigmaToVertex[i]      ->Write();
1048
1049     fhPt[i]                  ->Write();
1050     fhEta[i]                 ->Write();
1051     
1052     gDirectory->cd("../");
1053   }
1054
1055   gDirectory->cd("../");
1056 }
1057
1058 //____________________________________________________________________
1059 void AliESDtrackCuts::DrawHistograms()
1060 {
1061   // draws some histograms
1062
1063   TCanvas* canvas1 = new TCanvas(Form("%s_1", GetName()), "Track Quality Results1", 800, 800);
1064   canvas1->Divide(2, 2);
1065
1066   canvas1->cd(1);
1067   fhNClustersTPC[0]->SetStats(kFALSE);
1068   fhNClustersTPC[0]->Draw();
1069
1070   canvas1->cd(2);
1071   fhChi2PerClusterTPC[0]->SetStats(kFALSE);
1072   fhChi2PerClusterTPC[0]->Draw();
1073
1074   canvas1->cd(3);
1075   fhNSigmaToVertex[0]->SetStats(kFALSE);
1076   fhNSigmaToVertex[0]->GetXaxis()->SetRangeUser(0, 10);
1077   fhNSigmaToVertex[0]->Draw();
1078
1079   canvas1->SaveAs(Form("%s_%s.gif", GetName(), canvas1->GetName()));
1080
1081   TCanvas* canvas2 = new TCanvas(Form("%s_2", GetName()), "Track Quality Results2", 1200, 800);
1082   canvas2->Divide(3, 2);
1083
1084   canvas2->cd(1);
1085   fhC11[0]->SetStats(kFALSE);
1086   gPad->SetLogy();
1087   fhC11[0]->Draw();
1088
1089   canvas2->cd(2);
1090   fhC22[0]->SetStats(kFALSE);
1091   gPad->SetLogy();
1092   fhC22[0]->Draw();
1093
1094   canvas2->cd(3);
1095   fhC33[0]->SetStats(kFALSE);
1096   gPad->SetLogy();
1097   fhC33[0]->Draw();
1098
1099   canvas2->cd(4);
1100   fhC44[0]->SetStats(kFALSE);
1101   gPad->SetLogy();
1102   fhC44[0]->Draw();
1103
1104   canvas2->cd(5);
1105   fhC55[0]->SetStats(kFALSE);
1106   gPad->SetLogy();
1107   fhC55[0]->Draw();
1108
1109   canvas2->SaveAs(Form("%s_%s.gif", GetName(), canvas2->GetName()));
1110
1111   TCanvas* canvas3 = new TCanvas(Form("%s_3", GetName()), "Track Quality Results3", 1200, 800);
1112   canvas3->Divide(3, 2);
1113
1114   canvas3->cd(1);
1115   fhDXY[0]->SetStats(kFALSE);
1116   gPad->SetLogy();
1117   fhDXY[0]->Draw();
1118
1119   canvas3->cd(2);
1120   fhDZ[0]->SetStats(kFALSE);
1121   gPad->SetLogy();
1122   fhDZ[0]->Draw();
1123
1124   canvas3->cd(3);
1125   fhDXYvsDZ[0]->SetStats(kFALSE);
1126   gPad->SetLogz();
1127   gPad->SetRightMargin(0.15);
1128   fhDXYvsDZ[0]->Draw("COLZ");
1129
1130   canvas3->cd(4);
1131   fhDXYNormalized[0]->SetStats(kFALSE);
1132   gPad->SetLogy();
1133   fhDXYNormalized[0]->Draw();
1134
1135   canvas3->cd(5);
1136   fhDZNormalized[0]->SetStats(kFALSE);
1137   gPad->SetLogy();
1138   fhDZNormalized[0]->Draw();
1139
1140   canvas3->cd(6);
1141   fhDXYvsDZNormalized[0]->SetStats(kFALSE);
1142   gPad->SetLogz();
1143   gPad->SetRightMargin(0.15);
1144   fhDXYvsDZNormalized[0]->Draw("COLZ");
1145
1146   canvas3->SaveAs(Form("%s_%s.gif", GetName(), canvas3->GetName()));
1147
1148   TCanvas* canvas4 = new TCanvas(Form("%s_4", GetName()), "Track Quality Results4", 800, 500);
1149   canvas4->Divide(2, 1);
1150
1151   canvas4->cd(1);
1152   fhCutStatistics->SetStats(kFALSE);
1153   fhCutStatistics->LabelsOption("v");
1154   gPad->SetBottomMargin(0.3);
1155   fhCutStatistics->Draw();
1156
1157   canvas4->cd(2);
1158   fhCutCorrelation->SetStats(kFALSE);
1159   fhCutCorrelation->LabelsOption("v");
1160   gPad->SetBottomMargin(0.3);
1161   gPad->SetLeftMargin(0.3);
1162   fhCutCorrelation->Draw("COLZ");
1163
1164   canvas4->SaveAs(Form("%s_%s.gif", GetName(), canvas4->GetName()));
1165
1166   /*canvas->cd(1);
1167   fhDXYvsDZNormalized[0]->SetStats(kFALSE);
1168   fhDXYvsDZNormalized[0]->DrawCopy("COLZ");
1169
1170   canvas->cd(2);
1171   fhNClustersTPC[0]->SetStats(kFALSE);
1172   fhNClustersTPC[0]->DrawCopy();
1173
1174   canvas->cd(3);
1175   fhChi2PerClusterITS[0]->SetStats(kFALSE);
1176   fhChi2PerClusterITS[0]->DrawCopy();
1177   fhChi2PerClusterITS[1]->SetLineColor(2);
1178   fhChi2PerClusterITS[1]->DrawCopy("SAME");
1179
1180   canvas->cd(4);
1181   fhChi2PerClusterTPC[0]->SetStats(kFALSE);
1182   fhChi2PerClusterTPC[0]->DrawCopy();
1183   fhChi2PerClusterTPC[1]->SetLineColor(2);
1184   fhChi2PerClusterTPC[1]->DrawCopy("SAME");*/
1185 }
1186