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