]> git.uio.no Git - u/mrichter/AliRoot.git/blob - VZERO/AliVZERODigitizer.cxx
Updated materials in geometry (M. Sitta)
[u/mrichter/AliRoot.git] / VZERO / AliVZERODigitizer.cxx
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 /* $Id$ */
17
18 ///_________________________________________________________________________
19 ///
20 /// This class constructs Digits out of Hits
21 ///
22 ///
23
24 // --- Standard library ---
25
26 // --- ROOT system ---
27 #include <TMath.h>
28 #include <TTree.h>
29 #include <TMap.h>
30 #include <TGeoManager.h>
31 #include <TGeoPhysicalNode.h>
32 #include <AliGeomManager.h>
33 #include <TRandom.h>
34 #include <TF1.h>
35 #include <TH1F.h>
36
37 // --- AliRoot header files ---
38 #include "AliRun.h"
39 #include "AliVZERO.h"
40 #include "AliVZEROhit.h"
41 #include "AliRunLoader.h"
42 #include "AliLoader.h"
43 #include "AliGRPObject.h"
44 #include "AliRunDigitizer.h"
45 #include "AliCDBManager.h"
46 #include "AliCDBStorage.h"
47 #include "AliCDBEntry.h"
48 #include "AliVZEROCalibData.h"
49 #include "AliCTPTimeParams.h"
50 #include "AliLHCClockPhase.h"
51 #include "AliVZEROdigit.h"
52 #include "AliVZERODigitizer.h"
53 #include "AliVZEROSDigit.h"
54
55 ClassImp(AliVZERODigitizer)
56
57  AliVZERODigitizer::AliVZERODigitizer()
58                    :AliDigitizer(),
59                     fCalibData(GetCalibData()),
60                     fPhotoCathodeEfficiency(0.18),
61                     fNdigits(0),
62                     fDigits(0),
63                     fSignalShape(NULL),
64                     fPMResponse(NULL),
65                     fSinglePhESpectrum(NULL),
66                     fEvenOrOdd(kFALSE),
67                     fTask(kHits2Digits),
68                     fVZERO(NULL)
69 {
70   // default constructor
71   // Initialize OCDB and containers used in the digitization
72
73   Init();
74 }
75
76 //____________________________________________________________________________ 
77   AliVZERODigitizer::AliVZERODigitizer(AliVZERO *vzero, DigiTask_t task)
78                     :AliDigitizer(),
79                      fCalibData(GetCalibData()),
80                      fPhotoCathodeEfficiency(0.18),
81                      fNdigits(0),
82                      fDigits(0),
83                      fSignalShape(NULL),
84                      fPMResponse(NULL),
85                      fSinglePhESpectrum(NULL),
86                      fEvenOrOdd(kFALSE),
87                      fTask(task),
88                      fVZERO(vzero)
89 {
90   // constructor
91   // Initialize OCDB and containers used in the digitization
92
93   Init();
94 }
95            
96 //____________________________________________________________________________ 
97   AliVZERODigitizer::AliVZERODigitizer(AliRunDigitizer* manager)
98                     :AliDigitizer(manager),
99                      fCalibData(GetCalibData()),
100                      fPhotoCathodeEfficiency(0.18),
101                      fNdigits(0),
102                      fDigits(0),
103                      fSignalShape(NULL),
104                      fPMResponse(NULL),
105                      fSinglePhESpectrum(NULL),
106                      fEvenOrOdd(kFALSE),
107                      fTask(kHits2Digits),
108                      fVZERO(NULL)
109 {
110   // constructor
111   // Initialize OCDB and containers used in the digitization
112
113   Init();
114 }
115            
116 //____________________________________________________________________________ 
117   AliVZERODigitizer::~AliVZERODigitizer()
118 {
119   // destructor
120   
121   if (fDigits) {
122     fDigits->Delete();
123     delete fDigits;
124     fDigits=0; 
125   }
126
127   if (fSignalShape) {
128     delete fSignalShape;
129     fSignalShape = NULL;
130   }
131   if (fPMResponse) {
132     delete fPMResponse;
133     fPMResponse = NULL;
134   }
135   if (fSinglePhESpectrum) {
136     delete fSinglePhESpectrum;
137     fSinglePhESpectrum = NULL;
138   }
139
140   for(Int_t i = 0 ; i < 64; ++i) {
141     if (fTime[i]) delete [] fTime[i];
142   }
143 }
144
145 //_____________________________________________________________________________
146 Bool_t AliVZERODigitizer::Init()
147 {
148   // Initialises the digitizer
149   // Initialize OCDB and containers used in the digitization
150
151   // check if the digitizer was already initialized
152   if (fSignalShape) return kTRUE;
153
154   fSignalShape = new TF1("VZEROSignalShape",this,&AliVZERODigitizer::SignalShape,0,200,6,"AliVZERODigitizer","SignalShape");
155   //  fSignalShape->SetParameters(0,1.57345e1,-4.25603e-1,2.9,6.40982,3.69339e-01);
156   //  fSignalShape->SetParameters(1.34330e+00,1.13007e+02,-4.95705e-01,
157   //                          3.68911e+00,1.01040e+00, 3.94675e-01);
158   fSignalShape->SetParameters(-1.07335e+00,2.16002e+01,-1.26133e-01,
159                               1.41619e+00,5.50334e-01,3.86111e-01);
160   fPMResponse = new TF1("VZEROPMResponse",this,&AliVZERODigitizer::PMResponse,-kPMRespTime,2.*kPMRespTime,0,"AliVZERODigitizer","PMResponse");
161   fSinglePhESpectrum = new TF1("VZEROSinglePhESpectrum",this,&AliVZERODigitizer::SinglePhESpectrum,0,20,0,"AliVZERODigitizer","SinglePhESpectrum");
162   
163   // Now get the CTP L0->L1 delay
164   AliCDBEntry *entry = AliCDBManager::Instance()->Get("GRP/CTP/CTPtiming");
165   if (!entry) AliFatal("CTP timing parameters are not found in OCDB !");
166   AliCTPTimeParams *ctpParams = (AliCTPTimeParams*)entry->GetObject();
167   Float_t l1Delay = (Float_t)ctpParams->GetDelayL1L0()*25.0;
168
169   AliCDBEntry *entry1 = AliCDBManager::Instance()->Get("GRP/CTP/TimeAlign");
170   if (!entry1) AliFatal("CTP time-alignment is not found in OCDB !");
171   AliCTPTimeParams *ctpTimeAlign = (AliCTPTimeParams*)entry1->GetObject();
172   l1Delay += ((Float_t)ctpTimeAlign->GetDelayL1L0()*25.0);
173
174   AliCDBEntry *entry2 = AliCDBManager::Instance()->Get("VZERO/Calib/TimeDelays");
175   if (!entry2) AliFatal("VZERO time delays are not found in OCDB !");
176   TH1F *delays = (TH1F*)entry2->GetObject();
177
178   AliCDBEntry *entry3 = AliCDBManager::Instance()->Get("GRP/Calib/LHCClockPhase");
179   if (!entry3) AliFatal("LHC clock-phase shift is not found in OCDB !");
180   AliLHCClockPhase *phase = (AliLHCClockPhase*)entry3->GetObject();
181
182   for(Int_t i = 0 ; i < 64; ++i) {
183
184     for(Int_t j = 0; j < kNClocks; ++j) fAdc[i][j] = 0;
185     fLeadingTime[i] = fTimeWidth[i] = 0;
186
187     fPmGain[i] = fCalibData->GetGain(i);
188
189     fAdcPedestal[i][0] = fCalibData->GetPedestal(i);
190     fAdcSigma[i][0]    = fCalibData->GetSigma(i); 
191     fAdcPedestal[i][1] = fCalibData->GetPedestal(i+64);
192     fAdcSigma[i][1]    = fCalibData->GetSigma(i+64); 
193
194     Int_t board = AliVZEROCalibData::GetBoardNumber(i);
195     fNBins[i] = TMath::Nint(((Float_t)(fCalibData->GetMatchWindow(board)+1)*25.0+
196                              (Float_t)kMaxTDCWidth*fCalibData->GetWidthResolution(board))/
197                             fCalibData->GetTimeResolution(board));
198     fNBinsLT[i] = TMath::Nint(((Float_t)(fCalibData->GetMatchWindow(board)+1)*25.0)/
199                               fCalibData->GetTimeResolution(board));
200     fBinSize[i] = fCalibData->GetTimeResolution(board);
201     fHptdcOffset[i] = (((Float_t)fCalibData->GetTriggerCountOffset(board)-
202                         (Float_t)fCalibData->GetRollOver(board))*25.0+
203                        fCalibData->GetTimeOffset(i)-
204                        l1Delay-
205                        phase->GetMeanPhase()+
206                        delays->GetBinContent(i+1)+
207                        kV0Offset);
208
209     fTime[i] = new Float_t[fNBins[i]];
210     memset(fTime[i],0,fNBins[i]*sizeof(Float_t));
211   }
212
213   return kTRUE;
214 }
215
216 //____________________________________________________________________________
217 void AliVZERODigitizer::Exec(Option_t* /*option*/) 
218 {   
219   // Creates digits from hits
220   fNdigits = 0;  
221
222   if (fVZERO && !fManager) {
223     AliLoader *loader = fVZERO->GetLoader();
224     if (!loader) {
225       AliError("Can not get VZERO Loader via AliVZERO object!");
226       return;
227     }
228     AliRunLoader* runLoader = AliRunLoader::Instance();
229     for (Int_t iEvent = 0; iEvent < runLoader->GetNumberOfEvents(); ++iEvent) {
230       runLoader->GetEvent(iEvent);
231       if (fTask == kHits2Digits) {
232         DigitizeHits();
233         DigitizeSDigits();
234         WriteDigits(loader);
235       }
236       else {
237         DigitizeHits();
238         WriteSDigits(loader);
239       }
240     }
241   }
242   else if (fManager) {
243       ReadSDigits();
244       DigitizeSDigits();
245       AliRunLoader *currentLoader = AliRunLoader::GetRunLoader(fManager->GetOutputFolderName());
246       AliLoader *loader = currentLoader->GetLoader("VZEROLoader");
247       if (!loader) { 
248         AliError("Cannot get VZERO Loader via RunDigitizer!");
249         return;
250       }
251       WriteDigits(loader);
252   }
253   else {
254     AliFatal("Invalid digitization task! Exiting!");
255   }
256 }
257
258 //____________________________________________________________________________
259 void AliVZERODigitizer::AddDigit(Int_t pmnumber, Float_t time, Float_t width, Bool_t integrator, Short_t *chargeADC, Int_t *labels) 
260  { 
261  
262 // Adds Digit 
263  
264   TClonesArray &ldigits = *fDigits;  
265          
266   new(ldigits[fNdigits++]) AliVZEROdigit(pmnumber,time,width,integrator,chargeADC,labels);
267          
268 }
269 //____________________________________________________________________________
270 void AliVZERODigitizer::AddSDigit(Int_t pmnumber, Int_t nbins, Float_t *charges, Int_t *labels) 
271  { 
272  
273 // Adds SDigit 
274  
275   TClonesArray &ldigits = *fDigits;  
276          
277   new(ldigits[fNdigits++]) AliVZEROSDigit(pmnumber,nbins,charges,labels);
278          
279 }
280 //____________________________________________________________________________
281 void AliVZERODigitizer::ResetDigits()
282 {
283
284 // Clears Digits
285
286   fNdigits = 0;
287   if (fDigits) fDigits->Clear();
288 }
289
290 //____________________________________________________________________________
291 AliVZEROCalibData* AliVZERODigitizer::GetCalibData() const
292
293 {
294   AliCDBManager *man = AliCDBManager::Instance();
295
296   AliCDBEntry *entry=0;
297
298   entry = man->Get("VZERO/Calib/Data");
299
300 //   if(!entry){
301 //     AliWarning("Load of calibration data from default storage failed!");
302 //     AliWarning("Calibration data will be loaded from local storage ($ALICE_ROOT)");
303 //     Int_t runNumber = man->GetRun();
304 //     entry = man->GetStorage("local://$ALICE_ROOT/OCDB")
305 //       ->Get("VZERO/Calib/Data",runNumber);
306 //      
307 //   }
308
309   // Retrieval of data in directory VZERO/Calib/Data:
310
311
312   AliVZEROCalibData *calibdata = 0;
313
314   if (entry) calibdata = (AliVZEROCalibData*) entry->GetObject();
315   if (!calibdata)  AliFatal("No calibration data from calibration database !");
316
317   return calibdata;
318
319 }
320
321 double AliVZERODigitizer::SignalShape(double *x, double *par)
322 {
323   // this function simulates the time
324   // of arrival of the photons at the
325   // photocathode
326   Double_t xx = x[0];
327   if (xx <= par[0]) return 0;
328   Double_t a = 1./TMath::Power((xx-par[0])/par[1],1./par[2]);
329   if (xx <= par[3]) return a;
330   Double_t b = 1./TMath::Power((xx-par[3])/par[4],1./par[5]);
331   Double_t f = a*b/(a+b);
332   AliDebug(100,Form("x=%f func=%f",xx,f));
333   return f;
334 }
335
336 double AliVZERODigitizer::PMResponse(double *x, double * /* par */)
337 {
338   // this function describes the
339   // PM time response to a single
340   // photoelectron
341   Double_t xx = x[0]+kPMRespTime;
342   return xx*xx*TMath::Exp(-xx*xx/(kPMRespTime*kPMRespTime));
343 }
344
345 double AliVZERODigitizer::SinglePhESpectrum(double *x, double * /* par */)
346 {
347   // this function describes the
348   // PM amplitude response to a single
349   // photoelectron
350   Double_t xx = x[0];
351   if (xx < 0) return 0;
352   return (TMath::Poisson(xx,kPMNbOfSecElec)+kPMTransparency*TMath::Poisson(xx,1.0));
353 }
354
355 Int_t AliVZERODigitizer::Cell2Pmt(Int_t cell) const
356 {
357   // The method maps the scintillator
358   // indexes to the PM ones
359   if (cell < 0 || cell >= 80) {
360     AliError(Form("Wrong VZERO cell index %d",cell));
361     return -1;
362   }
363   if (cell < 16) return cell;
364   if (cell < 48) return 8 + cell/2;
365   return cell - 16;
366 }
367
368 void AliVZERODigitizer::DigitizeHits()
369 {
370   // Digitize the hits to the level of
371   // SDigits (fTime arrays)
372
373   for(Int_t i = 0 ; i < 64; ++i) {
374     memset(fTime[i],0,fNBins[i]*sizeof(Float_t));
375     fLabels[i][0] = fLabels[i][1] = fLabels[i][2] = -1;
376   }
377   Float_t integral = fPMResponse->Integral(-kPMRespTime,2.*kPMRespTime);
378   Float_t meansPhE = fSinglePhESpectrum->Mean(0,20);
379
380      AliLoader* loader = fVZERO->GetLoader();
381      if (!loader) {
382        AliError("Can not get VZERO Loader!");
383        return;
384      }
385      loader->LoadHits();
386      TTree* treeH = loader->TreeH();
387      if (!treeH) {
388        AliError("Cannot get TreeH!");
389        return;
390      }
391      TClonesArray* hits = fVZERO->Hits();
392
393 //  Now makes Digits from hits
394      Int_t nTracks = (Int_t) treeH->GetEntries();
395      for (Int_t iTrack = 0; iTrack < nTracks; iTrack++) {
396          fVZERO->ResetHits();
397          treeH->GetEvent(iTrack);
398          Int_t nHits = hits->GetEntriesFast();
399          for (Int_t iHit = 0; iHit < nHits; iHit++) {
400            AliVZEROhit* hit = (AliVZEROhit *)hits->UncheckedAt(iHit);
401            Int_t nPhot = hit->Nphot();
402            Int_t cell  = hit->Cell();                          
403            Int_t pmt = Cell2Pmt(cell);
404            if (pmt < 0) continue;
405            Int_t trackLabel = hit->GetTrack();
406            for(Int_t l = 0; l < 3; ++l) {
407              if (fLabels[pmt][l] < 0) {
408                fLabels[pmt][l] = trackLabel;
409                break;
410              }
411            }
412            Float_t dt_scintillator = gRandom->Gaus(0,kIntTimeRes);
413            Float_t t = dt_scintillator + 1e9*hit->Tof();
414            if (pmt < 32) t += kV0CDelayCables;
415            t += fHptdcOffset[pmt];
416            Int_t nPhE;
417            Float_t prob = fCalibData->GetLightYields(pmt)*fPhotoCathodeEfficiency; // Optical losses included!
418            if (nPhot > 100)
419              nPhE = (Int_t)gRandom->Gaus(prob*Float_t(nPhot)+0.5,
420                                          sqrt(Float_t(nPhot)*prob*(1.-prob)));
421            else
422              nPhE = gRandom->Binomial(nPhot,prob);
423            Float_t charge = TMath::Qe()*fPmGain[pmt]*fBinSize[pmt]/integral;
424            for (Int_t iPhE = 0; iPhE < nPhE; ++iPhE) {
425              Float_t tPhE = t + fSignalShape->GetRandom(0,fBinSize[pmt]*Float_t(fNBins[pmt]));
426              Float_t gainVar = fSinglePhESpectrum->GetRandom(0,20)/meansPhE;
427              Int_t firstBin = TMath::Max(0,(Int_t)((tPhE-kPMRespTime)/fBinSize[pmt]));
428              Int_t lastBin = TMath::Min(fNBins[pmt]-1,(Int_t)((tPhE+2.*kPMRespTime)/fBinSize[pmt]));
429              for(Int_t iBin = firstBin; iBin <= lastBin; ++iBin) {
430                Float_t tempT = fBinSize[pmt]*(0.5+iBin)-tPhE;
431                fTime[pmt][iBin] += gainVar*charge*fPMResponse->Eval(tempT);
432              }
433            }         // ph.e. loop
434          }           // hit loop
435      }               // track loop
436      loader->UnloadHits();
437 }
438
439
440 void AliVZERODigitizer::DigitizeSDigits()
441 {
442   // Digitize the fTime arrays (SDigits) to the level of
443   // Digits (fAdc arrays)
444   for(Int_t i = 0 ; i < 64; ++i) {
445     for(Int_t j = 0; j < kNClocks; ++j) fAdc[i][j] = 0;
446     fLeadingTime[i] = fTimeWidth[i] = 0;
447   }
448
449   Float_t maximum = 0.9*fSignalShape->GetMaximum(0,200); // Not exact, one needs to do this on the convoluted
450   Float_t integral2 = fSignalShape->Integral(0,200); // function. Anyway the effect is small <10% on the 2.5 ADC thr
451   for (Int_t ipmt = 0; ipmt < 64; ++ipmt) {
452     Float_t thr = fCalibData->GetDiscriThr(ipmt)*kChargePerADC*maximum*fBinSize[ipmt]/integral2;
453     Bool_t ltFound = kFALSE, ttFound = kFALSE;
454     for (Int_t iBin = 0; iBin < fNBins[ipmt]; ++iBin) {
455       Float_t t = fBinSize[ipmt]*Float_t(iBin);
456       if (fTime[ipmt][iBin] > thr) {
457         if (!ltFound && (iBin < fNBinsLT[ipmt])) {
458           ltFound = kTRUE;
459           fLeadingTime[ipmt] = t;
460         }
461       }
462       else {
463         if (ltFound) {
464           if (!ttFound) {
465             ttFound = kTRUE;
466             fTimeWidth[ipmt] = t - fLeadingTime[ipmt];
467           }
468         }
469       }
470       Float_t tadc = t - kClockOffset - fCalibData->GetTimeOffset(ipmt);
471       Int_t clock = kNClocks/2 - Int_t(tadc/25.0);
472       if (clock >= 0 && clock < kNClocks)
473         fAdc[ipmt][clock] += fTime[ipmt][iBin]/kChargePerADC;
474     }
475     Int_t board = AliVZEROCalibData::GetBoardNumber(ipmt);
476     if (ltFound && ttFound) {
477       fTimeWidth[ipmt] = fCalibData->GetWidthResolution(board)*
478         Float_t(Int_t(fTimeWidth[ipmt]/fCalibData->GetWidthResolution(board)));
479       if (fTimeWidth[ipmt] < Float_t(kMinTDCWidth)*fCalibData->GetWidthResolution(board))
480         fTimeWidth[ipmt] = Float_t(kMinTDCWidth)*fCalibData->GetWidthResolution(board);
481       if (fTimeWidth[ipmt] > Float_t(kMaxTDCWidth)*fCalibData->GetWidthResolution(board))
482         fTimeWidth[ipmt] = Float_t(kMaxTDCWidth)*fCalibData->GetWidthResolution(board);
483     }
484   }
485
486   fEvenOrOdd = gRandom->Integer(2);
487   for (Int_t j=0; j<64; ++j){
488     for (Int_t iClock = 0; iClock < kNClocks; ++iClock) {
489       Int_t integrator = (iClock + fEvenOrOdd) % 2;
490       fAdc[j][iClock]  += gRandom->Gaus(fAdcPedestal[j][integrator], fAdcSigma[j][integrator]);
491     }
492   }
493         
494 }
495
496 void AliVZERODigitizer::WriteDigits(AliLoader *loader)
497 {
498   // Take fAdc arrays filled by the previous
499   // method and produce and add digits to the digit Tree
500
501   loader->LoadDigits("UPDATE");
502
503   if (!loader->TreeD()) loader->MakeTree("D");
504   loader->MakeDigitsContainer();
505   TTree* treeD  = loader->TreeD();
506   DigitsArray();
507   treeD->Branch("VZERODigit", &fDigits); 
508   
509   Short_t *chargeADC = new Short_t[kNClocks];
510   for (Int_t i=0; i<64; i++) {      
511     for (Int_t j = 0; j < kNClocks; ++j) {
512       Int_t tempadc = Int_t(fAdc[i][j]);
513       if (tempadc > 1023) tempadc = 1023;
514       chargeADC[j] = tempadc;
515     }
516     AddDigit(i, fLeadingTime[i], fTimeWidth[i], Bool_t((10+fEvenOrOdd)%2), chargeADC, fLabels[i]);
517   }
518   delete [] chargeADC;
519
520   treeD->Fill();
521   loader->WriteDigits("OVERWRITE");  
522   loader->UnloadDigits();     
523   ResetDigits();
524 }
525
526 void AliVZERODigitizer::WriteSDigits(AliLoader *loader)
527 {
528   // Take fTime arrays filled by the previous
529   // method and produce and add sdigits to the sdigit Tree
530
531   loader->LoadSDigits("UPDATE");
532
533   if (!loader->TreeS()) loader->MakeTree("S");
534   loader->MakeSDigitsContainer();
535   TTree* treeS  = loader->TreeS();
536   SDigitsArray();
537   treeS->Branch("VZEROSDigit", &fDigits); 
538   
539   for (Int_t ipmt = 0; ipmt < 64; ++ipmt) {
540     AddSDigit(ipmt,fNBins[ipmt],fTime[ipmt],fLabels[ipmt]);
541   }
542
543   treeS->Fill();
544   loader->WriteSDigits("OVERWRITE");  
545   loader->UnloadSDigits();     
546   ResetDigits();
547 }
548
549 void AliVZERODigitizer::ReadSDigits()
550 {
551   // Read SDigits which are then to precessed
552   // in the following method
553   for(Int_t i = 0 ; i < 64; ++i) {
554     memset(fTime[i],0,fNBins[i]*sizeof(Float_t));
555     fLabels[i][0] = fLabels[i][1] = fLabels[i][2] = -1;
556   }
557
558   // Loop over input files
559   Int_t nFiles= fManager->GetNinputs();
560   for (Int_t inputFile = 0; inputFile < nFiles; inputFile++) {
561     // Get the current loader 
562     AliRunLoader* currentLoader = 
563       AliRunLoader::GetRunLoader(fManager->GetInputFolderName(inputFile));
564
565     AliLoader *loader = currentLoader->GetLoader("VZEROLoader");
566     loader->LoadSDigits("READ");
567   
568     // Get the tree of summable digits
569     TTree* sdigitsTree = loader->TreeS();
570     if (!sdigitsTree)  {
571       AliError("No sdigit tree from manager");
572       continue;
573     }
574
575     // Get the branch 
576     TBranch* sdigitsBranch = sdigitsTree->GetBranch("VZEROSDigit");
577     if (!sdigitsBranch) {
578       AliError("Failed to get sdigit branch");
579       return;
580     }
581
582     // Set the branch address
583     TClonesArray *sdigitsArray = NULL;
584     sdigitsBranch->SetAddress(&sdigitsArray);
585
586     // Sum contributions from the sdigits
587     // Get number of entries in the tree 
588     Int_t nentries  = Int_t(sdigitsBranch->GetEntries());
589     for (Int_t entry = 0; entry < nentries; ++entry)  {
590       sdigitsBranch->GetEntry(entry);
591       // Get the number of sdigits 
592       Int_t nsdigits = sdigitsArray->GetEntries();
593       for (Int_t sdigit = 0; sdigit < nsdigits; sdigit++) {
594         AliVZEROSDigit* sDigit = static_cast<AliVZEROSDigit*>(sdigitsArray->UncheckedAt(sdigit));
595         Int_t pmNumber = sDigit->PMNumber();
596         Int_t nbins = sDigit->GetNBins();
597         if (nbins != fNBins[pmNumber]) {
598           AliError(Form("Incompatible number of bins between digitizer (%d) and sdigit (%d) for PM %d! Skipping sdigit!",
599                         fNBins[pmNumber],nbins,pmNumber));
600           continue;
601         }
602         // Sum the charges
603         Float_t *charges = sDigit->GetCharges();
604         for(Int_t iBin = 0; iBin < nbins; ++iBin) fTime[pmNumber][iBin] += charges[iBin];
605         // and the labels
606         Int_t *labels = sDigit->GetTracks();
607         Int_t j = 0;
608         for(Int_t i = 0; i < 3; ++i) {
609           if (fLabels[pmNumber][i] < 0) {
610             if (labels[j] < 0) break;
611             fLabels[pmNumber][i] = labels[j];
612             j++;
613           }
614         }
615       }
616     }
617     loader->UnloadSDigits();
618   }
619 }
620
621 //____________________________________________________________________
622 TClonesArray*
623 AliVZERODigitizer::DigitsArray() 
624 {
625   // Initialize digit array if not already and
626   // return pointer to it. 
627   if (!fDigits) { 
628     fDigits = new TClonesArray("AliVZEROdigit", 64);
629     fNdigits = 0;
630   }
631   return fDigits;
632 }
633
634 //____________________________________________________________________
635 TClonesArray*
636 AliVZERODigitizer::SDigitsArray() 
637 {
638   // Initialize sdigit array if not already and
639   // return pointer to it. 
640   if (!fDigits) { 
641     fDigits = new TClonesArray("AliVZEROSDigit", 64);
642     fNdigits = 0;
643   }
644   return fDigits;
645 }