]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ANALYSIS/AliAnalysisEventCuts.cxx
Change in the signature of base class method Use to avoid compilation warnings (icc)
[u/mrichter/AliRoot.git] / ANALYSIS / AliAnalysisEventCuts.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 //           AliAnalysisEventCuts 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 #include "AliESD.h"
32
33 #include "AliAnalysisEventCuts.h"
34
35 ClassImp(AliAnalysisEventCuts)
36
37 //----------------------------------------//
38 AliAnalysisEventCuts::AliAnalysisEventCuts()
39 {
40   //Default constructor which calls the Reset method.
41   Reset();
42 }
43
44 //----------------------------------------//
45 AliAnalysisEventCuts::~AliAnalysisEventCuts()
46 {
47   //Defaut destructor.
48 }
49
50 //----------------------------------------//
51 void AliAnalysisEventCuts::Reset()
52 {
53   //Sets dummy values to every private member.
54   fVxMin = -1000.0;
55   fVxMax = 1000.0; 
56   fVyMin = -1000.0;
57   fVyMax = 1000.0;  
58   fVzMin = -1000.0;
59   fVzMax = 1000.0;
60   fMultMin = 0;
61   fMultMax = 100000;
62   
63   fMult = 0;  
64   fVx = 0;  
65   fVy = 0; 
66   fVz = 0; 
67   fTotalEvents = 0; 
68   fAcceptedEvents = 0; 
69
70   fFlagMult = 0;  
71   fFlagVx = 0;  
72   fFlagVy = 0; 
73   fFlagVz = 0; 
74 }
75
76 //----------------------------------------//
77 void AliAnalysisEventCuts::SetPrimaryVertexXRange(Float_t r1, Float_t r2)
78 {
79   //Sets the primary vertex x range.
80   fVxMin = r1;
81   fVxMax = r2; 
82   fFlagVx = 1;  
83 }
84
85 //----------------------------------------//
86 void AliAnalysisEventCuts::SetPrimaryVertexYRange(Float_t r1, Float_t r2)
87 {
88   //Sets the primary vertex y range.
89   fVyMin = r1;
90   fVyMax = r2; 
91   fFlagVy = 1;
92 }
93
94 //----------------------------------------//
95 void AliAnalysisEventCuts::SetPrimaryVertexZRange(Float_t r1, Float_t r2)
96 {
97   //Sets the primary vertex z range.
98   fVzMin = r1;
99   fVzMax = r2; 
100   fFlagVz = 1;
101 }
102
103 //----------------------------------------//
104 void AliAnalysisEventCuts::SetMultiplicityRange(Int_t n1, Int_t n2)
105 {
106   //Sets the multiplicity range.
107   fMultMin = n1;
108   fMultMax = n2; 
109   fFlagMult = 1;
110 }
111
112
113 //----------------------------------------//
114 Bool_t AliAnalysisEventCuts::IsAccepted(AliESD *esd)
115 {
116   //Returns true if the events is accepted otherwise false.
117   fTotalEvents++;
118   if((esd->GetNumberOfTracks() < fMultMin) || (esd->GetNumberOfTracks() > fMultMax))
119     {
120       fMult++;
121       AliInfo(Form("Event rejected due to multiplicity cut"));
122       return kFALSE;
123     }
124   const AliESDVertex *esdvertex = esd->GetVertex();
125   if((esdvertex->GetXv() < fVxMin) || (esdvertex->GetXv() > fVxMax))
126     {
127       fVx++;
128       AliInfo(Form("Event rejected due to Vx cut"));
129       return kFALSE;
130     }
131   if((esdvertex->GetYv() < fVyMin) || (esdvertex->GetYv() > fVyMax))
132     {
133       fVy++;
134       AliInfo(Form("Event rejected due to Vy cut"));
135       return kFALSE;
136     }
137  if((esdvertex->GetZv() < fVzMin) || (esdvertex->GetZv() > fVzMax))
138     {
139       fVz++;
140       AliInfo(Form("Event rejected due to Vz cut"));
141       return kFALSE;
142     }
143  fAcceptedEvents++;
144
145  return kTRUE;
146 }
147
148
149 //----------------------------------------//
150 TPaveText *AliAnalysisEventCuts::GetEventCuts()
151 {
152   //Shows a TPaveText with all the event cut stats.
153   TCanvas *ccuts = new TCanvas("ccuts","Event cuts",10,10,400,400);
154   ccuts->SetFillColor(10);
155   ccuts->SetHighLightColor(10);
156
157   TPaveText *pave = new TPaveText(0.01,0.01,0.98,0.98);
158   pave->SetFillColor(5);
159   Char_t cutName[256];
160  
161   TLine *l1 = pave->AddLine(0,0.78,1,0.78);
162   l1->SetLineWidth(2);
163   TLine *l2 = pave->AddLine(0,0.58,1,0.58);
164   l2->SetLineWidth(2);
165   TLine *l3 = pave->AddLine(0,0.38,1,0.38);
166   l3->SetLineWidth(2);
167   TLine *l4 = pave->AddLine(0,0.18,1,0.18);
168   l4->SetLineWidth(2);
169  
170   sprintf(cutName,"Total number of events: %d",fTotalEvents);
171   TText *t1 = pave->AddText(cutName);
172   t1->SetTextColor(4);
173   t1->SetTextSize(0.04);
174   t1->SetTextAlign(11);
175  
176   sprintf(cutName,"Total number of accepted events: %d",fAcceptedEvents);
177   t1 = pave->AddText(cutName);
178   t1->SetTextColor(4);
179   t1->SetTextSize(0.04);
180   t1->SetTextAlign(11);
181  
182   sprintf(cutName,"Multiplicity range: [%d,%d]",fMultMin,fMultMax);
183   t1 = pave->AddText(cutName);
184   t1->SetTextColor(4);
185   t1->SetTextSize(0.04);
186   t1->SetTextAlign(11);
187   sprintf(cutName,"Events rejected: %d",fMult);
188   t1 = pave->AddText(cutName);
189   t1->SetTextColor(4);
190   t1->SetTextSize(0.04);
191   t1->SetTextAlign(11);
192  
193   sprintf(cutName,"Vx range: [%f,%f]",fVxMin,fVxMax);
194   t1 = pave->AddText(cutName);
195   t1->SetTextColor(4);
196   t1->SetTextSize(0.04);
197   t1->SetTextAlign(11);
198   sprintf(cutName,"Events rejected: %d",fVx);
199   t1 = pave->AddText(cutName);
200   t1->SetTextColor(4);
201   t1->SetTextSize(0.04);
202   t1->SetTextAlign(11);
203  
204   sprintf(cutName,"Vy range: [%f,%f]",fVyMin,fVyMax);
205   t1 = pave->AddText(cutName);
206   t1->SetTextColor(4);
207   t1->SetTextSize(0.04);
208   t1->SetTextAlign(11);
209   sprintf(cutName,"Events rejected: %d",fVy);
210   t1 = pave->AddText(cutName);
211   t1->SetTextColor(4);
212   t1->SetTextSize(0.04);
213   t1->SetTextAlign(11);
214  
215   sprintf(cutName,"Vz range: [%f,%f]",fVzMin,fVzMax);
216   t1 = pave->AddText(cutName);
217   t1->SetTextColor(4);
218   t1->SetTextSize(0.04);
219   t1->SetTextAlign(11);
220   sprintf(cutName,"Events rejected: %d",fVz);
221   t1 = pave->AddText(cutName);
222   t1->SetTextColor(4);
223   t1->SetTextSize(0.04);
224   t1->SetTextAlign(11);
225  
226   return pave;
227 }
228
229 //----------------------------------------//
230 void AliAnalysisEventCuts::GetEventStats()
231 {
232   //Returns the total event stats.
233   //fTotalEvents is the total number of events.
234   //fAcceptedEvents is the total number of accepted events.  
235   AliInfo(Form("Total number of events: %d",fTotalEvents));
236   AliInfo(Form("Total number of accepted events: %d",fAcceptedEvents)); 
237 }
238
239 //----------------------------------------//
240 void AliAnalysisEventCuts::GetMultStats()
241 {
242   //Gets the multiplicity statistics.
243   //Prints the percentage of events rejected due to this cut. 
244   AliInfo(Form("Multiplicity range: [%d,%d]",fMultMin,fMultMax));
245   AliInfo(Form("Events rejected: %d",fMult));
246 }
247
248 //----------------------------------------//
249 void AliAnalysisEventCuts::GetVxStats()
250 {
251   //Gets the Vx statistics.
252   //Prints the percentage of events rejected due to this cut. 
253   AliInfo(Form("Vx range: [%f,%f]",fVxMin,fVxMax));
254   AliInfo(Form("Events rejected: %d",fVx));
255 }
256
257 //----------------------------------------//
258 void AliAnalysisEventCuts::GetVyStats()
259 {
260   //Gets the Vy statistics.
261   //Prints the percentage of events rejected due to this cut. 
262   AliInfo(Form("Vy range: [%f,%f]",fVyMin,fVyMax));
263   AliInfo(Form("Events rejected: %d",fVy));
264 }
265
266 //----------------------------------------//
267 void AliAnalysisEventCuts::GetVzStats()
268 {
269   //Gets the Vz statistics.
270   //Prints the percentage of events rejected due to this cut. 
271   AliInfo(Form("Vz range: [%f,%f]",fVzMin,fVzMax));
272   AliInfo(Form("Events rejected: %d",fVz));
273 }
274
275 //----------------------------------------//
276 void AliAnalysisEventCuts::PrintEventCuts()
277 {
278   //Prints the event stats.
279   //GetEventCuts()->Draw();  
280
281   AliInfo(Form("**************EVENT CUTS**************"));
282   GetEventStats();
283   if(fFlagMult)
284     GetMultStats();
285   if(fFlagVx)
286     GetVxStats();
287   if(fFlagVy)
288     GetVyStats();
289   if(fFlagVz)
290     GetVzStats();
291   AliInfo(Form("**************************************"));
292 }
293
294
295