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