]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG0/esdTrackCuts/AliESDtrackCuts.cxx
98ead8e1349c800487d0af42dabdb4851fbb2958
[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 Bool_t 
246 AliESDtrackCuts::AcceptTrack(AliESDtrack* esdTrack) {
247   // 
248   // figure out if the tracks survives all the track cuts defined
249   //
250   // the different quality parameter and kinematic values are first
251   // retrieved from the track. then it is found out what cuts the
252   // track did not survive and finally the cuts are imposed.
253
254
255
256   UInt_t status = esdTrack->GetStatus();
257
258   // dummy array
259   Int_t  fIdxInt[200];
260
261   // getting quality parameters from the ESD track
262   Int_t nClustersITS = esdTrack->GetITSclusters(fIdxInt);
263   Int_t nClustersTPC = esdTrack->GetTPCclusters(fIdxInt);
264   
265
266
267   Float_t chi2PerClusterITS = -1;
268   Float_t chi2PerClusterTPC = -1;
269   if (nClustersITS!=0)
270     chi2PerClusterITS = esdTrack->GetITSchi2()/Float_t(nClustersITS);
271   if (nClustersTPC!=0)
272     chi2PerClusterTPC = esdTrack->GetTPCchi2()/Float_t(nClustersTPC);  
273
274   Double_t extCov[15];
275   esdTrack->GetExternalCovariance(extCov);
276
277   // getting the track to vertex parameters
278   Float_t b[2];
279   Float_t bRes[2];
280   Float_t bCov[3];
281   esdTrack->GetImpactParameters(b,bCov);    
282   if (bCov[0]<=0 || bCov[2]<=0) {
283     AliDebug(1, "Estimated b resolution lower or equal zero!");
284     bCov[0]=0; bCov[2]=0;
285   } 
286   bRes[0] = TMath::Sqrt(bCov[0]);
287   bRes[1] = TMath::Sqrt(bCov[2]);
288
289   // -----------------------------------
290   // How to get to a n-sigma cut?
291   //
292   // The accumulated statistics from 0 to d is
293   //
294   // ->  Erf(d/Sqrt(2)) for a 1-dim gauss (d = n_sigma)
295   // ->  1 - Exp(-d**2) for a 2-dim gauss (d*d = dx*dx + dy*dy != n_sigma)
296   //
297   // It means that for a 2-dim gauss: n_sigma(d) = Sqrt(2)*ErfInv(1 - Exp((-x**2)/2)
298   // Can this be expressed in a different way?
299   //
300   //
301   // FIX: I don't think this is correct!!! Keeping d as n_sigma for now...
302
303   Float_t nSigmaToVertex = -1;
304   if (bRes[0]!=0 && bRes[1]!=0) {
305     Float_t d = TMath::Sqrt(TMath::Power(b[0]/bRes[0],2) + TMath::Power(b[1]/bRes[1],2));
306     nSigmaToVertex = d;//TMath::Sqrt(2)*(TMath::ErfInverse(1 - TMath::Exp(0.5*(-d*d))));
307     // JF solution: nSigmaToVertex = TMath::ErfInverse(TMath::Sqrt(2.0/TMath::Pi()) * TMath::Erf(d / TMath::Sqrt(2))) * TMath::Sqrt(2);
308   }
309
310   // getting the kinematic variables of the track 
311   // (assuming the mass is known)
312   Double_t p[3];
313   esdTrack->GetPxPyPz(p);
314   Float_t momentum = TMath::Sqrt(TMath::Power(p[0],2) + TMath::Power(p[1],2) + TMath::Power(p[2],2));
315   Float_t pt       = TMath::Sqrt(TMath::Power(p[0],2) + TMath::Power(p[1],2));
316   Float_t energy   = TMath::Sqrt(TMath::Power(esdTrack->GetMass(),2) + TMath::Power(momentum,2));
317
318
319   //y-eta related calculations
320   Float_t eta = -100.;
321   Float_t y   = -100.;
322   if((momentum != TMath::Abs(p[2]))&&(momentum != 0))
323     eta = 0.5*TMath::Log((momentum + p[2])/(momentum - p[2]));
324   if((energy != TMath::Abs(p[2]))&&(momentum != 0))
325     y = 0.5*TMath::Log((energy + p[2])/(energy - p[2]));
326
327   
328   //########################################################################
329   // cut the track?
330   
331   Bool_t cuts[kNCuts];
332   for (Int_t i=0; i<kNCuts; i++) cuts[i]=kFALSE;
333   
334   // track quality cuts
335   if (fCutRequireTPCRefit && (status&AliESDtrack::kTPCrefit)==0)
336     cuts[0]=kTRUE;
337   if (fCutRequireITSRefit && (status&AliESDtrack::kITSrefit)==0)
338     cuts[1]=kTRUE;
339   if (nClustersTPC<fCutMinNClusterTPC) 
340     cuts[2]=kTRUE;
341   if (nClustersITS<fCutMinNClusterITS) 
342     cuts[3]=kTRUE;
343   if (chi2PerClusterTPC>fCutMaxChi2PerClusterTPC) 
344     cuts[4]=kTRUE; 
345   if (chi2PerClusterITS>fCutMaxChi2PerClusterITS) 
346     cuts[5]=kTRUE;
347   if (extCov[0]  > fCutMaxC11) 
348     cuts[6]=kTRUE;  
349   if (extCov[2]  > fCutMaxC22) 
350     cuts[7]=kTRUE;  
351   if (extCov[5]  > fCutMaxC33) 
352     cuts[8]=kTRUE;  
353   if (extCov[9]  > fCutMaxC44) 
354     cuts[9]=kTRUE;  
355   if (extCov[14]  > fCutMaxC55) 
356     cuts[10]=kTRUE;  
357   if (nSigmaToVertex > fCutNsigmaToVertex) 
358     cuts[11] = kTRUE;
359   // if n sigma could not be calculated
360   if (nSigmaToVertex<0 && fCutSigmaToVertexRequired)   
361     cuts[12]=kTRUE;
362   if (!fCutAcceptKinkDaughters && esdTrack->GetKinkIndex(0)>0) 
363     cuts[13]=kTRUE;
364   // track kinematics cut
365   if((momentum < fPMin) || (momentum > fPMax)) 
366     cuts[14]=kTRUE;
367   if((pt < fPtMin) || (pt > fPtMax)) 
368     cuts[15] = kTRUE;
369   if((p[0] < fPxMin) || (p[0] > fPxMax)) 
370     cuts[16] = kTRUE;
371   if((p[1] < fPyMin) || (p[1] > fPyMax)) 
372     cuts[17] = kTRUE;
373   if((p[2] < fPzMin) || (p[2] > fPzMax)) 
374     cuts[18] = kTRUE;
375   if((eta < fEtaMin) || (eta > fEtaMax)) 
376     cuts[19] = kTRUE;
377   if((y < fRapMin) || (y > fRapMax)) 
378     cuts[20] = kTRUE;
379
380   Bool_t cut=kFALSE;
381   for (Int_t i=0; i<kNCuts; i++) 
382     if (cuts[i]) cut = kTRUE;
383   
384   //########################################################################
385   // filling histograms
386   if (fHistogramsOn) {
387     fhCutStatistics->Fill(fhCutStatistics->GetBinCenter(fhCutStatistics->GetXaxis()->FindBin("n tracks")));
388     
389     if (cut)
390       fhCutStatistics->Fill(fhCutStatistics->GetBinCenter(fhCutStatistics->GetXaxis()->FindBin("n cut tracks")));
391     
392     for (Int_t i=0; i<kNCuts; i++) {
393       if (cuts[i])
394         fhCutStatistics->Fill(fhCutStatistics->GetBinCenter(fhCutStatistics->GetXaxis()->FindBin(fgkCutNames[i])));
395       
396       for (Int_t j=i; j<kNCuts; j++) {
397         if (cuts[i] && cuts[j]) {
398           Float_t x = fhCutCorrelation->GetXaxis()->GetBinCenter(fhCutCorrelation->GetXaxis()->FindBin(fgkCutNames[i]));
399           Float_t y = fhCutCorrelation->GetYaxis()->GetBinCenter(fhCutCorrelation->GetYaxis()->FindBin(fgkCutNames[j]));
400           fhCutCorrelation->Fill(x,y);
401         }
402       }
403     }
404     
405
406     fhNClustersITS[0]->Fill(nClustersITS);        
407     fhNClustersTPC[0]->Fill(nClustersTPC);        
408     fhChi2PerClusterITS[0]->Fill(chi2PerClusterITS);
409     fhChi2PerClusterTPC[0]->Fill(chi2PerClusterTPC);   
410     
411     fhC11[0]->Fill(extCov[0]);                 
412     fhC22[0]->Fill(extCov[2]);                 
413     fhC33[0]->Fill(extCov[5]);                 
414     fhC44[0]->Fill(extCov[9]);                                  
415     fhC55[0]->Fill(extCov[14]);                                  
416     
417     fhDZ[0]->Fill(b[1]);     
418     fhDXY[0]->Fill(b[0]);    
419     fhDXYvsDZ[0]->Fill(b[1],b[0]);
420
421     if (bRes[0]!=0 && bRes[1]!=0) {
422       fhDZNormalized[0]->Fill(b[1]/bRes[1]);     
423       fhDXYNormalized[0]->Fill(b[0]/bRes[0]);    
424       fhDXYvsDZNormalized[0]->Fill(b[1]/bRes[1], b[0]/bRes[0]);
425     }
426   }
427
428   //########################################################################  
429   // cut the track!
430   if (cut) return kFALSE;
431
432   //########################################################################  
433   // filling histograms after cut
434   if (fHistogramsOn) {
435     fhNClustersITS[1]->Fill(nClustersITS);        
436     fhNClustersTPC[1]->Fill(nClustersTPC);        
437     fhChi2PerClusterITS[1]->Fill(chi2PerClusterITS);
438     fhChi2PerClusterTPC[1]->Fill(chi2PerClusterTPC);   
439     
440     fhC11[1]->Fill(extCov[0]);                 
441     fhC22[1]->Fill(extCov[2]);                 
442     fhC33[1]->Fill(extCov[5]);                 
443     fhC44[1]->Fill(extCov[9]);                                  
444     fhC55[1]->Fill(extCov[14]);                                  
445     
446     fhDZ[1]->Fill(b[1]);     
447     fhDXY[1]->Fill(b[0]);    
448     fhDXYvsDZ[1]->Fill(b[1],b[0]);
449
450     fhDZNormalized[1]->Fill(b[1]/bRes[1]);     
451     fhDXYNormalized[1]->Fill(b[0]/bRes[0]);    
452     fhDXYvsDZNormalized[1]->Fill(b[1]/bRes[1], b[0]/bRes[0]);
453   }
454   
455   return kTRUE;
456 }
457
458 //____________________________________________________________________
459 TObjArray*
460 AliESDtrackCuts::GetAcceptedTracks(AliESD* esd)
461 {
462   //
463   // returns an array of all tracks that pass the cuts
464   //
465
466   TObjArray* acceptedTracks = new TObjArray();
467
468   // loop over esd tracks
469   for (Int_t iTrack = 0; iTrack < esd->GetNumberOfTracks(); iTrack++) {
470     AliESDtrack* track = esd->GetTrack(iTrack);
471
472     if (AcceptTrack(track))
473       acceptedTracks->Add(track);
474   }
475
476   return acceptedTracks;
477 }
478
479 //____________________________________________________________________
480 Int_t
481 AliESDtrackCuts::CountAcceptedTracks(AliESD* esd)
482 {
483   //
484   // returns an the number of tracks that pass the cuts
485   //
486
487   Int_t count = 0;
488
489   // loop over esd tracks
490   for (Int_t iTrack = 0; iTrack < esd->GetNumberOfTracks(); iTrack++) {
491     AliESDtrack* track = esd->GetTrack(iTrack);
492
493     if (AcceptTrack(track))
494       count++;
495   }
496
497   return count;
498 }
499
500 //____________________________________________________________________
501  void AliESDtrackCuts::DefineHistograms(Int_t color) {
502    // 
503    // diagnostics histograms are defined
504    // 
505
506    fHistogramsOn=kTRUE;
507
508    //###################################################################################
509    // defining histograms
510
511    fhCutStatistics = new TH1F("cut_statistics","cut statistics",kNCuts+4,-0.5,kNCuts+3.5);
512
513    fhCutStatistics->GetXaxis()->SetBinLabel(1,"n tracks");
514    fhCutStatistics->GetXaxis()->SetBinLabel(2,"n cut tracks");
515
516    fhCutCorrelation = new TH2F("cut_correlation","cut correlation",kNCuts,-0.5,kNCuts-0.5,kNCuts,-0.5,kNCuts-0.5);;
517   
518    for (Int_t i=0; i<kNCuts; i++) {
519      fhCutStatistics->GetXaxis()->SetBinLabel(i+4,fgkCutNames[i]);
520      fhCutCorrelation->GetXaxis()->SetBinLabel(i+1,fgkCutNames[i]);
521      fhCutCorrelation->GetYaxis()->SetBinLabel(i+1,fgkCutNames[i]);
522    } 
523
524   fhCutStatistics  ->SetLineColor(color);
525   fhCutCorrelation ->SetLineColor(color);
526   fhCutStatistics  ->SetLineWidth(2);
527   fhCutCorrelation ->SetLineWidth(2);
528
529   Char_t str[256];
530   for (Int_t i=0; i<2; i++) {
531     if (i==0) sprintf(str," ");
532     else sprintf(str,"_cut");
533
534     fhNClustersITS[i]        = new TH1F(Form("nClustersITS%s",str),"",8,-0.5,7.5);
535     fhNClustersTPC[i]        = new TH1F(Form("nClustersTPC%s",str),"",165,-0.5,164.5);
536     fhChi2PerClusterITS[i]   = new TH1F(Form("chi2PerClusterITS%s",str),"",500,0,10);
537     fhChi2PerClusterTPC[i]   = new TH1F(Form("chi2PerClusterTPC%s",str),"",500,0,10);
538
539     fhC11[i]                 = new TH1F(Form("covMatrixDiagonal11%s",str),"",1000,0,5);
540     fhC22[i]                 = new  TH1F(Form("covMatrixDiagonal22%s",str),"",1000,0,5);
541     fhC33[i]                 = new  TH1F(Form("covMatrixDiagonal33%s",str),"",1000,0,0.5);
542     fhC44[i]                 = new  TH1F(Form("covMatrixDiagonal44%s",str),"",1000,0,5);
543     fhC55[i]                 = new  TH1F(Form("covMatrixDiagonal55%s",str),"",1000,0,5);
544
545     fhDXY[i]                 = new  TH1F(Form("dXY%s",str),"",500,-10,10);
546     fhDZ[i]                  = new  TH1F(Form("dZ%s",str),"",500,-10,10);
547     fhDXYvsDZ[i]             = new  TH2F(Form("dXYvsDZ%s",str),"",200,-10,10,200,-10,10);
548
549     fhDXYNormalized[i]       = new  TH1F(Form("dXYNormalized%s",str),"",500,-10,10);
550     fhDZNormalized[i]        = new  TH1F(Form("dZNormalized%s",str),"",500,-10,10);
551     fhDXYvsDZNormalized[i]   = new  TH2F(Form("dXYvsDZNormalized%s",str),"",200,-10,10,200,-10,10);
552
553
554     fhNClustersITS[i]->SetXTitle("n ITS clusters");
555     fhNClustersTPC[i]->SetXTitle("n TPC clusters");
556     fhChi2PerClusterITS[i]->SetXTitle("#Chi^{2} per ITS cluster");
557     fhChi2PerClusterTPC[i]->SetXTitle("#Chi^{2} per TPC cluster");
558
559     fhC11[i]->SetXTitle("cov 11 : #sigma_{y}^{2} [cm^{2}]");
560     fhC22[i]->SetXTitle("cov 22 : #sigma_{z}^{2} [cm^{2}]");
561     fhC33[i]->SetXTitle("cov 33 : #sigma_{sin(#phi)}^{2}");
562     fhC44[i]->SetXTitle("cov 44 : #sigma_{tan(#theta_{dip})}^{2}");
563     fhC55[i]->SetXTitle("cov 55 : #sigma_{1/p_{T}}^{2} [(c/GeV)^2]");
564
565     fhDXY[i]->SetXTitle("transverse impact parameter");
566     fhDZ[i]->SetXTitle("longitudinal impact parameter");
567     fhDXYvsDZ[i]->SetXTitle("longitudinal impact parameter");
568     fhDXYvsDZ[i]->SetYTitle("transverse impact parameter");
569
570     fhDXYNormalized[i]->SetXTitle("normalized trans impact par");
571     fhDZNormalized[i]->SetXTitle("normalized long impact par");
572     fhDXYvsDZNormalized[i]->SetXTitle("normalized long impact par"); 
573     fhDXYvsDZNormalized[i]->SetYTitle("normalized trans impact par");
574
575     fhNClustersITS[i]->SetLineColor(color);   fhNClustersITS[i]->SetLineWidth(2);
576     fhNClustersTPC[i]->SetLineColor(color);   fhNClustersTPC[i]->SetLineWidth(2);
577     fhChi2PerClusterITS[i]->SetLineColor(color);   fhChi2PerClusterITS[i]->SetLineWidth(2);
578     fhChi2PerClusterTPC[i]->SetLineColor(color);   fhChi2PerClusterTPC[i]->SetLineWidth(2);
579
580     fhC11[i]->SetLineColor(color);   fhC11[i]->SetLineWidth(2);
581     fhC22[i]->SetLineColor(color);   fhC22[i]->SetLineWidth(2);
582     fhC33[i]->SetLineColor(color);   fhC33[i]->SetLineWidth(2);
583     fhC44[i]->SetLineColor(color);   fhC44[i]->SetLineWidth(2);
584     fhC55[i]->SetLineColor(color);   fhC55[i]->SetLineWidth(2);
585
586     fhDXY[i]->SetLineColor(color);   fhDXY[i]->SetLineWidth(2);
587     fhDZ[i]->SetLineColor(color);   fhDZ[i]->SetLineWidth(2);
588
589     fhDXYNormalized[i]->SetLineColor(color);   fhDXYNormalized[i]->SetLineWidth(2);
590     fhDZNormalized[i]->SetLineColor(color);   fhDZNormalized[i]->SetLineWidth(2);
591   }
592 }
593
594 //____________________________________________________________________
595 void 
596 AliESDtrackCuts::Print(const Option_t*) const {
597   //
598   // print method - still to be implemented
599   //
600
601   AliInfo("AliESDtrackCuts...");
602 }
603
604
605 //____________________________________________________________________
606 void AliESDtrackCuts::SaveHistograms(Char_t* dir) {
607   // 
608   // saves the histograms in a directory (dir)
609   // 
610
611   
612   if (!fHistogramsOn) {
613     AliDebug(0, "Histograms not on - cannot save histograms!!!");
614     return;
615   }
616
617   gDirectory->mkdir(dir);
618   gDirectory->cd(dir);
619
620   gDirectory->mkdir("before_cuts");
621   gDirectory->mkdir("after_cuts");
622  
623   fhCutStatistics->Write();
624   fhCutCorrelation->Write();
625
626   for (Int_t i=0; i<2; i++) {
627     if (i==0)
628       gDirectory->cd("before_cuts");
629     else
630       gDirectory->cd("after_cuts");
631     
632     fhNClustersITS[i]        ->Write();
633     fhNClustersTPC[i]        ->Write();
634     fhChi2PerClusterITS[i]   ->Write();
635     fhChi2PerClusterTPC[i]   ->Write();
636     
637     fhC11[i]                 ->Write();
638     fhC22[i]                 ->Write();
639     fhC33[i]                 ->Write();
640     fhC44[i]                 ->Write();
641     fhC55[i]                 ->Write();
642
643     fhDXY[i]                 ->Write();
644     fhDZ[i]                  ->Write();
645     fhDXYvsDZ[i]             ->Write();
646     
647     fhDXYNormalized[i]       ->Write();
648     fhDZNormalized[i]        ->Write();
649     fhDXYvsDZNormalized[i]   ->Write();
650
651     gDirectory->cd("../");
652   }
653
654   gDirectory->cd("../");
655 }
656
657
658