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