3 //**************************************************************************
4 //* This file is property of and copyright by the ALICE HLT Project *
5 //* ALICE Experiment at CERN, All rights reserved. *
7 //* Primary Authors: Zhong-Bao Yin <Zhong-Bao.Yin@cern.ch> *
8 //* Matthias Richter <Matthias.Richter@ift.uib.no> *
9 //* for The ALICE HLT Project. *
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 //**************************************************************************
20 /** @file AliHLTTPCQADataMaker.cxx
21 @author Zhongbao Yin, Matthias Richter
23 @brief Container for the HLT offline QA
25 #include "AliHLTTPCQADataMaker.h"
26 #include "AliESDEvent.h"
32 #include "AliESDtrack.h"
36 /** ROOT macro for the implementation of ROOT specific class methods */
37 ClassImp(AliHLTTPCQADataMaker)
39 AliHLTTPCQADataMaker::AliHLTTPCQADataMaker()
40 : AliHLTQADataMakerBase()
42 // see header file for class documentation
44 // refer to README to build package
46 // visit http://web.ift.uib.no/~kjeks/doc/alice-hlt
49 AliHLTTPCQADataMaker::~AliHLTTPCQADataMaker()
51 // see header file for class documentation
54 void AliHLTTPCQADataMaker::StartOfDetectorCycle()
56 // see header file for class documentation
59 void AliHLTTPCQADataMaker::EndOfDetectorCycle(AliQAv1::TASKINDEX_t, TObjArray** /*list*/)
61 // see header file for class documentation
62 for (int itc=-1;itc<GetNTrigClasses();itc++) {
63 TH1* hp = GetESDsData(kPHLT,itc);
64 TH1* hr = GetESDsData(kPRatio,itc);
65 TH1* ho = GetESDsData(kPOffline,itc);
67 if (hp && hr && ho && hp->GetEntries()) {
73 hp = GetESDsData(kPHLTFired,itc);
74 hr = GetESDsData(kPRatioFired,itc);
75 ho = GetESDsData(kPOfflineFired,itc);
77 if (hp && hr && ho && hp->GetEntries()) {
85 void AliHLTTPCQADataMaker::MakeRaws(AliRawReader * rawReader)
87 // see header file for class documentation
88 if (!rawReader) return;
91 void AliHLTTPCQADataMaker::InitESDs(){
93 //create ESDs histograms in ESDs subdir
94 const Bool_t expert = kTRUE ;
95 const Bool_t image = kTRUE ;
97 TH2F * histESDMultiplicity =
98 new TH2F("hESDMultiplicity",
99 "Number of tracks from all events; Number of offline tracks; Number of HLT tracks",
100 300, 0, 300, 300, 0, 300);
101 histESDMultiplicity->Sumw2();
102 Add2ESDsList(histESDMultiplicity, kMultiplicity, !expert, image);
104 TH2F * histESDMultiplicityFired =
105 new TH2F("hESDMultiplicityFired",
106 "Number of tracks from HLT triggered events; Number of offline tracks; Number of HLT tracks",
107 300, 0, 300, 300, 0, 300);
108 histESDMultiplicityFired->Sumw2();
109 Add2ESDsList(histESDMultiplicityFired, kMultiplicityFired, !expert, image);
112 new TH2F("hESDNCls", "Mean number of TPC clusters from all events; Offline; HLT",
113 200, 0, 200, 200, 0, 200);
114 histESDNCls->Sumw2();
115 Add2ESDsList(histESDNCls, kNCls, !expert, image);
117 TH2F * histESDNClsFired =
118 new TH2F("hESDNClsFired", "Mean number of TPC clusters from triggered events; Offline; HLT",
119 200, 0, 200, 200, 0, 200);
120 histESDNClsFired->Sumw2();
121 Add2ESDsList(histESDNClsFired, kNClsFired, !expert, image);
124 = new TH1F("hPHLT", "P distribution for all events from HLT; P [GeV/c]",
127 Add2ESDsList(histPHLT, kPHLT, !expert, image);
130 = new TH1F("hPOffline",
131 "P distribution for all events from offline; P [GeV/c]",
133 histPOffline->Sumw2();
134 Add2ESDsList(histPOffline, kPOffline, !expert, image);
137 = new TH1F("hPHLTFired",
138 "P distribution for fired events from HLT; P [GeV/c]",
140 histPHLTFired->Sumw2();
141 Add2ESDsList(histPHLTFired, kPHLTFired, !expert, image);
143 TH1F * histPOfflineFired
144 = new TH1F("hPOfflineFired",
145 "P distribution for fired events from offline; P [GeV/c]",
147 histPOfflineFired->Sumw2();
148 Add2ESDsList(histPOfflineFired, kPOfflineFired, !expert, image);
151 = new TH1F("hPRatio",
152 "Ratio of P distribution for all events; P [GeV/c]",
155 Add2ESDsList(histPRatio, kPRatio, !expert, image);
157 TH1F * histPRatioFired
158 = new TH1F("hPRatioFired",
159 "Ratio of P distribution for fired events; P [GeV/c]",
161 histPRatioFired->Sumw2();
162 Add2ESDsList(histPRatioFired, kPRatioFired, !expert, image);
166 "P_{T} distribution for all events from HLT; P_{T} [GeV/c]",
169 Add2ESDsList(histPtHLT, kPtHLT, !expert, image);
172 = new TH1F("hPtOffline",
173 "P_{T} distribution for all events from offline; P_{T} [GeV/c]",
175 histPtOffline->Sumw2();
176 Add2ESDsList(histPtOffline, kPtOffline, !expert, image);
178 TH1F * histPtHLTFired
179 = new TH1F("hPtHLTFired",
180 "P_{T} distribution for fired events from HLT; P_{T} [GeV/c]",
182 histPtHLTFired->Sumw2();
183 Add2ESDsList(histPtHLTFired, kPtHLTFired, !expert, image);
185 TH1F * histPtOfflineFired
186 = new TH1F("hPtOfflineFired",
187 "P_{T} distribution for fired events from offline; P_{T} [GeV/c]",
189 histPtOfflineFired->Sumw2();
190 Add2ESDsList(histPtOfflineFired, kPtOfflineFired, !expert, image);
192 TH1F * histNClsPerTrkHLT
193 = new TH1F("hNClsPerTrkHLT", "Clusters per HLT track; N cluster; Counts",
195 histNClsPerTrkHLT->Sumw2();
196 Add2ESDsList(histNClsPerTrkHLT, kNClsPerTrkHLT, !expert, image);
198 TH1F * histNClsPerTrkOffline
199 = new TH1F("hNClsPerTrkOffline",
200 "Clusters per offline track; N cluster; Counts",
202 histNClsPerTrkOffline->Sumw2();
203 Add2ESDsList(histNClsPerTrkOffline, kNClsPerTrkOffline, !expert, image);
205 TH1F * histNClsPerTrkHLTFired
206 = new TH1F("hNClsPerTrkHLTFired",
207 "Clusters per HLT track from HLT fired events; N cluster; Counts",
209 histNClsPerTrkHLTFired->Sumw2();
210 Add2ESDsList(histNClsPerTrkHLTFired, kNClsPerTrkHLTFired, !expert, image);
212 TH1F * histNClsPerTrkOfflineFired
213 = new TH1F("hNClsPerTrkOfflineFired",
214 "Clusters per offline track from HLT fired events; N cluster; Counts",
216 histNClsPerTrkOfflineFired->Sumw2();
217 Add2ESDsList(histNClsPerTrkOfflineFired, kNClsPerTrkOfflineFired, !expert, image);
220 new TH1F("hPhiHLT", "Phi distribution of HLT tracks; Phi; Counts",
223 Add2ESDsList(histPhiHLT, kPhiHLT, !expert, image);
225 TH1F * histPhiOffline =
226 new TH1F("hPhiOffline",
227 "Phi distribution of offline tracks; Phi; Counts",
229 histPhiOffline->Sumw2();
230 Add2ESDsList(histPhiOffline, kPhiOffline, !expert, image);
232 TH1F * histPhiHLTFired =
233 new TH1F("hPhiHLTFired", "Phi distribution of HLT tracks from HLT fired event ; Phi; Counts",
235 histPhiHLTFired->Sumw2();
236 Add2ESDsList(histPhiHLTFired, kPhiHLTFired, !expert, image);
238 TH1F * histPhiOfflineFired =
239 new TH1F("hPhiOfflineFired",
240 "Phi distribution of offline tracks from HLT fired events; Phi; Counts",
242 histPhiOfflineFired->Sumw2();
243 Add2ESDsList(histPhiOfflineFired, kPhiOfflineFired, !expert, image);
246 new TH1F("hEtaHLT", "Eta distribution of HLT tracks; Eta; Counts",
249 Add2ESDsList(histEtaHLT, kEtaHLT, !expert, image);
251 TH1F * histEtaOffline =
252 new TH1F("hEtaOffline", "Eta distribution of offline tracks; Eta; Counts",
255 Add2ESDsList(histEtaOffline, kEtaOffline, !expert, image);
257 TH1F * histEtaHLTFired =
258 new TH1F("hEtaHLTFired",
259 "Eta distribution of HLT tracks from HLT fired events; Eta; Counts",
261 histEtaHLTFired->Sumw2();
262 Add2ESDsList(histEtaHLTFired, kEtaHLTFired, !expert, image);
264 TH1F * histEtaOfflineFired =
265 new TH1F("hEtaOfflineFired",
266 "Eta distribution of offline tracks from HLT fired events; Eta; Counts",
268 histEtaHLTFired->Sumw2();
269 Add2ESDsList(histEtaOfflineFired, kEtaOfflineFired, !expert, image);
273 void AliHLTTPCQADataMaker::MakeESDs(AliESDEvent * esd, AliESDEvent* hltesd)
276 if (!esd || !hltesd) {
277 AliError("invalid parameter: missing ESDs");
281 // make QA data from ESDs
283 const Int_t nESDTracks = esd->GetNumberOfTracks();
284 const Int_t nHLTesdTracks = hltesd->GetNumberOfTracks();
285 FillESDsData(kMultiplicity,nESDTracks, nHLTesdTracks);
288 Int_t nClsOffline = 0;
290 for(Int_t i = 0; i < nESDTracks; i++){
291 AliESDtrack * esdTrk = esd->GetTrack(i);
292 FillESDsData(kPOffline,esdTrk->P());
293 FillESDsData(kPtOffline,esdTrk->GetSignedPt());
294 FillESDsData(kNClsPerTrkOffline,esdTrk->GetTPCNcls());
295 FillESDsData(kPhiOffline,esdTrk->Phi()*TMath::RadToDeg());
296 FillESDsData(kEtaOffline,esdTrk->Eta());
297 nClsOffline+=esdTrk->GetTPCNcls();
300 for(Int_t i = 0; i < nHLTesdTracks; i++){
301 AliESDtrack * hltEsdTrk = hltesd->GetTrack(i);
302 FillESDsData(kPHLT,hltEsdTrk->P());
303 FillESDsData(kPtHLT,hltEsdTrk->GetSignedPt());
304 FillESDsData(kNClsPerTrkHLT,hltEsdTrk->GetTPCNcls());
305 FillESDsData(kPhiHLT,hltEsdTrk->Phi()*TMath::RadToDeg());
306 FillESDsData(kEtaHLT,hltEsdTrk->Eta());
307 nClsHLT += hltEsdTrk->GetTPCNcls();
311 nClsOffline /= nESDTracks;
313 nClsHLT /= nHLTesdTracks;
315 FillESDsData(kNCls,nClsOffline, nClsHLT);
318 if(hltesd->IsHLTTriggerFired()){
319 FillESDsData(kMultiplicityFired,nESDTracks, nHLTesdTracks);
320 FillESDsData(kNClsFired,nClsOffline, nClsHLT);
322 for(Int_t i = 0; i < nESDTracks; i++){
323 AliESDtrack * esdTrk = esd->GetTrack(i);
324 FillESDsData(kPOfflineFired,esdTrk->P());
325 FillESDsData(kPtOfflineFired,esdTrk->GetSignedPt());
326 FillESDsData(kNClsPerTrkOfflineFired,esdTrk->GetTPCNcls());
327 FillESDsData(kPhiOfflineFired,esdTrk->Phi()*TMath::RadToDeg());
328 FillESDsData(kEtaOfflineFired,esdTrk->Eta());
331 for(Int_t i = 0; i < nHLTesdTracks; i++){
332 AliESDtrack * hltEsdTrk = hltesd->GetTrack(i);
333 FillESDsData(kPHLTFired,hltEsdTrk->P());
334 FillESDsData(kPtHLTFired,hltEsdTrk->GetSignedPt());
335 FillESDsData(kNClsPerTrkHLTFired,hltEsdTrk->GetTPCNcls());
336 FillESDsData(kPhiHLTFired,hltEsdTrk->Phi()*TMath::RadToDeg());
337 FillESDsData(kEtaHLTFired,hltEsdTrk->Eta());