]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/QA/AliHLTTPCQADataMaker.cxx
updated jet response taks (B. Bathen)
[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   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);
66     //
67     if (hp && hr && ho && hp->GetEntries()) {
68       hp->Sumw2();
69       hr->Add(ho);
70       hr->Divide(hp);
71     }
72     //
73     hp = GetESDsData(kPHLTFired,itc);
74     hr = GetESDsData(kPRatioFired,itc);
75     ho = GetESDsData(kPOfflineFired,itc);
76     //
77     if (hp && hr && ho && hp->GetEntries()) {
78       hp->Sumw2();
79       hr->Add(ho);
80       hr->Divide(hp);
81     }
82   } 
83 }
84
85 void AliHLTTPCQADataMaker::MakeRaws(AliRawReader * rawReader)
86 {
87   // see header file for class documentation
88   if (!rawReader) return;
89 }
90
91 void AliHLTTPCQADataMaker::InitESDs(){
92   
93   //create ESDs histograms in ESDs subdir
94   const Bool_t expert   = kTRUE ; 
95   const Bool_t image    = kTRUE ; 
96
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);
103
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);
110
111   TH2F * histESDNCls = 
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);
116   
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);
122
123   TH1F * histPHLT 
124     = new TH1F("hPHLT", "P distribution for all events from HLT; P [GeV/c]", 
125                100, -0.5, 99.5);
126   histPHLT->Sumw2();
127   Add2ESDsList(histPHLT, kPHLT, !expert, image);
128
129   TH1F * histPOffline 
130     = new TH1F("hPOffline", 
131                "P distribution for all events from offline; P [GeV/c]",  
132                100, -0.5, 99.5);
133   histPOffline->Sumw2();
134   Add2ESDsList(histPOffline, kPOffline, !expert, image);
135   
136   TH1F * histPHLTFired 
137     = new TH1F("hPHLTFired", 
138                "P distribution for fired events from HLT; P [GeV/c]",  
139                100, -0.5, 99.5);
140   histPHLTFired->Sumw2();
141   Add2ESDsList(histPHLTFired, kPHLTFired, !expert, image); 
142  
143   TH1F * histPOfflineFired 
144     = new TH1F("hPOfflineFired",  
145                "P distribution for fired events from offline; P [GeV/c]",   
146                100, -0.5, 99.5);
147   histPOfflineFired->Sumw2();
148   Add2ESDsList(histPOfflineFired, kPOfflineFired, !expert, image);
149
150   TH1F * histPRatio 
151     = new TH1F("hPRatio", 
152                "Ratio of P distribution for all events; P [GeV/c]",
153                100, -0.5, 99.5);
154   histPRatio->Sumw2();
155   Add2ESDsList(histPRatio, kPRatio, !expert, image);
156   
157   TH1F * histPRatioFired  
158     = new TH1F("hPRatioFired",  
159                "Ratio of P distribution for fired events; P [GeV/c]",
160                100, -0.5, 99.5);
161   histPRatioFired->Sumw2();
162   Add2ESDsList(histPRatioFired, kPRatioFired, !expert, image);
163
164   TH1F * histPtHLT  
165     = new TH1F("hPtHLT", 
166                "P_{T} distribution for all events from HLT; P_{T} [GeV/c]",  
167                200, -100.5, 99.5);
168   histPtHLT->Sumw2();
169   Add2ESDsList(histPtHLT, kPtHLT, !expert, image); 
170  
171   TH1F * histPtOffline  
172     = new TH1F("hPtOffline",  
173                "P_{T} distribution for all events from offline; P_{T} [GeV/c]",   
174                200, -100.5, 99.5);
175   histPtOffline->Sumw2();
176   Add2ESDsList(histPtOffline, kPtOffline, !expert, image);
177
178   TH1F * histPtHLTFired   
179     = new TH1F("hPtHLTFired",  
180                "P_{T} distribution for fired events from HLT; P_{T} [GeV/c]",   
181                200, -100.5, 99.5);
182   histPtHLTFired->Sumw2();
183   Add2ESDsList(histPtHLTFired, kPtHLTFired, !expert, image);  
184   
185   TH1F * histPtOfflineFired   
186     = new TH1F("hPtOfflineFired",   
187                "P_{T} distribution for fired events from offline; P_{T} [GeV/c]", 
188                200, -100.5, 99.5);
189   histPtOfflineFired->Sumw2();
190   Add2ESDsList(histPtOfflineFired, kPtOfflineFired, !expert, image); 
191
192   TH1F * histNClsPerTrkHLT
193     = new TH1F("hNClsPerTrkHLT", "Clusters per HLT track; N cluster; Counts",
194                200, 0, 200);
195   histNClsPerTrkHLT->Sumw2();
196   Add2ESDsList(histNClsPerTrkHLT, kNClsPerTrkHLT, !expert, image);
197
198   TH1F * histNClsPerTrkOffline 
199     = new TH1F("hNClsPerTrkOffline", 
200                "Clusters per offline track; N cluster; Counts", 
201                200, 0, 200); 
202   histNClsPerTrkOffline->Sumw2();
203   Add2ESDsList(histNClsPerTrkOffline, kNClsPerTrkOffline, !expert, image);
204
205 TH1F * histNClsPerTrkHLTFired 
206   = new TH1F("hNClsPerTrkHLTFired", 
207              "Clusters per HLT track from HLT fired events; N cluster; Counts", 
208              200, 0, 200); 
209  histNClsPerTrkHLTFired->Sumw2();
210  Add2ESDsList(histNClsPerTrkHLTFired, kNClsPerTrkHLTFired, !expert, image); 
211  
212   TH1F * histNClsPerTrkOfflineFired  
213     = new TH1F("hNClsPerTrkOfflineFired",  
214                "Clusters per offline track from HLT fired events; N cluster; Counts",  
215                200, 0, 200);  
216   histNClsPerTrkOfflineFired->Sumw2();
217   Add2ESDsList(histNClsPerTrkOfflineFired, kNClsPerTrkOfflineFired, !expert, image);
218
219   TH1F * histPhiHLT = 
220     new TH1F("hPhiHLT", "Phi distribution of HLT tracks; Phi; Counts", 
221              360, 0, 360);
222   histPhiHLT->Sumw2();
223   Add2ESDsList(histPhiHLT, kPhiHLT, !expert, image);
224
225   TH1F * histPhiOffline =  
226     new TH1F("hPhiOffline", 
227              "Phi distribution of offline tracks; Phi; Counts",  
228              360, 0, 360); 
229   histPhiOffline->Sumw2(); 
230   Add2ESDsList(histPhiOffline, kPhiOffline, !expert, image);
231
232   TH1F * histPhiHLTFired =  
233     new TH1F("hPhiHLTFired", "Phi distribution of HLT tracks from HLT fired event ; Phi; Counts",  
234              360, 0, 360); 
235   histPhiHLTFired->Sumw2(); 
236   Add2ESDsList(histPhiHLTFired, kPhiHLTFired, !expert, image); 
237  
238   TH1F * histPhiOfflineFired =   
239     new TH1F("hPhiOfflineFired",  
240              "Phi distribution of offline tracks from HLT fired events; Phi; Counts",   
241              360, 0, 360);  
242   histPhiOfflineFired->Sumw2();  
243   Add2ESDsList(histPhiOfflineFired, kPhiOfflineFired, !expert, image); 
244
245   TH1F * histEtaHLT =
246     new TH1F("hEtaHLT", "Eta distribution of HLT tracks; Eta; Counts",
247              200, -1, 1);
248   histEtaHLT->Sumw2();
249   Add2ESDsList(histEtaHLT, kEtaHLT, !expert, image);
250
251   TH1F * histEtaOffline = 
252     new TH1F("hEtaOffline", "Eta distribution of offline tracks; Eta; Counts", 
253              200, -1, 1); 
254   histEtaHLT->Sumw2(); 
255   Add2ESDsList(histEtaOffline, kEtaOffline, !expert, image);
256
257   TH1F * histEtaHLTFired = 
258     new TH1F("hEtaHLTFired", 
259              "Eta distribution of HLT tracks from HLT fired events; Eta; Counts", 
260              200, -1, 1); 
261   histEtaHLTFired->Sumw2(); 
262   Add2ESDsList(histEtaHLTFired, kEtaHLTFired, !expert, image); 
263  
264   TH1F * histEtaOfflineFired =  
265     new TH1F("hEtaOfflineFired", 
266              "Eta distribution of offline tracks from HLT fired events; Eta; Counts",  
267              200, -1, 1);  
268   histEtaHLTFired->Sumw2();  
269   Add2ESDsList(histEtaOfflineFired, kEtaOfflineFired, !expert, image);
270
271 }
272
273 void AliHLTTPCQADataMaker::MakeESDs(AliESDEvent * esd, AliESDEvent* hltesd)
274 {
275   // HLT QA on ESDs
276   if (!esd || !hltesd) {
277     AliError("invalid parameter: missing ESDs");
278     return;
279   }
280
281   // make QA data from ESDs
282
283   const Int_t nESDTracks = esd->GetNumberOfTracks();
284   const Int_t nHLTesdTracks = hltesd->GetNumberOfTracks();
285   FillESDsData(kMultiplicity,nESDTracks, nHLTesdTracks);
286
287   Int_t nClsHLT = 0;
288   Int_t nClsOffline = 0;
289   
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();
298   }
299
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();
308   } 
309   
310   if(nESDTracks)
311     nClsOffline /= nESDTracks;
312   if(nHLTesdTracks)
313     nClsHLT /= nHLTesdTracks;
314
315   FillESDsData(kNCls,nClsOffline, nClsHLT);
316
317
318   if(hltesd->IsHLTTriggerFired()){
319     FillESDsData(kMultiplicityFired,nESDTracks, nHLTesdTracks);
320     FillESDsData(kNClsFired,nClsOffline, nClsHLT);
321     
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());
329     } 
330    
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());
338     }  
339
340     
341   }
342
343 }