]> git.uio.no Git - u/mrichter/AliRoot.git/blame - HLT/QA/AliHLTQADataMakerRec.cxx
Fixing coding violations (Livio)
[u/mrichter/AliRoot.git] / HLT / QA / AliHLTQADataMakerRec.cxx
CommitLineData
341aab10 1// $Id$
ed292d06 2
3//**************************************************************************
4//* This file is property of and copyright by the ALICE HLT Project *
5//* ALICE Experiment at CERN, All rights reserved. *
6//* *
019f778f 7//* Primary Authors: Zhong-Bao Yin <Zhong-Bao.Yin@cern.ch> *
8//* Matthias Richter <Matthias.Richter@ift.uib.no> *
ed292d06 9//* for The ALICE HLT Project. *
10//* *
11//* Permission to use, copy, modify and distribute this software and its *
12//* documentation strictly for non-commercial purposes is hereby granted *
13//* without fee, provided that the above copyright notice appears in all *
14//* copies and that both the copyright notice and this permission notice *
15//* appear in the supporting documentation. The authors make no claims *
16//* about the suitability of this software for any purpose. It is *
17//* provided "as is" without express or implied warranty. *
18//**************************************************************************
19
20/** @file AliHLTQADataMakerRec.cxx
019f778f 21 @author Zhongbao Yin, Matthias Richter
ed292d06 22 @date 2009-05-14
23 @brief Container for the HLT offline QA
24*/
25#include "AliHLTQADataMakerRec.h"
4c080b28 26#include "AliESDEvent.h"
27#include <iostream>
28
019f778f 29#include "TH1F.h"
30#include "TH2F.h"
31
32#include "AliESDtrack.h"
33
4c080b28 34using namespace std;
ed292d06 35
36/** ROOT macro for the implementation of ROOT specific class methods */
37ClassImp(AliHLTQADataMakerRec)
38
39AliHLTQADataMakerRec::AliHLTQADataMakerRec()
40{
41 // see header file for class documentation
42 // or
43 // refer to README to build package
44 // or
45 // visit http://web.ift.uib.no/~kjeks/doc/alice-hlt
46}
47
48AliHLTQADataMakerRec::~AliHLTQADataMakerRec()
49{
50 // see header file for class documentation
51}
52
4c080b28 53void AliHLTQADataMakerRec::Exec(AliQAv1::TASKINDEX_t task, TObject * data)
54{
55 // special handling for esds
56 if ( task == AliQAv1::kESDS ) {
57 AliESDEvent * esd = NULL;
58 AliESDEvent * hltesd = NULL;
59 if (data->IsA() == AliESDEvent::Class()) {
60 // silently skip this. Currently HLT QA is still called as
61 // part of AliQAManager::RunOneEvent with the esd
62 return;
63 }
64 if (data->InheritsFrom("TObjArray")) {
65 TObjArray* array=dynamic_cast<TObjArray*>(data);
66 if (array && array->GetEntriesFast()>0) {
67 esd = dynamic_cast<AliESDEvent *>(array->At(0)) ;
68 }
69 if (array && array->GetEntriesFast()>1) {
70 hltesd = dynamic_cast<AliESDEvent *>(array->At(1)) ;
71 }
72 } else {
73 esd = static_cast<AliESDEvent *>(data) ;
74 }
75
76 if (esd && strcmp(esd->ClassName(), "AliESDEvent") == 0) {
77 if (hltesd) {
78 MakeESDs(esd, hltesd);
79 } else {
80 AliError(Form("HLT ESD missing or wrong class type (%p), skipping HLT QA task kESDs", hltesd));
81 }
82 } else {
83 AliError(Form("ESD missing or wrong class type (%p), skipping HLT QA task kESDSs", esd));
84 }
85 } else {
86 // forward for all other types
87 AliQADataMakerRec::Exec(task, data);
88 }
89}
90
ed292d06 91void AliHLTQADataMakerRec::StartOfDetectorCycle()
92{
93 // see header file for class documentation
94}
95
96void AliHLTQADataMakerRec::EndOfDetectorCycle(AliQAv1::TASKINDEX_t, TObjArray** /*list*/)
97{
98 // see header file for class documentation
019f778f 99
100 if(GetESDsData(kPHLT)->GetEntries()){
101 GetESDsData(kPRatio)->Sumw2();
102 GetESDsData(kPRatio)->Add(GetESDsData(kPOffline));
103 GetESDsData(kPRatio)->Divide(GetESDsData(kPHLT));
104 }
105
106 if(GetESDsData(kPHLTFired)->GetEntries()){
107 GetESDsData(kPRatioFired)->Sumw2();
108 GetESDsData(kPRatioFired)->Add(GetESDsData(kPOfflineFired));
109 GetESDsData(kPRatioFired)->Divide(GetESDsData(kPHLTFired));
110 }
111
112
ed292d06 113}
4c080b28 114
115void AliHLTQADataMakerRec::MakeRaws(AliRawReader * rawReader)
116{
117 // see header file for class documentation
118 if (!rawReader) return;
119}
120
019f778f 121void AliHLTQADataMakerRec::InitESDs(){
122
123 //create ESDs histograms in ESDs subdir
124 const Bool_t expert = kTRUE ;
125 const Bool_t image = kTRUE ;
126
127 TH2F * histESDMultiplicity =
128 new TH2F("hESDMultiplicity",
129 "Number of tracks from all events; Number of offline tracks; Number of HLT tracks",
130 300, 0, 300, 300, 0, 300);
131 histESDMultiplicity->Sumw2();
132 Add2ESDsList(histESDMultiplicity, kMultiplicity, !expert, image);
133
134 TH2F * histESDMultiplicityFired =
135 new TH2F("hESDMultiplicityFired",
136 "Number of tracks from HLT triggered events; Number of offline tracks; Number of HLT tracks",
137 300, 0, 300, 300, 0, 300);
138 histESDMultiplicityFired->Sumw2();
139 Add2ESDsList(histESDMultiplicityFired, kMultiplicityFired, !expert, image);
140
141 TH2F * histESDNCls =
142 new TH2F("hESDNCls", "Mean number of TPC clusters from all events; Offline; HLT",
143 200, 0, 200, 200, 0, 200);
144 histESDNCls->Sumw2();
145 Add2ESDsList(histESDNCls, kNCls, !expert, image);
146
147 TH2F * histESDNClsFired =
148 new TH2F("hESDNClsFired", "Mean number of TPC clusters from triggered events; Offline; HLT",
149 200, 0, 200, 200, 0, 200);
150 histESDNClsFired->Sumw2();
151 Add2ESDsList(histESDNClsFired, kNClsFired, !expert, image);
152
153 TH1F * histPHLT
154 = new TH1F("hPHLT", "P distribution for all events from HLT; P [GeV/c]",
155 100, -0.5, 99.5);
156 histPHLT->Sumw2();
157 Add2ESDsList(histPHLT, kPHLT, !expert, image);
158
159 TH1F * histPOffline
160 = new TH1F("hPOffline",
161 "P distribution for all events from offline; P [GeV/c]",
162 100, -0.5, 99.5);
163 histPOffline->Sumw2();
164 Add2ESDsList(histPOffline, kPOffline, !expert, image);
165
166 TH1F * histPHLTFired
167 = new TH1F("hPHLTFired",
168 "P distribution for fired events from HLT; P [GeV/c]",
169 100, -0.5, 99.5);
170 histPHLTFired->Sumw2();
171 Add2ESDsList(histPHLTFired, kPHLTFired, !expert, image);
172
173 TH1F * histPOfflineFired
174 = new TH1F("hPOfflineFired",
175 "P distribution for fired events from offline; P [GeV/c]",
176 100, -0.5, 99.5);
177 histPOfflineFired->Sumw2();
178 Add2ESDsList(histPOfflineFired, kPOfflineFired, !expert, image);
179
180 TH1F * histPRatio
181 = new TH1F("hPRatio",
182 "Ratio of P distribution for all events; P [GeV/c]",
183 100, -0.5, 99.5);
184 histPRatio->Sumw2();
185 Add2ESDsList(histPRatio, kPRatio, !expert, image);
186
187 TH1F * histPRatioFired
188 = new TH1F("hPRatioFired",
189 "Ratio of P distribution for fired events; P [GeV/c]",
190 100, -0.5, 99.5);
191 histPRatioFired->Sumw2();
192 Add2ESDsList(histPRatioFired, kPRatioFired, !expert, image);
193
194 TH1F * histPtHLT
195 = new TH1F("hPtHLT",
196 "P_{T} distribution for all events from HLT; P_{T} [GeV/c]",
197 200, -100.5, 99.5);
198 histPtHLT->Sumw2();
199 Add2ESDsList(histPtHLT, kPtHLT, !expert, image);
200
201 TH1F * histPtOffline
202 = new TH1F("hPtOffline",
203 "P_{T} distribution for all events from offline; P_{T} [GeV/c]",
204 200, -100.5, 99.5);
205 histPtOffline->Sumw2();
206 Add2ESDsList(histPtOffline, kPtOffline, !expert, image);
207
208 TH1F * histPtHLTFired
209 = new TH1F("hPtHLTFired",
210 "P_{T} distribution for fired events from HLT; P_{T} [GeV/c]",
211 200, -100.5, 99.5);
212 histPtHLTFired->Sumw2();
213 Add2ESDsList(histPtHLTFired, kPtHLTFired, !expert, image);
214
215 TH1F * histPtOfflineFired
216 = new TH1F("hPtOfflineFired",
217 "P_{T} distribution for fired events from offline; P_{T} [GeV/c]",
218 200, -100.5, 99.5);
219 histPtOfflineFired->Sumw2();
220 Add2ESDsList(histPtOfflineFired, kPtOfflineFired, !expert, image);
221
222 TH1F * histNClsPerTrkHLT
223 = new TH1F("hNClsPerTrkHLT", "Clusters per HLT track; N cluster; Counts",
224 200, 0, 200);
225 histNClsPerTrkHLT->Sumw2();
226 Add2ESDsList(histNClsPerTrkHLT, kNClsPerTrkHLT, !expert, image);
227
228 TH1F * histNClsPerTrkOffline
229 = new TH1F("hNClsPerTrkOffline",
230 "Clusters per offline track; N cluster; Counts",
231 200, 0, 200);
232 histNClsPerTrkOffline->Sumw2();
233 Add2ESDsList(histNClsPerTrkOffline, kNClsPerTrkOffline, !expert, image);
234
235TH1F * histNClsPerTrkHLTFired
236 = new TH1F("hNClsPerTrkHLTFired",
237 "Clusters per HLT track from HLT fired events; N cluster; Counts",
238 200, 0, 200);
239 histNClsPerTrkHLTFired->Sumw2();
240 Add2ESDsList(histNClsPerTrkHLTFired, kNClsPerTrkHLTFired, !expert, image);
241
242 TH1F * histNClsPerTrkOfflineFired
243 = new TH1F("hNClsPerTrkOfflineFired",
244 "Clusters per offline track from HLT fired events; N cluster; Counts",
245 200, 0, 200);
246 histNClsPerTrkOfflineFired->Sumw2();
247 Add2ESDsList(histNClsPerTrkOfflineFired, kNClsPerTrkOfflineFired, !expert, image);
248
249 TH1F * histPhiHLT =
250 new TH1F("hPhiHLT", "Phi distribution of HLT tracks; Phi; Counts",
251 360, 0, 360);
252 histPhiHLT->Sumw2();
253 Add2ESDsList(histPhiHLT, kPhiHLT, !expert, image);
254
255 TH1F * histPhiOffline =
256 new TH1F("hPhiOffline",
257 "Phi distribution of offline tracks; Phi; Counts",
258 360, 0, 360);
259 histPhiOffline->Sumw2();
260 Add2ESDsList(histPhiOffline, kPhiOffline, !expert, image);
261
262 TH1F * histPhiHLTFired =
263 new TH1F("hPhiHLTFired", "Phi distribution of HLT tracks from HLT fired event ; Phi; Counts",
264 360, 0, 360);
265 histPhiHLTFired->Sumw2();
266 Add2ESDsList(histPhiHLTFired, kPhiHLTFired, !expert, image);
267
268 TH1F * histPhiOfflineFired =
269 new TH1F("hPhiOfflineFired",
270 "Phi distribution of offline tracks from HLT fired events; Phi; Counts",
271 360, 0, 360);
272 histPhiOfflineFired->Sumw2();
273 Add2ESDsList(histPhiOfflineFired, kPhiOfflineFired, !expert, image);
274
275 TH1F * histEtaHLT =
276 new TH1F("hEtaHLT", "Eta distribution of HLT tracks; Eta; Counts",
277 200, -1, 1);
278 histEtaHLT->Sumw2();
279 Add2ESDsList(histEtaHLT, kEtaHLT, !expert, image);
280
281 TH1F * histEtaOffline =
282 new TH1F("hEtaOffline", "Eta distribution of offline tracks; Eta; Counts",
283 200, -1, 1);
284 histEtaHLT->Sumw2();
285 Add2ESDsList(histEtaOffline, kEtaOffline, !expert, image);
286
287 TH1F * histEtaHLTFired =
288 new TH1F("hEtaHLTFired",
289 "Eta distribution of HLT tracks from HLT fired events; Eta; Counts",
290 200, -1, 1);
291 histEtaHLTFired->Sumw2();
292 Add2ESDsList(histEtaHLTFired, kEtaHLTFired, !expert, image);
293
294 TH1F * histEtaOfflineFired =
295 new TH1F("hEtaOfflineFired",
296 "Eta distribution of offline tracks from HLT fired events; Eta; Counts",
297 200, -1, 1);
298 histEtaHLTFired->Sumw2();
299 Add2ESDsList(histEtaOfflineFired, kEtaOfflineFired, !expert, image);
300
301}
302
4c080b28 303void AliHLTQADataMakerRec::MakeESDs(AliESDEvent * esd)
304{
305 // see header file for class documentation
306
307 // as an extension in the QA interface also the hlt esd can be sent
308 // in order to preserve backward compatibility, a new function has been
309 // introduced.
310 //
311 // NOTE: This function is not the place for HLT QA
312 if (!esd) return;
313}
314
315void AliHLTQADataMakerRec::MakeESDs(AliESDEvent * esd, AliESDEvent* hltesd)
316{
317 // HLT QA on ESDs
318 if (!esd || !hltesd) {
319 AliError("invalid parameter: missing ESDs");
320 return;
321 }
019f778f 322
323 // make QA data from ESDs
324
325 const Int_t nESDTracks = esd->GetNumberOfTracks();
326 const Int_t nHLTesdTracks = hltesd->GetNumberOfTracks();
327 GetESDsData(kMultiplicity)->Fill(nESDTracks, nHLTesdTracks);
328
329 Int_t nClsHLT = 0;
330 Int_t nClsOffline = 0;
331
332 for(Int_t i = 0; i < nESDTracks; i++){
333 AliESDtrack * esdTrk = esd->GetTrack(i);
334 GetESDsData(kPOffline)->Fill(esdTrk->P());
335 GetESDsData(kPtOffline)->Fill(esdTrk->GetSignedPt());
336 GetESDsData(kNClsPerTrkOffline)->Fill(esdTrk->GetTPCNcls());
337 GetESDsData(kPhiOffline)->Fill(esdTrk->Phi()*TMath::RadToDeg());
338 GetESDsData(kEtaOffline)->Fill(esdTrk->Eta());
339 nClsOffline+=esdTrk->GetTPCNcls();
340 }
341
342 for(Int_t i = 0; i < nHLTesdTracks; i++){
343 AliESDtrack * hltEsdTrk = hltesd->GetTrack(i);
344 GetESDsData(kPHLT)->Fill(hltEsdTrk->P());
345 GetESDsData(kPtHLT)->Fill(hltEsdTrk->GetSignedPt());
346 GetESDsData(kNClsPerTrkHLT)->Fill(hltEsdTrk->GetTPCNcls());
347 GetESDsData(kPhiHLT)->Fill(hltEsdTrk->Phi()*TMath::RadToDeg());
348 GetESDsData(kEtaHLT)->Fill(hltEsdTrk->Eta());
349 nClsHLT += hltEsdTrk->GetTPCNcls();
350 }
351
352 if(nESDTracks)
353 nClsOffline /= Float_t(nESDTracks);
354 if(nHLTesdTracks)
355 nClsHLT /= Float_t(nHLTesdTracks);
356
357 GetESDsData(kNCls)->Fill(nClsOffline, nClsHLT);
358
359
360 if(hltesd->IsHLTTriggerFired()){
361 GetESDsData(kMultiplicityFired)->Fill(nESDTracks, nHLTesdTracks);
362 GetESDsData(kNClsFired)->Fill(nClsOffline, nClsHLT);
363
364 for(Int_t i = 0; i < nESDTracks; i++){
365 AliESDtrack * esdTrk = esd->GetTrack(i);
366 GetESDsData(kPOfflineFired)->Fill(esdTrk->P());
367 GetESDsData(kPtOfflineFired)->Fill(esdTrk->GetSignedPt());
368 GetESDsData(kNClsPerTrkOfflineFired)->Fill(esdTrk->GetTPCNcls());
369 GetESDsData(kPhiOfflineFired)->Fill(esdTrk->Phi()*TMath::RadToDeg());
370 GetESDsData(kEtaOfflineFired)->Fill(esdTrk->Eta());
371 }
372
373 for(Int_t i = 0; i < nHLTesdTracks; i++){
374 AliESDtrack * hltEsdTrk = hltesd->GetTrack(i);
375 GetESDsData(kPHLTFired)->Fill(hltEsdTrk->P());
376 GetESDsData(kPtHLTFired)->Fill(hltEsdTrk->GetSignedPt());
377 GetESDsData(kNClsPerTrkHLTFired)->Fill(hltEsdTrk->GetTPCNcls());
378 GetESDsData(kPhiHLTFired)->Fill(hltEsdTrk->Phi()*TMath::RadToDeg());
379 GetESDsData(kEtaHLTFired)->Fill(hltEsdTrk->Eta());
380 }
381
382
383 }
384
4c080b28 385}