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