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