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