fixing (finally) the nsigma to vertex cut
[u/mrichter/AliRoot.git] / PWG0 / esdTrackCuts / AliESDtrackCuts.cxx
1 #include "AliESDtrackCuts.h"
2
3
4 #include <AliESDtrack.h>
5 #include <AliESD.h>
6 #include <AliLog.h>
7
8 //____________________________________________________________________
9 ClassImp(AliESDtrackCuts)
10
11 // Cut names
12 const Char_t* AliESDtrackCuts::fgkCutNames[kNCuts] = {
13  "require TPC refit",
14  "require ITS refit",
15  "n clusters TPC",
16  "n clusters ITS",
17  "#Chi^{2}/clusters TPC",
18  "#Chi^{2}/clusters ITS",
19  "cov 11",
20  "cov 22",
21  "cov 33",
22  "cov 44",
23  "cov 55",
24  "trk-to-vtx",
25  "trk-to-vtx failed",
26  "kink daughters",
27
28  "p",
29  "p_{T}",
30  "p_{x}",
31  "p_{y}",
32  "p_{z}",
33  "y",
34  "eta"
35 };
36
37 //____________________________________________________________________
38 AliESDtrackCuts::AliESDtrackCuts()
39 {
40   //
41   // constructor
42   //
43
44   Init();
45
46   //##############################################################################
47   // setting default cuts
48   SetMinNClustersTPC();
49   SetMinNClustersITS();
50   SetMaxChi2PerClusterTPC();
51   SetMaxChi2PerClusterITS();                                
52   SetMaxCovDiagonalElements();                                      
53   SetRequireTPCRefit();
54   SetRequireITSRefit();
55   SetAcceptKingDaughters();
56   SetMinNsigmaToVertex();
57   SetRequireSigmaToVertex();
58   SetPRange();
59   SetPtRange();
60   SetPxRange();
61   SetPyRange();
62   SetPzRange();
63   SetEtaRange();
64   SetRapRange();
65
66   SetHistogramsOn();
67 }
68
69 //_____________________________________________________________________________
70 AliESDtrackCuts::AliESDtrackCuts(const AliESDtrackCuts &c) : TObject(c)
71 {
72   //
73   // copy constructor
74   //
75
76   ((AliESDtrackCuts &) c).Copy(*this);
77 }
78
79 AliESDtrackCuts::~AliESDtrackCuts()
80 {
81   //
82   // destructor
83   //
84
85   // ## TODO to be implemented
86 }
87
88 void AliESDtrackCuts::Init()
89 {
90   //
91   // sets everything to zero
92   //
93
94   fCutMinNClusterTPC = 0;
95   fCutMinNClusterITS = 0;
96
97   fCutMaxChi2PerClusterTPC = 0;
98   fCutMaxChi2PerClusterITS = 0;
99
100   fCutMaxC11 = 0;
101   fCutMaxC22 = 0;
102   fCutMaxC33 = 0;
103   fCutMaxC44 = 0;
104   fCutMaxC55 = 0;
105
106   fCutAcceptKinkDaughters = 0;
107   fCutRequireTPCRefit = 0;
108   fCutRequireITSRefit = 0;
109
110   fCutNsigmaToVertex = 0;
111   fCutSigmaToVertexRequired = 0;
112
113   fPMin = 0;
114   fPMax = 0;
115   fPtMin = 0;
116   fPtMax = 0;
117   fPxMin = 0;
118   fPxMax = 0;
119   fPyMin = 0;
120   fPyMax = 0;
121   fPzMin = 0;
122   fPzMax = 0;
123   fEtaMin = 0;
124   fEtaMax = 0;
125   fRapMin = 0;
126   fRapMax = 0;
127
128   fHistogramsOn = kFALSE;
129
130   for (Int_t i=0; i<2; ++i)
131   {
132     fhNClustersITS[i] = 0;
133     fhNClustersTPC[i] = 0;
134
135     fhChi2PerClusterITS[i] = 0;
136     fhChi2PerClusterTPC[i] = 0;
137
138     fhC11[i] = 0;
139     fhC22[i] = 0;
140     fhC33[i] = 0;
141     fhC44[i] = 0;
142     fhC55[i] = 0;
143
144     fhDXY[i] = 0;
145     fhDZ[i] = 0;
146     fhDXYvsDZ[i] = 0;
147
148     fhDXYNormalized[i] = 0;
149     fhDZNormalized[i] = 0;
150     fhDXYvsDZNormalized[i] = 0;
151   }
152
153   fhCutStatistics = 0;
154   fhCutCorrelation = 0;
155 }
156
157 //_____________________________________________________________________________
158 AliESDtrackCuts &AliESDtrackCuts::operator=(const AliESDtrackCuts &c)
159 {
160   //
161   // Assignment operator
162   //
163
164   if (this != &c) ((AliESDtrackCuts &) c).Copy(*this);
165   return *this;
166 }
167
168 //_____________________________________________________________________________
169 void AliESDtrackCuts::Copy(TObject &c) const
170 {
171   //
172   // Copy function
173   //
174
175   AliESDtrackCuts& target = (AliESDtrackCuts &) c;
176
177   target.Init();
178
179   target.fCutMinNClusterTPC = fCutMinNClusterTPC;
180   target.fCutMinNClusterITS = fCutMinNClusterITS;
181
182   target.fCutMaxChi2PerClusterTPC = fCutMaxChi2PerClusterTPC;
183   target.fCutMaxChi2PerClusterITS = fCutMaxChi2PerClusterITS;
184
185   target.fCutMaxC11 = fCutMaxC11;
186   target.fCutMaxC22 = fCutMaxC22;
187   target.fCutMaxC33 = fCutMaxC33;
188   target.fCutMaxC44 = fCutMaxC44;
189   target.fCutMaxC55 = fCutMaxC55;
190
191   target.fCutAcceptKinkDaughters = fCutAcceptKinkDaughters;
192   target.fCutRequireTPCRefit = fCutRequireTPCRefit;
193   target.fCutRequireITSRefit = fCutRequireITSRefit;
194
195   target.fCutNsigmaToVertex = fCutNsigmaToVertex;
196   target.fCutSigmaToVertexRequired = fCutSigmaToVertexRequired;
197
198   target.fPMin = fPMin;
199   target.fPMax = fPMax;
200   target.fPtMin = fPtMin;
201   target.fPtMax = fPtMax;
202   target.fPxMin = fPxMin;
203   target.fPxMax = fPxMax;
204   target.fPyMin = fPyMin;
205   target.fPyMax = fPyMax;
206   target.fPzMin = fPzMin;
207   target.fPzMax = fPzMax;
208   target.fEtaMin = fEtaMin;
209   target.fEtaMax = fEtaMax;
210   target.fRapMin = fRapMin;
211   target.fRapMax = fRapMax;
212
213   target.fHistogramsOn = fHistogramsOn;
214
215   for (Int_t i=0; i<2; ++i)
216   {
217     if (fhNClustersITS[i]) target.fhNClustersITS[i] = (TH1F*) fhNClustersITS[i]->Clone();
218     if (fhNClustersTPC[i]) target.fhNClustersTPC[i] = (TH1F*) fhNClustersTPC[i]->Clone();
219
220     if (fhChi2PerClusterITS[i]) target.fhChi2PerClusterITS[i] = (TH1F*) fhChi2PerClusterITS[i]->Clone();
221     if (fhChi2PerClusterTPC[i]) target.fhChi2PerClusterTPC[i] = (TH1F*) fhChi2PerClusterTPC[i]->Clone();
222
223     if (fhC11[i]) target.fhC11[i] = (TH1F*) fhC11[i]->Clone();
224     if (fhC22[i]) target.fhC22[i] = (TH1F*) fhC22[i]->Clone();
225     if (fhC33[i]) target.fhC33[i] = (TH1F*) fhC33[i]->Clone();
226     if (fhC44[i]) target.fhC44[i] = (TH1F*) fhC44[i]->Clone();
227     if (fhC55[i]) target.fhC55[i] = (TH1F*) fhC55[i]->Clone();
228
229     if (fhDXY[i]) target.fhDXY[i] = (TH1F*) fhDXY[i]->Clone();
230     if (fhDZ[i]) target.fhDZ[i] = (TH1F*) fhDZ[i]->Clone();
231     if (fhDXYvsDZ[i]) target.fhDXYvsDZ[i] = (TH2F*) fhDXYvsDZ[i]->Clone();
232
233     if (fhDXYNormalized[i]) target.fhDXYNormalized[i] = (TH1F*) fhDXYNormalized[i]->Clone();
234     if (fhDZNormalized[i]) target.fhDZNormalized[i] = (TH1F*) fhDZNormalized[i]->Clone();
235     if (fhDXYvsDZNormalized[i]) target.fhDXYvsDZNormalized[i] = (TH2F*) fhDXYvsDZNormalized[i]->Clone();
236   }
237
238   if (fhCutStatistics) target.fhCutStatistics = (TH1F*) fhCutStatistics->Clone();
239   if (fhCutCorrelation) target.fhCutCorrelation = (TH2F*) fhCutCorrelation->Clone();
240
241   TObject::Copy(c);
242 }
243
244 //____________________________________________________________________
245 Float_t AliESDtrackCuts::GetSigmaToVertex(AliESDtrack* esdTrack)
246 {
247   //
248
249   Float_t b[2];
250   Float_t bRes[2];
251   Float_t bCov[3];
252   esdTrack->GetImpactParameters(b,bCov);
253   if (bCov[0]<=0 || bCov[2]<=0) {
254     AliDebug(1, "Estimated b resolution lower or equal zero!");
255     bCov[0]=0; bCov[2]=0;
256   }
257   bRes[0] = TMath::Sqrt(bCov[0]);
258   bRes[1] = TMath::Sqrt(bCov[2]);
259
260   // -----------------------------------
261   // How to get to a n-sigma cut?
262   //
263   // The accumulated statistics from 0 to d is
264   //
265   // ->  Erf(d/Sqrt(2)) for a 1-dim gauss (d = n_sigma)
266   // ->  1 - Exp(-d**2) for a 2-dim gauss (d*d = dx*dx + dy*dy != n_sigma)
267   //
268   // It means that for a 2-dim gauss: n_sigma(d) = Sqrt(2)*ErfInv(1 - Exp((-x**2)/2)
269   // Can this be expressed in a different way?
270   //
271   //
272   // FIX: I don't think this is correct!!! Keeping d as n_sigma for now...
273
274   if (bRes[0] == 0 || bRes[1] ==0)
275     return -1;
276
277   Float_t d = TMath::Sqrt(TMath::Power(b[0]/bRes[0],2) + TMath::Power(b[1]/bRes[1],2));
278
279   // stupid rounding problem screws up everything:
280   // if d is too big, TMath::Exp(...) gets 0, and TMath::ErfInverse(1) that should be infinite, gets 0 :(
281   if (TMath::Exp(-d * d / 2) < 1e-10)
282     return 1000;
283
284   d = TMath::ErfInverse(1 - TMath::Exp(-d * d / 2)) * TMath::Sqrt(2);
285   return d;
286 }
287
288 //____________________________________________________________________
289 Bool_t
290 AliESDtrackCuts::AcceptTrack(AliESDtrack* esdTrack) {
291   // 
292   // figure out if the tracks survives all the track cuts defined
293   //
294   // the different quality parameter and kinematic values are first
295   // retrieved from the track. then it is found out what cuts the
296   // track did not survive and finally the cuts are imposed.
297
298
299
300   UInt_t status = esdTrack->GetStatus();
301
302   // dummy array
303   Int_t  fIdxInt[200];
304
305   // getting quality parameters from the ESD track
306   Int_t nClustersITS = esdTrack->GetITSclusters(fIdxInt);
307   Int_t nClustersTPC = esdTrack->GetTPCclusters(fIdxInt);
308   
309
310
311   Float_t chi2PerClusterITS = -1;
312   Float_t chi2PerClusterTPC = -1;
313   if (nClustersITS!=0)
314     chi2PerClusterITS = esdTrack->GetITSchi2()/Float_t(nClustersITS);
315   if (nClustersTPC!=0)
316     chi2PerClusterTPC = esdTrack->GetTPCchi2()/Float_t(nClustersTPC);  
317
318   Double_t extCov[15];
319   esdTrack->GetExternalCovariance(extCov);
320
321   // getting the track to vertex parameters
322   Float_t nSigmaToVertex = GetSigmaToVertex(esdTrack);
323
324   // getting the kinematic variables of the track
325   // (assuming the mass is known)
326   Double_t p[3];
327   esdTrack->GetPxPyPz(p);
328   Float_t momentum = TMath::Sqrt(TMath::Power(p[0],2) + TMath::Power(p[1],2) + TMath::Power(p[2],2));
329   Float_t pt       = TMath::Sqrt(TMath::Power(p[0],2) + TMath::Power(p[1],2));
330   Float_t energy   = TMath::Sqrt(TMath::Power(esdTrack->GetMass(),2) + TMath::Power(momentum,2));
331
332
333   //y-eta related calculations
334   Float_t eta = -100.;
335   Float_t y   = -100.;
336   if((momentum != TMath::Abs(p[2]))&&(momentum != 0))
337     eta = 0.5*TMath::Log((momentum + p[2])/(momentum - p[2]));
338   if((energy != TMath::Abs(p[2]))&&(momentum != 0))
339     y = 0.5*TMath::Log((energy + p[2])/(energy - p[2]));
340
341   
342   //########################################################################
343   // cut the track?
344   
345   Bool_t cuts[kNCuts];
346   for (Int_t i=0; i<kNCuts; i++) cuts[i]=kFALSE;
347   
348   // track quality cuts
349   if (fCutRequireTPCRefit && (status&AliESDtrack::kTPCrefit)==0)
350     cuts[0]=kTRUE;
351   if (fCutRequireITSRefit && (status&AliESDtrack::kITSrefit)==0)
352     cuts[1]=kTRUE;
353   if (nClustersTPC<fCutMinNClusterTPC)
354     cuts[2]=kTRUE;
355   if (nClustersITS<fCutMinNClusterITS) 
356     cuts[3]=kTRUE;
357   if (chi2PerClusterTPC>fCutMaxChi2PerClusterTPC) 
358     cuts[4]=kTRUE; 
359   if (chi2PerClusterITS>fCutMaxChi2PerClusterITS) 
360     cuts[5]=kTRUE;
361   if (extCov[0]  > fCutMaxC11) 
362     cuts[6]=kTRUE;  
363   if (extCov[2]  > fCutMaxC22) 
364     cuts[7]=kTRUE;  
365   if (extCov[5]  > fCutMaxC33) 
366     cuts[8]=kTRUE;  
367   if (extCov[9]  > fCutMaxC44) 
368     cuts[9]=kTRUE;  
369   if (extCov[14]  > fCutMaxC55) 
370     cuts[10]=kTRUE;  
371   if (nSigmaToVertex > fCutNsigmaToVertex)
372     cuts[11] = kTRUE;
373   // if n sigma could not be calculated
374   if (nSigmaToVertex<0 && fCutSigmaToVertexRequired)
375     cuts[12]=kTRUE;
376   if (!fCutAcceptKinkDaughters && esdTrack->GetKinkIndex(0)>0) 
377     cuts[13]=kTRUE;
378   // track kinematics cut
379   if((momentum < fPMin) || (momentum > fPMax)) 
380     cuts[14]=kTRUE;
381   if((pt < fPtMin) || (pt > fPtMax)) 
382     cuts[15] = kTRUE;
383   if((p[0] < fPxMin) || (p[0] > fPxMax)) 
384     cuts[16] = kTRUE;
385   if((p[1] < fPyMin) || (p[1] > fPyMax)) 
386     cuts[17] = kTRUE;
387   if((p[2] < fPzMin) || (p[2] > fPzMax))
388     cuts[18] = kTRUE;
389   if((eta < fEtaMin) || (eta > fEtaMax)) 
390     cuts[19] = kTRUE;
391   if((y < fRapMin) || (y > fRapMax)) 
392     cuts[20] = kTRUE;
393
394   Bool_t cut=kFALSE;
395   for (Int_t i=0; i<kNCuts; i++) 
396     if (cuts[i]) cut = kTRUE;
397   
398   //########################################################################
399   // filling histograms
400   if (fHistogramsOn) {
401     fhCutStatistics->Fill(fhCutStatistics->GetBinCenter(fhCutStatistics->GetXaxis()->FindBin("n tracks")));
402     
403     if (cut)
404       fhCutStatistics->Fill(fhCutStatistics->GetBinCenter(fhCutStatistics->GetXaxis()->FindBin("n cut tracks")));
405     
406     for (Int_t i=0; i<kNCuts; i++) {
407       if (cuts[i])
408         fhCutStatistics->Fill(fhCutStatistics->GetBinCenter(fhCutStatistics->GetXaxis()->FindBin(fgkCutNames[i])));
409       
410       for (Int_t j=i; j<kNCuts; j++) {
411         if (cuts[i] && cuts[j]) {
412           Float_t x = fhCutCorrelation->GetXaxis()->GetBinCenter(fhCutCorrelation->GetXaxis()->FindBin(fgkCutNames[i]));
413           Float_t y = fhCutCorrelation->GetYaxis()->GetBinCenter(fhCutCorrelation->GetYaxis()->FindBin(fgkCutNames[j]));
414           fhCutCorrelation->Fill(x,y);
415         }
416       }
417     }
418     
419
420     fhNClustersITS[0]->Fill(nClustersITS);
421     fhNClustersTPC[0]->Fill(nClustersTPC);
422     fhChi2PerClusterITS[0]->Fill(chi2PerClusterITS);
423     fhChi2PerClusterTPC[0]->Fill(chi2PerClusterTPC);
424
425     fhC11[0]->Fill(extCov[0]);
426     fhC22[0]->Fill(extCov[2]);
427     fhC33[0]->Fill(extCov[5]);
428     fhC44[0]->Fill(extCov[9]);
429     fhC55[0]->Fill(extCov[14]);
430
431     Float_t b[2];
432     Float_t bRes[2];
433     Float_t bCov[3];
434     esdTrack->GetImpactParameters(b,bCov);
435     if (bCov[0]<=0 || bCov[2]<=0) {
436       AliDebug(1, "Estimated b resolution lower or equal zero!");
437       bCov[0]=0; bCov[2]=0;
438     }
439     bRes[0] = TMath::Sqrt(bCov[0]);
440     bRes[1] = TMath::Sqrt(bCov[2]);
441
442     fhDZ[0]->Fill(b[1]);
443     fhDXY[0]->Fill(b[0]);
444     fhDXYvsDZ[0]->Fill(b[1],b[0]);
445
446     if (bRes[0]!=0 && bRes[1]!=0) {
447       fhDZNormalized[0]->Fill(b[1]/bRes[1]);
448       fhDXYNormalized[0]->Fill(b[0]/bRes[0]);
449       fhDXYvsDZNormalized[0]->Fill(b[1]/bRes[1], b[0]/bRes[0]);
450     }
451   }
452
453   //########################################################################
454   // cut the track!
455   if (cut) return kFALSE;
456
457   //########################################################################
458   // filling histograms after cut
459   if (fHistogramsOn) {
460     fhNClustersITS[1]->Fill(nClustersITS);
461     fhNClustersTPC[1]->Fill(nClustersTPC);
462     fhChi2PerClusterITS[1]->Fill(chi2PerClusterITS);
463     fhChi2PerClusterTPC[1]->Fill(chi2PerClusterTPC);
464
465     fhC11[1]->Fill(extCov[0]);
466     fhC22[1]->Fill(extCov[2]);
467     fhC33[1]->Fill(extCov[5]);
468     fhC44[1]->Fill(extCov[9]);
469     fhC55[1]->Fill(extCov[14]);
470
471     Float_t b[2];
472     Float_t bRes[2];
473     Float_t bCov[3];
474     esdTrack->GetImpactParameters(b,bCov);
475     if (bCov[0]<=0 || bCov[2]<=0) {
476       AliDebug(1, "Estimated b resolution lower or equal zero!");
477       bCov[0]=0; bCov[2]=0;
478     }
479     bRes[0] = TMath::Sqrt(bCov[0]);
480     bRes[1] = TMath::Sqrt(bCov[2]);
481
482     fhDZ[1]->Fill(b[1]);
483     fhDXY[1]->Fill(b[0]);
484     fhDXYvsDZ[1]->Fill(b[1],b[0]);
485
486     if (bRes[0]!=0 && bRes[1]!=0)
487     {
488       fhDZNormalized[1]->Fill(b[1]/bRes[1]);
489       fhDXYNormalized[1]->Fill(b[0]/bRes[0]);
490       fhDXYvsDZNormalized[1]->Fill(b[1]/bRes[1], b[0]/bRes[0]);
491     }
492   }
493
494   return kTRUE;
495 }
496
497 //____________________________________________________________________
498 TObjArray*
499 AliESDtrackCuts::GetAcceptedTracks(AliESD* esd)
500 {
501   //
502   // returns an array of all tracks that pass the cuts
503   //
504
505   TObjArray* acceptedTracks = new TObjArray();
506
507   // loop over esd tracks
508   for (Int_t iTrack = 0; iTrack < esd->GetNumberOfTracks(); iTrack++) {
509     AliESDtrack* track = esd->GetTrack(iTrack);
510
511     if (AcceptTrack(track))
512       acceptedTracks->Add(track);
513   }
514
515   return acceptedTracks;
516 }
517
518 //____________________________________________________________________
519 Int_t
520 AliESDtrackCuts::CountAcceptedTracks(AliESD* esd)
521 {
522   //
523   // returns an the number of tracks that pass the cuts
524   //
525
526   Int_t count = 0;
527
528   // loop over esd tracks
529   for (Int_t iTrack = 0; iTrack < esd->GetNumberOfTracks(); iTrack++) {
530     AliESDtrack* track = esd->GetTrack(iTrack);
531
532     if (AcceptTrack(track))
533       count++;
534   }
535
536   return count;
537 }
538
539 //____________________________________________________________________
540  void AliESDtrackCuts::DefineHistograms(Int_t color) {
541    // 
542    // diagnostics histograms are defined
543    // 
544
545    fHistogramsOn=kTRUE;
546
547    //###################################################################################
548    // defining histograms
549
550    fhCutStatistics = new TH1F("cut_statistics","cut statistics",kNCuts+4,-0.5,kNCuts+3.5);
551
552    fhCutStatistics->GetXaxis()->SetBinLabel(1,"n tracks");
553    fhCutStatistics->GetXaxis()->SetBinLabel(2,"n cut tracks");
554
555    fhCutCorrelation = new TH2F("cut_correlation","cut correlation",kNCuts,-0.5,kNCuts-0.5,kNCuts,-0.5,kNCuts-0.5);;
556   
557    for (Int_t i=0; i<kNCuts; i++) {
558      fhCutStatistics->GetXaxis()->SetBinLabel(i+4,fgkCutNames[i]);
559      fhCutCorrelation->GetXaxis()->SetBinLabel(i+1,fgkCutNames[i]);
560      fhCutCorrelation->GetYaxis()->SetBinLabel(i+1,fgkCutNames[i]);
561    } 
562
563   fhCutStatistics  ->SetLineColor(color);
564   fhCutCorrelation ->SetLineColor(color);
565   fhCutStatistics  ->SetLineWidth(2);
566   fhCutCorrelation ->SetLineWidth(2);
567
568   Char_t str[256];
569   for (Int_t i=0; i<2; i++) {
570     if (i==0) sprintf(str," ");
571     else sprintf(str,"_cut");
572
573     fhNClustersITS[i]        = new TH1F(Form("nClustersITS%s",str),"",8,-0.5,7.5);
574     fhNClustersTPC[i]        = new TH1F(Form("nClustersTPC%s",str),"",165,-0.5,164.5);
575     fhChi2PerClusterITS[i]   = new TH1F(Form("chi2PerClusterITS%s",str),"",500,0,10);
576     fhChi2PerClusterTPC[i]   = new TH1F(Form("chi2PerClusterTPC%s",str),"",500,0,10);
577
578     fhC11[i]                 = new TH1F(Form("covMatrixDiagonal11%s",str),"",1000,0,5);
579     fhC22[i]                 = new  TH1F(Form("covMatrixDiagonal22%s",str),"",1000,0,5);
580     fhC33[i]                 = new  TH1F(Form("covMatrixDiagonal33%s",str),"",1000,0,0.5);
581     fhC44[i]                 = new  TH1F(Form("covMatrixDiagonal44%s",str),"",1000,0,5);
582     fhC55[i]                 = new  TH1F(Form("covMatrixDiagonal55%s",str),"",1000,0,5);
583
584     fhDXY[i]                 = new  TH1F(Form("dXY%s",str),"",500,-10,10);
585     fhDZ[i]                  = new  TH1F(Form("dZ%s",str),"",500,-10,10);
586     fhDXYvsDZ[i]             = new  TH2F(Form("dXYvsDZ%s",str),"",200,-10,10,200,-10,10);
587
588     fhDXYNormalized[i]       = new  TH1F(Form("dXYNormalized%s",str),"",500,-10,10);
589     fhDZNormalized[i]        = new  TH1F(Form("dZNormalized%s",str),"",500,-10,10);
590     fhDXYvsDZNormalized[i]   = new  TH2F(Form("dXYvsDZNormalized%s",str),"",200,-10,10,200,-10,10);
591
592
593     fhNClustersITS[i]->SetXTitle("n ITS clusters");
594     fhNClustersTPC[i]->SetXTitle("n TPC clusters");
595     fhChi2PerClusterITS[i]->SetXTitle("#Chi^{2} per ITS cluster");
596     fhChi2PerClusterTPC[i]->SetXTitle("#Chi^{2} per TPC cluster");
597
598     fhC11[i]->SetXTitle("cov 11 : #sigma_{y}^{2} [cm^{2}]");
599     fhC22[i]->SetXTitle("cov 22 : #sigma_{z}^{2} [cm^{2}]");
600     fhC33[i]->SetXTitle("cov 33 : #sigma_{sin(#phi)}^{2}");
601     fhC44[i]->SetXTitle("cov 44 : #sigma_{tan(#theta_{dip})}^{2}");
602     fhC55[i]->SetXTitle("cov 55 : #sigma_{1/p_{T}}^{2} [(c/GeV)^2]");
603
604     fhDXY[i]->SetXTitle("transverse impact parameter");
605     fhDZ[i]->SetXTitle("longitudinal impact parameter");
606     fhDXYvsDZ[i]->SetXTitle("longitudinal impact parameter");
607     fhDXYvsDZ[i]->SetYTitle("transverse impact parameter");
608
609     fhDXYNormalized[i]->SetXTitle("normalized trans impact par");
610     fhDZNormalized[i]->SetXTitle("normalized long impact par");
611     fhDXYvsDZNormalized[i]->SetXTitle("normalized long impact par"); 
612     fhDXYvsDZNormalized[i]->SetYTitle("normalized trans impact par");
613
614     fhNClustersITS[i]->SetLineColor(color);   fhNClustersITS[i]->SetLineWidth(2);
615     fhNClustersTPC[i]->SetLineColor(color);   fhNClustersTPC[i]->SetLineWidth(2);
616     fhChi2PerClusterITS[i]->SetLineColor(color);   fhChi2PerClusterITS[i]->SetLineWidth(2);
617     fhChi2PerClusterTPC[i]->SetLineColor(color);   fhChi2PerClusterTPC[i]->SetLineWidth(2);
618
619     fhC11[i]->SetLineColor(color);   fhC11[i]->SetLineWidth(2);
620     fhC22[i]->SetLineColor(color);   fhC22[i]->SetLineWidth(2);
621     fhC33[i]->SetLineColor(color);   fhC33[i]->SetLineWidth(2);
622     fhC44[i]->SetLineColor(color);   fhC44[i]->SetLineWidth(2);
623     fhC55[i]->SetLineColor(color);   fhC55[i]->SetLineWidth(2);
624
625     fhDXY[i]->SetLineColor(color);   fhDXY[i]->SetLineWidth(2);
626     fhDZ[i]->SetLineColor(color);   fhDZ[i]->SetLineWidth(2);
627
628     fhDXYNormalized[i]->SetLineColor(color);   fhDXYNormalized[i]->SetLineWidth(2);
629     fhDZNormalized[i]->SetLineColor(color);   fhDZNormalized[i]->SetLineWidth(2);
630   }
631 }
632
633 //____________________________________________________________________
634 void 
635 AliESDtrackCuts::Print(const Option_t*) const {
636   //
637   // print method - still to be implemented
638   //
639
640   AliInfo("AliESDtrackCuts...");
641 }
642
643
644 //____________________________________________________________________
645 void AliESDtrackCuts::SaveHistograms(Char_t* dir) {
646   // 
647   // saves the histograms in a directory (dir)
648   // 
649
650   
651   if (!fHistogramsOn) {
652     AliDebug(0, "Histograms not on - cannot save histograms!!!");
653     return;
654   }
655
656   gDirectory->mkdir(dir);
657   gDirectory->cd(dir);
658
659   gDirectory->mkdir("before_cuts");
660   gDirectory->mkdir("after_cuts");
661  
662   fhCutStatistics->Write();
663   fhCutCorrelation->Write();
664
665   for (Int_t i=0; i<2; i++) {
666     if (i==0)
667       gDirectory->cd("before_cuts");
668     else
669       gDirectory->cd("after_cuts");
670     
671     fhNClustersITS[i]        ->Write();
672     fhNClustersTPC[i]        ->Write();
673     fhChi2PerClusterITS[i]   ->Write();
674     fhChi2PerClusterTPC[i]   ->Write();
675     
676     fhC11[i]                 ->Write();
677     fhC22[i]                 ->Write();
678     fhC33[i]                 ->Write();
679     fhC44[i]                 ->Write();
680     fhC55[i]                 ->Write();
681
682     fhDXY[i]                 ->Write();
683     fhDZ[i]                  ->Write();
684     fhDXYvsDZ[i]             ->Write();
685     
686     fhDXYNormalized[i]       ->Write();
687     fhDZNormalized[i]        ->Write();
688     fhDXYvsDZNormalized[i]   ->Write();
689
690     gDirectory->cd("../");
691   }
692
693   gDirectory->cd("../");
694 }
695
696
697