]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ANALYSIS/AliAnalysisTrackCuts.cxx
Updated flags
[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 <TROOT.h>
26 #include <TObject.h>
27 #include <TSystem.h>
28 #include <TObject.h>
29 #include <TPaveText.h>
30 #include <TText.h>
31 #include <TLine.h>
32 #include <TCanvas.h>
33
34 #include "AliLog.h"
35 #include "AliESDtrack.h"
36 #include "AliESD.h"
37
38 #include "AliAnalysisTrackCuts.h"
39
40 ClassImp(AliAnalysisTrackCuts)
41
42 //----------------------------------------//
43 AliAnalysisTrackCuts::AliAnalysisTrackCuts()
44 {
45   Reset();
46 }
47
48 //----------------------------------------//
49 AliAnalysisTrackCuts::~AliAnalysisTrackCuts()
50 {
51 }
52
53 //----------------------------------------//
54 void AliAnalysisTrackCuts::Reset()
55 {
56   fPxMin = -1000.0;
57   fPxMax = 1000.0; 
58   fPyMin = -1000.0;
59   fPyMax = 1000.0;  
60   fPzMin = -1000.0;
61   fPzMax = 1000.0;
62   fPtMin = 0.0;
63   fPtMax = 1000.0;
64   fPMin = 0.0;
65   fPMax = 1000.0;
66   fBrMin = 0.0;
67   fBrMax = 1000.0;
68   fBzMin = 0.0;
69   fBzMax = 1000.0;
70   fEtaMin = -100.0;
71   fEtaMax = 100.0;
72   fRapMin = -100.0;
73   fRapMax = 100.0;
74   
75   fP = 0;  
76   fPt = 0; 
77   fPx = 0;  
78   fPy = 0; 
79   fPz = 0; 
80   fbr = 0; 
81   fbz = 0; 
82   fEta = 0;  
83   fRap = 0; 
84   fTotalTracks = 0; 
85   fAcceptedTracks = 0; 
86
87   fFlagP = 0;  
88   fFlagPt = 0;  
89   fFlagPx = 0;  
90   fFlagPy = 0;  
91   fFlagPz = 0;  
92   fFlagEta = 0;  
93   fFlagRap = 0;  
94   fFlagbr = 0;  
95   fFlagbz = 0;  
96 }
97
98 //----------------------------------------//
99 void AliAnalysisTrackCuts::SetPxRange(Float_t r1, Float_t r2)
100 {
101   fPxMin = r1;
102   fPxMax = r2;
103   fFlagPx = 1;
104 }
105
106 //----------------------------------------//
107 void AliAnalysisTrackCuts::SetPyRange(Float_t r1, Float_t r2)
108 {
109   fPyMin = r1;
110   fPyMax = r2; 
111   fFlagPy = 1;
112 }
113
114 //----------------------------------------//
115 void AliAnalysisTrackCuts::SetPzRange(Float_t r1, Float_t r2)
116 {
117   fPzMin = r1;
118   fPzMax = r2; 
119   fFlagPy = 1;
120 }
121
122 //----------------------------------------//
123 void AliAnalysisTrackCuts::SetPRange(Float_t r1, Float_t r2)
124 {
125   fPMin = r1;
126   fPMax = r2; 
127   fFlagPz = 1;
128 }
129
130 //----------------------------------------//
131 void AliAnalysisTrackCuts::SetPtRange(Float_t r1, Float_t r2)
132 {
133   fPtMin = r1;
134   fPtMax = r2; 
135   fFlagPt = 1;
136 }
137
138 //----------------------------------------//
139 void AliAnalysisTrackCuts::SetBrRange(Float_t r1, Float_t r2)
140 {
141   fBrMin = r1;
142   fBrMax = r2; 
143   fFlagbr = 1;
144 }
145
146 //----------------------------------------//
147 void AliAnalysisTrackCuts::SetBzRange(Float_t r1, Float_t r2)
148 {
149   fBzMin = r1;
150   fBzMax = r2; 
151   fFlagbz = 1;
152 }
153
154 //----------------------------------------//
155 void AliAnalysisTrackCuts::SetEtaRange(Float_t r1, Float_t r2)
156 {
157   fEtaMin = r1;
158   fEtaMax = r2; 
159   fFlagEta = 1;
160 }
161
162 //----------------------------------------//
163 void AliAnalysisTrackCuts::SetRapRange(Float_t r1, Float_t r2)
164 {
165   fRapMin = r1;
166   fRapMax = r2; 
167   fFlagRap = 1;
168 }
169
170 //----------------------------------------//
171 void AliAnalysisTrackCuts::GetTrackStats()
172 {
173   AliInfo(Form("Total number of tracks: %d",fTotalTracks));
174   AliInfo(Form("Total number of accepted tracks: %d",fAcceptedTracks)); 
175 }
176
177 //----------------------------------------//
178 void AliAnalysisTrackCuts::GetPStats()
179 {
180   AliInfo(Form("P range: [%f,%f]",fPMin,fPMax));
181   if(fTotalTracks != 0)
182     AliInfo(Form("Tracks rejected: %f",100.0*fP/fTotalTracks)); 
183 }
184
185 //----------------------------------------//
186 void AliAnalysisTrackCuts::GetPtStats()
187 {
188   AliInfo(Form("Pt range: [%f,%f]",fPtMin,fPtMax));
189   if(fTotalTracks != 0) 
190     AliInfo(Form("Tracks rejected: %f",100.0*fPt/fTotalTracks)); 
191 }
192
193 //----------------------------------------//
194 void AliAnalysisTrackCuts::GetPxStats()
195 {
196   AliInfo(Form("Px range: [%f,%f]",fPxMin,fPxMax));
197   if(fTotalTracks != 0) 
198     AliInfo(Form("Tracks rejected: %f",100.0*fPx/fTotalTracks)); 
199 }
200
201 //----------------------------------------//
202 void AliAnalysisTrackCuts::GetPyStats()
203 {
204   AliInfo(Form("Py range: [%f,%f]",fPyMin,fPyMax));
205   if(fTotalTracks != 0) 
206     AliInfo(Form("Tracks rejected: %f",100.0*fPy/fTotalTracks)); 
207 }
208
209 //----------------------------------------//
210 void AliAnalysisTrackCuts::GetPzStats()
211 {
212   AliInfo(Form("Pz range: [%f,%f]",fPzMin,fPzMax));
213   if(fTotalTracks != 0) 
214     AliInfo(Form("Tracks rejected: %f",100.0*fPz/fTotalTracks)); 
215 }
216
217 //----------------------------------------//
218 void AliAnalysisTrackCuts::GetEtaStats()
219 {
220   AliInfo(Form("eta range: [%f,%f]",fEtaMin,fEtaMax));
221   if(fTotalTracks != 0)
222     AliInfo(Form("Tracks rejected: %f",100.0*fEta/fTotalTracks)); 
223 }
224
225 //----------------------------------------//
226 void AliAnalysisTrackCuts::GetRapStats()
227 {
228   AliInfo(Form("y range: [%f,%f]",fRapMin,fRapMax));
229   if(fTotalTracks != 0)
230     AliInfo(Form("Tracks rejected: %f",100.0*fRap/fTotalTracks)); 
231 }
232
233 //----------------------------------------//
234 void AliAnalysisTrackCuts::GetBrStats()
235 {
236   AliInfo(Form("br range: [%f,%f]",fBrMin,fBrMax));
237   if(fTotalTracks != 0) 
238     AliInfo(Form("Tracks rejected: %f",100.0*fbr/fTotalTracks)); 
239 }
240
241 //----------------------------------------//
242 void AliAnalysisTrackCuts::GetBzStats()
243 {
244   AliInfo(Form("bz range: [%f,%f]",fBzMin,fBzMax));
245   if(fTotalTracks != 0)
246     AliInfo(Form("Tracks rejected: %f",100.0*fbz/fTotalTracks)); 
247 }
248
249
250 //----------------------------------------//
251 Bool_t AliAnalysisTrackCuts::IsAccepted(AliESD *esd ,AliESDtrack *esdtrack)
252 {
253   fTotalTracks++;
254   
255   //momentum related calculations
256   Double_t p[3];
257   esdtrack->GetPxPyPz(p);
258   Float_t P = sqrt(pow(p[0],2) + pow(p[1],2) + pow(p[2],2));
259   Float_t Pt = sqrt(pow(p[0],2) + pow(p[1],2));
260   Float_t E = sqrt(pow(esdtrack->GetMass(),2) + pow(P,2));
261
262   //y-eta related calculations
263   Float_t eta = 0.5*log((P - p[2])/(P + p[2]));
264   Float_t y = 0.5*log((E - p[2])/(E + p[2]));
265  
266   //impact parameter related calculations
267   Double_t TrackPosition[3];
268   esdtrack->GetXYZ(TrackPosition);
269   const AliESDVertex * VertexIn = esd->GetVertex();
270   Double_t VertexPosition[3];
271   VertexIn->GetXYZ(VertexPosition);                 
272   for (Int_t ii=0; ii<3; ii++) TrackPosition[ii] -= VertexPosition[ii];
273                     
274   Float_t br = Float_t(TMath::Sqrt(pow(TrackPosition[0],2) + pow(TrackPosition[1],2)));
275   Float_t bz = Float_t(TMath::Abs(TrackPosition[2]));
276  
277   if((P < fPMin) || (P > fPMax))
278     {
279       fP++;
280       return kFALSE;
281     }
282   if((Pt < fPtMin) || (Pt > fPtMax))
283     {
284       fPt++;
285       return kFALSE;
286     }
287   if((p[0] < fPxMin) || (p[0] > fPxMax))
288     {
289       fPx++;
290       return kFALSE;
291     }
292   if((p[1] < fPyMin) || (p[1] > fPyMax))
293     {
294       fPy++;
295       return kFALSE;
296     }
297   if((p[2] < fPzMin) || (p[2] > fPzMax))
298     {
299       fPz++;
300       return kFALSE;
301     } 
302   if((br < fBrMin) || (br > fBrMax))
303     {
304       fbr++;
305       return kFALSE;
306     }
307   if((bz < fBzMin) || (bz > fBzMax))
308     {
309       fbz++;
310       return kFALSE;
311     }
312   if((eta < fEtaMin) || (eta > fEtaMax))
313     {
314       fEta++;
315       return kFALSE;
316     }
317   if((y < fRapMin) || (y > fRapMax))
318     {
319       fRap++;
320       return kFALSE;
321     }
322   
323   fAcceptedTracks++;
324   
325   return kTRUE;
326 }
327
328
329 //----------------------------------------//
330 TPaveText *AliAnalysisTrackCuts::GetTrackCuts()
331 {
332   TCanvas *ccuts2 = new TCanvas("ccuts2","Track cuts",410,10,400,400);
333   ccuts2->SetFillColor(10);
334   ccuts2->SetHighLightColor(10);
335
336   TPaveText *pave = new TPaveText(0.01,0.01,0.98,0.98);
337   pave->SetFillColor(3);
338   Char_t CutName[256];
339  
340   TLine *l1 = pave->AddLine(0,0.89,1,0.89);
341   l1->SetLineWidth(2);
342   TLine *l2 = pave->AddLine(0,0.79,1,0.79);
343   l2->SetLineWidth(2);
344   TLine *l3 = pave->AddLine(0,0.69,1,0.69);
345   l3->SetLineWidth(2);
346   TLine *l4 = pave->AddLine(0,0.59,1,0.59);
347   l4->SetLineWidth(2);
348   TLine *l5 = pave->AddLine(0,0.49,1,0.49);
349   l5->SetLineWidth(2);
350   TLine *l6 = pave->AddLine(0,0.39,1,0.39);
351   l6->SetLineWidth(2);
352   TLine *l7 = pave->AddLine(0,0.29,1,0.29);
353   l7->SetLineWidth(2);
354   TLine *l8 = pave->AddLine(0,0.19,1,0.19);
355   l8->SetLineWidth(2);
356   TLine *l9 = pave->AddLine(0,0.09,1,0.09);
357   l9->SetLineWidth(2);
358
359   sprintf(CutName,"Total number of tracks: %d",fTotalTracks);
360   TText *t1 = pave->AddText(CutName);
361   t1->SetTextColor(1);
362   t1->SetTextSize(0.04);
363   t1->SetTextAlign(11);
364  
365   sprintf(CutName,"Total number of accepted tracks: %d",fAcceptedTracks);
366   t1 = pave->AddText(CutName);
367   t1->SetTextColor(1);
368   t1->SetTextSize(0.04);
369   t1->SetTextAlign(11);
370  
371   sprintf(CutName,"P range: [%f,%f]",fPMin,fPMax);
372   t1 = pave->AddText(CutName);
373   t1->SetTextColor(1);
374   t1->SetTextSize(0.04);
375   t1->SetTextAlign(11);
376   sprintf(CutName,"Tracks rejected: %f",100.0*fP/fTotalTracks);
377   t1 = pave->AddText(CutName);
378   t1->SetTextColor(1);
379   t1->SetTextSize(0.04);
380   t1->SetTextAlign(11);
381  
382   sprintf(CutName,"Pt range: [%f,%f]",fPtMin,fPtMax);
383   t1 = pave->AddText(CutName);
384   t1->SetTextColor(1);
385   t1->SetTextSize(0.04);
386   t1->SetTextAlign(11);
387   sprintf(CutName,"Tracks rejected: %f",100.0*fPt/fTotalTracks);
388   t1 = pave->AddText(CutName);
389   t1->SetTextColor(1);
390   t1->SetTextSize(0.04);
391   t1->SetTextAlign(11);
392
393   sprintf(CutName,"Px range: [%f,%f]",fPxMin,fPxMax);
394   t1 = pave->AddText(CutName);
395   t1->SetTextColor(1);
396   t1->SetTextSize(0.04);
397   t1->SetTextAlign(11);
398   sprintf(CutName,"Tracks rejected: %f",100.0*fPx/fTotalTracks);
399   t1 = pave->AddText(CutName);
400   t1->SetTextColor(1);
401   t1->SetTextSize(0.04);
402   t1->SetTextAlign(11);
403  
404   sprintf(CutName,"Py range: [%f,%f]",fPyMin,fPyMax);
405   t1 = pave->AddText(CutName);
406   t1->SetTextColor(1);
407   t1->SetTextSize(0.04);
408   t1->SetTextAlign(11);
409   sprintf(CutName,"Tracks rejected: %f",100.0*fPy/fTotalTracks);
410   t1 = pave->AddText(CutName);
411   t1->SetTextColor(1);
412   t1->SetTextSize(0.04);
413   t1->SetTextAlign(11);
414  
415   sprintf(CutName,"Pz range: [%f,%f]",fPzMin,fPzMax);
416   t1 = pave->AddText(CutName);
417   t1->SetTextColor(1);
418   t1->SetTextSize(0.04);
419   t1->SetTextAlign(11);
420   sprintf(CutName,"Tracks rejected: %f",100.0*fPz/fTotalTracks);
421   t1 = pave->AddText(CutName);
422   t1->SetTextColor(1);
423   t1->SetTextSize(0.04);
424   t1->SetTextAlign(11);
425   
426   sprintf(CutName,"br range: [%f,%f]",fBrMin,fBrMax);
427   t1 = pave->AddText(CutName);
428   t1->SetTextColor(1);
429   t1->SetTextSize(0.04);
430   t1->SetTextAlign(11);
431   sprintf(CutName,"Tracks rejected: %f",100.0*fbr/fTotalTracks);
432   t1 = pave->AddText(CutName);
433   t1->SetTextColor(1);
434   t1->SetTextSize(0.04);
435   t1->SetTextAlign(11);
436  
437   sprintf(CutName,"bz range: [%f,%f]",fBzMin,fBzMax);
438   t1 = pave->AddText(CutName);
439   t1->SetTextColor(1);
440   t1->SetTextSize(0.04);
441   t1->SetTextAlign(11);
442   sprintf(CutName,"Tracks rejected: %f",100.0*fbz/fTotalTracks);
443   t1 = pave->AddText(CutName);
444   t1->SetTextColor(1);
445   t1->SetTextSize(0.04);
446   t1->SetTextAlign(11);
447
448   sprintf(CutName,"eta range: [%f,%f]",fEtaMin,fEtaMax);
449   t1 = pave->AddText(CutName);
450   t1->SetTextColor(1);
451   t1->SetTextSize(0.04);
452   t1->SetTextAlign(11);
453   sprintf(CutName,"Tracks rejected: %f",100.0*fEta/fTotalTracks);
454   t1 = pave->AddText(CutName);
455   t1->SetTextColor(1);
456   t1->SetTextSize(0.04);
457   t1->SetTextAlign(11);
458
459   sprintf(CutName,"y range: [%f,%f]",fRapMin,fRapMax);
460   t1 = pave->AddText(CutName);
461   t1->SetTextColor(1);
462   t1->SetTextSize(0.04);
463   t1->SetTextAlign(11);
464   sprintf(CutName,"Tracks rejected: %f",100.0*fRap/fTotalTracks);
465   t1 = pave->AddText(CutName);
466   t1->SetTextColor(1);
467   t1->SetTextSize(0.04);
468   t1->SetTextAlign(11);
469
470   return pave;
471 }
472
473 //----------------------------------------//
474 void AliAnalysisTrackCuts::PrintTrackCuts()
475 {
476   //GetTrackCuts()->Draw();
477
478   AliInfo(Form("**************TRACK CUTS**************"));
479   GetTrackStats();
480   if(fFlagP)
481     GetPStats();
482   if(fFlagPt)
483     GetPtStats();
484   if(fFlagPx)
485     GetPxStats();
486   if(fFlagPy)
487     GetPyStats();
488   if(fFlagPz)
489     GetPzStats();
490   if(fFlagEta)
491     GetEtaStats();
492   if(fFlagRap)
493     GetRapStats();
494   if(fFlagbr)
495     GetBrStats();
496   if(fFlagbz)
497     GetBzStats(); 
498   AliInfo(Form("**************************************"));
499 }