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