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