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