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