]> git.uio.no Git - u/mrichter/AliRoot.git/blob - VZERO/AliVZERODigitizer.cxx
Fix for coverity (AdC)
[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 "AliDigitizationInput.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(AliDigitizationInput* digInput)
98                     :AliDigitizer(digInput),
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->GetRollOver(board)-
202                         (Float_t)fCalibData->GetTriggerCountOffset(board))*25.0+
203                        fCalibData->GetTimeOffset(i)-
204                        l1Delay-
205                        phase->GetMeanPhase()+
206                        delays->GetBinContent(i+1)+
207                        kV0Offset);
208     fClockOffset[i] = (((Float_t)fCalibData->GetRollOver(board)-
209                         (Float_t)fCalibData->GetTriggerCountOffset(board))*25.0+
210                        fCalibData->GetTimeOffset(i)-
211                        l1Delay+
212                        kV0Offset);
213
214     fTime[i] = new Float_t[fNBins[i]];
215     memset(fTime[i],0,fNBins[i]*sizeof(Float_t));
216   }
217
218   return kTRUE;
219 }
220
221 //____________________________________________________________________________
222 void AliVZERODigitizer::Digitize(Option_t* /*option*/) 
223 {   
224   // Creates digits from hits
225   fNdigits = 0;  
226
227   if (fVZERO && !fDigInput) {
228     AliLoader *loader = fVZERO->GetLoader();
229     if (!loader) {
230       AliError("Can not get VZERO Loader via AliVZERO object!");
231       return;
232     }
233     AliRunLoader* runLoader = AliRunLoader::Instance();
234     for (Int_t iEvent = 0; iEvent < runLoader->GetNumberOfEvents(); ++iEvent) {
235       runLoader->GetEvent(iEvent);
236       if (fTask == kHits2Digits) {
237         DigitizeHits();
238         DigitizeSDigits();
239         WriteDigits(loader);
240       }
241       else {
242         DigitizeHits();
243         WriteSDigits(loader);
244       }
245     }
246   }
247   else if (fDigInput) {
248       ReadSDigits();
249       DigitizeSDigits();
250       AliRunLoader *currentLoader = AliRunLoader::GetRunLoader(fDigInput->GetOutputFolderName());
251       AliLoader *loader = currentLoader->GetLoader("VZEROLoader");
252       if (!loader) { 
253         AliError("Cannot get VZERO Loader via RunDigitizer!");
254         return;
255       }
256       WriteDigits(loader);
257   }
258   else {
259     AliFatal("Invalid digitization task! Exiting!");
260   }
261 }
262
263 //____________________________________________________________________________
264 void AliVZERODigitizer::AddDigit(Int_t pmnumber, Float_t time, Float_t width, Bool_t integrator, Short_t *chargeADC, Int_t *labels) 
265  { 
266  
267 // Adds Digit 
268  
269   TClonesArray &ldigits = *fDigits;  
270          
271   new(ldigits[fNdigits++]) AliVZEROdigit(pmnumber,time,width,integrator,chargeADC,labels);
272          
273 }
274 //____________________________________________________________________________
275 void AliVZERODigitizer::AddSDigit(Int_t pmnumber, Int_t nbins, Float_t *charges, Int_t *labels) 
276  { 
277  
278 // Adds SDigit 
279  
280   TClonesArray &ldigits = *fDigits;  
281          
282   new(ldigits[fNdigits++]) AliVZEROSDigit(pmnumber,nbins,charges,labels);
283          
284 }
285 //____________________________________________________________________________
286 void AliVZERODigitizer::ResetDigits()
287 {
288
289 // Clears Digits
290
291   fNdigits = 0;
292   if (fDigits) fDigits->Clear();
293 }
294
295 //____________________________________________________________________________
296 AliVZEROCalibData* AliVZERODigitizer::GetCalibData() const
297
298 {
299   AliCDBManager *man = AliCDBManager::Instance();
300
301   AliCDBEntry *entry=0;
302
303   entry = man->Get("VZERO/Calib/Data");
304
305 //   if(!entry){
306 //     AliWarning("Load of calibration data from default storage failed!");
307 //     AliWarning("Calibration data will be loaded from local storage ($ALICE_ROOT)");
308 //     Int_t runNumber = man->GetRun();
309 //     entry = man->GetStorage("local://$ALICE_ROOT/OCDB")
310 //       ->Get("VZERO/Calib/Data",runNumber);
311 //      
312 //   }
313
314   // Retrieval of data in directory VZERO/Calib/Data:
315
316
317   AliVZEROCalibData *calibdata = 0;
318
319   if (entry) calibdata = (AliVZEROCalibData*) entry->GetObject();
320   if (!calibdata)  AliFatal("No calibration data from calibration database !");
321
322   return calibdata;
323
324 }
325
326 double AliVZERODigitizer::SignalShape(double *x, double *par)
327 {
328   // this function simulates the time
329   // of arrival of the photons at the
330   // photocathode
331   Double_t xx = x[0];
332   if (xx <= par[0]) return 0;
333   Double_t a = 1./TMath::Power((xx-par[0])/par[1],1./par[2]);
334   if (xx <= par[3]) return a;
335   Double_t b = 1./TMath::Power((xx-par[3])/par[4],1./par[5]);
336   Double_t f = a*b/(a+b);
337   AliDebug(100,Form("x=%f func=%f",xx,f));
338   return f;
339 }
340
341 double AliVZERODigitizer::PMResponse(double *x, double * /* par */)
342 {
343   // this function describes the
344   // PM time response to a single
345   // photoelectron
346   Double_t xx = x[0]+kPMRespTime;
347   return xx*xx*TMath::Exp(-xx*xx/(kPMRespTime*kPMRespTime));
348 }
349
350 double AliVZERODigitizer::SinglePhESpectrum(double *x, double * /* par */)
351 {
352   // this function describes the
353   // PM amplitude response to a single
354   // photoelectron
355   Double_t xx = x[0];
356   if (xx < 0) return 0;
357   return (TMath::Poisson(xx,kPMNbOfSecElec)+kPMTransparency*TMath::Poisson(xx,1.0));
358 }
359
360 Int_t AliVZERODigitizer::Cell2Pmt(Int_t cell) const
361 {
362   // The method maps the scintillator
363   // indexes to the PM ones
364   if (cell < 0 || cell >= 80) {
365     AliError(Form("Wrong VZERO cell index %d",cell));
366     return -1;
367   }
368   if (cell < 16) return cell;
369   if (cell < 48) return 8 + cell/2;
370   return cell - 16;
371 }
372
373 void AliVZERODigitizer::DigitizeHits()
374 {
375   // Digitize the hits to the level of
376   // SDigits (fTime arrays)
377
378   for(Int_t i = 0 ; i < 64; ++i) {
379     memset(fTime[i],0,fNBins[i]*sizeof(Float_t));
380     fLabels[i][0] = fLabels[i][1] = fLabels[i][2] = -1;
381   }
382   Float_t integral = fPMResponse->Integral(-kPMRespTime,2.*kPMRespTime);
383   Float_t meansPhE = fSinglePhESpectrum->Mean(0,20);
384
385      AliLoader* loader = fVZERO->GetLoader();
386      if (!loader) {
387        AliError("Can not get VZERO Loader!");
388        return;
389      }
390      loader->LoadHits();
391      TTree* treeH = loader->TreeH();
392      if (!treeH) {
393        AliError("Cannot get TreeH!");
394        return;
395      }
396      TClonesArray* hits = fVZERO->Hits();
397
398 //  Now makes Digits from hits
399      Int_t nTracks = (Int_t) treeH->GetEntries();
400      for (Int_t iTrack = 0; iTrack < nTracks; iTrack++) {
401          fVZERO->ResetHits();
402          treeH->GetEvent(iTrack);
403          Int_t nHits = hits->GetEntriesFast();
404          for (Int_t iHit = 0; iHit < nHits; iHit++) {
405            AliVZEROhit* hit = (AliVZEROhit *)hits->UncheckedAt(iHit);
406            Int_t nPhot = hit->Nphot();
407            Int_t cell  = hit->Cell();                          
408            Int_t pmt = Cell2Pmt(cell);
409            if (pmt < 0) continue;
410            Int_t trackLabel = hit->GetTrack();
411            for(Int_t l = 0; l < 3; ++l) {
412              if (fLabels[pmt][l] < 0) {
413                fLabels[pmt][l] = trackLabel;
414                break;
415              }
416            }
417            Float_t dt_scintillator = gRandom->Gaus(0,kIntTimeRes);
418            Float_t t = dt_scintillator + 1e9*hit->Tof();
419            if (pmt < 32) t += kV0CDelayCables;
420            t += fHptdcOffset[pmt];
421            Int_t nPhE;
422            Float_t prob = fCalibData->GetLightYields(pmt)*fPhotoCathodeEfficiency; // Optical losses included!
423            if (nPhot > 100)
424              nPhE = (Int_t)gRandom->Gaus(prob*Float_t(nPhot)+0.5,
425                                          sqrt(Float_t(nPhot)*prob*(1.-prob)));
426            else
427              nPhE = gRandom->Binomial(nPhot,prob);
428            Float_t charge = TMath::Qe()*fPmGain[pmt]*fBinSize[pmt]/integral;
429            for (Int_t iPhE = 0; iPhE < nPhE; ++iPhE) {
430              Float_t tPhE = t + fSignalShape->GetRandom(0,fBinSize[pmt]*Float_t(fNBins[pmt]));
431              Float_t gainVar = fSinglePhESpectrum->GetRandom(0,20)/meansPhE;
432              Int_t firstBin = TMath::Max(0,(Int_t)((tPhE-kPMRespTime)/fBinSize[pmt]));
433              Int_t lastBin = TMath::Min(fNBins[pmt]-1,(Int_t)((tPhE+2.*kPMRespTime)/fBinSize[pmt]));
434              for(Int_t iBin = firstBin; iBin <= lastBin; ++iBin) {
435                Float_t tempT = fBinSize[pmt]*(0.5+iBin)-tPhE;
436                fTime[pmt][iBin] += gainVar*charge*fPMResponse->Eval(tempT);
437              }
438            }         // ph.e. loop
439          }           // hit loop
440      }               // track loop
441      loader->UnloadHits();
442 }
443
444
445 void AliVZERODigitizer::DigitizeSDigits()
446 {
447   // Digitize the fTime arrays (SDigits) to the level of
448   // Digits (fAdc arrays)
449   for(Int_t i = 0 ; i < 64; ++i) {
450     for(Int_t j = 0; j < kNClocks; ++j) fAdc[i][j] = 0;
451     fLeadingTime[i] = fTimeWidth[i] = 0;
452   }
453
454   Float_t maximum = 0.9*fSignalShape->GetMaximum(0,200); // Not exact, one needs to do this on the convoluted
455   Float_t integral2 = fSignalShape->Integral(0,200); // function. Anyway the effect is small <10% on the 2.5 ADC thr
456   for (Int_t ipmt = 0; ipmt < 64; ++ipmt) {
457     Float_t thr = fCalibData->GetCalibDiscriThr(ipmt,kFALSE)*kChargePerADC*maximum*fBinSize[ipmt]/integral2;
458     Bool_t ltFound = kFALSE, ttFound = kFALSE;
459     for (Int_t iBin = 0; iBin < fNBins[ipmt]; ++iBin) {
460       Float_t t = fBinSize[ipmt]*Float_t(iBin);
461       if (fTime[ipmt][iBin] > thr) {
462         if (!ltFound && (iBin < fNBinsLT[ipmt])) {
463           ltFound = kTRUE;
464           fLeadingTime[ipmt] = t;
465         }
466       }
467       else {
468         if (ltFound) {
469           if (!ttFound) {
470             ttFound = kTRUE;
471             fTimeWidth[ipmt] = t - fLeadingTime[ipmt];
472           }
473         }
474       }
475       Float_t tadc = t - fClockOffset[ipmt];
476       Int_t clock = kNClocks/2 - Int_t(tadc/25.0);
477       if (clock >= 0 && clock < kNClocks)
478         fAdc[ipmt][clock] += fTime[ipmt][iBin]/kChargePerADC;
479     }
480     AliDebug(1,Form("Channel %d Offset %f Time %f",ipmt,fClockOffset[ipmt],fLeadingTime[ipmt]));
481     Int_t board = AliVZEROCalibData::GetBoardNumber(ipmt);
482     if (ltFound && ttFound) {
483       fTimeWidth[ipmt] = fCalibData->GetWidthResolution(board)*
484         Float_t(Int_t(fTimeWidth[ipmt]/fCalibData->GetWidthResolution(board)));
485       if (fTimeWidth[ipmt] < Float_t(kMinTDCWidth)*fCalibData->GetWidthResolution(board))
486         fTimeWidth[ipmt] = Float_t(kMinTDCWidth)*fCalibData->GetWidthResolution(board);
487       if (fTimeWidth[ipmt] > Float_t(kMaxTDCWidth)*fCalibData->GetWidthResolution(board))
488         fTimeWidth[ipmt] = Float_t(kMaxTDCWidth)*fCalibData->GetWidthResolution(board);
489     }
490   }
491
492   fEvenOrOdd = gRandom->Integer(2);
493   for (Int_t j=0; j<64; ++j){
494     for (Int_t iClock = 0; iClock < kNClocks; ++iClock) {
495       Int_t integrator = (iClock + fEvenOrOdd) % 2;
496       AliDebug(1,Form("ADC %d %d %f",j,iClock,fAdc[j][iClock]));
497       fAdc[j][iClock]  += gRandom->Gaus(fAdcPedestal[j][integrator], fAdcSigma[j][integrator]);
498     }
499   }
500         
501 }
502
503 void AliVZERODigitizer::WriteDigits(AliLoader *loader)
504 {
505   // Take fAdc arrays filled by the previous
506   // method and produce and add digits to the digit Tree
507
508   loader->LoadDigits("UPDATE");
509
510   if (!loader->TreeD()) loader->MakeTree("D");
511   loader->MakeDigitsContainer();
512   TTree* treeD  = loader->TreeD();
513   DigitsArray();
514   treeD->Branch("VZERODigit", &fDigits); 
515   
516   Short_t *chargeADC = new Short_t[kNClocks];
517   for (Int_t i=0; i<64; i++) {      
518     for (Int_t j = 0; j < kNClocks; ++j) {
519       Int_t tempadc = Int_t(fAdc[i][j]);
520       if (tempadc > 1023) tempadc = 1023;
521       chargeADC[j] = tempadc;
522     }
523     AddDigit(i, fLeadingTime[i], fTimeWidth[i], Bool_t((10+fEvenOrOdd)%2), chargeADC, fLabels[i]);
524   }
525   delete [] chargeADC;
526
527   treeD->Fill();
528   loader->WriteDigits("OVERWRITE");  
529   loader->UnloadDigits();     
530   ResetDigits();
531 }
532
533 void AliVZERODigitizer::WriteSDigits(AliLoader *loader)
534 {
535   // Take fTime arrays filled by the previous
536   // method and produce and add sdigits to the sdigit Tree
537
538   loader->LoadSDigits("UPDATE");
539
540   if (!loader->TreeS()) loader->MakeTree("S");
541   loader->MakeSDigitsContainer();
542   TTree* treeS  = loader->TreeS();
543   SDigitsArray();
544   treeS->Branch("VZEROSDigit", &fDigits); 
545   
546   for (Int_t ipmt = 0; ipmt < 64; ++ipmt) {
547     AddSDigit(ipmt,fNBins[ipmt],fTime[ipmt],fLabels[ipmt]);
548   }
549
550   treeS->Fill();
551   loader->WriteSDigits("OVERWRITE");  
552   loader->UnloadSDigits();     
553   ResetDigits();
554 }
555
556 void AliVZERODigitizer::ReadSDigits()
557 {
558   // Read SDigits which are then to precessed
559   // in the following method
560   for(Int_t i = 0 ; i < 64; ++i) {
561     memset(fTime[i],0,fNBins[i]*sizeof(Float_t));
562     fLabels[i][0] = fLabels[i][1] = fLabels[i][2] = -1;
563   }
564
565   // Loop over input files
566   Int_t nFiles= fDigInput->GetNinputs();
567   for (Int_t inputFile = 0; inputFile < nFiles; inputFile++) {
568     // Get the current loader 
569     AliRunLoader* currentLoader = 
570       AliRunLoader::GetRunLoader(fDigInput->GetInputFolderName(inputFile));
571
572     AliLoader *loader = currentLoader->GetLoader("VZEROLoader");
573     loader->LoadSDigits("READ");
574   
575     // Get the tree of summable digits
576     TTree* sdigitsTree = loader->TreeS();
577     if (!sdigitsTree)  {
578       AliError("No sdigit tree from digInput");
579       continue;
580     }
581
582     // Get the branch 
583     TBranch* sdigitsBranch = sdigitsTree->GetBranch("VZEROSDigit");
584     if (!sdigitsBranch) {
585       AliError("Failed to get sdigit branch");
586       return;
587     }
588
589     // Set the branch address
590     TClonesArray *sdigitsArray = NULL;
591     sdigitsBranch->SetAddress(&sdigitsArray);
592
593     // Sum contributions from the sdigits
594     // Get number of entries in the tree 
595     Int_t nentries  = Int_t(sdigitsBranch->GetEntries());
596     for (Int_t entry = 0; entry < nentries; ++entry)  {
597       sdigitsBranch->GetEntry(entry);
598       // Get the number of sdigits 
599       Int_t nsdigits = sdigitsArray->GetEntries();
600       for (Int_t sdigit = 0; sdigit < nsdigits; sdigit++) {
601         AliVZEROSDigit* sDigit = static_cast<AliVZEROSDigit*>(sdigitsArray->UncheckedAt(sdigit));
602         Int_t pmNumber = sDigit->PMNumber();
603         Int_t nbins = sDigit->GetNBins();
604         if (nbins != fNBins[pmNumber]) {
605           AliError(Form("Incompatible number of bins between digitizer (%d) and sdigit (%d) for PM %d! Skipping sdigit!",
606                         fNBins[pmNumber],nbins,pmNumber));
607           continue;
608         }
609         // Sum the charges
610         Float_t *charges = sDigit->GetCharges();
611         for(Int_t iBin = 0; iBin < nbins; ++iBin) fTime[pmNumber][iBin] += charges[iBin];
612         // and the labels
613         Int_t *labels = sDigit->GetTracks();
614         Int_t j = 0;
615         for(Int_t i = 0; i < 3; ++i) {
616           if (fLabels[pmNumber][i] < 0) {
617             if (labels[j] < 0) break;
618             fLabels[pmNumber][i] = labels[j];
619             j++;
620           }
621         }
622       }
623     }
624     loader->UnloadSDigits();
625   }
626 }
627
628 //____________________________________________________________________
629 TClonesArray*
630 AliVZERODigitizer::DigitsArray() 
631 {
632   // Initialize digit array if not already and
633   // return pointer to it. 
634   if (!fDigits) { 
635     fDigits = new TClonesArray("AliVZEROdigit", 64);
636     fNdigits = 0;
637   }
638   return fDigits;
639 }
640
641 //____________________________________________________________________
642 TClonesArray*
643 AliVZERODigitizer::SDigitsArray() 
644 {
645   // Initialize sdigit array if not already and
646   // return pointer to it. 
647   if (!fDigits) { 
648     fDigits = new TClonesArray("AliVZEROSDigit", 64);
649     fNdigits = 0;
650   }
651   return fDigits;
652 }