]> git.uio.no Git - u/mrichter/AliRoot.git/blame - PWG2/FORWARD/analysis2/AliFMDEventInspector.cxx
Added debug flags
[u/mrichter/AliRoot.git] / PWG2 / FORWARD / analysis2 / AliFMDEventInspector.cxx
CommitLineData
8565b10b 1#include "AliFMDEventInspector.h"
2#include "AliLog.h"
3#include "AliESDEvent.h"
4#include "AliMultiplicity.h"
5#include "AliAnalysisManager.h"
6#include "AliInputEventHandler.h"
7#include "AliTriggerAnalysis.h"
8#include "AliPhysicsSelection.h"
9#include "AliAODForwardMult.h"
10#include <TH1.h>
11#include <TList.h>
12#include <TDirectory.h>
13
14//====================================================================
15AliFMDEventInspector::AliFMDEventInspector()
16 : TNamed(),
17 fHEventsTr(0),
18 fHEventsTrVtx(0),
19 fHTriggers(0),
20 fLowFluxCut(1000),
21 fMaxVzErr(0.1),
22 fList(0),
23 fDebug(0)
24{
25}
26
27//____________________________________________________________________
28AliFMDEventInspector::AliFMDEventInspector(const char* name)
29 : TNamed("fmdEventInspector", name),
30 fHEventsTr(0),
31 fHEventsTrVtx(0),
32 fHTriggers(0),
33 fLowFluxCut(1000),
34 fMaxVzErr(0.1),
35 fList(0),
36 fDebug(0)
37{
38}
39
40//____________________________________________________________________
41AliFMDEventInspector::AliFMDEventInspector(const AliFMDEventInspector& o)
42 : TNamed(o),
43 fHEventsTr(o.fHEventsTr),
44 fHEventsTrVtx(o.fHEventsTrVtx),
45 fHTriggers(o.fHTriggers),
46 fLowFluxCut(o.fMaxVzErr),
47 fMaxVzErr(o.fMaxVzErr),
48 fList(o.fList),
49 fDebug(0)
50{
51}
52
53//____________________________________________________________________
54AliFMDEventInspector::~AliFMDEventInspector()
55{
56 if (fHEventsTr) delete fHEventsTr;
57 if (fHEventsTrVtx) delete fHEventsTrVtx;
58 if (fHTriggers) delete fHTriggers;
59 if (fList) delete fList;
60}
61//____________________________________________________________________
62AliFMDEventInspector&
63AliFMDEventInspector::operator=(const AliFMDEventInspector& o)
64{
65 TNamed::operator=(o);
66 fHEventsTr = o.fHEventsTr;
67 fHEventsTrVtx = o.fHEventsTrVtx;
68 fHTriggers = o.fHTriggers;
69 fLowFluxCut = o.fLowFluxCut;
70 fMaxVzErr = o.fMaxVzErr;
71 fDebug = o.fDebug;
72 fList = (o.fList ? new TList : 0);
73 if (fList) {
74 fList->SetName(GetName());
75 if (fHEventsTr) fList->Add(fHEventsTr);
76 if (fHEventsTrVtx) fList->Add(fHEventsTrVtx);
77 if (fHTriggers) fList->Add(fHTriggers);
78 }
79 return *this;
80}
81
82//____________________________________________________________________
83Bool_t
84AliFMDEventInspector::FetchHistograms(TList* d,
85 TH1I*& hEventsTr,
86 TH1I*& hEventsTrVtx,
87 TH1I*& hTriggers) const
88{
89 hEventsTr = 0;
90 hEventsTrVtx = 0;
91 hTriggers = 0;
92 TList* dd = dynamic_cast<TList*>(d->FindObject(GetName()));
93 if (!dd) return kFALSE;
94
95 hEventsTr = dynamic_cast<TH1I*>(dd->FindObject("nEventsTr"));
96 hEventsTrVtx = dynamic_cast<TH1I*>(dd->FindObject("nEventsTrVtx"));
97 hTriggers = dynamic_cast<TH1I*>(dd->FindObject("triggers"));
98
99 if (!hEventsTr || !hEventsTrVtx || !hTriggers) return kFALSE;
100 return kTRUE;
101}
102//____________________________________________________________________
103void
104AliFMDEventInspector::Init(const TAxis& vtxAxis)
105{
106 fHEventsTr = new TH1I("nEventsTr", "Number of events w/trigger",
107 vtxAxis.GetNbins(),
108 vtxAxis.GetXmin(),
109 vtxAxis.GetXmax());
110 fHEventsTr->SetXTitle("v_{z} [cm]");
111 fHEventsTr->SetYTitle("# of events");
112 fHEventsTr->SetFillColor(kRed+1);
113 fHEventsTr->SetFillStyle(3001);
114 fHEventsTr->SetDirectory(0);
115 // fHEventsTr->Sumw2();
116 fList->Add(fHEventsTr);
117
118 fHEventsTrVtx = new TH1I("nEventsTrVtx",
119 "Number of events w/trigger and vertex",
120 vtxAxis.GetNbins(),
121 vtxAxis.GetXmin(),
122 vtxAxis.GetXmax());
123 fHEventsTrVtx->SetXTitle("v_{z} [cm]");
124 fHEventsTrVtx->SetYTitle("# of events");
125 fHEventsTrVtx->SetFillColor(kBlue+1);
126 fHEventsTrVtx->SetFillStyle(3001);
127 fHEventsTrVtx->SetDirectory(0);
128 // fHEventsTrVtx->Sumw2();
129 fList->Add(fHEventsTrVtx);
130
131
132 fHTriggers = new TH1I("triggers", "Triggers", 10, 0, 10);
133 fHTriggers->SetFillColor(kRed+1);
134 fHTriggers->SetFillStyle(3001);
135 fHTriggers->SetStats(0);
136 fHTriggers->SetDirectory(0);
137 fHTriggers->GetXaxis()->SetBinLabel(kInel +1,"INEL");
138 fHTriggers->GetXaxis()->SetBinLabel(kInelGt0+1,"INEL>0");
139 fHTriggers->GetXaxis()->SetBinLabel(kNSD +1,"NSD");
140 fHTriggers->GetXaxis()->SetBinLabel(kEmpty +1,"Empty");
141 fHTriggers->GetXaxis()->SetBinLabel(kA +1,"A");
142 fHTriggers->GetXaxis()->SetBinLabel(kB +1,"B");
143 fHTriggers->GetXaxis()->SetBinLabel(kC +1,"C");
144 fHTriggers->GetXaxis()->SetBinLabel(kE +1,"E");
145 fHTriggers->GetXaxis()->SetBinLabel(9, "spare1");
146 fHTriggers->GetXaxis()->SetBinLabel(10, "spare2");
147 fList->Add(fHTriggers);
148}
149
150//____________________________________________________________________
151void
152AliFMDEventInspector::DefineOutput(TList* dir)
153{
154 fList = new TList;
155 fList->SetName(GetName());
156 dir->Add(fList);
157}
158
159//____________________________________________________________________
160UInt_t
161AliFMDEventInspector::Process(const AliESDEvent* event,
162 UInt_t& triggers,
163 Bool_t& lowFlux,
164 Int_t& ivz,
165 Double_t& vz)
166{
167 // Check that we have an event
168 if (!event) {
169 AliWarning("No ESD event found for input event");
170 return kNoEvent;
171 }
172
173 // Read trigger information from the ESD and store in AOD object
174 if (!ReadTriggers(event, triggers)) {
175 if (fDebug > 2)
176 AliWarning("Failed to read triggers from ESD");
177 return kNoTriggers;
178 }
179
180 // Check if this is a high-flux event
181 const AliMultiplicity* testmult = event->GetMultiplicity();
182 if (!testmult) {
183 if (fDebug > 3)
184 AliWarning("No central multiplicity object found");
185 return kNoSPD;
186 }
187 lowFlux = testmult->GetNumberOfTracklets() < fLowFluxCut;
188
189 // Check the FMD ESD data
190 if (!event->GetFMDData()) {
191 if (fDebug > 3)
192 AliWarning("No FMD data found in ESD");
193 return kNoFMD;
194 }
195
196 // Get the vertex information
197 vz = 0;
198 Bool_t vzOk = ReadVertex(event, vz);
199
200 fHEventsTr->Fill(vz);
201 if (!vzOk) {
202 if (fDebug > 3)
203 AliWarning("Failed to read vertex from ESD");
204 return kNoVertex;
205 }
206 fHEventsTrVtx->Fill(vz);
207
208 // Get the vertex bin
209 ivz = fHEventsTr->GetXaxis()->FindBin(vz)-1;
210 if (ivz < 0 || ivz >= fHEventsTr->GetXaxis()->GetNbins()) {
211 if (fDebug > 3)
212 AliWarning(Form("Vertex @ %f outside of range [%f,%f]",
213 vz, fHEventsTr->GetXaxis()->GetXmin(),
214 fHEventsTr->GetXaxis()->GetXmax()));
215 ivz = -1;
216 return kBadVertex;
217 }
218 return kOk;
219}
220
221//____________________________________________________________________
222Bool_t
223AliFMDEventInspector::ReadTriggers(const AliESDEvent* esd, UInt_t& triggers)
224{
225 triggers = 0;
226
227 // Get the analysis manager - should always be there
228 AliAnalysisManager* am = AliAnalysisManager::GetAnalysisManager();
229 if (!am) {
230 AliWarning("No analysis manager defined!");
231 return kFALSE;
232 }
233
234 // Get the input handler - should always be there
235 AliInputEventHandler* ih =
236 static_cast<AliInputEventHandler*>(am->GetInputEventHandler());
237 if (!ih) {
238 AliWarning("No input handler");
239 return kFALSE;
240 }
241
242 // Get the physics selection - add that by using the macro
243 // AddTaskPhysicsSelection.C
244 AliPhysicsSelection* ps =
245 static_cast<AliPhysicsSelection*>(ih->GetEventSelection());
246 if (!ps) {
247 AliWarning("No physics selection");
248 return kFALSE;
249 }
250
251 // Check if this is a collision candidate (INEL)
252 Bool_t inel = ps->IsCollisionCandidate(esd);
253 if (inel) {
254 triggers |= AliAODForwardMult::kInel;
255 fHTriggers->Fill(kInel+0.5);
256 }
257
258 // IF this is inel, see if we have a tracklet
259 if (inel) {
260 const AliMultiplicity* spdmult = esd->GetMultiplicity();
261 if (!spdmult) {
262 AliWarning("No SPD multiplicity");
263 }
264 else {
265 Int_t n = spdmult->GetNumberOfTracklets();
266 for (Int_t j = 0; j < n; j++) {
267 if(TMath::Abs(spdmult->GetEta(j)) < 1) {
268 triggers |= AliAODForwardMult::kInelGt0;
269 fHTriggers->Fill(kInelGt0+.5);
270 break;
271 }
272 }
273 }
274 }
275
276 // Analyse some trigger stuff
277 AliTriggerAnalysis ta;
278 if (ta.IsOfflineTriggerFired(esd, AliTriggerAnalysis::kNSD1)) {
279 triggers |= AliAODForwardMult::kNSD;
280 fHTriggers->Fill(kNSD+.5);
281 }
282
283 // Get trigger stuff
284 TString trigStr = esd->GetFiredTriggerClasses();
285 if (trigStr.Contains("CBEAMB-ABCE-NOPF-ALL")) {
286 triggers |= AliAODForwardMult::kEmpty;
287 fHTriggers->Fill(kEmpty+.5);
288 }
289
290 if (trigStr.Contains("CINT1A-ABCE-NOPF-ALL")) {
291 triggers |= AliAODForwardMult::kA;
292 fHTriggers->Fill(kA+.5);
293 }
294
295 if (trigStr.Contains("CINT1B-ABCE-NOPF-ALL")) {
296 triggers |= AliAODForwardMult::kB;
297 fHTriggers->Fill(kB+.5);
298 }
299
300
301 if (trigStr.Contains("CINT1C-ABCE-NOPF-ALL")) {
302 triggers |= AliAODForwardMult::kC;
303 fHTriggers->Fill(kC+.5);
304 }
305
306 if (trigStr.Contains("CINT1-E-NOPF-ALL")) {
307 triggers |= AliAODForwardMult::kE;
308 fHTriggers->Fill(kE+.5);
309 }
310
311 return kTRUE;
312}
313//____________________________________________________________________
314Bool_t
315AliFMDEventInspector::ReadVertex(const AliESDEvent* esd, Double_t& vz)
316{
317 vz = 0;
318 // Get the vertex
319 const AliESDVertex* vertex = esd->GetPrimaryVertexSPD();
320 if (!vertex) {
321 if (fDebug > 2)
322 AliWarning("No SPD vertex found in ESD");
323 return kFALSE;
324 }
325
326 // Check that enough tracklets contributed
327 if(vertex->GetNContributors() <= 0) {
328 if (fDebug > 2)
329 AliWarning(Form("Number of contributors to vertex is %d<=0",
330 vertex->GetNContributors()));
331 vz = 0;
332 return kFALSE;
333 }
334
335 // Check that the uncertainty isn't too large
336 if (vertex->GetZRes() > fMaxVzErr) {
337 if (fDebug > 2)
338 AliWarning(Form("Uncertaintity in Z of vertex is too large %f > %f",
339 vertex->GetZRes(), fMaxVzErr));
340 return kFALSE;
341 }
342
343 // Get the z coordiante
344 vz = vertex->GetZ();
345 return kTRUE;
346}
347//
348// EOF
349//
350