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