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