]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ZDC/AliZDCDigitizer.cxx
Index corrected in AliZDCReconstructor.cxx + getting rid of unused calib object in...
[u/mrichter/AliRoot.git] / ZDC / AliZDCDigitizer.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 //                      ZDC digitizer class                                  //
21 //                                                                           //
22 ///////////////////////////////////////////////////////////////////////////////
23
24 #include <stdlib.h>
25
26 // --- ROOT system
27 #include <TTree.h>
28 #include <TFile.h>
29 #include <TNtuple.h>
30 #include <TRandom.h>
31
32 // --- AliRoot header files
33 #include "AliLog.h"
34 #include "AliRun.h"
35 #include "AliHeader.h"
36 #include "AliGenHijingEventHeader.h"
37 #include "AliRunDigitizer.h"
38 #include "AliRunLoader.h"
39 #include "AliGRPObject.h"
40 #include "AliCDBManager.h"
41 #include "AliCDBEntry.h"
42 #include "AliZDCSDigit.h"
43 #include "AliZDCDigit.h"
44 #include "AliZDCFragment.h"
45 #include "AliZDCv3.h"
46 #include "AliZDCDigitizer.h"
47
48 class AliCDBStorage;
49 class AliZDCPedestals;
50
51 ClassImp(AliZDCDigitizer)
52
53
54 //____________________________________________________________________________
55 AliZDCDigitizer::AliZDCDigitizer() :
56   fIsCalibration(0), 
57   fIsSignalInADCGate(kFALSE),
58   fFracLostSignal(0.),
59   fPedData(0), 
60   fSpectators2Track(kFALSE)
61 {
62   // Default constructor    
63
64 }
65
66 //____________________________________________________________________________
67 AliZDCDigitizer::AliZDCDigitizer(AliRunDigitizer* manager):
68   AliDigitizer(manager),
69   fIsCalibration(0), //By default the simulation doesn't create calib. data
70   fIsSignalInADCGate(kFALSE),
71   fFracLostSignal(0.),
72   fPedData(GetPedData()), 
73   fSpectators2Track(kFALSE)
74 {
75   // Get calibration data
76   if(fIsCalibration!=0) printf("\n\t AliZDCDigitizer -> Creating calibration data (pedestals)\n");
77   for(Int_t i=0; i<6; i++){
78     for(Int_t j=0; j<5; j++)
79        fPMGain[i][j] = 0.;
80   }
81 }
82
83 //____________________________________________________________________________
84 AliZDCDigitizer::~AliZDCDigitizer()
85 {
86 // Destructor
87 // Not implemented
88 }
89
90
91 //____________________________________________________________________________
92 AliZDCDigitizer::AliZDCDigitizer(const AliZDCDigitizer &digitizer):
93   AliDigitizer(),
94   fIsCalibration(digitizer.fIsCalibration),
95   fIsSignalInADCGate(digitizer.fIsSignalInADCGate),
96   fFracLostSignal(digitizer.fFracLostSignal),
97   fPedData(digitizer.fPedData),
98   fSpectators2Track(digitizer.fSpectators2Track)
99 {
100   // Copy constructor
101
102   for(Int_t i=0; i<6; i++){
103      for(Int_t j=0; j<5; j++){
104         fPMGain[i][j]   = digitizer.fPMGain[i][j];           
105      }
106   }
107   for(Int_t i=0; i<2; i++) fADCRes[i] = digitizer.fADCRes[i];
108
109
110 }
111
112 //____________________________________________________________________________
113 Bool_t AliZDCDigitizer::Init()
114 {
115   // Initialize the digitizer
116   
117   AliCDBEntry*  entry = AliCDBManager::Instance()->Get("GRP/GRP/Data");
118   AliGRPObject* grpData = 0x0;
119   if(entry){
120     TMap* m = dynamic_cast<TMap*>(entry->GetObject());  // old GRP entry
121     if(m){
122       //m->Print();
123       grpData = new AliGRPObject();
124       grpData->ReadValuesFromMap(m);
125     }
126     else{
127       grpData = dynamic_cast<AliGRPObject*>(entry->GetObject());  // new GRP entry
128     }
129     entry->SetOwner(0);
130     AliCDBManager::Instance()->UnloadFromCache("GRP/GRP/Data");
131   }
132   if(!grpData) AliError("No GRP entry found in OCDB!");
133   
134   TString beamType = grpData->GetBeamType();
135   if(beamType==AliGRPObject::GetInvalidString()){
136     AliError("\t UNKNOWN beam type from GRP obj -> PMT gains not set in ZDC digitizer!!!\n");
137   }
138   
139   Float_t beamEnergy = grpData->GetBeamEnergy();
140   if(beamEnergy==AliGRPObject::GetInvalidFloat()){
141     AliWarning("GRP/GRP/Data entry:  missing value for the beam energy ! Using 0.");
142     AliError("\t UNKNOWN beam type from GRP obj -> PMT gains not set in ZDC digitizer!!!\n");
143     beamEnergy = 0.;
144   }
145
146   if((beamType.CompareTo("P-P")) == 0){
147     //PTM gains rescaled to beam energy for p-p
148     if(beamEnergy != 0){
149       for(Int_t j = 0; j < 5; j++){
150         fPMGain[0][j] = (661.444/beamEnergy+0.000740671)*10000000;
151         fPMGain[1][j] = (864.350/beamEnergy+0.002344)*10000000;
152         fPMGain[2][j] = (1.32312-0.000101515*beamEnergy)*10000000;
153         fPMGain[3][j] = fPMGain[0][j];
154         fPMGain[4][j] = fPMGain[1][j] ;
155       }
156       AliInfo(Form("    PMT gains for p-p @ %1.0f+%1.0f: ZN(%1.0f), ZP(%1.0f), ZEM(%1.0f)\n",
157         beamEnergy, beamEnergy, fPMGain[0][0], fPMGain[1][0], fPMGain[2][0]));
158     }
159   }
160   else{
161     // PTM gains for Pb-Pb @ 2.7_2.7 A TeV ***************
162     for(Int_t j = 0; j < 5; j++){
163       fPMGain[0][j] = 50000.;
164       fPMGain[1][j] = 100000.;
165       fPMGain[2][j] = 100000.;
166       fPMGain[3][j] = 50000.;
167       fPMGain[4][j] = 100000.;
168       fPMGain[5][j] = 100000.;
169     }
170     AliInfo(Form("    PMT gains for Pb-Pb @ %1.0f+%1.0f: ZN(%1.0f), ZP(%1.0f), ZEM(%1.0f)\n",
171         beamEnergy, beamEnergy, fPMGain[0][0], fPMGain[1][0], fPMGain[2][0]));
172   }
173     
174   // ADC Caen V965
175   fADCRes[0] = 0.0000008; // ADC Resolution high gain: 200 fC/adcCh
176   fADCRes[1] = 0.0000064; // ADC Resolution low gain:  25  fC/adcCh
177
178   return kTRUE;
179 }
180
181 //____________________________________________________________________________
182 void AliZDCDigitizer::Exec(Option_t* /*option*/)
183 {
184   // Execute digitization
185
186   // ------------------------------------------------------------
187   // --- pm[0][...] = light in ZN side C  [C, Q1, Q2, Q3, Q4]
188   // --- pm[1][...] = light in ZP side C [C, Q1, Q2, Q3, Q4]
189   // --- pm[2][...] = light in ZEM [x, 1, 2, x, x]
190   // --- pm[3][...] = light in ZN side A [C, Q1, Q2, Q3, Q4] ->NEW!
191   // --- pm[4][...] = light in ZP side A [C, Q1, Q2, Q3, Q4] ->NEW!
192   // ------------------------------------------------------------
193   Float_t pm[5][5]; 
194   for(Int_t iSector1=0; iSector1<5; iSector1++) 
195     for(Int_t iSector2=0; iSector2<5; iSector2++){
196       pm[iSector1][iSector2] = 0;
197     }
198     
199   // ------------------------------------------------------------
200   // ### Out of time ADC added (22 channels)
201   // --- same codification as for signal PTMs (see above)
202   // ------------------------------------------------------------
203   Float_t pmoot[5][5];
204   for(Int_t iSector1=0; iSector1<5; iSector1++) 
205     for(Int_t iSector2=0; iSector2<5; iSector2++){
206       pmoot[iSector1][iSector2] = 0;
207     }
208
209   // impact parameter and number of spectators
210   Float_t impPar = -1;
211   Int_t specNTarg = 0, specPTarg = 0;
212   Int_t specNProj = 0, specPProj = 0;
213   Float_t signalTime0 = 0.;
214
215   // loop over input streams
216   for(Int_t iInput = 0; iInput<fManager->GetNinputs(); iInput++){
217
218     // get run loader and ZDC loader
219     AliRunLoader* runLoader = 
220       AliRunLoader::GetRunLoader(fManager->GetInputFolderName(iInput));
221     AliLoader* loader = runLoader->GetLoader("ZDCLoader");
222     if(!loader) continue;
223
224     // load sdigits
225     loader->LoadSDigits();
226     TTree* treeS = loader->TreeS();
227     if(!treeS) continue;
228     AliZDCSDigit sdigit;
229     AliZDCSDigit* psdigit = &sdigit;
230     treeS->SetBranchAddress("ZDC", &psdigit);
231
232     // loop over sdigits
233     for(Int_t iSDigit=0; iSDigit<treeS->GetEntries(); iSDigit++){
234       treeS->GetEntry(iSDigit);
235       //
236       if(!psdigit) continue;
237       if((sdigit.GetSector(1) < 0) || (sdigit.GetSector(1) > 4)){
238         AliError(Form("\nsector[0] = %d, sector[1] = %d\n", 
239                       sdigit.GetSector(0), sdigit.GetSector(1)));
240         continue;
241       }
242       // Checking if signal is inside ADC gate
243       if(iSDigit==0) signalTime0 = sdigit.GetTrackTime();
244       else{
245         // Assuming a signal lenght of 20 ns, signal is in gate if 
246         // signal ENDS in signalTime0+50. -> sdigit.GetTrackTime()+20<=signalTime0+50.
247         if(sdigit.GetTrackTime()<=signalTime0+30.) fIsSignalInADCGate = kTRUE;
248         if(sdigit.GetTrackTime()>signalTime0+30.){
249           fIsSignalInADCGate = kFALSE;
250           // Vedi quaderno per spiegazione approx. usata 
251           // nel calcolo della fraz. di segnale perso
252           fFracLostSignal = (sdigit.GetTrackTime()-30)*(sdigit.GetTrackTime()-30)/280.;
253         }
254       }
255       
256       Float_t sdSignal = sdigit.GetLightPM();
257       if(fIsSignalInADCGate == kFALSE){
258         AliDebug(2,Form("\t Signal time %f -> fraction %f of ZDC signal for det.(%d, %d) out of ADC gate\n",
259                 sdigit.GetTrackTime(),fFracLostSignal,sdigit.GetSector(0),sdigit.GetSector(1)));
260         sdSignal = (1-fFracLostSignal)*sdSignal;
261       }
262       
263       pm[(sdigit.GetSector(0))-1][sdigit.GetSector(1)] += sdigit.GetLightPM();
264       //Ch. debug
265       /*printf("\t Detector %d, Tower %d -> pm[%d][%d] = %.0f \n",
266           sdigit.GetSector(0), sdigit.GetSector(1),sdigit.GetSector(0)-1,
267           sdigit.GetSector(1), pm[sdigit.GetSector(0)-1][sdigit.GetSector(1)]); 
268       */
269       
270     }
271
272     loader->UnloadSDigits();
273
274     // get the impact parameter and the number of spectators in case of hijing
275     if(!runLoader->GetAliRun()) runLoader->LoadgAlice();
276     AliHeader* header = runLoader->GetHeader();
277     if(!header) continue;
278     AliGenEventHeader* genHeader = header->GenEventHeader();
279     if(!genHeader) continue;
280     if(!genHeader->InheritsFrom(AliGenHijingEventHeader::Class())) continue;
281     
282     if(fSpectators2Track==kTRUE){
283       impPar = ((AliGenHijingEventHeader*) genHeader)->ImpactParameter(); 
284       specNProj = ((AliGenHijingEventHeader*) genHeader)->ProjSpectatorsn();
285       specPProj = ((AliGenHijingEventHeader*) genHeader)->ProjSpectatorsp();
286       specNTarg = ((AliGenHijingEventHeader*) genHeader)->TargSpectatorsn();
287       specPTarg = ((AliGenHijingEventHeader*) genHeader)->TargSpectatorsp();
288       printf("\n\t AliZDCDigitizer: b = %1.2f fm\n"
289       " \t    PROJ.:  #spectator n %d, #spectator p %d\n"
290       " \t    TARG.:  #spectator n %d, #spectator p %d\n\n", 
291       impPar, specNProj, specPProj, specNTarg, specPTarg);
292     }
293     
294   }
295
296   // Applying fragmentation algorithm and adding spectator signal
297   if(fSpectators2Track==kTRUE && impPar) {
298     Int_t freeSpecNProj, freeSpecPProj;
299     Fragmentation(impPar, specNProj, specPProj, freeSpecNProj, freeSpecPProj);
300     Int_t freeSpecNTarg, freeSpecPTarg;
301     Fragmentation(impPar, specNTarg, specPTarg, freeSpecNTarg, freeSpecPTarg);
302     SpectatorSignal(1, freeSpecNProj, pm);
303     //printf("    AliZDCDigitizer -> Signal for %d PROJ free spectator n",freeSpecNProj);
304     SpectatorSignal(2, freeSpecPProj, pm);
305     //printf(" and %d free spectator p added\n",freeSpecPProj);
306     SpectatorSignal(3, freeSpecNTarg, pm);
307     //printf("    AliZDCDigitizer -> Signal for %d TARG free spectator n",freeSpecNTarg);
308     SpectatorSignal(4, freeSpecPTarg, pm);
309     //printf("and %d free spectator p added\n",freeSpecPTarg);
310   }
311
312
313   // get the output run loader and loader
314   AliRunLoader* runLoader = 
315     AliRunLoader::GetRunLoader(fManager->GetOutputFolderName());
316   AliLoader* loader = runLoader->GetLoader("ZDCLoader");
317   if(!loader) {
318     AliError("no ZDC loader found");
319     return;
320   }
321
322   // create the output tree
323   const char* mode = "update";
324   if(runLoader->GetEventNumber() == 0) mode = "recreate";
325   loader->LoadDigits(mode);
326   loader->MakeTree("D");
327   TTree* treeD = loader->TreeD();
328   AliZDCDigit digit;
329   AliZDCDigit* pdigit = &digit;
330   const Int_t kBufferSize = 4000;
331   treeD->Branch("ZDC", "AliZDCDigit", &pdigit, kBufferSize);
332
333   // Create digits
334   Int_t sector[2];
335   Int_t digi[2], digioot[2];
336   for(sector[0]=1; sector[0]<6; sector[0]++){
337     for(sector[1]=0; sector[1]<5; sector[1]++){
338         if((sector[0]==3) && ((sector[1]<1) || (sector[1]>2))) continue;
339         for(Int_t res=0; res<2; res++){
340            digi[res] = Phe2ADCch(sector[0], sector[1], pm[sector[0]-1][sector[1]], res) 
341                     + Pedestal(sector[0], sector[1], res);
342         }
343         new(pdigit) AliZDCDigit(sector, digi);
344         treeD->Fill();
345
346         //Ch. debug
347         //printf("\t DIGIT added -> det %d quad %d - digi[0,1] = [%d, %d]\n",
348         //     sector[0], sector[1], digi[0], digi[1]); // Chiara debugging!
349         
350     }
351   } // Loop over detector
352   // Adding in-time digits for 2 reference PTM signals (after signal ch.)
353   // (for the moment the ref. signal is completely invented assuming a PMgain of 5*10^4!)
354   Int_t sectorRef[2];
355   sectorRef[1] = 5;
356   Int_t sigRef[2];
357   // Reference signal are set to 100 (high gain chain) and 800 (low gain chain)
358   if(fIsCalibration==0) {sigRef[0]=100;  sigRef[1]=800;}
359   else {sigRef[0]=0;  sigRef[1]=0;} // calibration -> simulation of pedestal values
360   //
361   for(Int_t iref=0; iref<2; iref++){
362      sectorRef[0] = 3*iref+1;
363      for(Int_t res=0; res<2; res++){
364        sigRef[res] += Pedestal(sectorRef[0], sectorRef[1], res);
365      }
366      new(pdigit) AliZDCDigit(sectorRef, sigRef);
367      treeD->Fill();     
368      
369      //Ch. debug
370      //printf("\t RefDigit added -> det = %d, quad = %d - digi[0,1] = [%d, %d]\n",
371      //    sectorRef[0], sectorRef[1], sigRef[0], sigRef[1]); // Chiara debugging!     
372   }
373   //
374   // --- Adding digits for out-of-time channels after signal digits
375   for(sector[0]=1; sector[0]<6; sector[0]++){
376     for(sector[1]=0; sector[1]<5; sector[1]++){
377         if((sector[0]==3) && ((sector[1]<1) || (sector[1]>2))) continue;
378         for(Int_t res=0; res<2; res++){
379            digioot[res] = Pedestal(sector[0], sector[1], res); // out-of-time ADCs
380         }
381         new(pdigit) AliZDCDigit(sector, digioot);
382         treeD->Fill();
383
384         //Ch. debug
385         //printf("\t DIGIToot added -> det = %d, quad = %d - digi[0,1] = [%d, %d]\n",
386         //     sector[0], sector[1], digioot[0], digioot[1]); // Chiara debugging!      
387     }
388   }
389   // Adding out-of-time digits for 2 reference PTM signals (after out-of-time ch.)
390   Int_t sigRefoot[2];
391   for(Int_t iref=0; iref<2; iref++){
392      sectorRef[0] = 3*iref+1;
393      for(Int_t res=0; res<2; res++){
394        sigRefoot[res] = Pedestal(sectorRef[0], sectorRef[1], res);
395      }
396      new(pdigit) AliZDCDigit(sectorRef, sigRefoot);
397      treeD->Fill();
398      //Ch. debug
399      //printf("\t RefDigitoot added -> det = %d, quad = %d - digi[0,1] = [%d, %d]\n",
400      //    sectorRef[0], sectorRef[1], sigRefoot[0], sigRefoot[1]); // Chiara debugging!
401      
402   }
403   //printf("\t AliZDCDigitizer -> TreeD has %d entries\n",(Int_t) treeD->GetEntries());
404
405   // write the output tree
406   loader->WriteDigits("OVERWRITE");
407   loader->UnloadDigits();
408 }
409
410
411 //_____________________________________________________________________________
412 void AliZDCDigitizer::Fragmentation(Float_t impPar, Int_t specN, Int_t specP,
413                                     Int_t &freeSpecN, Int_t &freeSpecP) const
414 {
415 // simulate fragmentation of spectators
416
417   AliZDCFragment frag(impPar);
418
419   // Fragments generation
420   frag.GenerateIMF();
421   Int_t nAlpha = frag.GetNalpha();
422
423   // Attach neutrons
424   Int_t ztot = frag.GetZtot();
425   Int_t ntot = frag.GetNtot();
426   frag.AttachNeutrons();
427   freeSpecN = specN-ntot-2*nAlpha;
428   freeSpecP = specP-ztot-2*nAlpha;
429   // Removing deuterons
430   Int_t ndeu = (Int_t) (freeSpecN*frag.DeuteronNumber());
431   freeSpecN -= ndeu;
432   //
433   if(freeSpecN<0) freeSpecN=0;
434   if(freeSpecP<0) freeSpecP=0;
435   AliDebug(2, Form("FreeSpn = %d, FreeSpp = %d", freeSpecN, freeSpecP));
436 }
437
438 //_____________________________________________________________________________
439 void AliZDCDigitizer::SpectatorSignal(Int_t SpecType, Int_t numEvents, 
440                                       Float_t pm[5][5]) const
441 {
442 // add signal of the spectators
443   
444   TString hfn; 
445   if(SpecType == 1) {      // --- Signal for projectile spectator neutrons
446     hfn = "$ALICE_ROOT/ZDC/ZNCSignal.root";
447   } 
448   else if(SpecType == 2) { // --- Signal for projectile spectator protons
449     hfn = "$ALICE_ROOT/ZDC/ZPCSignal.root";
450   }
451   else if(SpecType == 3) { // --- Signal for target spectator neutrons
452     hfn = "$ALICE_ROOT/ZDC/ZNASignal.root";
453   }
454   else if(SpecType == 4) { // --- Signal for target spectator protons
455     hfn = "$ALICE_ROOT/ZDC/ZPASignal.root";
456   }
457   
458   TFile* file = TFile::Open(hfn);
459   if(!file || !file->IsOpen()) {
460     AliError((Form(" Opening file %s failed\n",hfn.Data())));
461     return;
462   }
463
464   TNtuple* zdcSignal = (TNtuple*) file->Get("ZDCSignal");
465   Int_t nentries = (Int_t) zdcSignal->GetEntries();
466   
467   Float_t *entry;
468   Int_t pl, i, k, iev=0, rnd[125], volume[2];
469   for(pl=0;pl<125;pl++) rnd[pl] = 0;
470   if(numEvents > 125) {
471     AliDebug(2,Form("numEvents (%d) is larger than 125", numEvents));
472     numEvents = 125;
473   }
474   for(pl=0;pl<numEvents;pl++){
475      rnd[pl] = (Int_t) (9999*gRandom->Rndm());
476      if(rnd[pl] >= 9999) rnd[pl] = 9998;
477      //printf(" rnd[%d] = %d\n",pl,rnd[pl]);     
478   }
479   // Sorting vector in ascending order with C function QSORT 
480   qsort((void*)rnd,numEvents,sizeof(Int_t),comp);
481   do{
482      for(i=0; i<nentries; i++){  
483         zdcSignal->GetEvent(i);
484         entry = zdcSignal->GetArgs();
485         if(entry[0] == rnd[iev]){
486           for(k=0; k<2; k++) volume[k] = (Int_t) entry[k+1];
487           //
488           Float_t lightQ = entry[7];
489           Float_t lightC = entry[8];
490           //
491           if(volume[0] != 3) {  // ZN or ZP
492             pm[volume[0]-1][0] += lightC;
493             pm[volume[0]-1][volume[1]] += lightQ;
494             //printf("\n   pm[%d][0] = %.0f, pm[%d][%d] = %.0f\n",(volume[0]-1),pm[volume[0]-1][0],
495             //  (volume[0]-1),volume[1],pm[volume[0]-1][volume[1]]);
496           } 
497           else { 
498             if(volume[1] == 1) pm[2][1] += lightC; // ZEM 1
499             else               pm[2][2] += lightQ; // ZEM 2
500             //printf("\n   pm[2][1] = %.0f, pm[2][2] = %.0f\n",pm[2][1],pm[2][2]);
501           }
502         }
503         else if(entry[0] > rnd[iev]){
504           iev++;
505           continue;
506         }
507      }
508   }while(iev<numEvents);
509   
510   file->Close();
511   delete file;
512 }
513
514
515 //_____________________________________________________________________________
516 Int_t AliZDCDigitizer::Phe2ADCch(Int_t Det, Int_t Quad, Float_t Light, 
517                                  Int_t Res) const
518 {
519   // Evaluation of the ADC channel corresponding to the light yield Light
520   Int_t vADCch = (Int_t) (Light * fPMGain[Det-1][Quad] * fADCRes[Res]);
521   // Ch. debug
522   //printf("\t Phe2ADCch -> det %d quad %d - PMGain[%d][%d] %1.0f phe %1.2f  ADC %d\n", 
523   //    Det,Quad,Det-1,Quad,fPMGain[Det-1][Quad],Light,vADCch);
524
525   return vADCch;
526 }
527
528 //_____________________________________________________________________________
529 Int_t AliZDCDigitizer::Pedestal(Int_t Det, Int_t Quad, Int_t Res) const
530 {
531   // Returns a pedestal for detector det, PM quad, channel with res.
532   //
533   Float_t pedValue;
534   // Normal run
535   if(fIsCalibration == 0){
536     Int_t index=0, kNch=24;
537     if(Quad!=5){
538       if(Det==1)        index = Quad+kNch*Res;    // ZNC
539       else if(Det==2)   index = (Quad+5)+kNch*Res;  // ZPC
540       else if(Det==3)   index = (Quad+9)+kNch*Res;  // ZEM
541       else if(Det==4)   index = (Quad+12)+kNch*Res; // ZNA
542       else if(Det==5)   index = (Quad+17)+kNch*Res; // ZPA
543     }
544     else index = (Det-1)/3+22+kNch*Res; // Reference PMs
545     //
546     Float_t meanPed = fPedData->GetMeanPed(index);
547     Float_t pedWidth = fPedData->GetMeanPedWidth(index);
548     pedValue = gRandom->Gaus(meanPed,pedWidth);
549     //
550     /*printf("\t  AliZDCDigitizer::Pedestal -> det %d quad %d res %d - Ped[%d] = %d\n",
551         Det, Quad, Res, index,(Int_t) pedValue); // Chiara debugging!
552     */
553   }
554   // To create calibration object
555   else{
556     if(Res == 0) pedValue = gRandom->Gaus((35.+10.*gRandom->Rndm()),(0.5+0.2*gRandom->Rndm())); //High gain
557     else  pedValue = gRandom->Gaus((250.+100.*gRandom->Rndm()),(3.5+2.*gRandom->Rndm())); //Low gain
558   }
559
560   return (Int_t) pedValue;
561 }
562
563 //_____________________________________________________________________________
564 AliCDBStorage* AliZDCDigitizer::SetStorage(const char *uri) 
565 {
566
567   Bool_t deleteManager = kFALSE;
568   
569   AliCDBManager *manager = AliCDBManager::Instance();
570   AliCDBStorage *defstorage = manager->GetDefaultStorage();
571   
572   if(!defstorage || !(defstorage->Contains("ZDC"))){ 
573      AliWarning("No default storage set or default storage doesn't contain ZDC!");
574      manager->SetDefaultStorage(uri);
575      deleteManager = kTRUE;
576   }
577  
578   AliCDBStorage *storage = manager->GetDefaultStorage();
579
580   if(deleteManager){
581     AliCDBManager::Instance()->UnsetDefaultStorage();
582     defstorage = 0;   // the storage is killed by AliCDBManager::Instance()->Destroy()
583   }
584
585   return storage; 
586 }
587
588 //_____________________________________________________________________________
589 AliZDCPedestals* AliZDCDigitizer::GetPedData() const
590 {
591
592   // Getting pedestal calibration object for ZDC set
593
594   AliCDBEntry  *entry = AliCDBManager::Instance()->Get("ZDC/Calib/Pedestals");
595   if(!entry) AliFatal("No calibration data loaded!");  
596
597   AliZDCPedestals *calibdata = dynamic_cast<AliZDCPedestals*>  (entry->GetObject());
598   if(!calibdata)  AliFatal("Wrong calibration object in calibration  file!");
599
600   return calibdata;
601 }