]> git.uio.no Git - u/mrichter/AliRoot.git/blame - PWG1/AliAnaVZEROQA.cxx
constant bining for q/pt
[u/mrichter/AliRoot.git] / PWG1 / AliAnaVZEROQA.cxx
CommitLineData
ecb38463 1/**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3 * *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
6 * *
7 * Permission to use, copy, modify and distribute this software and its *
8 * documentation strictly for non-commercial purposes is hereby granted *
9 * without fee, provided that the above copyright notice appears in all *
10 * copies and that both the copyright notice and this permission notice *
11 * appear in the supporting documentation. The authors make no claims *
12 * about the suitability of this software for any purpose. It is *
13 * provided "as is" without express or implied warranty. *
14 **************************************************************************/
15
16//------------------------------
17// Analysis task for quality-assurance
18// of VZERO ESD
19//
20// 05/12/2009 cvetan.cheshkov@cern.ch
21//------------------------------
22
23#include "TChain.h"
24#include "TROOT.h"
25#include "TFile.h"
26#include "TH1F.h"
27#include "TH2F.h"
ecb38463 28#include "TCanvas.h"
29#include "AliLog.h"
30#include "AliESDEvent.h"
31#include "AliESDVZERO.h"
32
33#include "AliAnaVZEROQA.h"
34
35ClassImp(AliAnaVZEROQA)
36
37AliAnaVZEROQA::AliAnaVZEROQA():
38AliAnalysisTaskSE("AliAnaVZEROQA"),
39 fListOfHistos(0),
40
41 fhAdcNoTimeA(0),
42 fhAdcWithTimeA(0),
43 fhAdcNoTimeC(0),
44 fhAdcWithTimeC(0),
45
46 fhAdcPMTNoTime(0),
47 fhAdcPMTWithTime(0),
48
49 fhTimeA(0),
50 fhTimeC(0),
51
52 fhWidthA(0),
53 fhWidthC(0),
54
55 fhTimePMT(0),
56 fhWidthPMT(0),
57
58 fhAdcWidthA(0),
59 fhAdcWidthC(0),
60
61 fhTimeCorr(0),
62
63 fhAdcTimeA(0),
64 fhAdcTimeC(0),
65
66 fV0a(0),
67 fV0c(0),
68 fV0multA(0),
69 fV0multC(0),
892effee 70 fV0ampl(0),
71
72 fhTimePMTCorr(0),
73 fhEvents(0),
74
75 fhVtxXYBB(0),
76 fhVtxZBB(0),
77 fhVtxXYBGA(0),
78 fhVtxZBGA(0),
79 fhVtxXYBGC(0),
80 fhVtxZBGC(0)
ecb38463 81{
82 // Default constructor
83 // Define input and output slots here
84 // Input slot #0 works with a TChain
85 DefineInput(0, TChain::Class());
86 // Output slot #1 TList
87 DefineOutput(1, TList::Class());
88}
89
90AliAnaVZEROQA::AliAnaVZEROQA(const char* name):
91AliAnalysisTaskSE(name),
92 fListOfHistos(0),
93 fhAdcNoTimeA(0),
94 fhAdcWithTimeA(0),
95 fhAdcNoTimeC(0),
96 fhAdcWithTimeC(0),
97
98 fhAdcPMTNoTime(0),
99 fhAdcPMTWithTime(0),
100
101 fhTimeA(0),
102 fhTimeC(0),
103
104 fhWidthA(0),
105 fhWidthC(0),
106
107 fhTimePMT(0),
108 fhWidthPMT(0),
109
110 fhAdcWidthA(0),
111 fhAdcWidthC(0),
112
113 fhTimeCorr(0),
114
115 fhAdcTimeA(0),
116 fhAdcTimeC(0),
117
118 fV0a(0),
119 fV0c(0),
120 fV0multA(0),
121 fV0multC(0),
892effee 122 fV0ampl(0),
123
124 fhTimePMTCorr(0),
125 fhEvents(0),
126
127 fhVtxXYBB(0),
128 fhVtxZBB(0),
129 fhVtxXYBGA(0),
130 fhVtxZBGA(0),
131 fhVtxXYBGC(0),
132 fhVtxZBGC(0)
ecb38463 133{
134 // Constructor
135 AliInfo("Constructor AliAnaVZEROQA");
136 // Define input and output slots here
137 // Input slot #0 works with a TChain
138 DefineInput(0, TChain::Class());
139 // Output slot #1 TList
140 DefineOutput(1, TList::Class());
141}
142
143TH1F * AliAnaVZEROQA::CreateHisto1D(const char* name, const char* title,Int_t nBins,
144 Double_t xMin, Double_t xMax,
145 const char* xLabel, const char* yLabel)
146{
147 // create a histogram
148 TH1F* result = new TH1F(name, title, nBins, xMin, xMax);
149 result->SetOption("E");
150 if (xLabel) result->GetXaxis()->SetTitle(xLabel);
151 if (yLabel) result->GetYaxis()->SetTitle(yLabel);
152 result->SetMarkerStyle(kFullCircle);
153 return result;
154}
155
156TH2F * AliAnaVZEROQA::CreateHisto2D(const char* name, const char* title,Int_t nBinsX,
157 Double_t xMin, Double_t xMax,
158 Int_t nBinsY,
159 Double_t yMin, Double_t yMax,
160 const char* xLabel, const char* yLabel)
161{
162 // create a histogram
163 TH2F* result = new TH2F(name, title, nBinsX, xMin, xMax, nBinsY, yMin, yMax);
164 if (xLabel) result->GetXaxis()->SetTitle(xLabel);
165 if (yLabel) result->GetYaxis()->SetTitle(yLabel);
166 return result;
167}
168
169void AliAnaVZEROQA::UserCreateOutputObjects()
170{
171 // Create histograms
172 AliInfo("AliAnaVZEROQA::UserCreateOutputObjects");
173 // Create output container
174 fListOfHistos = new TList();
13317536 175 fListOfHistos->SetOwner();
ecb38463 176
177 fhAdcNoTimeA = CreateHisto1D("hAdcNoTimeA","ADC (no Leading Time) V0A",200,0,200,"ADC charge","Entries");
178 fhAdcWithTimeA = CreateHisto1D("hAdcWithTimeA","ADC ( with Leading Time) V0A",200,0,200,"ADC charge","Entries");
179 fhAdcNoTimeC = CreateHisto1D("hAdcNoTimeC","ADC (no Leading Time) V0C",200,0,200,"ADC charge","Entries");
180 fhAdcWithTimeC = CreateHisto1D("hAdcWithTimeC","ADC ( with Leading Time) V0C",200,0,200,"ADC charge","Entries");
181
182 fhAdcPMTNoTime = CreateHisto2D("hadcpmtnotime","ADC vs PMT index (no leading time)",64,-0.5,63.5,200,0,200,"PMT index","ADC charge");
183 fhAdcPMTWithTime = CreateHisto2D("hadcpmtwithtime","ADC vs PMT index (with leading time)",64,-0.5,63.5,200,0,200,"PMT index","ADC charge");
184
185 fhTimeA = CreateHisto1D("htimepmtA","Time measured by TDC V0A",400,-100,100,"Leading time (ns)","Entries");
186 fhTimeC = CreateHisto1D("htimepmtC","Time measured by TDC V0C",400,-100,100,"Leading time (ns)","Entries");
187
188 fhWidthA = CreateHisto1D("hwidthA","Signal width measured by TDC V0A",200,0,100,"Signal width (ns)","Entries");
189 fhWidthC = CreateHisto1D("hwidthC","Signal width measured by TDC V0C",200,0,100,"Signal width (ns)","Entries");
190
191 fhTimePMT = CreateHisto2D("htimepmt","Time measured by TDC vs PMT index",64,-0.5,63.5,200,0,100,"PMT Index","Leading time (ns)");
192 fhWidthPMT = CreateHisto2D("hwidthpmt","Time width vs PMT index",64,-0.5,63.5,200,0,100,"PMT Index","Signal width (ns)");
193
194 fhAdcWidthA = CreateHisto2D("hadcwidthA","Time width vs ADC V0A",200,0,200,200,0,100,"ADC charge","Width (ns)");
195 fhAdcWidthC = CreateHisto2D("hadcwidthC","Time width vs ADC V0C",200,0,200,200,0,100,"ADC charge","Width (ns)");
196
197 fhTimeCorr = CreateHisto2D("htimecorr","Average time C side vs. A side",200,0,100,200,0,100,"Time V0A (ns)","Time V0C (ns");
198
199 fhAdcTimeA = CreateHisto2D("hAdcTimeA","ADC vs Time V0A",1000,-100,100,200,0,200,"Time (ns)","ADC charge");
200 fhAdcTimeC = CreateHisto2D("hAdcTimeC","ADC vs Time V0C",1000,-100,100,200,0,200,"Time (ns)","ADC charge");
201
202 fV0a = CreateHisto1D("hV0a","Number of fired PMTs (V0A)",65,-0.5,64.5);
203 fV0c = CreateHisto1D("hV0c","Number of fired PMTs (V0C)",65,-0.5,64.5);
204 fV0multA = CreateHisto1D("hV0multA","Total reconstructed multiplicity (V0A)",100,0.,1000.);
205 fV0multC = CreateHisto1D("hV0multC","Total reconstructed multiplicity (V0C)",100,0.,1000.);
206 fV0ampl = CreateHisto1D("hV0ampl","V0 multiplicity in single channel (all V0 channels)",400,-0.5,99.5);
207
892effee 208 fhTimePMTCorr = CreateHisto2D("htimepmtcorr","Time measured by TDC (corrected for slewing, channels aligned) vs PMT index",64,-0.5,63.5,200,0,100,"PMT Index","Leading time (ns)");
209 fhEvents = CreateHisto2D("hEvents","V0C vs V0A (empty,bb,bg)",3,-0.5,2.5,3,-0.5,2.5);
210
211 fhVtxXYBB = CreateHisto2D("fhVtxXYBB","XY SPD vertex (bb)",200,-2,2,200,-2,2);
212 fhVtxZBB = CreateHisto1D("fhVtxZBB","Z SPD vertex (bb)",400,-50,50);
213 fhVtxXYBGA = CreateHisto2D("fhVtxXYBGA","XY SPD vertex (bga)",200,-2,2,200,-2,2);
214 fhVtxZBGA = CreateHisto1D("fhVtxZBGA","Z SPD vertex (bga)",400,-50,50);
215 fhVtxXYBGC = CreateHisto2D("fhVtxXYBGC","XY SPD vertex (bgc)",200,-2,2,200,-2,2);
216 fhVtxZBGC = CreateHisto1D("fhVtxZBGC","Z SPD vertex (bgc)",400,-50,50);
217
ecb38463 218 fListOfHistos->Add(fhAdcNoTimeA);
219 fListOfHistos->Add(fhAdcWithTimeA);
220 fListOfHistos->Add(fhAdcNoTimeC);
221 fListOfHistos->Add(fhAdcWithTimeC);
222
223 fListOfHistos->Add(fhAdcPMTNoTime);
224 fListOfHistos->Add(fhAdcPMTWithTime);
225
226 fListOfHistos->Add(fhTimeA);
227 fListOfHistos->Add(fhTimeC);
228
229 fListOfHistos->Add(fhWidthA);
230 fListOfHistos->Add(fhWidthC);
231
232 fListOfHistos->Add(fhTimePMT);
233 fListOfHistos->Add(fhWidthPMT);
234
235 fListOfHistos->Add(fhAdcWidthA);
236 fListOfHistos->Add(fhAdcWidthC);
237
238 fListOfHistos->Add(fhTimeCorr);
239
240 fListOfHistos->Add(fhAdcTimeA);
241 fListOfHistos->Add(fhAdcTimeC);
242
243 fListOfHistos->Add(fV0a);
244 fListOfHistos->Add(fV0c);
245 fListOfHistos->Add(fV0multA);
246 fListOfHistos->Add(fV0multC);
247 fListOfHistos->Add(fV0ampl);
892effee 248
249 fListOfHistos->Add(fhTimePMTCorr);
250 fListOfHistos->Add(fhEvents);
251
252 fListOfHistos->Add(fhVtxXYBB);
253 fListOfHistos->Add(fhVtxZBB);
254 fListOfHistos->Add(fhVtxXYBGA);
255 fListOfHistos->Add(fhVtxZBGA);
256 fListOfHistos->Add(fhVtxXYBGC);
257 fListOfHistos->Add(fhVtxZBGC);
3b5284bd 258
259 PostData(1, fListOfHistos);
ecb38463 260}
261
262void AliAnaVZEROQA::UserExec(Option_t */*option*/)
263{
f0a651f0 264 // Fill the QA histograms
265 // using ESD data
ecb38463 266 AliVEvent* event = InputEvent();
267 if (!event) {
268 Printf("ERROR: Could not retrieve event");
269 return;
270 }
271
272 AliESDEvent* esd = dynamic_cast<AliESDEvent*>(event);
713a03ae 273 if (!esd) {
274 Printf("ERROR: No ESD event");
275 return;
276 }
ecb38463 277 AliESDVZERO* esdV0 = esd->GetVZEROData();
713a03ae 278 if (!esdV0) {
279 Printf("ERROR: No ESD VZERO");
280 return;
281 }
ecb38463 282
283 Float_t timeA = 0,timeC = 0;
284 Int_t ntimeA = 0, ntimeC = 0;
285 for (Int_t i=0; i<64; ++i) {
286 if (esdV0->GetTime(i) < 1e-6) {
287 if (i >= 32) {
288 fhAdcNoTimeA->Fill(esdV0->GetAdc(i));
289 }
290 else {
291 fhAdcNoTimeC->Fill(esdV0->GetAdc(i));
292 }
293 fhAdcPMTNoTime->Fill(i,esdV0->GetAdc(i));
294 }
295 else {
296 if (i >= 32) {
297 fhAdcWithTimeA->Fill(esdV0->GetAdc(i));
298 }
299 else {
300 fhAdcWithTimeC->Fill(esdV0->GetAdc(i));
301 }
302 fhAdcPMTWithTime->Fill(i,esdV0->GetAdc(i));
303 }
304
305 if (i >= 32) {
306 fhTimeA->Fill(esdV0->GetTime(i));
307 fhWidthA->Fill(esdV0->GetWidth(i));
308 fhAdcWidthA->Fill(esdV0->GetAdc(i),esdV0->GetTime(i));
309 fhAdcTimeA->Fill(esdV0->GetTime(i),esdV0->GetAdc(i));
310 }
311 else {
312 fhTimeC->Fill(esdV0->GetTime(i));
313 fhWidthC->Fill(esdV0->GetWidth(i));
314 fhAdcWidthC->Fill(esdV0->GetAdc(i),esdV0->GetTime(i));
315 fhAdcTimeC->Fill(esdV0->GetTime(i),esdV0->GetAdc(i));
316 }
317 fhTimePMT->Fill(i,esdV0->GetTime(i));
318 fhWidthPMT->Fill(i,esdV0->GetWidth(i));
319
892effee 320 Float_t correctedTime = CorrectLeadingTime(i,esdV0->GetTime(i),esdV0->GetAdc(i));
321 fhTimePMTCorr->Fill(i,correctedTime);
322
ecb38463 323 if (esdV0->GetTime(i) > 1e-6) {
324 if (i >= 32) {
892effee 325 timeA += correctedTime;
ecb38463 326 ntimeA++;
327 }
328 else {
892effee 329 timeC += correctedTime;
ecb38463 330 ntimeC++;
331 }
332 }
333 }
334
335 if (ntimeA > 0) timeA = timeA/ntimeA;
336 if (ntimeC > 0) timeC = timeC/ntimeC;
337
338 fhTimeCorr->Fill(timeA,timeC);
339
892effee 340 Int_t flaga = 0, flagc = 0;
341 TString stra("emptyA"),strc("emptyC");
342 if (timeA > 48 && timeA < 62) { flaga = 1; stra = "BBA"; }
343 if (timeA > 26 && timeA < 33) { flaga = 2; stra = "BGA"; }
344 if (timeC > 49 && timeC < 60) { flagc = 1; strc = "BBC"; }
345 if (timeC > 43 && timeC < 48.5) { flagc = 2; strc = "BGC"; }
346
347 fhEvents->Fill(flaga,flagc);
348
349 const AliESDVertex *vtx = esd->GetPrimaryVertexSPD();
350
351 if (flaga <= 1 && flagc <=1) {
352 fhVtxXYBB->Fill(vtx->GetXv(),vtx->GetYv());
353 fhVtxZBB->Fill(vtx->GetZv());
354 }
355 else {
356 if (flaga == 2) {
357 fhVtxXYBGA->Fill(vtx->GetXv(),vtx->GetYv());
358 fhVtxZBGA->Fill(vtx->GetZv());
359 }
360 if (flagc == 2) {
361 fhVtxXYBGC->Fill(vtx->GetXv(),vtx->GetYv());
362 fhVtxZBGC->Fill(vtx->GetZv());
363 }
364 }
365
ecb38463 366 fV0a->Fill(esdV0->GetNbPMV0A());
367 fV0c->Fill(esdV0->GetNbPMV0C());
368 fV0multA->Fill(esdV0->GetMTotV0A());
369 fV0multC->Fill(esdV0->GetMTotV0C());
370 for(Int_t i = 0; i < 64; i++) {
371 fV0ampl->Fill(esdV0->GetMultiplicity(i));
372 }
373
374 // Post output data.
375 PostData(1, fListOfHistos);
376}
377
378void AliAnaVZEROQA::Terminate(Option_t *)
379{
f0a651f0 380 // Store the output histograms
381 // from the list
ecb38463 382 fListOfHistos = dynamic_cast<TList*>(GetOutputData(1));
383 if (!fListOfHistos) {
384 Printf("ERROR: fListOfHistos not available");
385 return;
386 }
387
ecb38463 388 Info("AliAnaVZEROQA", "Successfully finished");
389}
390
f0a651f0 391Float_t AliAnaVZEROQA::CorrectLeadingTime(Int_t i, Float_t time, Float_t adc) const
892effee 392{
393 // Correct for slewing and align the channels
394
f0a651f0 395 if (time < 1e-6) return 0;
892effee 396
397 // Time offsets between channels
398 Float_t timeShift[64] = {30.2914 , 30.0019 , 30.7429 , 30.1997 , 30.1511 , 29.6437 , 30.0609 , 29.5452 , 30.1437 , 30.745 , 30.7537 , 30.446 , 30.2771 , 30.838 , 30.3748 , 30.0635 , 30.1786 , 30.282 , 31.0992 , 30.7491 , 30.624 , 30.9268 , 30.6585 , 30.4895 , 31.5815 , 31.3871 , 31.2032 , 31.5778 , 31.0838 , 31.2259 , 31.2122 , 31.5989 , 28.3792 , 28.8325 , 27.8719 , 28.3475 , 26.9925 , 27.9300 , 28.4223 , 28.4996 , 28.2934 , 28.1281 , 27.209 , 28.5327 , 28.1181 , 28.0888 , 29.5111 , 28.6601 , 29.7705 , 29.6531 , 30.3373 , 30.2345 , 30.5935 , 29.8164 , 30.2235 , 29.6505 , 30.1225 , 31.2045 , 30.8399 , 30.6789 , 30.2784 , 31.7028 , 31.4239 , 30.1814};
399 time -= timeShift[i];
400
401 // Slewing correction
f0a651f0 402 if (adc < 1e-6) return time;
892effee 403
404 time += 30.;
405 if (adc > 300.) adc = 300.;
406 if (adc > 70.) {
407 return (time -
408 2.93028e+01 +
409 adc*1.25188e-02 -
410 adc*adc*2.68348e-05);
411 }
412 else {
413 return (time -
414 3.52314e+01 +
415 adc*5.99289e-01 -
416 adc*adc*2.74668e-02 +
417 adc*adc*adc*6.61224e-04 -
418 adc*adc*adc*adc*7.77105e-06 +
419 adc*adc*adc*adc*adc*3.51229e-08);
420 }
421}
ecb38463 422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482