]> git.uio.no Git - u/mrichter/AliRoot.git/blob - VZERO/AliVZERODigitizer.cxx
8d1c8a3f751d1411af3a95ccdb69f9614069e042
[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 "AliVZEROdigit.h"
51 #include "AliVZERODigitizer.h"
52
53 ClassImp(AliVZERODigitizer)
54
55  AliVZERODigitizer::AliVZERODigitizer()
56                    :AliDigitizer(),
57                     fCalibData(GetCalibData()),
58                     fPhotoCathodeEfficiency(0.18),
59                     fNdigits(0),
60                     fDigits(0),
61                     fSignalShape(NULL),
62                     fPMResponse(NULL),
63                     fSinglePhESpectrum(NULL)
64    
65 {
66   // default constructor
67
68 //    fNdigits = 0;
69 //    fDigits  = 0;
70 //   
71 //    fPhotoCathodeEfficiency =   0.18;
72 //    fPMVoltage              =  768.0;
73 //    fPMGain = TMath::Power((fPMVoltage / 112.5) ,7.04277); 
74    
75 //   fCalibData = GetCalibData();
76
77   fSignalShape = new TF1("VZEROSignalShape",this,&AliVZERODigitizer::SignalShape,0,200,6,"AliVZERODigitizer","SignalShape");
78   //  fSignalShape->SetParameters(0,1.57345e1,-4.25603e-1,2.9,6.40982,3.69339e-01);
79   fSignalShape->SetParameters(1.34330e+00,1.13007e+02,-4.95705e-01,
80                               3.68911e+00,1.01040e+00, 3.94675e-01);
81   fPMResponse = new TF1("VZEROPMResponse",this,&AliVZERODigitizer::PMResponse,-kPMRespTime,2.*kPMRespTime,0,"AliVZERODigitizer","PMResponse");
82   fSinglePhESpectrum = new TF1("VZEROSinglePhESpectrum",this,&AliVZERODigitizer::SinglePhESpectrum,0,20,0,"AliVZERODigitizer","SinglePhESpectrum");
83
84   // Now get the CTP L0->L1 delay
85   AliCDBEntry *entry = AliCDBManager::Instance()->Get("GRP/CTP/CTPtiming");
86   if (!entry) AliFatal("CTP timing parameters are not found in OCDB !");
87   AliCTPTimeParams *ctpParams = (AliCTPTimeParams*)entry->GetObject();
88   Float_t l1Delay = (Float_t)ctpParams->GetDelayL1L0()*25.0;
89
90   AliCDBEntry *entry2 = AliCDBManager::Instance()->Get("VZERO/Calib/TimeDelays");
91   if (!entry2) AliFatal("VZERO time delays are not found in OCDB !");
92   TH1F *delays = (TH1F*)entry2->GetObject();
93
94   for(Int_t i = 0 ; i < 64; ++i) {
95
96     for(Int_t j = 0; j < kNClocks; ++j) fAdc[i][j] = 0;
97     fLeadingTime[i] = fTimeWidth[i] = 0;
98
99     fPmGain[i] = fCalibData->GetGain(i);
100
101     fAdcPedestal[i][0] = fCalibData->GetPedestal(i);
102     fAdcSigma[i][0]    = fCalibData->GetSigma(i); 
103     fAdcPedestal[i][1] = fCalibData->GetPedestal(i+64);
104     fAdcSigma[i][1]    = fCalibData->GetSigma(i+64); 
105
106     Int_t board = i / 8;
107     fNBins[i] = TMath::Nint(((Float_t)(fCalibData->GetMatchWindow(board)+1)*25.0+
108                              (Float_t)kMaxTDCWidth*fCalibData->GetWidthResolution(board))/
109                             fCalibData->GetTimeResolution(board));
110     fNBinsLT[i] = TMath::Nint(((Float_t)(fCalibData->GetMatchWindow(board)+1)*25.0)/
111                               fCalibData->GetTimeResolution(board));
112     fBinSize[i] = fCalibData->GetTimeResolution(board);
113     fHptdcOffset[i] = (((Float_t)fCalibData->GetTriggerCountOffset(board)-
114                         (Float_t)fCalibData->GetRollOver(board))*25.0+
115                        fCalibData->GetTimeOffset(i)+
116                        l1Delay+
117                        delays->GetBinContent(i+1)+
118                        kV0Offset);
119
120     fTime[i] = new Float_t[fNBins[i]];
121     memset(fTime[i],0,fNBins[i]*sizeof(Float_t));
122   }
123
124 }
125
126 //____________________________________________________________________________ 
127   AliVZERODigitizer::AliVZERODigitizer(AliRunDigitizer* manager)
128                     :AliDigitizer(manager),
129                      fCalibData(GetCalibData()),
130                      fPhotoCathodeEfficiency(0.18),
131                      fNdigits(0),
132                      fDigits(0),
133                      fSignalShape(NULL),
134                      fPMResponse(NULL),
135                      fSinglePhESpectrum(NULL)
136                                         
137 {
138   // constructor
139   // Initialize OCDB and containers used in the digitization
140   
141   fSignalShape = new TF1("VZEROSignalShape",this,&AliVZERODigitizer::SignalShape,0,200,6,"AliVZERODigitizer","SignalShape");
142   //  fSignalShape->SetParameters(0,1.57345e1,-4.25603e-1,2.9,6.40982,3.69339e-01);
143   fSignalShape->SetParameters(1.34330e+00,1.13007e+02,-4.95705e-01,
144                               3.68911e+00,1.01040e+00, 3.94675e-01);
145   fPMResponse = new TF1("VZEROPMResponse",this,&AliVZERODigitizer::PMResponse,-kPMRespTime,2.*kPMRespTime,0,"AliVZERODigitizer","PMResponse");
146   fSinglePhESpectrum = new TF1("VZEROSinglePhESpectrum",this,&AliVZERODigitizer::SinglePhESpectrum,0,20,0,"AliVZERODigitizer","SinglePhESpectrum");
147   
148   // Now get the CTP L0->L1 delay
149   AliCDBEntry *entry = AliCDBManager::Instance()->Get("GRP/CTP/CTPtiming");
150   if (!entry) AliFatal("CTP timing parameters are not found in OCDB !");
151   AliCTPTimeParams *ctpParams = (AliCTPTimeParams*)entry->GetObject();
152   Float_t l1Delay = (Float_t)ctpParams->GetDelayL1L0()*25.0;
153
154   AliCDBEntry *entry2 = AliCDBManager::Instance()->Get("VZERO/Calib/TimeDelays");
155   if (!entry2) AliFatal("VZERO time delays are not found in OCDB !");
156   TH1F *delays = (TH1F*)entry2->GetObject();
157
158   for(Int_t i = 0 ; i < 64; ++i) {
159
160     for(Int_t j = 0; j < kNClocks; ++j) fAdc[i][j] = 0;
161     fLeadingTime[i] = fTimeWidth[i] = 0;
162
163     fPmGain[i] = fCalibData->GetGain(i);
164
165     fAdcPedestal[i][0] = fCalibData->GetPedestal(i);
166     fAdcSigma[i][0]    = fCalibData->GetSigma(i); 
167     fAdcPedestal[i][1] = fCalibData->GetPedestal(i+64);
168     fAdcSigma[i][1]    = fCalibData->GetSigma(i+64); 
169
170     Int_t board = i / 8;
171     fNBins[i] = TMath::Nint(((Float_t)(fCalibData->GetMatchWindow(board)+1)*25.0+
172                              (Float_t)kMaxTDCWidth*fCalibData->GetWidthResolution(board))/
173                             fCalibData->GetTimeResolution(board));
174     fNBinsLT[i] = TMath::Nint(((Float_t)(fCalibData->GetMatchWindow(board)+1)*25.0)/
175                               fCalibData->GetTimeResolution(board));
176     fBinSize[i] = fCalibData->GetTimeResolution(board);
177     fHptdcOffset[i] = (((Float_t)fCalibData->GetTriggerCountOffset(board)-
178                         (Float_t)fCalibData->GetRollOver(board))*25.0+
179                        fCalibData->GetTimeOffset(i)+
180                        l1Delay+
181                        delays->GetBinContent(i+1)+
182                        kV0Offset);
183
184     fTime[i] = new Float_t[fNBins[i]];
185     memset(fTime[i],0,fNBins[i]*sizeof(Float_t));
186   }
187 }
188            
189 //____________________________________________________________________________ 
190   AliVZERODigitizer::~AliVZERODigitizer()
191 {
192   // destructor
193   
194   if (fDigits) {
195     fDigits->Delete();
196     delete fDigits;
197     fDigits=0; 
198   }
199
200   if (fSignalShape) {
201     delete fSignalShape;
202     fSignalShape = NULL;
203   }
204   if (fPMResponse) {
205     delete fPMResponse;
206     fPMResponse = NULL;
207   }
208   if (fSinglePhESpectrum) {
209     delete fSinglePhESpectrum;
210     fSinglePhESpectrum = NULL;
211   }
212
213   for(Int_t i = 0 ; i < 64; ++i) {
214     if (fTime[i]) delete [] fTime[i];
215   }
216 }
217
218 //_____________________________________________________________________________
219 Bool_t AliVZERODigitizer::Init()
220 {
221   // Initialises the digitizer
222
223   // Initialises the Digit array
224   fDigits = new TClonesArray ("AliVZEROdigit", 1000);
225   
226   return kTRUE;
227 }
228
229 //____________________________________________________________________________
230 void AliVZERODigitizer::Exec(Option_t* /*option*/) 
231 {   
232   // Creates digits from hits
233   fNdigits     =    0;  
234
235   for(Int_t i = 0 ; i < 64; ++i) {
236     memset(fTime[i],0,fNBins[i]*sizeof(Float_t));
237     for(Int_t j = 0; j < kNClocks; ++j) fAdc[i][j] = 0;
238     fLeadingTime[i] = fTimeWidth[i] = 0;
239   }
240   Float_t integral = fPMResponse->Integral(-kPMRespTime,2.*kPMRespTime);
241   Float_t meansPhE = fSinglePhESpectrum->Mean(0,20);
242
243   AliRunLoader* outRunLoader =  AliRunLoader::GetRunLoader(fManager->GetOutputFolderName());    
244   if (!outRunLoader) {
245     Error("Exec", "Can not get output Run Loader");
246     return;
247   }
248     
249   AliLoader* outLoader = outRunLoader->GetLoader("VZEROLoader");
250
251   if (!outLoader) {
252     Error("Exec", "Can not get output VZERO Loader");
253     return;
254   }
255
256   const char* mode = "update";
257   if(outRunLoader->GetEventNumber() == 0) mode = "recreate";
258   outLoader->LoadDigits(mode);
259
260   if (!outLoader->TreeD()) outLoader->MakeTree("D");
261   outLoader->MakeDigitsContainer();
262   TTree* treeD  = outLoader->TreeD();
263   Int_t bufsize = 16000;
264   treeD->Branch("VZERODigit", &fDigits, bufsize); 
265   
266   for (Int_t iInput = 0; iInput < fManager->GetNinputs(); iInput++) {
267      AliRunLoader* runLoader = AliRunLoader::GetRunLoader(fManager->GetInputFolderName(iInput));
268      AliLoader* loader = runLoader->GetLoader("VZEROLoader");
269      if (!loader) {
270        Error("Exec", "Can not get VZERO Loader for input %d", iInput);
271        continue;
272          }
273       
274      if (!runLoader->GetAliRun()) runLoader->LoadgAlice();
275
276      AliVZERO* vzero = (AliVZERO*) runLoader->GetAliRun()->GetDetector("VZERO");
277      if (!vzero) {
278        Error("Exec", "No VZERO detector for input %d", iInput);
279        continue;
280          }
281       
282      loader->LoadHits();
283      TTree* treeH = loader->TreeH();
284      if (!treeH) {
285        Error("Exec", "Cannot get TreeH for input %d", iInput);
286        continue; 
287          }
288        
289      TClonesArray* hits = vzero->Hits();
290
291      //     Float_t lightYieldCorr[64] = {0.00707,0.00517,0.00520,0.00537,0.00735,0.00537,0.00733,0.00605,0.00778,0.00749,0.00701,0.00755,0.00732,0.00617,0.00669,0.00525,0.00752,0.00820,0.00797,0.01107,0.01080,0.00889,0.00880,0.01712,0.00866,0.00701,0.00811,0.00602,0.00879,0.00821,0.00861,0.01433,0.00061,0.00032,0.00099,0.00061,0.00034,0.00046,0.00031,0.00122,0.00155,0.00091,0.00032,0.00096,0.00120,0.00067,0.00113,0.00060,0.00158,0.00136,0.00340,0.00066,0.00076,0.00119,0.00129,0.00147,0.00137,0.00117,0.00088,0.00164,0.00128,0.00081,0.00121,0.00250};
292      Float_t lightYieldCorr[64] = {0.01173,0.00874,0.00878,0.00886,0.01151,0.00925,0.01167,0.00983,0.01181,0.01243,0.01115,0.01220,0.01228,0.01053,0.01021,0.00930,0.01270,0.01411,0.01316,0.01894,0.01923,0.01860,0.01738,0.00305,0.01584,0.01251,0.01344,0.00310,0.01302,0.01266,0.01407,0.00338,0.00089,0.00100,0.00130,0.00081,0.00052,0.01230,0.00059,0.02452,0.02980,0.00137,0.01421,0.00116,0.00141,0.00092,0.02480,0.00096,0.00182,0.00174,0.00218,0.00106,0.00116,0.00160,0.00162,0.03097,0.00194,0.00171,0.00132,0.00239,0.00173,0.00118,0.00163,0.00262};
293 //  Now makes Digits from hits
294      Int_t nTracks = (Int_t) treeH->GetEntries();
295      for (Int_t iTrack = 0; iTrack < nTracks; iTrack++) {
296          vzero->ResetHits();
297          treeH->GetEvent(iTrack);
298          Int_t nHits = hits->GetEntriesFast();
299          for (Int_t iHit = 0; iHit < nHits; iHit++) {
300            AliVZEROhit* hit = (AliVZEROhit *)hits->UncheckedAt(iHit);
301            Int_t nPhot = hit->Nphot();
302            Int_t cell  = hit->Cell();                          
303            Int_t pmt = Cell2Pmt(cell);
304            Float_t dt_scintillator = gRandom->Gaus(0,kIntTimeRes);
305            Float_t t = dt_scintillator + 1e9*hit->Tof();
306            if (pmt < 32) t += kV0CDelayCables;
307            t += fHptdcOffset[pmt];
308            Int_t nPhE;
309            Float_t prob = lightYieldCorr[pmt]*fPhotoCathodeEfficiency; // Optical losses included!
310            if (nPhot > 100)
311              nPhE =gRandom->Gaus(prob*Float_t(nPhot)+0.5,
312                                  sqrt(Float_t(nPhot)*prob*(1.-prob)));
313            else
314              nPhE = gRandom->Binomial(nPhot,prob);
315            Float_t charge = TMath::Qe()*fPmGain[pmt]*fBinSize[pmt]/integral;
316            for (Int_t iPhE = 0; iPhE < nPhE; ++iPhE) {
317              Float_t tPhE = t + fSignalShape->GetRandom(0,fBinSize[pmt]*Float_t(fNBins[pmt]));
318              Float_t gainVar = fSinglePhESpectrum->GetRandom(0,20)/meansPhE;
319              Int_t firstBin = TMath::Max(0,(Int_t)((tPhE-kPMRespTime)/fBinSize[pmt]));
320              Int_t lastBin = TMath::Min(fNBins[pmt]-1,(Int_t)((tPhE+2.*kPMRespTime)/fBinSize[pmt]));
321              for(Int_t iBin = firstBin; iBin <= lastBin; ++iBin) {
322                Float_t tempT = fBinSize[pmt]*(0.5+iBin)-tPhE;
323                fTime[pmt][iBin] += gainVar*charge*fPMResponse->Eval(tempT);
324              }
325            }         // ph.e. loop
326          }           // hit loop
327      }               // track loop
328      loader->UnloadHits();
329   }                  // input loop
330
331   Float_t maximum = 0.9*fSignalShape->GetMaximum(0,200); // Not exact, one needs to do this on the convoluted
332   Float_t integral2 = fSignalShape->Integral(0,200); // function. Anyway the effect is small <10% on the 2.5 ADC thr
333   for (Int_t ipmt = 0; ipmt < 64; ++ipmt) {
334     Float_t thr = fCalibData->GetDiscriThr(ipmt)*kChargePerADC*maximum*fBinSize[ipmt]/integral2;
335     Bool_t ltFound = kFALSE, ttFound = kFALSE;
336     for (Int_t iBin = 0; iBin < fNBins[ipmt]; ++iBin) {
337       Float_t t = fBinSize[ipmt]*Float_t(iBin+0.5);
338       if (fTime[ipmt][iBin] > thr) {
339         if (!ltFound && (iBin < fNBinsLT[ipmt])) {
340           ltFound = kTRUE;
341           fLeadingTime[ipmt] = t;
342         }
343       }
344       else {
345         if (ltFound) {
346           if (!ttFound) {
347             ttFound = kTRUE;
348             fTimeWidth[ipmt] = t - fLeadingTime[ipmt];
349           }
350         }
351       }
352       Float_t tadc = t - kADCTimeOffset;
353       Int_t clock = kNClocks - Int_t(tadc/25.0) - 1;
354       if (clock >= 0 && clock < kNClocks)
355         fAdc[ipmt][clock] += fTime[ipmt][iBin]/kChargePerADC;
356     }
357     Int_t board = ipmt / 8;
358     if (ltFound && ttFound) {
359       fTimeWidth[ipmt] = fCalibData->GetWidthResolution(board)*
360         Float_t(Int_t(fTimeWidth[ipmt]/fCalibData->GetWidthResolution(board)));
361       if (fTimeWidth[ipmt] < Float_t(kMinTDCWidth)*fCalibData->GetWidthResolution(board))
362         fTimeWidth[ipmt] = Float_t(kMinTDCWidth)*fCalibData->GetWidthResolution(board);
363       if (fTimeWidth[ipmt] > Float_t(kMaxTDCWidth)*fCalibData->GetWidthResolution(board))
364         fTimeWidth[ipmt] = Float_t(kMaxTDCWidth)*fCalibData->GetWidthResolution(board);
365     }
366   }
367
368   Int_t evenOrOdd = gRandom->Integer(2);
369   for (Int_t j=0; j<64; ++j){
370     for (Int_t iClock = 0; iClock < kNClocks; ++iClock) {
371       Int_t integrator = (iClock + evenOrOdd) % 2;
372       fAdc[j][iClock]  += gRandom->Gaus(fAdcPedestal[j][integrator], fAdcSigma[j][integrator]);
373     }
374   }
375         
376   // Now add digits to the digit Tree 
377
378   for (Int_t i=0; i<64; i++) {      
379     //    printf(" cell %d adc ",i);
380     //    for (Int_t j = 0; j < kNClocks; ++j)
381     //      printf("%d ", Int_t(fAdc[i][j]));
382     //    printf("\n");
383     Float_t totADC = 0;
384     for (Int_t j = 8; j <= 11; ++j) {
385       Int_t tempadc = Int_t(fAdc[i][j]);
386       if (tempadc > 1023) tempadc = 1023;
387       Int_t integrator = (j + evenOrOdd) % 2;
388       if ((Float_t(tempadc) - fAdcPedestal[i][integrator]) > (2.*fAdcSigma[i][integrator]))
389         totADC += (Float_t(tempadc) - fAdcPedestal[i][integrator]);
390     }
391     totADC += fAdcPedestal[i][(10+evenOrOdd)%2];
392     AddDigit(i, totADC, fLeadingTime[i], fTimeWidth[i], Bool_t((10+evenOrOdd)%2));
393   }
394
395   treeD->Fill();
396   outLoader->WriteDigits("OVERWRITE");  
397   outLoader->UnloadDigits();     
398   ResetDigit();
399 }
400
401 //____________________________________________________________________________
402 void AliVZERODigitizer::AddDigit(Int_t PMnumber, Float_t adc, Float_t time, Float_t width, Bool_t integrator) 
403  { 
404  
405 // Adds Digit 
406  
407   TClonesArray &ldigits = *fDigits;  
408          
409   new(ldigits[fNdigits++]) AliVZEROdigit(PMnumber,adc,time,width,kFALSE,kFALSE,integrator);
410          
411 }
412 //____________________________________________________________________________
413 void AliVZERODigitizer::ResetDigit()
414 {
415
416 // Clears Digits
417
418   fNdigits = 0;
419   if (fDigits) fDigits->Delete();
420 }
421
422 //____________________________________________________________________________
423 AliVZEROCalibData* AliVZERODigitizer::GetCalibData() const
424
425 {
426   AliCDBManager *man = AliCDBManager::Instance();
427
428   AliCDBEntry *entry=0;
429
430   entry = man->Get("VZERO/Calib/Data");
431
432 //   if(!entry){
433 //     AliWarning("Load of calibration data from default storage failed!");
434 //     AliWarning("Calibration data will be loaded from local storage ($ALICE_ROOT)");
435 //     Int_t runNumber = man->GetRun();
436 //     entry = man->GetStorage("local://$ALICE_ROOT/OCDB")
437 //       ->Get("VZERO/Calib/Data",runNumber);
438 //      
439 //   }
440
441   // Retrieval of data in directory VZERO/Calib/Data:
442
443
444   AliVZEROCalibData *calibdata = 0;
445
446   if (entry) calibdata = (AliVZEROCalibData*) entry->GetObject();
447   if (!calibdata)  AliFatal("No calibration data from calibration database !");
448
449   return calibdata;
450
451 }
452
453 double AliVZERODigitizer::SignalShape(double *x, double *par)
454 {
455   // this function simulates the time
456   // of arrival of the photons at the
457   // photocathode
458   Double_t xx = x[0];
459   if (xx <= par[0]) return 0;
460   Double_t a = 1./TMath::Power((xx-par[0])/par[1],1./par[2]);
461   if (xx <= par[3]) return a;
462   Double_t b = 1./TMath::Power((xx-par[3])/par[4],1./par[5]);
463   Double_t f = a*b/(a+b);
464   AliDebug(100,Form("x=%f func=%f",xx,f));
465   return f;
466 }
467
468 double AliVZERODigitizer::PMResponse(double *x, double * /* par */)
469 {
470   // this function describes the
471   // PM time response to a single
472   // photoelectron
473   Double_t xx = x[0]+kPMRespTime;
474   return xx*xx*TMath::Exp(-xx*xx/(kPMRespTime*kPMRespTime));
475 }
476
477 double AliVZERODigitizer::SinglePhESpectrum(double *x, double * /* par */)
478 {
479   // this function describes the
480   // PM amplitude response to a single
481   // photoelectron
482   Double_t xx = x[0];
483   if (xx < 0) return 0;
484   return (TMath::Poisson(xx,kPMNbOfSecElec)+kPMTransparency*TMath::Poisson(xx,1.0));
485 }
486
487 Int_t AliVZERODigitizer::Cell2Pmt(Int_t cell) const
488 {
489   // The method maps the scintillator
490   // indexes to the PM ones
491   if (cell < 0 || cell >= 80) {
492     AliError(Form("Wrong VZERO cell index %d",cell));
493     return -1;
494   }
495   if (cell < 16) return cell;
496   if (cell < 48) return 8 + cell/2;
497   return cell - 16;
498 }