]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGPP/VZERO/AliAnaVZEROPbPb.cxx
Merge branch 'master' of https://git.cern.ch/reps/AliRoot
[u/mrichter/AliRoot.git] / PWGPP / VZERO / AliAnaVZEROPbPb.cxx
1 #include "TChain.h"
2 #include "TTree.h"
3 #include "TH1F.h"
4 #include "TH2F.h"
5
6 #include "AliAnalysisTask.h"
7 #include "AliAnalysisManager.h"
8
9 #include "AliESDEvent.h"
10 #include "AliESDInputHandler.h"
11 #include "AliAnaVZEROPbPb.h"
12 #include "AliMultiplicity.h"
13 #include "AliESDUtils.h"
14 #include "AliCentrality.h"
15
16 // VZERO includes
17 #include "AliESDVZERO.h"
18
19 ClassImp(AliAnaVZEROPbPb)
20
21 AliAnaVZEROPbPb::AliAnaVZEROPbPb() 
22   : AliAnalysisTaskSE(), fESD(0), fEsdV0(0), fOutputList(0), fNClasses(0), fClassesNames(0),
23   fNFlags(0),
24   fhAdcPMTNoTime(0),
25   fhAdcPMTWithTime(0),
26   fhTimePMT(0),
27   fhWidthPMT(0),
28   fhTimeCorr(0),
29   fhPmtMult(0),
30   fhV0ampl(0),
31   fhEvents(0),
32   fhVtxXYBB(0),
33   fhVtxZBB(0),
34   fhVtxXYBGA(0),
35   fhVtxZBGA(0),
36   fhVtxXYBGC(0),
37   fhVtxZBGC(0),
38   fhL2Triggers(0),
39 fhOnlineCharge(0),
40 fhRecoMult(0),
41 fhRecoMultPMT(0),
42 fhV0vsSPDCentrality(0),
43 fhTriggerBits(0),
44 fhTotRecoMult(0),
45 fhCentrality(0),
46 fhEqualizedMult(0),
47 fhEqualizedMultSum(0),
48 fNBinTotMult(100),
49 fTotMultMax(30000),
50 fNBinMult(100),
51 fV0AMultMax(20000),
52 fV0CMultMax(30000),
53 fNBinOnlineCharge(100),
54 fV0AOnlineChargeMax(20000),
55 fV0COnlineChargeMax(30000),
56 fNBinEquaMult(400),
57 fEquaMultMax(1000),
58 fNBinSumsEqMult(100),
59 fV0AEqMultMax(4000.),
60 fV0CEqMultMax(4000.),
61 fZvtxCut(10.)
62
63 {
64   // Constructor
65   // Init arrays
66   for(Int_t i = 0; i < 2; ++i) {
67         fhAdcNoTime[i] = fhAdcWithTime[i] =  fhWidth[i] =  fhTime[i] = NULL;
68         fhAdcTime[i] =  fhAdcWidth[i] = NULL;
69   }
70
71 }
72
73 //________________________________________________________________________
74 AliAnaVZEROPbPb::AliAnaVZEROPbPb(const char *name) 
75   : AliAnalysisTaskSE(name), fESD(0), fEsdV0(0), fOutputList(0),fNClasses(0),fClassesNames(0),
76   fNFlags(0),
77   fhAdcPMTNoTime(0),
78   fhAdcPMTWithTime(0),
79   fhTimePMT(0),
80   fhWidthPMT(0),
81   fhTimeCorr(0),
82   fhPmtMult(0),
83   fhV0ampl(0),
84   fhEvents(0),
85   fhVtxXYBB(0),
86   fhVtxZBB(0),
87   fhVtxXYBGA(0),
88   fhVtxZBGA(0),
89   fhVtxXYBGC(0),
90   fhVtxZBGC(0),
91   fhL2Triggers(0),
92 fhOnlineCharge(0),
93 fhRecoMult(0),
94 fhRecoMultPMT(0),
95 fhV0vsSPDCentrality(0),
96 fhTriggerBits(0),
97 fhTotRecoMult(0),
98 fhCentrality(0),
99 fhEqualizedMult(0),
100 fhEqualizedMultSum(0),
101 fNBinTotMult(100),
102 fTotMultMax(30000),
103 fNBinMult(100),
104 fV0AMultMax(20000),
105 fV0CMultMax(30000),
106 fNBinOnlineCharge(100),
107 fV0AOnlineChargeMax(20000),
108 fV0COnlineChargeMax(30000),
109 fNBinEquaMult(400),
110 fEquaMultMax(1000),
111 fNBinSumsEqMult(100),
112 fV0AEqMultMax(4000.),
113 fV0CEqMultMax(4000.),
114 fZvtxCut(10.)
115 {
116   // Constructor
117   // Init arrays
118   for(Int_t i = 0; i < 2; ++i) {
119         fhAdcNoTime[i] = fhAdcWithTime[i] =  fhWidth[i] =  fhTime[i] = NULL;
120         fhAdcTime[i] =  fhAdcWidth[i] = NULL;
121   }
122   // Define input and output slots here
123   // Input slot #0 works with a TChain
124   DefineInput(0, TChain::Class());
125   // Output slot #0 id reserved by the base class for AOD
126   // Output slot #1 writes into a TH1 container
127   DefineOutput(1, TList::Class());
128 }
129 //________________________________________________________________________
130 void AliAnaVZEROPbPb::SetClassesNames(const Char_t * nameList){
131   // Initialize the class names
132   // which are used in the trigger split
133         TString  names("AllClasses,");
134         names += nameList;
135         fClassesNames = names.Tokenize(",");
136         fNClasses = fClassesNames->GetEntriesFast();
137 }
138 //________________________________________________________________________
139 TH1F * AliAnaVZEROPbPb::CreateHisto1D(const char* name, const char* title,Int_t nBins, 
140                                     Double_t xMin, Double_t xMax,
141                                     const char* xLabel, const char* yLabel)
142 {
143   // create a histogram
144   // and set the axis labels and the style
145   TH1F* result = new TH1F(name, title, nBins, xMin, xMax);
146   result->SetOption("hist");
147   if (xLabel) result->GetXaxis()->SetTitle(xLabel);
148   if (yLabel) result->GetYaxis()->SetTitle(yLabel);
149   result->SetMarkerStyle(kFullCircle);
150   return result;
151 }
152
153 //________________________________________________________________________
154 TH2F * AliAnaVZEROPbPb::CreateHisto2D(const char* name, const char* title,Int_t nBinsX, 
155                                     Double_t xMin, Double_t xMax,
156                                     Int_t nBinsY,
157                                     Double_t yMin, Double_t yMax,
158                                     const char* xLabel, const char* yLabel)
159 {
160   // create a histogram
161   // and set the axis labels and the style
162   TH2F* result = new TH2F(name, title, nBinsX, xMin, xMax, nBinsY, yMin, yMax);
163   if (xLabel) result->GetXaxis()->SetTitle(xLabel);
164   if (yLabel) result->GetYaxis()->SetTitle(yLabel);
165   return result;
166 }
167
168
169 //________________________________________________________________________
170 void AliAnaVZEROPbPb::UserCreateOutputObjects()
171 {
172   // Create histograms
173   // Called once
174         if(fNClasses==0) {
175                 AliFatal("No Classes Defined");
176                 return;
177         }
178
179   fOutputList = new TList();
180   fOutputList->SetOwner(kTRUE);
181
182
183         CreateHistosPerL2Trigger();
184     CreateQAHistos();
185
186   PostData(1, fOutputList);
187  }
188 //________________________________________________________________________
189 void AliAnaVZEROPbPb::Init()
190 {
191   // Nothing here
192   // ...
193 }
194
195 //________________________________________________________________________
196 void AliAnaVZEROPbPb::UserExec(Option_t *) 
197 {
198   // Main loop
199   // Called for each event
200
201   fESD = dynamic_cast<AliESDEvent*>(InputEvent());
202   if (!fESD) {
203     printf("ERROR: fESD not available\n");
204     return;
205   }
206   
207   fEsdV0 = fESD->GetVZEROData();
208   if (!fEsdV0) {
209     Printf("ERROR: esd V0  not available");
210     return;
211   }
212
213   FillQAHistos();
214
215   Bool_t isSelected = kFALSE;
216
217   const AliESDVertex *vtx = fESD->GetPrimaryVertexSPD();
218
219   if ((vtx != NULL) && (vtx->GetStatus()) && (TMath::Abs(vtx->GetZv())<fZvtxCut) && (fEsdV0->GetV0ADecision()==1) && (fEsdV0->GetV0CDecision()==1)){
220         isSelected = kTRUE;
221   }
222
223   if (isSelected) FillPerL2TriggerHistos();
224
225         
226
227   PostData(1, fOutputList);
228 }      
229 //________________________________________________________________________
230 void AliAnaVZEROPbPb::CreateHistosPerL2Trigger(){
231   // Create the histograms
232   // for all the required L2 trigger classes    
233         fhOnlineCharge= new TH2F*[fNClasses];
234         fhCentrality= new TH1F*[fNClasses];
235         fhV0vsSPDCentrality= new TH2F*[fNClasses];
236         fhRecoMult= new TH2F*[fNClasses];
237         fhRecoMultPMT= new TH2F*[fNClasses];
238         fhTotRecoMult= new TH1F*[fNClasses];
239         fhTriggerBits= new TH1F*[fNClasses];
240         fhEqualizedMult= new TH2F*[fNClasses];
241         fhEqualizedMultSum = new TH2F*[fNClasses];
242         
243         fhL2Triggers = CreateHisto1D("hL2Triggers","L2 Triggers",fNClasses,0,fNClasses);
244         fOutputList->Add(fhL2Triggers);   
245         
246         TIter iter(fClassesNames);
247         TObjString* name;
248         Int_t iClass =0;
249         while((name = (TObjString*) iter.Next())){
250                 fhL2Triggers->GetXaxis()->SetBinLabel(iClass+1,name->String().Data());
251                 
252                 fhOnlineCharge[iClass] = CreateHisto2D(Form("hOnlineCharge_%s",name->String().Data()),Form("Online Charge for %s",name->String().Data()),fNBinOnlineCharge,0.,fV0AOnlineChargeMax,fNBinOnlineCharge,0.,fV0COnlineChargeMax,"V0A","V0C");
253                 fOutputList->Add(fhOnlineCharge[iClass]);         
254
255                 fhCentrality[iClass] = CreateHisto1D(Form("hV0Centrality_%s",name->String().Data()),Form("V0 centrality for %s",name->String().Data()),100,0.,100.);
256                 fOutputList->Add(fhCentrality[iClass]);   
257
258                 fhV0vsSPDCentrality[iClass] = CreateHisto2D(Form("hV0vsSPDCentrality_%s",name->String().Data()),Form("Centrality for %s",name->String().Data()),100,0.,100.,100,0.,100,"SPD Centrality (CL1)","V0 Centrality (V0M)");
259                 fOutputList->Add(fhV0vsSPDCentrality[iClass]);    
260
261                 fhRecoMult[iClass] = CreateHisto2D(Form("hRecoMult_%s",name->String().Data()),Form("Reco Multiplicity for %s",name->String().Data()),fNBinMult,0.,fV0AMultMax,fNBinMult,0.,fV0CMultMax,"V0A Offline Mult","V0C Offline Mult");
262                 fOutputList->Add(fhRecoMult[iClass]);     
263
264                 fhRecoMultPMT[iClass] = CreateHisto2D(Form("hRecoMultPMT_%s",name->String().Data()),Form("Reco Multiplicity per PMT for %s",name->String().Data()),64,-0.5,63.5,200,0.,2000.,"PMT channel","Offline Mult");
265                 fOutputList->Add(fhRecoMultPMT[iClass]);          
266
267                 fhTotRecoMult[iClass] = CreateHisto1D(Form("hTotRecoMult_%s",name->String().Data()),Form("Total Reco Multiplicity for %s",name->String().Data()),fNBinTotMult,0.,fTotMultMax,"V0A + V0C Offline Mult");
268                 fOutputList->Add(fhTotRecoMult[iClass]);          
269
270                 fhTriggerBits[iClass] = CreateHisto1D(Form("hTriggerBits_%s",name->String().Data()),Form("Trigger Bits for %s",name->String().Data()),16,-0.5,15.5);
271                 fhTriggerBits[iClass]->GetXaxis()->SetBinLabel(1,"BBA_AND_BBC");
272                 fhTriggerBits[iClass]->GetXaxis()->SetBinLabel(2,"BBA_OR_BBC");
273                 fhTriggerBits[iClass]->GetXaxis()->SetBinLabel(3,"BGA_AND_BBC");
274                 fhTriggerBits[iClass]->GetXaxis()->SetBinLabel(4,"BGA");
275                 fhTriggerBits[iClass]->GetXaxis()->SetBinLabel(5,"BBA_AND_BGC");
276                 fhTriggerBits[iClass]->GetXaxis()->SetBinLabel(6,"BGC");
277                 fhTriggerBits[iClass]->GetXaxis()->SetBinLabel(7,"CTA1_AND_CTC1");
278                 fhTriggerBits[iClass]->GetXaxis()->SetBinLabel(8,"CTA1_OR_CTC1");
279                 fhTriggerBits[iClass]->GetXaxis()->SetBinLabel(9,"CTA2_AND_CTC2");
280                 fhTriggerBits[iClass]->GetXaxis()->SetBinLabel(10,"CTA2_OR_CTC2");
281                 fhTriggerBits[iClass]->GetXaxis()->SetBinLabel(11,"MTA_AND_MTC");
282                 fhTriggerBits[iClass]->GetXaxis()->SetBinLabel(12,"MTA_OR_MTC");
283                 fhTriggerBits[iClass]->GetXaxis()->SetBinLabel(13,"BBA");
284                 fhTriggerBits[iClass]->GetXaxis()->SetBinLabel(14,"BBC");
285                 fhTriggerBits[iClass]->GetXaxis()->SetBinLabel(15,"BGA_OR_BGC");
286                 fhTriggerBits[iClass]->GetXaxis()->SetBinLabel(16,"All True BG");
287                 fOutputList->Add(fhTriggerBits[iClass]);
288         
289                 fhEqualizedMult[iClass] = CreateHisto2D(Form("hEqualizedMult_%s",name->String().Data()),Form("Equalized Multiplicity for %s",name->String().Data()),64,-0.5,63.5,fNBinEquaMult,0.,fEquaMultMax,"PMT channel","Equalized Multiplicity");
290                 fOutputList->Add(fhEqualizedMult[iClass]);
291         
292                 fhEqualizedMultSum[iClass] = CreateHisto2D(Form("hEqualizedMultSum_%s",name->String().Data()),Form("Summed Equalized Multiplicity for %s",name->String().Data()),fNBinSumsEqMult,0.,fV0AEqMultMax,fNBinSumsEqMult,0.,fV0CEqMultMax,"V0A","V0C");
293                 fOutputList->Add(fhEqualizedMultSum[iClass]);
294
295                 
296                 iClass++;
297         }
298         
299 }
300 //________________________________________________________________________
301 void AliAnaVZEROPbPb::CreateQAHistos(){
302   // Create the main
303   // QA histos
304   fNFlags    = CreateHisto2D("hNFlags","BB Flags",33,-0.5,32.5,33,-0.5,32.5,"V0A","V0C");
305   fOutputList->Add(fNFlags);
306
307         for(int iSide = 0; iSide < 2; ++iSide){
308                 TString side;
309                 if(iSide) side = "V0A";
310                 else side = "V0C";
311                 
312                 fhAdcNoTime[iSide] = CreateHisto1D(Form("hAdcNoTime%s",side.Data()),Form("ADC (no Leading Time) %s",side.Data()),200,0,200,"ADC charge","Entries");
313                 fOutputList->Add(fhAdcNoTime[iSide]);     
314                 fhAdcWithTime[iSide] = CreateHisto1D(Form("hAdcWithTime%s",side.Data()),Form("ADC (with Leading Time) %s",side.Data()),200,0,200,"ADC charge","Entries");
315                 fOutputList->Add(fhAdcWithTime[iSide]);
316                 fhTime[iSide] = CreateHisto1D(Form("htimepmt%s",side.Data()),Form("Time measured by TDC %s",side.Data()),1000,-100,100,"Leading time (ns)","Entries");
317                 fOutputList->Add(fhTime[iSide]);
318                 fhWidth[iSide] = CreateHisto1D(Form("hwidth%s",side.Data()),Form("Signal width measured by TDC %s",side.Data()),128,0,800,"Signal width (ns)","Entries");
319                 fOutputList->Add(fhWidth[iSide]);
320                 fhAdcWidth[iSide] = CreateHisto2D(Form("hadcwidth%s",side.Data()),Form("Time width vs ADC %s",side.Data()),200,0,1200,128,0,800,"ADC charge","Width (ns)");
321                 fOutputList->Add(fhAdcWidth[iSide]);
322                 fhAdcTime[iSide] = CreateHisto2D(Form("hAdcTime%s",side.Data()),Form("ADC vs Time %s",side.Data()),200,-100,100,200,0,200,"Time (ns)","ADC charge");
323                 fOutputList->Add(fhAdcTime[iSide]);
324         }
325
326   fhAdcPMTNoTime = CreateHisto2D("hadcpmtnotime","ADC vs PMT index (no leading time)",64,-0.5,63.5,200,0,200,"PMT index","ADC charge");
327   fhAdcPMTWithTime = CreateHisto2D("hadcpmtwithtime","ADC vs PMT index (with leading time)",64,-0.5,63.5,200,0,2000,"PMT index","ADC charge");
328
329   fhTimePMT = CreateHisto2D("htimepmt","Time measured by TDC vs PMT index",64,-0.5,63.5,200,-50,50,"PMT Index","Leading time (ns)");
330   fhWidthPMT = CreateHisto2D("hwidthpmt","Time width vs PMT index",64,-0.5,63.5,128,0,800,"PMT Index","Signal width (ns)");
331
332   fhTimeCorr = CreateHisto2D("htimecorr","Average time C side vs. A side",200,-50,50,200,-50,50,"Time V0A (ns)","Time V0C (ns");
333
334   fhV0ampl  = CreateHisto1D("hV0ampl","V0 multiplicity in single channel (all V0 channels)",500,-0.5,499.5);
335
336   fhEvents = CreateHisto2D("hTriggerDecision","Trigger Decision",3,-0.5,2.5,3,-0.5,2.5,"V0A Decision","V0C Decision");
337   fhEvents->GetXaxis()->SetBinLabel(1,"Fake");
338   fhEvents->GetXaxis()->SetBinLabel(2,"BB");
339   fhEvents->GetXaxis()->SetBinLabel(3,"BG");
340   fhEvents->GetYaxis()->SetBinLabel(1,"Fake");
341   fhEvents->GetYaxis()->SetBinLabel(2,"BB");
342   fhEvents->GetYaxis()->SetBinLabel(3,"BG");
343
344   fhVtxXYBB = CreateHisto2D("fhVtxXYBB","XY SPD vertex (bb)",200,-2,2,200,-2,2);
345   fhVtxZBB = CreateHisto1D("fhVtxZBB","Z SPD vertex (bb)",400,-50,50);
346   fhVtxXYBGA = CreateHisto2D("fhVtxXYBGA","XY SPD vertex (bga)",200,-2,2,200,-2,2);
347   fhVtxZBGA = CreateHisto1D("fhVtxZBGA","Z SPD vertex (bga)",400,-50,50);
348   fhVtxXYBGC = CreateHisto2D("fhVtxXYBGC","XY SPD vertex (bgc)",200,-2,2,200,-2,2);
349   fhVtxZBGC = CreateHisto1D("fhVtxZBGC","Z SPD vertex (bgc)",400,-50,50);
350
351   fhPmtMult = CreateHisto2D("hV0CellMult","Number of fired PMTs (V0C vs V0A)",33,-0.5,32.5,33,-0.5,32.5,"# Cell (V0A)","# Cell (V0C)");
352         fOutputList->Add(fhPmtMult);
353
354
355   fOutputList->Add(fhAdcPMTNoTime);
356   fOutputList->Add(fhAdcPMTWithTime);
357
358   fOutputList->Add(fhTimePMT);
359   fOutputList->Add(fhWidthPMT);
360
361   fOutputList->Add(fhTimeCorr);
362   fOutputList->Add(fhV0ampl);
363
364   fOutputList->Add(fhEvents);
365
366   fOutputList->Add(fhVtxXYBB);
367   fOutputList->Add(fhVtxZBB);
368   fOutputList->Add(fhVtxXYBGA);
369   fOutputList->Add(fhVtxZBGA);
370   fOutputList->Add(fhVtxXYBGC);
371   fOutputList->Add(fhVtxZBGC);
372   
373 }
374 //________________________________________________________________________
375 void AliAnaVZEROPbPb::FillPerL2TriggerHistos(){
376   // Fill the histos which are split
377   // by L2 trigger class
378    TString trigStr(fESD->GetFiredTriggerClasses());
379   
380         TIter iter(fClassesNames);
381         TObjString* name;
382         Int_t iClass =0;
383         
384    
385         TObjArray * tokens = trigStr.Tokenize(" ");
386     Int_t ntokens = tokens->GetEntriesFast();
387                 
388   while((name = (TObjString*) iter.Next())){
389         
390         Bool_t goodTrig = kFALSE;
391         if(iClass>0){
392         for (Int_t itoken = 0; itoken < ntokens; ++itoken) {
393                         if (!((((TObjString*)tokens->At(itoken))->String()).Contains("-B-"))) continue;
394                         if ((((TObjString*)tokens->At(itoken))->String()).BeginsWith(name->String().Data())) {
395                                 goodTrig = kTRUE;
396                                 break;
397                         }
398                 }
399         } else goodTrig = kTRUE;
400         
401         if (!goodTrig) {
402                 iClass++;
403                 continue;
404         }
405 //      if (!trigStr.Contains(name->String().Data())) continue;
406         
407         fhOnlineCharge[iClass]->Fill(fEsdV0->GetTriggerChargeA(),fEsdV0->GetTriggerChargeC());
408         
409         fhL2Triggers->Fill(iClass);
410                 
411         fhRecoMult[iClass]->Fill(fEsdV0->GetMTotV0A(),fEsdV0->GetMTotV0C());
412         fhTotRecoMult[iClass]->Fill(fEsdV0->GetMTotV0A()+fEsdV0->GetMTotV0C());
413
414         for(int iTrig = 0; iTrig < 16; ++iTrig){
415                 if(fEsdV0->GetTriggerBits() & (1<<iTrig)) fhTriggerBits[iClass]->Fill(iTrig);   
416         }
417         
418     AliCentrality *centrality = fESD->GetCentrality();
419         Float_t percentile = centrality->GetCentralityPercentile("V0M");
420         Float_t spdPercentile = centrality->GetCentralityPercentile("CL1");
421         //      if (spdPercentile < 0) spdPercentile = 0;
422         //      if (percentile < 0) percentile = 0;
423         fhCentrality[iClass]->Fill(percentile);
424         fhV0vsSPDCentrality[iClass]->Fill(spdPercentile,percentile);
425         
426         Float_t sumEqMult[2] = {0.,0.};
427         for(int iCh = 0; iCh < 64; ++iCh){
428                 if(fEsdV0->GetTime(iCh) < -1024.+ 1e-6) continue;
429                 fhRecoMultPMT[iClass]->Fill(iCh,fEsdV0->GetMultiplicity(iCh));
430                 Int_t side = iCh / 32 ;
431                 sumEqMult[side] += fESD->GetVZEROEqMultiplicity(iCh);
432                 fhEqualizedMult[iClass]->Fill(iCh,fESD->GetVZEROEqMultiplicity(iCh));
433         }
434         fhEqualizedMultSum[iClass]->Fill(sumEqMult[1],sumEqMult[0]);
435         
436         iClass++;
437   }
438         delete tokens;
439
440 }
441 //________________________________________________________________________
442 void AliAnaVZEROPbPb::FillQAHistos(){
443   // Fill the main QA histos
444   // Count V0 flags
445   Int_t nV0A = 0, nV0C = 0;
446   for(Int_t i = 0; i < 32; ++i) {
447     if (fEsdV0->GetBBFlag(i)) nV0C++;
448     if (fEsdV0->GetBBFlag(i+32)) nV0A++;
449   }
450   fNFlags->Fill((Float_t)nV0A,(Float_t)nV0C);
451
452   for (Int_t i=0; i<64; ++i) {
453         Int_t side = i/32;
454         if (fEsdV0->GetTime(i) < 1e-6) {
455                 fhAdcNoTime[side]->Fill(fEsdV0->GetAdc(i));
456                 fhAdcPMTNoTime->Fill(i,fEsdV0->GetAdc(i));
457     } else {
458                 fhAdcWithTime[side]->Fill(fEsdV0->GetAdc(i));
459                 fhAdcPMTWithTime->Fill(i,fEsdV0->GetAdc(i));
460     }
461
462         fhTime[side]->Fill(fEsdV0->GetTime(i));
463         fhWidth[side]->Fill(fEsdV0->GetWidth(i));
464         fhAdcWidth[side]->Fill(fEsdV0->GetAdc(i),fEsdV0->GetWidth(i));
465         fhAdcTime[side]->Fill(fEsdV0->GetTime(i),fEsdV0->GetAdc(i));
466
467     fhTimePMT->Fill(i,fEsdV0->GetTime(i));
468     fhWidthPMT->Fill(i,fEsdV0->GetWidth(i));
469
470   }
471
472   fhTimeCorr->Fill(fEsdV0->GetV0ATime(),fEsdV0->GetV0CTime());
473
474   AliESDVZERO::Decision flaga = fEsdV0->GetV0ADecision();
475   AliESDVZERO::Decision flagc = fEsdV0->GetV0CDecision();
476
477   fhEvents->Fill(flaga,flagc);
478
479   const AliESDVertex *vtx = fESD->GetPrimaryVertexSPD();
480
481   if (flaga <= 1 && flagc <=1) {
482     fhVtxXYBB->Fill(vtx->GetXv(),vtx->GetYv());
483     fhVtxZBB->Fill(vtx->GetZv());
484   }
485   else {
486     if (flaga == 2) {
487       fhVtxXYBGA->Fill(vtx->GetXv(),vtx->GetYv());
488       fhVtxZBGA->Fill(vtx->GetZv());
489     }
490     if (flagc == 2) {
491       fhVtxXYBGC->Fill(vtx->GetXv(),vtx->GetYv());
492       fhVtxZBGC->Fill(vtx->GetZv());
493     }
494   }
495
496   fhPmtMult->Fill(fEsdV0->GetNbPMV0A(),fEsdV0->GetNbPMV0C());
497   for(Int_t i = 0; i < 64; i++) {
498     fhV0ampl->Fill(fEsdV0->GetMultiplicity(i));
499   }
500   
501 }
502 //________________________________________________________________________
503 void AliAnaVZEROPbPb::Terminate(Option_t *) 
504 {
505   // Check if the output list is there
506   // Called once at the end of the query
507
508   fOutputList = dynamic_cast<TList*> (GetOutputData(1));
509   if (!fOutputList) {
510     printf("ERROR: Output list not available\n");
511     return;
512   }
513
514 }
515 //________________________________________________________________________
516 void AliAnaVZEROPbPb::SetOnlineChargeRange(Int_t nbins, Float_t maxA, Float_t maxC){
517   // Set Trigger charge
518   // range used for histogramming
519         fNBinOnlineCharge = nbins;
520         fV0AOnlineChargeMax = maxA;
521         fV0COnlineChargeMax = maxC;
522 }
523 //________________________________________________________________________
524 void AliAnaVZEROPbPb::SetTotalMultiplicityRange(Int_t nbins, Float_t max){
525   // Set Multiplicity 
526   // range used for histogramming
527         fNBinTotMult = nbins;
528         fTotMultMax = max;
529 }
530 //________________________________________________________________________
531 void AliAnaVZEROPbPb::SetMultiplicityRange(Int_t nbins, Float_t maxA, Float_t maxC){
532   // Set Multiplicity 
533   // range used for histogramming
534         fNBinMult = nbins;
535         fV0AMultMax = maxA;
536         fV0CMultMax = maxC;
537 }
538 //________________________________________________________________________
539 void AliAnaVZEROPbPb::SetEquaMultRange(Int_t nbins, Float_t max){
540   // Set Equalized multiplicity
541   // range used for histogramming
542         fNBinEquaMult = nbins;
543         fEquaMultMax = max;
544 }
545 //________________________________________________________________________
546 void AliAnaVZEROPbPb::SetSumEquaMultRange(Int_t nbins, Float_t maxA, Float_t maxC){
547   // Set Equalized multiplicity
548   // range used for histogramming
549         fNBinSumsEqMult = nbins;
550         fV0AEqMultMax = maxA;
551         fV0CEqMultMax = maxC;
552 }
553