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