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