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