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