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