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