]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ANALYSIS/AliESDtrackCuts.cxx
Fixed signatures of some methods
[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   // stupid rounding problem screws up everything:
496   // if d is too big, TMath::Exp(...) gets 0, and TMath::ErfInverse(1) that should be infinite, gets 0 :(
497   if (TMath::Exp(-d * d / 2) < 1e-10)
498     return 1000;
499
500   d = TMath::ErfInverse(1 - TMath::Exp(-d * d / 2)) * TMath::Sqrt(2);
501   return d;
502 }
503
504 void AliESDtrackCuts::EnableNeededBranches(TTree* tree)
505 {
506   // enables the branches needed by AcceptTrack, for a list see comment of AcceptTrack
507
508   tree->SetBranchStatus("fTracks.fFlags", 1);
509   tree->SetBranchStatus("fTracks.fITSncls", 1);
510   tree->SetBranchStatus("fTracks.fTPCncls", 1);
511   tree->SetBranchStatus("fTracks.fITSchi2", 1);
512   tree->SetBranchStatus("fTracks.fTPCchi2", 1);
513   tree->SetBranchStatus("fTracks.fC*", 1);
514   tree->SetBranchStatus("fTracks.fD", 1);
515   tree->SetBranchStatus("fTracks.fZ", 1);
516   tree->SetBranchStatus("fTracks.fCdd", 1);
517   tree->SetBranchStatus("fTracks.fCdz", 1);
518   tree->SetBranchStatus("fTracks.fCzz", 1);
519   tree->SetBranchStatus("fTracks.fP*", 1);
520   tree->SetBranchStatus("fTracks.fR*", 1);
521   tree->SetBranchStatus("fTracks.fKinkIndexes*", 1);
522 }
523
524 //____________________________________________________________________
525 Bool_t
526 AliESDtrackCuts::AcceptTrack(AliESDtrack* esdTrack) {
527   // 
528   // figure out if the tracks survives all the track cuts defined
529   //
530   // the different quality parameter and kinematic values are first
531   // retrieved from the track. then it is found out what cuts the
532   // track did not survive and finally the cuts are imposed.
533
534   // this function needs the following branches:
535   // fTracks.fFlags
536   // fTracks.fITSncls
537   // fTracks.fTPCncls
538   // fTracks.fITSchi2
539   // fTracks.fTPCchi2
540   // fTracks.fC   //GetExternalCovariance
541   // fTracks.fD   //GetImpactParameters
542   // fTracks.fZ   //GetImpactParameters
543   // fTracks.fCdd //GetImpactParameters
544   // fTracks.fCdz //GetImpactParameters
545   // fTracks.fCzz //GetImpactParameters
546   // fTracks.fP   //GetPxPyPz
547   // fTracks.fR   //GetMass
548   // fTracks.fP   //GetMass
549   // fTracks.fKinkIndexes
550
551
552   UInt_t status = esdTrack->GetStatus();
553
554   // getting quality parameters from the ESD track
555   Int_t nClustersITS = esdTrack->GetITSclusters(0);
556   Int_t nClustersTPC = esdTrack->GetTPCclusters(0);
557
558   Float_t chi2PerClusterITS = -1;
559   Float_t chi2PerClusterTPC = -1;
560   if (nClustersITS!=0)
561     chi2PerClusterITS = esdTrack->GetITSchi2()/Float_t(nClustersITS);
562   if (nClustersTPC!=0)
563     chi2PerClusterTPC = esdTrack->GetTPCchi2()/Float_t(nClustersTPC);
564   Double_t extCov[15];
565   esdTrack->GetExternalCovariance(extCov);
566
567   // getting the track to vertex parameters
568   Float_t nSigmaToVertex = GetSigmaToVertex(esdTrack);
569       
570   Float_t b[2];
571   Float_t bCov[3];
572   esdTrack->GetImpactParameters(b,bCov);
573   if (bCov[0]<=0 || bCov[2]<=0) {
574     AliDebug(1, "Estimated b resolution lower or equal zero!");
575     bCov[0]=0; bCov[2]=0;
576   }
577   Float_t dcaToVertex   = TMath::Sqrt(b[0]*b[0] + b[1]*b[1]);
578  
579   Float_t dcaToVertexXY = b[0];
580   
581   // getting the kinematic variables of the track
582   // (assuming the mass is known)
583   Double_t p[3];
584   esdTrack->GetPxPyPz(p);
585
586   Float_t momentum = TMath::Sqrt(TMath::Power(p[0],2) + TMath::Power(p[1],2) + TMath::Power(p[2],2));
587   Float_t pt       = TMath::Sqrt(TMath::Power(p[0],2) + TMath::Power(p[1],2));
588   Float_t energy   = TMath::Sqrt(TMath::Power(esdTrack->GetMass(),2) + TMath::Power(momentum,2));
589
590
591   //y-eta related calculations
592   Float_t eta = -100.;
593   Float_t y   = -100.;
594   if((momentum != TMath::Abs(p[2]))&&(momentum != 0))
595     eta = 0.5*TMath::Log((momentum + p[2])/(momentum - p[2]));
596   if((energy != TMath::Abs(p[2]))&&(momentum != 0))
597     y = 0.5*TMath::Log((energy + p[2])/(energy - p[2]));
598
599   
600   //########################################################################
601   // cut the track?
602   
603   Bool_t cuts[kNCuts];
604   for (Int_t i=0; i<kNCuts; i++) cuts[i]=kFALSE;
605   
606   // track quality cuts
607   if (fCutRequireTPCRefit && (status&AliESDtrack::kTPCrefit)==0)
608     cuts[0]=kTRUE;
609   if (fCutRequireITSRefit && (status&AliESDtrack::kITSrefit)==0)
610     cuts[1]=kTRUE;
611   if (nClustersTPC<fCutMinNClusterTPC)
612     cuts[2]=kTRUE;
613   if (nClustersITS<fCutMinNClusterITS) 
614     cuts[3]=kTRUE;
615   if (chi2PerClusterTPC>fCutMaxChi2PerClusterTPC) 
616     cuts[4]=kTRUE; 
617   if (chi2PerClusterITS>fCutMaxChi2PerClusterITS) 
618     cuts[5]=kTRUE;
619   if (extCov[0]  > fCutMaxC11) 
620     cuts[6]=kTRUE;  
621   if (extCov[2]  > fCutMaxC22) 
622     cuts[7]=kTRUE;  
623   if (extCov[5]  > fCutMaxC33) 
624     cuts[8]=kTRUE;  
625   if (extCov[9]  > fCutMaxC44) 
626     cuts[9]=kTRUE;  
627   if (extCov[14]  > fCutMaxC55) 
628     cuts[10]=kTRUE;  
629   if (nSigmaToVertex > fCutNsigmaToVertex && fCutSigmaToVertexRequired)
630     cuts[11] = kTRUE;
631   // if n sigma could not be calculated
632   if (nSigmaToVertex<0 && fCutSigmaToVertexRequired)
633     cuts[12]=kTRUE;
634   if (!fCutAcceptKinkDaughters && esdTrack->GetKinkIndex(0)>0)
635     cuts[13]=kTRUE;
636   // track kinematics cut
637   if((momentum < fPMin) || (momentum > fPMax)) 
638     cuts[14]=kTRUE;
639   if((pt < fPtMin) || (pt > fPtMax)) 
640     cuts[15] = kTRUE;
641   if((p[0] < fPxMin) || (p[0] > fPxMax)) 
642     cuts[16] = kTRUE;
643   if((p[1] < fPyMin) || (p[1] > fPyMax)) 
644     cuts[17] = kTRUE;
645   if((p[2] < fPzMin) || (p[2] > fPzMax))
646     cuts[18] = kTRUE;
647   if((eta < fEtaMin) || (eta > fEtaMax)) 
648     cuts[19] = kTRUE;
649   if((y < fRapMin) || (y > fRapMax)) 
650     cuts[20] = kTRUE;
651   if (dcaToVertex > fCutDCAToVertex)
652     cuts[21] = kTRUE;
653   if (dcaToVertexXY > fCutDCAToVertexXY)
654     cuts[22] = kTRUE;
655
656   Bool_t cut=kFALSE;
657   for (Int_t i=0; i<kNCuts; i++) 
658     if (cuts[i]) cut = kTRUE;
659
660
661
662   //########################################################################
663   // filling histograms
664   if (fHistogramsOn) {
665     fhCutStatistics->Fill(fhCutStatistics->GetBinCenter(fhCutStatistics->GetXaxis()->FindBin("n tracks")));
666     if (cut)
667       fhCutStatistics->Fill(fhCutStatistics->GetBinCenter(fhCutStatistics->GetXaxis()->FindBin("n cut tracks")));
668     
669     for (Int_t i=0; i<kNCuts; i++) {
670       if (cuts[i])
671         fhCutStatistics->Fill(fhCutStatistics->GetBinCenter(fhCutStatistics->GetXaxis()->FindBin(fgkCutNames[i])));
672
673       for (Int_t j=i; j<kNCuts; j++) {
674         if (cuts[i] && cuts[j]) {
675           Float_t xC = fhCutCorrelation->GetXaxis()->GetBinCenter(fhCutCorrelation->GetXaxis()->FindBin(fgkCutNames[i]));
676           Float_t yC = fhCutCorrelation->GetYaxis()->GetBinCenter(fhCutCorrelation->GetYaxis()->FindBin(fgkCutNames[j]));
677           fhCutCorrelation->Fill(xC, yC);
678         }
679       }
680     }
681   }
682
683   // now we loop over the filling of the histograms twice: once "before" the cut, once "after"
684   // the code is not in a function due to too many local variables that would need to be passed
685
686   for (Int_t id = 0; id < 2; id++)
687   {
688     // id = 0 --> before cut
689     // id = 1 --> after cut
690
691     if (fHistogramsOn)
692     {
693       fhNClustersITS[id]->Fill(nClustersITS);
694       fhNClustersTPC[id]->Fill(nClustersTPC);
695       fhChi2PerClusterITS[id]->Fill(chi2PerClusterITS);
696       fhChi2PerClusterTPC[id]->Fill(chi2PerClusterTPC);
697
698       fhC11[id]->Fill(extCov[0]);
699       fhC22[id]->Fill(extCov[2]);
700       fhC33[id]->Fill(extCov[5]);
701       fhC44[id]->Fill(extCov[9]);
702       fhC55[id]->Fill(extCov[14]);
703
704       fhPt[id]->Fill(pt);
705       fhEta[id]->Fill(eta);
706
707       Float_t bRes[2];
708       bRes[0] = TMath::Sqrt(bCov[0]);
709       bRes[1] = TMath::Sqrt(bCov[2]);
710
711       fhDZ[id]->Fill(b[1]);
712       fhDXY[id]->Fill(b[0]);
713       fhDXYDZ[id]->Fill(dcaToVertex);
714       fhDXYvsDZ[id]->Fill(b[1],b[0]);
715
716       if (bRes[0]!=0 && bRes[1]!=0) {
717         fhDZNormalized[id]->Fill(b[1]/bRes[1]);
718         fhDXYNormalized[id]->Fill(b[0]/bRes[0]);
719         fhDXYvsDZNormalized[id]->Fill(b[1]/bRes[1], b[0]/bRes[0]);
720         fhNSigmaToVertex[id]->Fill(nSigmaToVertex);
721       }
722     }
723
724     // cut the track
725     if (cut)
726       return kFALSE;
727   }
728
729   return kTRUE;
730 }
731
732 //____________________________________________________________________
733 AliESDtrack* AliESDtrackCuts::GetTPCOnlyTrack(AliESDEvent* esd, Int_t iTrack)
734 {
735   // creates a TPC only track from the given esd track
736   // the track has to be deleted by the user
737   //
738   // NB. most of the functionality to get a TPC only track from an ESD track is in AliESDtrack, where it should be
739   // there are only missing propagations here that are needed for old data
740   // this function will therefore become obsolete
741   //
742   // adapted from code provided by CKB
743
744   if (!esd->GetPrimaryVertexTPC())
745     return 0; // No TPC vertex no TPC tracks
746
747   AliESDtrack* track = esd->GetTrack(iTrack);
748   if (!track)
749     return 0;
750
751   AliESDtrack *tpcTrack = new AliESDtrack();
752
753   // This should have been done during the reconstruction
754   // fixed by Juri in r26675
755   // but recalculate for older data CKB
756   Float_t p[2],cov[3];
757   track->GetImpactParametersTPC(p,cov);
758   if(p[0]==0&&p[1]==0)
759     track->RelateToVertexTPC(esd->GetPrimaryVertexTPC(),esd->GetMagneticField(),kVeryBig);
760   // BKC
761
762   // only true if we have a tpc track
763   if (!track->FillTPCOnlyTrack(*tpcTrack))
764   {
765     delete tpcTrack;
766     return 0;
767   }
768
769   // propagate to Vertex
770   // not needed for normal reconstructed ESDs...
771   // Double_t pTPC[2],covTPC[3];
772   // tpcTrack->PropagateToDCA(esd->GetPrimaryVertexTPC(), esd->GetMagneticField(), 10000,  pTPC, covTPC);
773
774   return tpcTrack;
775 }
776
777 //____________________________________________________________________
778 TObjArray* AliESDtrackCuts::GetAcceptedTracks(AliESDEvent* esd,Bool_t bTPC)
779 {
780   //
781   // returns an array of all tracks that pass the cuts
782   // or an array of TPC only tracks (propagated to the TPC vertex during reco)
783   // tracks that pass the cut
784
785   TObjArray* acceptedTracks = new TObjArray();
786
787   // loop over esd tracks
788   for (Int_t iTrack = 0; iTrack < esd->GetNumberOfTracks(); iTrack++) {
789     if(bTPC){
790       if(!esd->GetPrimaryVertexTPC())return acceptedTracks; // No TPC vertex no TPC tracks
791
792       AliESDtrack *tpcTrack = GetTPCOnlyTrack(esd, iTrack);
793       if (!tpcTrack)
794         continue;
795
796       if (AcceptTrack(tpcTrack)) {
797         acceptedTracks->Add(tpcTrack);
798       }
799       else
800         delete tpcTrack;
801     }
802     else
803     {
804       AliESDtrack* track = esd->GetTrack(iTrack);
805       if(AcceptTrack(track))
806         acceptedTracks->Add(track);
807     }
808   } 
809   if(bTPC)acceptedTracks->SetOwner(kTRUE);
810   return acceptedTracks;
811 }
812
813 //____________________________________________________________________
814 Int_t AliESDtrackCuts::CountAcceptedTracks(AliESDEvent* esd)
815 {
816   //
817   // returns an the number of tracks that pass the cuts
818   //
819
820   Int_t count = 0;
821
822   // loop over esd tracks
823   for (Int_t iTrack = 0; iTrack < esd->GetNumberOfTracks(); iTrack++) {
824     AliESDtrack* track = esd->GetTrack(iTrack);
825     if (AcceptTrack(track))
826       count++;
827   }
828
829   return count;
830 }
831
832 //____________________________________________________________________
833  void AliESDtrackCuts::DefineHistograms(Int_t color) {
834    // 
835    // diagnostics histograms are defined
836    // 
837
838    fHistogramsOn=kTRUE;
839    
840    Bool_t oldStatus = TH1::AddDirectoryStatus();
841    TH1::AddDirectory(kFALSE);
842    
843    //###################################################################################
844    // defining histograms
845
846    fhCutStatistics = new TH1F("cut_statistics","cut statistics",kNCuts+4,-0.5,kNCuts+3.5);
847
848    fhCutStatistics->GetXaxis()->SetBinLabel(1,"n tracks");
849    fhCutStatistics->GetXaxis()->SetBinLabel(2,"n cut tracks");
850
851    fhCutCorrelation = new TH2F("cut_correlation","cut correlation",kNCuts,-0.5,kNCuts-0.5,kNCuts,-0.5,kNCuts-0.5);;
852   
853    for (Int_t i=0; i<kNCuts; i++) {
854      fhCutStatistics->GetXaxis()->SetBinLabel(i+4,fgkCutNames[i]);
855      fhCutCorrelation->GetXaxis()->SetBinLabel(i+1,fgkCutNames[i]);
856      fhCutCorrelation->GetYaxis()->SetBinLabel(i+1,fgkCutNames[i]);
857    } 
858
859   fhCutStatistics  ->SetLineColor(color);
860   fhCutCorrelation ->SetLineColor(color);
861   fhCutStatistics  ->SetLineWidth(2);
862   fhCutCorrelation ->SetLineWidth(2);
863
864   for (Int_t i=0; i<2; i++) {
865     fhNClustersITS[i]        = new TH1F("nClustersITS"    ,"",8,-0.5,7.5);
866     fhNClustersTPC[i]        = new TH1F("nClustersTPC"     ,"",165,-0.5,164.5);
867     fhChi2PerClusterITS[i]   = new TH1F("chi2PerClusterITS","",500,0,10);
868     fhChi2PerClusterTPC[i]   = new TH1F("chi2PerClusterTPC","",500,0,10);
869
870     fhC11[i]                 = new TH1F("covMatrixDiagonal11","",2000,0,20);
871     fhC22[i]                 = new TH1F("covMatrixDiagonal22","",2000,0,20);
872     fhC33[i]                 = new TH1F("covMatrixDiagonal33","",1000,0,1);
873     fhC44[i]                 = new TH1F("covMatrixDiagonal44","",1000,0,5);
874     fhC55[i]                 = new TH1F("covMatrixDiagonal55","",1000,0,5);
875
876     fhDXY[i]                 = new TH1F("dXY"    ,"",500,-10,10);
877     fhDZ[i]                  = new TH1F("dZ"     ,"",500,-10,10);
878     fhDXYDZ[i]               = new TH1F("dXYDZ"  ,"",500,0,10);
879     fhDXYvsDZ[i]             = new TH2F("dXYvsDZ","",200,-10,10,200,-10,10);
880
881     fhDXYNormalized[i]       = new TH1F("dXYNormalized"    ,"",500,-10,10);
882     fhDZNormalized[i]        = new TH1F("dZNormalized"     ,"",500,-10,10);
883     fhDXYvsDZNormalized[i]   = new TH2F("dXYvsDZNormalized","",200,-10,10,200,-10,10);
884
885     fhNSigmaToVertex[i]      = new TH1F("nSigmaToVertex","",500,0,10);
886
887     fhPt[i]                  = new TH1F("pt"     ,"p_{T} distribution;p_{T} (GeV/c)",500,0.0,100.0);
888     fhEta[i]                 = new TH1F("eta"     ,"#eta distribution;#eta",40,-2.0,2.0);
889     
890     fhNClustersITS[i]->SetTitle("n ITS clusters");
891     fhNClustersTPC[i]->SetTitle("n TPC clusters");
892     fhChi2PerClusterITS[i]->SetTitle("#Chi^{2} per ITS cluster");
893     fhChi2PerClusterTPC[i]->SetTitle("#Chi^{2} per TPC cluster");
894
895     fhC11[i]->SetTitle("cov 11 : #sigma_{y}^{2} [cm^{2}]");
896     fhC22[i]->SetTitle("cov 22 : #sigma_{z}^{2} [cm^{2}]");
897     fhC33[i]->SetTitle("cov 33 : #sigma_{sin(#phi)}^{2}");
898     fhC44[i]->SetTitle("cov 44 : #sigma_{tan(#theta_{dip})}^{2}");
899     fhC55[i]->SetTitle("cov 55 : #sigma_{1/p_{T}}^{2} [(c/GeV)^2]");
900
901     fhDXY[i]->SetTitle("transverse impact parameter");
902     fhDZ[i]->SetTitle("longitudinal impact parameter");
903     fhDXYDZ[i]->SetTitle("absolute impact parameter;sqrt(dXY**2 + dZ**2) in cm");
904     fhDXYvsDZ[i]->SetXTitle("longitudinal impact parameter");
905     fhDXYvsDZ[i]->SetYTitle("transverse impact parameter");
906
907     fhDXYNormalized[i]->SetTitle("normalized trans impact par");
908     fhDZNormalized[i]->SetTitle("normalized long impact par");
909     fhDXYvsDZNormalized[i]->SetTitle("normalized long impact par");
910     fhDXYvsDZNormalized[i]->SetYTitle("normalized trans impact par");
911     fhNSigmaToVertex[i]->SetTitle("n #sigma to vertex");
912
913     fhNClustersITS[i]->SetLineColor(color);   fhNClustersITS[i]->SetLineWidth(2);
914     fhNClustersTPC[i]->SetLineColor(color);   fhNClustersTPC[i]->SetLineWidth(2);
915     fhChi2PerClusterITS[i]->SetLineColor(color);   fhChi2PerClusterITS[i]->SetLineWidth(2);
916     fhChi2PerClusterTPC[i]->SetLineColor(color);   fhChi2PerClusterTPC[i]->SetLineWidth(2);
917
918     fhC11[i]->SetLineColor(color);   fhC11[i]->SetLineWidth(2);
919     fhC22[i]->SetLineColor(color);   fhC22[i]->SetLineWidth(2);
920     fhC33[i]->SetLineColor(color);   fhC33[i]->SetLineWidth(2);
921     fhC44[i]->SetLineColor(color);   fhC44[i]->SetLineWidth(2);
922     fhC55[i]->SetLineColor(color);   fhC55[i]->SetLineWidth(2);
923
924     fhDXY[i]->SetLineColor(color);   fhDXY[i]->SetLineWidth(2);
925     fhDZ[i]->SetLineColor(color);    fhDZ[i]->SetLineWidth(2);
926     fhDXYDZ[i]->SetLineColor(color); fhDXYDZ[i]->SetLineWidth(2);
927
928     fhDXYNormalized[i]->SetLineColor(color);   fhDXYNormalized[i]->SetLineWidth(2);
929     fhDZNormalized[i]->SetLineColor(color);    fhDZNormalized[i]->SetLineWidth(2);
930     fhNSigmaToVertex[i]->SetLineColor(color);  fhNSigmaToVertex[i]->SetLineWidth(2); 
931   }
932
933   // The number of sigmas to the vertex is per definition gaussian
934   ffDTheoretical = new TF1("nSigmaToVertexTheoretical","([0]/2.506628274)*exp(-(x**2)/2)",0,50);
935   ffDTheoretical->SetParameter(0,1);
936
937   TH1::AddDirectory(oldStatus);
938 }
939
940 //____________________________________________________________________
941 Bool_t AliESDtrackCuts::LoadHistograms(const Char_t* dir)
942 {
943   //
944   // loads the histograms from a file
945   // if dir is empty a directory with the name of this object is taken (like in SaveHistogram)
946   //
947
948   if (!dir)
949     dir = GetName();
950
951   if (!gDirectory->cd(dir))
952     return kFALSE;
953
954   ffDTheoretical = dynamic_cast<TF1*> (gDirectory->Get("nSigmaToVertexTheory"));
955
956   fhCutStatistics = dynamic_cast<TH1F*> (gDirectory->Get("cut_statistics"));
957   fhCutCorrelation = dynamic_cast<TH2F*> (gDirectory->Get("cut_correlation"));
958
959   for (Int_t i=0; i<2; i++) {
960     if (i==0)
961     {
962       gDirectory->cd("before_cuts");
963     }
964     else
965       gDirectory->cd("after_cuts");
966
967     fhNClustersITS[i]      = dynamic_cast<TH1F*> (gDirectory->Get("nClustersITS"     ));
968     fhNClustersTPC[i]      = dynamic_cast<TH1F*> (gDirectory->Get("nClustersTPC"     ));
969     fhChi2PerClusterITS[i] = dynamic_cast<TH1F*> (gDirectory->Get("chi2PerClusterITS"));
970     fhChi2PerClusterTPC[i] = dynamic_cast<TH1F*> (gDirectory->Get("chi2PerClusterTPC"));
971
972     fhC11[i] = dynamic_cast<TH1F*> (gDirectory->Get("covMatrixDiagonal11"));
973     fhC22[i] = dynamic_cast<TH1F*> (gDirectory->Get("covMatrixDiagonal22"));
974     fhC33[i] = dynamic_cast<TH1F*> (gDirectory->Get("covMatrixDiagonal33"));
975     fhC44[i] = dynamic_cast<TH1F*> (gDirectory->Get("covMatrixDiagonal44"));
976     fhC55[i] = dynamic_cast<TH1F*> (gDirectory->Get("covMatrixDiagonal55"));
977
978     fhDXY[i] =     dynamic_cast<TH1F*> (gDirectory->Get("dXY"    ));
979     fhDZ[i] =      dynamic_cast<TH1F*> (gDirectory->Get("dZ"     ));
980     fhDXYDZ[i] =   dynamic_cast<TH1F*> (gDirectory->Get("dXYDZ"));
981     fhDXYvsDZ[i] = dynamic_cast<TH2F*> (gDirectory->Get("dXYvsDZ"));
982
983     fhDXYNormalized[i] =     dynamic_cast<TH1F*> (gDirectory->Get("dXYNormalized"    ));
984     fhDZNormalized[i] =      dynamic_cast<TH1F*> (gDirectory->Get("dZNormalized"     ));
985     fhDXYvsDZNormalized[i] = dynamic_cast<TH2F*> (gDirectory->Get("dXYvsDZNormalized"));
986     fhNSigmaToVertex[i] = dynamic_cast<TH1F*> (gDirectory->Get("nSigmaToVertex"));
987
988     fhPt[i] = dynamic_cast<TH1F*> (gDirectory->Get("pt"));
989     fhEta[i] = dynamic_cast<TH1F*> (gDirectory->Get("eta"));
990
991     gDirectory->cd("../");
992   }
993
994   gDirectory->cd("..");
995
996   return kTRUE;
997 }
998
999 //____________________________________________________________________
1000 void AliESDtrackCuts::SaveHistograms(const Char_t* dir) {
1001   //
1002   // saves the histograms in a directory (dir)
1003   //
1004
1005   if (!fHistogramsOn) {
1006     AliDebug(0, "Histograms not on - cannot save histograms!!!");
1007     return;
1008   }
1009
1010   if (!dir)
1011     dir = GetName();
1012
1013   gDirectory->mkdir(dir);
1014   gDirectory->cd(dir);
1015
1016   gDirectory->mkdir("before_cuts");
1017   gDirectory->mkdir("after_cuts");
1018  
1019   // a factor of 2 is needed since n sigma is positive
1020   ffDTheoretical->SetParameter(0,2*fhNSigmaToVertex[0]->Integral("width"));
1021   ffDTheoretical->Write("nSigmaToVertexTheory");
1022
1023   fhCutStatistics->Write();
1024   fhCutCorrelation->Write();
1025
1026   for (Int_t i=0; i<2; i++) {
1027     if (i==0)
1028       gDirectory->cd("before_cuts");
1029     else
1030       gDirectory->cd("after_cuts");
1031
1032     fhNClustersITS[i]        ->Write();
1033     fhNClustersTPC[i]        ->Write();
1034     fhChi2PerClusterITS[i]   ->Write();
1035     fhChi2PerClusterTPC[i]   ->Write();
1036
1037     fhC11[i]                 ->Write();
1038     fhC22[i]                 ->Write();
1039     fhC33[i]                 ->Write();
1040     fhC44[i]                 ->Write();
1041     fhC55[i]                 ->Write();
1042
1043     fhDXY[i]                 ->Write();
1044     fhDZ[i]                  ->Write();
1045     fhDXYDZ[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     fhPt[i]                  ->Write();
1054     fhEta[i]                 ->Write();
1055     
1056     gDirectory->cd("../");
1057   }
1058
1059   gDirectory->cd("../");
1060 }
1061
1062 //____________________________________________________________________
1063 void AliESDtrackCuts::DrawHistograms()
1064 {
1065   // draws some histograms
1066
1067   TCanvas* canvas1 = new TCanvas(Form("%s_1", GetName()), "Track Quality Results1", 800, 800);
1068   canvas1->Divide(2, 2);
1069
1070   canvas1->cd(1);
1071   fhNClustersTPC[0]->SetStats(kFALSE);
1072   fhNClustersTPC[0]->Draw();
1073
1074   canvas1->cd(2);
1075   fhChi2PerClusterTPC[0]->SetStats(kFALSE);
1076   fhChi2PerClusterTPC[0]->Draw();
1077
1078   canvas1->cd(3);
1079   fhNSigmaToVertex[0]->SetStats(kFALSE);
1080   fhNSigmaToVertex[0]->GetXaxis()->SetRangeUser(0, 10);
1081   fhNSigmaToVertex[0]->Draw();
1082
1083   canvas1->SaveAs(Form("%s_%s.gif", GetName(), canvas1->GetName()));
1084
1085   TCanvas* canvas2 = new TCanvas(Form("%s_2", GetName()), "Track Quality Results2", 1200, 800);
1086   canvas2->Divide(3, 2);
1087
1088   canvas2->cd(1);
1089   fhC11[0]->SetStats(kFALSE);
1090   gPad->SetLogy();
1091   fhC11[0]->Draw();
1092
1093   canvas2->cd(2);
1094   fhC22[0]->SetStats(kFALSE);
1095   gPad->SetLogy();
1096   fhC22[0]->Draw();
1097
1098   canvas2->cd(3);
1099   fhC33[0]->SetStats(kFALSE);
1100   gPad->SetLogy();
1101   fhC33[0]->Draw();
1102
1103   canvas2->cd(4);
1104   fhC44[0]->SetStats(kFALSE);
1105   gPad->SetLogy();
1106   fhC44[0]->Draw();
1107
1108   canvas2->cd(5);
1109   fhC55[0]->SetStats(kFALSE);
1110   gPad->SetLogy();
1111   fhC55[0]->Draw();
1112
1113   canvas2->SaveAs(Form("%s_%s.gif", GetName(), canvas2->GetName()));
1114
1115   TCanvas* canvas3 = new TCanvas(Form("%s_3", GetName()), "Track Quality Results3", 1200, 800);
1116   canvas3->Divide(3, 2);
1117
1118   canvas3->cd(1);
1119   fhDXY[0]->SetStats(kFALSE);
1120   gPad->SetLogy();
1121   fhDXY[0]->Draw();
1122
1123   canvas3->cd(2);
1124   fhDZ[0]->SetStats(kFALSE);
1125   gPad->SetLogy();
1126   fhDZ[0]->Draw();
1127
1128   canvas3->cd(3);
1129   fhDXYvsDZ[0]->SetStats(kFALSE);
1130   gPad->SetLogz();
1131   gPad->SetRightMargin(0.15);
1132   fhDXYvsDZ[0]->Draw("COLZ");
1133
1134   canvas3->cd(4);
1135   fhDXYNormalized[0]->SetStats(kFALSE);
1136   gPad->SetLogy();
1137   fhDXYNormalized[0]->Draw();
1138
1139   canvas3->cd(5);
1140   fhDZNormalized[0]->SetStats(kFALSE);
1141   gPad->SetLogy();
1142   fhDZNormalized[0]->Draw();
1143
1144   canvas3->cd(6);
1145   fhDXYvsDZNormalized[0]->SetStats(kFALSE);
1146   gPad->SetLogz();
1147   gPad->SetRightMargin(0.15);
1148   fhDXYvsDZNormalized[0]->Draw("COLZ");
1149
1150   canvas3->SaveAs(Form("%s_%s.gif", GetName(), canvas3->GetName()));
1151
1152   TCanvas* canvas4 = new TCanvas(Form("%s_4", GetName()), "Track Quality Results4", 800, 500);
1153   canvas4->Divide(2, 1);
1154
1155   canvas4->cd(1);
1156   fhCutStatistics->SetStats(kFALSE);
1157   fhCutStatistics->LabelsOption("v");
1158   gPad->SetBottomMargin(0.3);
1159   fhCutStatistics->Draw();
1160
1161   canvas4->cd(2);
1162   fhCutCorrelation->SetStats(kFALSE);
1163   fhCutCorrelation->LabelsOption("v");
1164   gPad->SetBottomMargin(0.3);
1165   gPad->SetLeftMargin(0.3);
1166   fhCutCorrelation->Draw("COLZ");
1167
1168   canvas4->SaveAs(Form("%s_%s.gif", GetName(), canvas4->GetName()));
1169
1170   /*canvas->cd(1);
1171   fhDXYvsDZNormalized[0]->SetStats(kFALSE);
1172   fhDXYvsDZNormalized[0]->DrawCopy("COLZ");
1173
1174   canvas->cd(2);
1175   fhNClustersTPC[0]->SetStats(kFALSE);
1176   fhNClustersTPC[0]->DrawCopy();
1177
1178   canvas->cd(3);
1179   fhChi2PerClusterITS[0]->SetStats(kFALSE);
1180   fhChi2PerClusterITS[0]->DrawCopy();
1181   fhChi2PerClusterITS[1]->SetLineColor(2);
1182   fhChi2PerClusterITS[1]->DrawCopy("SAME");
1183
1184   canvas->cd(4);
1185   fhChi2PerClusterTPC[0]->SetStats(kFALSE);
1186   fhChi2PerClusterTPC[0]->DrawCopy();
1187   fhChi2PerClusterTPC[1]->SetLineColor(2);
1188   fhChi2PerClusterTPC[1]->DrawCopy("SAME");*/
1189 }
1190