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