]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ANALYSIS/AliAnalysisTrackCuts.cxx
new digitization and reconstruction corresponded to new data format
[u/mrichter/AliRoot.git] / ANALYSIS / AliAnalysisTrackCuts.cxx
1 /**************************************************************************
2  * Author: Panos Christakoglou.                                           *
3  * Contributors are mentioned in the code where appropriate.              *
4  *                                                                        *
5  * Permission to use, copy, modify and distribute this software and its   *
6  * documentation strictly for non-commercial purposes is hereby granted   *
7  * without fee, provided that the above copyright notice appears in all   *
8  * copies and that both the copyright notice and this permission notice   *
9  * appear in the supporting documentation. The authors make no claims     *
10  * about the suitability of this software for any purpose. It is          *
11  * provided "as is" without express or implied warranty.                  *
12  **************************************************************************/
13
14 /* $Id$ */
15
16 //-----------------------------------------------------------------
17 //           AliAnalysisTrackCuts class
18 //   This is the class to deal with the event and track level cuts
19 //   Origin: Panos Christakoglou, UOA-CERN, Panos.Christakoglou@cern.ch
20 //-----------------------------------------------------------------
21
22
23
24 //ROOT
25 #include <TPaveText.h>
26 #include <TText.h>
27 #include <TLine.h>
28 #include <TCanvas.h>
29
30 #include "AliLog.h"
31
32 #include "AliESDtrack.h"
33 #include "AliESD.h"
34
35 #include "AliAnalysisTrackCuts.h"
36
37 ClassImp(AliAnalysisTrackCuts)
38
39 //----------------------------------------//
40 AliAnalysisTrackCuts::AliAnalysisTrackCuts()
41 {
42   //Default constructor.
43   //Calls the Reset method.
44   Reset();
45 }
46
47 //----------------------------------------//
48 AliAnalysisTrackCuts::~AliAnalysisTrackCuts()
49 {
50   //Destructor.
51 }
52
53 //----------------------------------------//
54 void AliAnalysisTrackCuts::Reset()
55 {
56   //Assigns dummy values to every private member.
57   fPxMin = -1000.0;
58   fPxMax = 1000.0; 
59   fPyMin = -1000.0;
60   fPyMax = 1000.0;  
61   fPzMin = -1000.0;
62   fPzMax = 1000.0;
63   fPtMin = 0.0;
64   fPtMax = 1000.0;
65   fPMin = 0.0;
66   fPMax = 1000.0;
67   fBrMin = 0.0;
68   fBrMax = 1000.0;
69   fBzMin = 0.0;
70   fBzMax = 1000.0;
71   fEtaMin = -100.0;
72   fEtaMax = 100.0;
73   fRapMin = -100.0;
74   fRapMax = 100.0;
75   
76   fP = 0;  
77   fPt = 0; 
78   fPx = 0;  
79   fPy = 0; 
80   fPz = 0; 
81   fbr = 0; 
82   fbz = 0; 
83   fEta = 0;  
84   fRap = 0; 
85   fTotalTracks = 0; 
86   fAcceptedTracks = 0; 
87
88   fFlagP = 0;  
89   fFlagPt = 0;  
90   fFlagPx = 0;  
91   fFlagPy = 0;  
92   fFlagPz = 0;  
93   fFlagEta = 0;  
94   fFlagRap = 0;  
95   fFlagbr = 0;  
96   fFlagbz = 0;  
97 }
98
99 //----------------------------------------//
100 void AliAnalysisTrackCuts::SetPxRange(Float_t r1, Float_t r2)
101 {
102   //Sets the range for the momentum x component. 
103   fPxMin = r1;
104   fPxMax = r2;
105   fFlagPx = 1;
106 }
107
108 //----------------------------------------//
109 void AliAnalysisTrackCuts::SetPyRange(Float_t r1, Float_t r2)
110 {
111   //Sets the range for the momentum y component. 
112   fPyMin = r1;
113   fPyMax = r2; 
114   fFlagPy = 1;
115 }
116
117 //----------------------------------------//
118 void AliAnalysisTrackCuts::SetPzRange(Float_t r1, Float_t r2)
119 {
120   //Sets the range for the momentum z component. 
121   fPzMin = r1;
122   fPzMax = r2; 
123   fFlagPy = 1;
124 }
125
126 //----------------------------------------//
127 void AliAnalysisTrackCuts::SetPRange(Float_t r1, Float_t r2)
128 {
129   //Sets the range for the momentum. 
130   fPMin = r1;
131   fPMax = r2; 
132   fFlagPz = 1;
133 }
134
135 //----------------------------------------//
136 void AliAnalysisTrackCuts::SetPtRange(Float_t r1, Float_t r2)
137 {
138   //Sets the range for the teransverse momentum. 
139   fPtMin = r1;
140   fPtMax = r2; 
141   fFlagPt = 1;
142 }
143
144 //----------------------------------------//
145 void AliAnalysisTrackCuts::SetBrRange(Float_t r1, Float_t r2)
146 {
147   //Sets the range of the closest approach of the track 
148   //to the primary vertex in the r-phi plane. 
149   fBrMin = r1;
150   fBrMax = r2; 
151   fFlagbr = 1;
152 }
153
154 //----------------------------------------//
155 void AliAnalysisTrackCuts::SetBzRange(Float_t r1, Float_t r2)
156 {
157   //Sets the range of the closest approach of the track 
158   //to the primary vertex in the beam axis. 
159   fBzMin = r1;
160   fBzMax = r2; 
161   fFlagbz = 1;
162 }
163
164 //----------------------------------------//
165 void AliAnalysisTrackCuts::SetEtaRange(Float_t r1, Float_t r2)
166 {
167   //Sets the range of the pseudo-rapidity. 
168   fEtaMin = r1;
169   fEtaMax = r2; 
170   fFlagEta = 1;
171 }
172
173 //----------------------------------------//
174 void AliAnalysisTrackCuts::SetRapRange(Float_t r1, Float_t r2)
175 {
176   //Sets the range of the rapidity. 
177   fRapMin = r1;
178   fRapMax = r2; 
179   fFlagRap = 1;
180 }
181
182 //----------------------------------------//
183 void AliAnalysisTrackCuts::GetTrackStats()
184 {
185   //Gets the statistics.
186   //fTotalTracks is the total number of tracks.
187   //fAcceptedTracks is the number of accepted tracks after the cuts. 
188   AliInfo(Form("Total number of tracks: %d",fTotalTracks));
189   AliInfo(Form("Total number of accepted tracks: %d",fAcceptedTracks)); 
190 }
191
192 //----------------------------------------//
193 void AliAnalysisTrackCuts::GetPStats()
194 {
195   //Gets the momentum statistics.
196   //Prints the percentage of tracks rejected due to this cut. 
197    AliInfo(Form("P range: [%f,%f]",fPMin,fPMax));
198   if(fTotalTracks != 0)
199     AliInfo(Form("Tracks rejected: %f",100.0*fP/fTotalTracks)); 
200 }
201
202 //----------------------------------------//
203 void AliAnalysisTrackCuts::GetPtStats()
204 {
205   //Gets the transverse momentum statistics.
206   //Prints the percentage of tracks rejected due to this cut. 
207   AliInfo(Form("Pt range: [%f,%f]",fPtMin,fPtMax));
208   if(fTotalTracks != 0) 
209     AliInfo(Form("Tracks rejected: %f",100.0*fPt/fTotalTracks)); 
210 }
211
212 //----------------------------------------//
213 void AliAnalysisTrackCuts::GetPxStats()
214 {
215   //Gets the x momentum statistics.
216   //Prints the percentage of tracks rejected due to this cut. 
217   AliInfo(Form("Px range: [%f,%f]",fPxMin,fPxMax));
218   if(fTotalTracks != 0) 
219     AliInfo(Form("Tracks rejected: %f",100.0*fPx/fTotalTracks)); 
220 }
221
222 //----------------------------------------//
223 void AliAnalysisTrackCuts::GetPyStats()
224 {
225   //Gets the y momentum statistics.
226   //Prints the percentage of tracks rejected due to this cut. 
227   AliInfo(Form("Py range: [%f,%f]",fPyMin,fPyMax));
228   if(fTotalTracks != 0) 
229     AliInfo(Form("Tracks rejected: %f",100.0*fPy/fTotalTracks)); 
230 }
231
232 //----------------------------------------//
233 void AliAnalysisTrackCuts::GetPzStats()
234 {
235   //Gets the z momentum statistics.
236   //Prints the percentage of tracks rejected due to this cut. 
237   AliInfo(Form("Pz range: [%f,%f]",fPzMin,fPzMax));
238   if(fTotalTracks != 0) 
239     AliInfo(Form("Tracks rejected: %f",100.0*fPz/fTotalTracks)); 
240 }
241
242 //----------------------------------------//
243 void AliAnalysisTrackCuts::GetEtaStats()
244 {
245   //Gets the pseudo-rapidity statistics.
246   //Prints the percentage of tracks rejected due to this cut. 
247   AliInfo(Form("eta range: [%f,%f]",fEtaMin,fEtaMax));
248   if(fTotalTracks != 0)
249     AliInfo(Form("Tracks rejected: %f",100.0*fEta/fTotalTracks)); 
250 }
251
252 //----------------------------------------//
253 void AliAnalysisTrackCuts::GetRapStats()
254 {
255   //Gets the rapidity statistics.
256   //Prints the percentage of tracks rejected due to this cut. 
257   AliInfo(Form("y range: [%f,%f]",fRapMin,fRapMax));
258   if(fTotalTracks != 0)
259     AliInfo(Form("Tracks rejected: %f",100.0*fRap/fTotalTracks)); 
260 }
261
262 //----------------------------------------//
263 void AliAnalysisTrackCuts::GetBrStats()
264 {
265   //Gets the statistics fro the closest distance of 
266   //the track to the primary vertex in the r-phi plane.
267   //Prints the percentage of tracks rejected due to this cut. 
268   AliInfo(Form("br range: [%f,%f]",fBrMin,fBrMax));
269   if(fTotalTracks != 0) 
270     AliInfo(Form("Tracks rejected: %f",100.0*fbr/fTotalTracks)); 
271 }
272
273 //----------------------------------------//
274 void AliAnalysisTrackCuts::GetBzStats()
275 {
276   //Gets the statistics fro the closest distance of 
277   //the track to the primary vertex in the beam axis.
278   //Prints the percentage of tracks rejected due to this cut. 
279   AliInfo(Form("bz range: [%f,%f]",fBzMin,fBzMax));
280   if(fTotalTracks != 0)
281     AliInfo(Form("Tracks rejected: %f",100.0*fbz/fTotalTracks)); 
282 }
283
284
285 //----------------------------------------//
286 Bool_t AliAnalysisTrackCuts::IsAccepted(AliESD *esd ,AliESDtrack *esdtrack)
287 {
288   //Returns true if the tracks is accepted otherwise false.
289   fTotalTracks++;
290   
291   //momentum related calculations
292   Double_t p[3];
293   esdtrack->GetPxPyPz(p);
294   Float_t momentum = sqrt(pow(p[0],2) + pow(p[1],2) + pow(p[2],2));
295   Float_t pt = sqrt(pow(p[0],2) + pow(p[1],2));
296   Float_t energy = sqrt(pow(esdtrack->GetMass(),2) + pow(momentum,2));
297
298   //y-eta related calculations
299   Float_t eta = 0.5*log((momentum - p[2])/(momentum + p[2]));
300   Float_t y = 0.5*log((energy - p[2])/(energy + p[2]));
301  
302   //impact parameter related calculations
303   Double_t trackPosition[3];
304   esdtrack->GetXYZ(trackPosition);
305   const AliESDVertex * vertexIn = esd->GetVertex();
306   Double_t vertexPosition[3];
307   vertexIn->GetXYZ(vertexPosition);                 
308   for (Int_t ii=0; ii<3; ii++) trackPosition[ii] -= vertexPosition[ii];
309                     
310   Float_t br = Float_t(TMath::Sqrt(pow(trackPosition[0],2) + pow(trackPosition[1],2)));
311   Float_t bz = Float_t(TMath::Abs(trackPosition[2]));
312  
313   if((momentum < fPMin) || (momentum > fPMax))
314     {
315       fP++;
316       return kFALSE;
317     }
318   if((pt < fPtMin) || (pt > fPtMax))
319     {
320       fPt++;
321       return kFALSE;
322     }
323   if((p[0] < fPxMin) || (p[0] > fPxMax))
324     {
325       fPx++;
326       return kFALSE;
327     }
328   if((p[1] < fPyMin) || (p[1] > fPyMax))
329     {
330       fPy++;
331       return kFALSE;
332     }
333   if((p[2] < fPzMin) || (p[2] > fPzMax))
334     {
335       fPz++;
336       return kFALSE;
337     } 
338   if((br < fBrMin) || (br > fBrMax))
339     {
340       fbr++;
341       return kFALSE;
342     }
343   if((bz < fBzMin) || (bz > fBzMax))
344     {
345       fbz++;
346       return kFALSE;
347     }
348   if((eta < fEtaMin) || (eta > fEtaMax))
349     {
350       fEta++;
351       return kFALSE;
352     }
353   if((y < fRapMin) || (y > fRapMax))
354     {
355       fRap++;
356       return kFALSE;
357     }
358   
359   fAcceptedTracks++;
360   
361   return kTRUE;
362 }
363
364
365 //----------------------------------------//
366 TPaveText *AliAnalysisTrackCuts::GetTrackCuts()
367 {
368   //Shows a TPaveText with all the track cuts stats.
369   TCanvas *ccuts2 = new TCanvas("ccuts2","Track cuts",410,10,400,400);
370   ccuts2->SetFillColor(10);
371   ccuts2->SetHighLightColor(10);
372
373   TPaveText *pave = new TPaveText(0.01,0.01,0.98,0.98);
374   pave->SetFillColor(3);
375   Char_t cutName[256];
376  
377   TLine *l1 = pave->AddLine(0,0.89,1,0.89);
378   l1->SetLineWidth(2);
379   TLine *l2 = pave->AddLine(0,0.79,1,0.79);
380   l2->SetLineWidth(2);
381   TLine *l3 = pave->AddLine(0,0.69,1,0.69);
382   l3->SetLineWidth(2);
383   TLine *l4 = pave->AddLine(0,0.59,1,0.59);
384   l4->SetLineWidth(2);
385   TLine *l5 = pave->AddLine(0,0.49,1,0.49);
386   l5->SetLineWidth(2);
387   TLine *l6 = pave->AddLine(0,0.39,1,0.39);
388   l6->SetLineWidth(2);
389   TLine *l7 = pave->AddLine(0,0.29,1,0.29);
390   l7->SetLineWidth(2);
391   TLine *l8 = pave->AddLine(0,0.19,1,0.19);
392   l8->SetLineWidth(2);
393   TLine *l9 = pave->AddLine(0,0.09,1,0.09);
394   l9->SetLineWidth(2);
395
396   sprintf(cutName,"Total number of tracks: %d",fTotalTracks);
397   TText *t1 = pave->AddText(cutName);
398   t1->SetTextColor(1);
399   t1->SetTextSize(0.04);
400   t1->SetTextAlign(11);
401  
402   sprintf(cutName,"Total number of accepted tracks: %d",fAcceptedTracks);
403   t1 = pave->AddText(cutName);
404   t1->SetTextColor(1);
405   t1->SetTextSize(0.04);
406   t1->SetTextAlign(11);
407  
408   sprintf(cutName,"P range: [%f,%f]",fPMin,fPMax);
409   t1 = pave->AddText(cutName);
410   t1->SetTextColor(1);
411   t1->SetTextSize(0.04);
412   t1->SetTextAlign(11);
413   sprintf(cutName,"Tracks rejected: %f",100.0*fP/fTotalTracks);
414   t1 = pave->AddText(cutName);
415   t1->SetTextColor(1);
416   t1->SetTextSize(0.04);
417   t1->SetTextAlign(11);
418  
419   sprintf(cutName,"Pt range: [%f,%f]",fPtMin,fPtMax);
420   t1 = pave->AddText(cutName);
421   t1->SetTextColor(1);
422   t1->SetTextSize(0.04);
423   t1->SetTextAlign(11);
424   sprintf(cutName,"Tracks rejected: %f",100.0*fPt/fTotalTracks);
425   t1 = pave->AddText(cutName);
426   t1->SetTextColor(1);
427   t1->SetTextSize(0.04);
428   t1->SetTextAlign(11);
429
430   sprintf(cutName,"Px range: [%f,%f]",fPxMin,fPxMax);
431   t1 = pave->AddText(cutName);
432   t1->SetTextColor(1);
433   t1->SetTextSize(0.04);
434   t1->SetTextAlign(11);
435   sprintf(cutName,"Tracks rejected: %f",100.0*fPx/fTotalTracks);
436   t1 = pave->AddText(cutName);
437   t1->SetTextColor(1);
438   t1->SetTextSize(0.04);
439   t1->SetTextAlign(11);
440  
441   sprintf(cutName,"Py range: [%f,%f]",fPyMin,fPyMax);
442   t1 = pave->AddText(cutName);
443   t1->SetTextColor(1);
444   t1->SetTextSize(0.04);
445   t1->SetTextAlign(11);
446   sprintf(cutName,"Tracks rejected: %f",100.0*fPy/fTotalTracks);
447   t1 = pave->AddText(cutName);
448   t1->SetTextColor(1);
449   t1->SetTextSize(0.04);
450   t1->SetTextAlign(11);
451  
452   sprintf(cutName,"Pz range: [%f,%f]",fPzMin,fPzMax);
453   t1 = pave->AddText(cutName);
454   t1->SetTextColor(1);
455   t1->SetTextSize(0.04);
456   t1->SetTextAlign(11);
457   sprintf(cutName,"Tracks rejected: %f",100.0*fPz/fTotalTracks);
458   t1 = pave->AddText(cutName);
459   t1->SetTextColor(1);
460   t1->SetTextSize(0.04);
461   t1->SetTextAlign(11);
462   
463   sprintf(cutName,"br range: [%f,%f]",fBrMin,fBrMax);
464   t1 = pave->AddText(cutName);
465   t1->SetTextColor(1);
466   t1->SetTextSize(0.04);
467   t1->SetTextAlign(11);
468   sprintf(cutName,"Tracks rejected: %f",100.0*fbr/fTotalTracks);
469   t1 = pave->AddText(cutName);
470   t1->SetTextColor(1);
471   t1->SetTextSize(0.04);
472   t1->SetTextAlign(11);
473  
474   sprintf(cutName,"bz range: [%f,%f]",fBzMin,fBzMax);
475   t1 = pave->AddText(cutName);
476   t1->SetTextColor(1);
477   t1->SetTextSize(0.04);
478   t1->SetTextAlign(11);
479   sprintf(cutName,"Tracks rejected: %f",100.0*fbz/fTotalTracks);
480   t1 = pave->AddText(cutName);
481   t1->SetTextColor(1);
482   t1->SetTextSize(0.04);
483   t1->SetTextAlign(11);
484
485   sprintf(cutName,"eta range: [%f,%f]",fEtaMin,fEtaMax);
486   t1 = pave->AddText(cutName);
487   t1->SetTextColor(1);
488   t1->SetTextSize(0.04);
489   t1->SetTextAlign(11);
490   sprintf(cutName,"Tracks rejected: %f",100.0*fEta/fTotalTracks);
491   t1 = pave->AddText(cutName);
492   t1->SetTextColor(1);
493   t1->SetTextSize(0.04);
494   t1->SetTextAlign(11);
495
496   sprintf(cutName,"y range: [%f,%f]",fRapMin,fRapMax);
497   t1 = pave->AddText(cutName);
498   t1->SetTextColor(1);
499   t1->SetTextSize(0.04);
500   t1->SetTextAlign(11);
501   sprintf(cutName,"Tracks rejected: %f",100.0*fRap/fTotalTracks);
502   t1 = pave->AddText(cutName);
503   t1->SetTextColor(1);
504   t1->SetTextSize(0.04);
505   t1->SetTextAlign(11);
506
507   return pave;
508 }
509
510 //----------------------------------------//
511 void AliAnalysisTrackCuts::PrintTrackCuts()
512 {
513   //Prints the track cut stats.
514   //GetTrackCuts()->Draw();
515
516   AliInfo(Form("**************TRACK CUTS**************"));
517   GetTrackStats();
518   if(fFlagP)
519     GetPStats();
520   if(fFlagPt)
521     GetPtStats();
522   if(fFlagPx)
523     GetPxStats();
524   if(fFlagPy)
525     GetPyStats();
526   if(fFlagPz)
527     GetPzStats();
528   if(fFlagEta)
529     GetEtaStats();
530   if(fFlagRap)
531     GetRapStats();
532   if(fFlagbr)
533     GetBrStats();
534   if(fFlagbz)
535     GetBzStats(); 
536   AliInfo(Form("**************************************"));
537 }