]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/QA/AliHLTTPCQADataMaker.cxx
bug https://savannah.cern.ch/bugs/?69885
[u/mrichter/AliRoot.git] / HLT / QA / AliHLTTPCQADataMaker.cxx
1 // $Id$
2
3 //**************************************************************************
4 //* This file is property of and copyright by the ALICE HLT Project        * 
5 //* ALICE Experiment at CERN, All rights reserved.                         *
6 //*                                                                        *
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.                            *
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   AliHLTTPCQADataMaker.cxx
21     @author Zhongbao Yin, Matthias Richter
22     @date   2009-05-14
23     @brief  Container for the HLT offline QA
24 */
25 #include "AliHLTTPCQADataMaker.h"
26 #include "AliESDEvent.h"
27 #include <iostream>
28
29 #include "TH1F.h"
30 #include "TH2F.h"
31
32 #include "AliESDtrack.h"
33
34 using namespace std;
35
36 /** ROOT macro for the implementation of ROOT specific class methods */
37 ClassImp(AliHLTTPCQADataMaker)
38
39 AliHLTTPCQADataMaker::AliHLTTPCQADataMaker()
40   : AliHLTQADataMakerBase()
41 {
42   // see header file for class documentation
43   // or
44   // refer to README to build package
45   // or
46   // visit http://web.ift.uib.no/~kjeks/doc/alice-hlt
47 }
48
49 AliHLTTPCQADataMaker::~AliHLTTPCQADataMaker()
50 {
51   // see header file for class documentation
52 }
53
54 void AliHLTTPCQADataMaker::StartOfDetectorCycle()
55 {
56   // see header file for class documentation
57 }
58
59 void AliHLTTPCQADataMaker::EndOfDetectorCycle(AliQAv1::TASKINDEX_t, TObjArray** /*list*/)
60 {
61   // see header file for class documentation
62   
63   if(GetESDsData(kPHLT)->GetEntries()){
64     GetESDsData(kPRatio)->Sumw2();
65     GetESDsData(kPRatio)->Add(GetESDsData(kPOffline));
66     GetESDsData(kPRatio)->Divide(GetESDsData(kPHLT));
67   }
68   
69   if(GetESDsData(kPHLTFired)->GetEntries()){ 
70     GetESDsData(kPRatioFired)->Sumw2(); 
71     GetESDsData(kPRatioFired)->Add(GetESDsData(kPOfflineFired)); 
72     GetESDsData(kPRatioFired)->Divide(GetESDsData(kPHLTFired));
73   } 
74   
75
76 }
77
78 void AliHLTTPCQADataMaker::MakeRaws(AliRawReader * rawReader)
79 {
80   // see header file for class documentation
81   if (!rawReader) return;
82 }
83
84 void AliHLTTPCQADataMaker::InitESDs(){
85   
86   //create ESDs histograms in ESDs subdir
87   const Bool_t expert   = kTRUE ; 
88   const Bool_t image    = kTRUE ; 
89
90   TH2F * histESDMultiplicity =
91     new TH2F("hESDMultiplicity", 
92              "Number of tracks from all events; Number of offline tracks; Number of HLT tracks", 
93              300, 0, 300, 300, 0, 300);
94   histESDMultiplicity->Sumw2();
95   Add2ESDsList(histESDMultiplicity, kMultiplicity, !expert, image);
96
97   TH2F * histESDMultiplicityFired = 
98     new TH2F("hESDMultiplicityFired",  
99              "Number of tracks from HLT triggered events; Number of offline tracks; Number of HLT tracks",  
100              300, 0, 300, 300, 0, 300);
101   histESDMultiplicityFired->Sumw2();
102   Add2ESDsList(histESDMultiplicityFired, kMultiplicityFired, !expert, image);
103
104   TH2F * histESDNCls = 
105     new TH2F("hESDNCls", "Mean number of TPC clusters from all events; Offline; HLT",
106              200, 0, 200, 200, 0, 200);
107   histESDNCls->Sumw2();
108   Add2ESDsList(histESDNCls, kNCls, !expert, image);
109   
110   TH2F * histESDNClsFired =  
111     new TH2F("hESDNClsFired", "Mean number of TPC clusters from triggered events; Offline; HLT", 
112              200, 0, 200, 200, 0, 200); 
113   histESDNClsFired->Sumw2();
114   Add2ESDsList(histESDNClsFired, kNClsFired, !expert, image);
115
116   TH1F * histPHLT 
117     = new TH1F("hPHLT", "P distribution for all events from HLT; P [GeV/c]", 
118                100, -0.5, 99.5);
119   histPHLT->Sumw2();
120   Add2ESDsList(histPHLT, kPHLT, !expert, image);
121
122   TH1F * histPOffline 
123     = new TH1F("hPOffline", 
124                "P distribution for all events from offline; P [GeV/c]",  
125                100, -0.5, 99.5);
126   histPOffline->Sumw2();
127   Add2ESDsList(histPOffline, kPOffline, !expert, image);
128   
129   TH1F * histPHLTFired 
130     = new TH1F("hPHLTFired", 
131                "P distribution for fired events from HLT; P [GeV/c]",  
132                100, -0.5, 99.5);
133   histPHLTFired->Sumw2();
134   Add2ESDsList(histPHLTFired, kPHLTFired, !expert, image); 
135  
136   TH1F * histPOfflineFired 
137     = new TH1F("hPOfflineFired",  
138                "P distribution for fired events from offline; P [GeV/c]",   
139                100, -0.5, 99.5);
140   histPOfflineFired->Sumw2();
141   Add2ESDsList(histPOfflineFired, kPOfflineFired, !expert, image);
142
143   TH1F * histPRatio 
144     = new TH1F("hPRatio", 
145                "Ratio of P distribution for all events; P [GeV/c]",
146                100, -0.5, 99.5);
147   histPRatio->Sumw2();
148   Add2ESDsList(histPRatio, kPRatio, !expert, image);
149   
150   TH1F * histPRatioFired  
151     = new TH1F("hPRatioFired",  
152                "Ratio of P distribution for fired events; P [GeV/c]",
153                100, -0.5, 99.5);
154   histPRatioFired->Sumw2();
155   Add2ESDsList(histPRatioFired, kPRatioFired, !expert, image);
156
157   TH1F * histPtHLT  
158     = new TH1F("hPtHLT", 
159                "P_{T} distribution for all events from HLT; P_{T} [GeV/c]",  
160                200, -100.5, 99.5);
161   histPtHLT->Sumw2();
162   Add2ESDsList(histPtHLT, kPtHLT, !expert, image); 
163  
164   TH1F * histPtOffline  
165     = new TH1F("hPtOffline",  
166                "P_{T} distribution for all events from offline; P_{T} [GeV/c]",   
167                200, -100.5, 99.5);
168   histPtOffline->Sumw2();
169   Add2ESDsList(histPtOffline, kPtOffline, !expert, image);
170
171   TH1F * histPtHLTFired   
172     = new TH1F("hPtHLTFired",  
173                "P_{T} distribution for fired events from HLT; P_{T} [GeV/c]",   
174                200, -100.5, 99.5);
175   histPtHLTFired->Sumw2();
176   Add2ESDsList(histPtHLTFired, kPtHLTFired, !expert, image);  
177   
178   TH1F * histPtOfflineFired   
179     = new TH1F("hPtOfflineFired",   
180                "P_{T} distribution for fired events from offline; P_{T} [GeV/c]", 
181                200, -100.5, 99.5);
182   histPtOfflineFired->Sumw2();
183   Add2ESDsList(histPtOfflineFired, kPtOfflineFired, !expert, image); 
184
185   TH1F * histNClsPerTrkHLT
186     = new TH1F("hNClsPerTrkHLT", "Clusters per HLT track; N cluster; Counts",
187                200, 0, 200);
188   histNClsPerTrkHLT->Sumw2();
189   Add2ESDsList(histNClsPerTrkHLT, kNClsPerTrkHLT, !expert, image);
190
191   TH1F * histNClsPerTrkOffline 
192     = new TH1F("hNClsPerTrkOffline", 
193                "Clusters per offline track; N cluster; Counts", 
194                200, 0, 200); 
195   histNClsPerTrkOffline->Sumw2();
196   Add2ESDsList(histNClsPerTrkOffline, kNClsPerTrkOffline, !expert, image);
197
198 TH1F * histNClsPerTrkHLTFired 
199   = new TH1F("hNClsPerTrkHLTFired", 
200              "Clusters per HLT track from HLT fired events; N cluster; Counts", 
201              200, 0, 200); 
202  histNClsPerTrkHLTFired->Sumw2();
203  Add2ESDsList(histNClsPerTrkHLTFired, kNClsPerTrkHLTFired, !expert, image); 
204  
205   TH1F * histNClsPerTrkOfflineFired  
206     = new TH1F("hNClsPerTrkOfflineFired",  
207                "Clusters per offline track from HLT fired events; N cluster; Counts",  
208                200, 0, 200);  
209   histNClsPerTrkOfflineFired->Sumw2();
210   Add2ESDsList(histNClsPerTrkOfflineFired, kNClsPerTrkOfflineFired, !expert, image);
211
212   TH1F * histPhiHLT = 
213     new TH1F("hPhiHLT", "Phi distribution of HLT tracks; Phi; Counts", 
214              360, 0, 360);
215   histPhiHLT->Sumw2();
216   Add2ESDsList(histPhiHLT, kPhiHLT, !expert, image);
217
218   TH1F * histPhiOffline =  
219     new TH1F("hPhiOffline", 
220              "Phi distribution of offline tracks; Phi; Counts",  
221              360, 0, 360); 
222   histPhiOffline->Sumw2(); 
223   Add2ESDsList(histPhiOffline, kPhiOffline, !expert, image);
224
225   TH1F * histPhiHLTFired =  
226     new TH1F("hPhiHLTFired", "Phi distribution of HLT tracks from HLT fired event ; Phi; Counts",  
227              360, 0, 360); 
228   histPhiHLTFired->Sumw2(); 
229   Add2ESDsList(histPhiHLTFired, kPhiHLTFired, !expert, image); 
230  
231   TH1F * histPhiOfflineFired =   
232     new TH1F("hPhiOfflineFired",  
233              "Phi distribution of offline tracks from HLT fired events; Phi; Counts",   
234              360, 0, 360);  
235   histPhiOfflineFired->Sumw2();  
236   Add2ESDsList(histPhiOfflineFired, kPhiOfflineFired, !expert, image); 
237
238   TH1F * histEtaHLT =
239     new TH1F("hEtaHLT", "Eta distribution of HLT tracks; Eta; Counts",
240              200, -1, 1);
241   histEtaHLT->Sumw2();
242   Add2ESDsList(histEtaHLT, kEtaHLT, !expert, image);
243
244   TH1F * histEtaOffline = 
245     new TH1F("hEtaOffline", "Eta distribution of offline tracks; Eta; Counts", 
246              200, -1, 1); 
247   histEtaHLT->Sumw2(); 
248   Add2ESDsList(histEtaOffline, kEtaOffline, !expert, image);
249
250   TH1F * histEtaHLTFired = 
251     new TH1F("hEtaHLTFired", 
252              "Eta distribution of HLT tracks from HLT fired events; Eta; Counts", 
253              200, -1, 1); 
254   histEtaHLTFired->Sumw2(); 
255   Add2ESDsList(histEtaHLTFired, kEtaHLTFired, !expert, image); 
256  
257   TH1F * histEtaOfflineFired =  
258     new TH1F("hEtaOfflineFired", 
259              "Eta distribution of offline tracks from HLT fired events; Eta; Counts",  
260              200, -1, 1);  
261   histEtaHLTFired->Sumw2();  
262   Add2ESDsList(histEtaOfflineFired, kEtaOfflineFired, !expert, image);
263
264 }
265
266 void AliHLTTPCQADataMaker::MakeESDs(AliESDEvent * esd, AliESDEvent* hltesd)
267 {
268   // HLT QA on ESDs
269   if (!esd || !hltesd) {
270     AliError("invalid parameter: missing ESDs");
271     return;
272   }
273
274   // make QA data from ESDs
275
276   const Int_t nESDTracks = esd->GetNumberOfTracks();
277   const Int_t nHLTesdTracks = hltesd->GetNumberOfTracks();
278   GetESDsData(kMultiplicity)->Fill(nESDTracks, nHLTesdTracks);
279
280   Int_t nClsHLT = 0;
281   Int_t nClsOffline = 0;
282   
283   for(Int_t i = 0; i < nESDTracks; i++){
284     AliESDtrack *  esdTrk = esd->GetTrack(i);
285     GetESDsData(kPOffline)->Fill(esdTrk->P());
286     GetESDsData(kPtOffline)->Fill(esdTrk->GetSignedPt());
287     GetESDsData(kNClsPerTrkOffline)->Fill(esdTrk->GetTPCNcls());
288     GetESDsData(kPhiOffline)->Fill(esdTrk->Phi()*TMath::RadToDeg());
289     GetESDsData(kEtaOffline)->Fill(esdTrk->Eta());
290     nClsOffline+=esdTrk->GetTPCNcls();
291   }
292
293   for(Int_t i = 0; i < nHLTesdTracks; i++){ 
294     AliESDtrack *  hltEsdTrk = hltesd->GetTrack(i); 
295     GetESDsData(kPHLT)->Fill(hltEsdTrk->P()); 
296     GetESDsData(kPtHLT)->Fill(hltEsdTrk->GetSignedPt());
297     GetESDsData(kNClsPerTrkHLT)->Fill(hltEsdTrk->GetTPCNcls());
298     GetESDsData(kPhiHLT)->Fill(hltEsdTrk->Phi()*TMath::RadToDeg());
299     GetESDsData(kEtaHLT)->Fill(hltEsdTrk->Eta());
300     nClsHLT += hltEsdTrk->GetTPCNcls();
301   } 
302   
303   if(nESDTracks)
304     nClsOffline /= Float_t(nESDTracks);
305   if(nHLTesdTracks)
306     nClsHLT /= Float_t(nHLTesdTracks);
307
308   GetESDsData(kNCls)->Fill(nClsOffline, nClsHLT);
309
310
311   if(hltesd->IsHLTTriggerFired()){
312     GetESDsData(kMultiplicityFired)->Fill(nESDTracks, nHLTesdTracks);
313     GetESDsData(kNClsFired)->Fill(nClsOffline, nClsHLT);
314     
315     for(Int_t i = 0; i < nESDTracks; i++){ 
316       AliESDtrack *  esdTrk = esd->GetTrack(i); 
317       GetESDsData(kPOfflineFired)->Fill(esdTrk->P()); 
318       GetESDsData(kPtOfflineFired)->Fill(esdTrk->GetSignedPt());
319       GetESDsData(kNClsPerTrkOfflineFired)->Fill(esdTrk->GetTPCNcls());
320       GetESDsData(kPhiOfflineFired)->Fill(esdTrk->Phi()*TMath::RadToDeg());
321       GetESDsData(kEtaOfflineFired)->Fill(esdTrk->Eta());
322     } 
323    
324     for(Int_t i = 0; i < nHLTesdTracks; i++){  
325       AliESDtrack *  hltEsdTrk = hltesd->GetTrack(i);  
326       GetESDsData(kPHLTFired)->Fill(hltEsdTrk->P());  
327       GetESDsData(kPtHLTFired)->Fill(hltEsdTrk->GetSignedPt());
328       GetESDsData(kNClsPerTrkHLTFired)->Fill(hltEsdTrk->GetTPCNcls());
329       GetESDsData(kPhiHLTFired)->Fill(hltEsdTrk->Phi()*TMath::RadToDeg());
330       GetESDsData(kEtaHLTFired)->Fill(hltEsdTrk->Eta());
331     }  
332
333     
334   }
335
336 }