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