1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
18 ///////////////////////////////////////////////////////////////////////////////
20 // ZDC digitizer class //
22 ///////////////////////////////////////////////////////////////////////////////
32 // --- AliRoot header files
35 #include "AliHeader.h"
36 #include "AliGenHijingEventHeader.h"
37 #include "AliRunDigitizer.h"
38 #include "AliRunLoader.h"
39 #include "AliGRPObject.h"
40 #include "AliCDBManager.h"
41 #include "AliCDBEntry.h"
42 #include "AliZDCSDigit.h"
43 #include "AliZDCDigit.h"
44 #include "AliZDCFragment.h"
46 #include "AliZDCDigitizer.h"
49 class AliZDCPedestals;
51 class AliZDCTowerCalib;
53 ClassImp(AliZDCDigitizer)
56 //____________________________________________________________________________
57 AliZDCDigitizer::AliZDCDigitizer() :
59 fIsSignalInADCGate(kFALSE),
64 fSpectators2Track(kFALSE)
66 // Default constructor
70 //____________________________________________________________________________
71 AliZDCDigitizer::AliZDCDigitizer(AliRunDigitizer* manager):
72 AliDigitizer(manager),
73 fIsCalibration(0), //By default the simulation doesn't create calib. data
74 fIsSignalInADCGate(kFALSE),
76 fPedData(GetPedData()),
77 fEnCalibData(GetEnCalibData()),
78 fTowCalibData(GetTowCalibData()),
79 fSpectators2Track(kFALSE)
81 // Get calibration data
82 if(fIsCalibration!=0) printf("\n\t AliZDCDigitizer -> Creating calibration data (pedestals)\n");
83 for(Int_t i=0; i<6; i++){
84 for(Int_t j=0; j<5; j++)
89 //____________________________________________________________________________
90 AliZDCDigitizer::~AliZDCDigitizer()
97 //____________________________________________________________________________
98 AliZDCDigitizer::AliZDCDigitizer(const AliZDCDigitizer &digitizer):
100 fIsCalibration(digitizer.fIsCalibration),
101 fIsSignalInADCGate(digitizer.fIsSignalInADCGate),
102 fFracLostSignal(digitizer.fFracLostSignal),
103 fPedData(digitizer.fPedData),
104 fEnCalibData(digitizer.fEnCalibData),
105 fTowCalibData(digitizer.fTowCalibData),
106 fSpectators2Track(digitizer.fSpectators2Track)
110 for(Int_t i=0; i<6; i++){
111 for(Int_t j=0; j<5; j++){
112 fPMGain[i][j] = digitizer.fPMGain[i][j];
115 for(Int_t i=0; i<2; i++) fADCRes[i] = digitizer.fADCRes[i];
120 //____________________________________________________________________________
121 Bool_t AliZDCDigitizer::Init()
123 // Initialize the digitizer
125 AliCDBEntry* entry = AliCDBManager::Instance()->Get("GRP/GRP/Data");
126 AliGRPObject* grpData = 0x0;
128 TMap* m = dynamic_cast<TMap*>(entry->GetObject()); // old GRP entry
131 grpData = new AliGRPObject();
132 grpData->ReadValuesFromMap(m);
135 grpData = dynamic_cast<AliGRPObject*>(entry->GetObject()); // new GRP entry
138 AliCDBManager::Instance()->UnloadFromCache("GRP/GRP/Data");
140 if(!grpData) AliError("No GRP entry found in OCDB!");
142 TString beamType = grpData->GetBeamType();
143 if(beamType==AliGRPObject::GetInvalidString()){
144 AliError("\t UNKNOWN beam type from GRP obj -> PMT gains not set in ZDC digitizer!!!\n");
147 Float_t beamEnergy = grpData->GetBeamEnergy()/2.;
148 if(beamEnergy==AliGRPObject::GetInvalidFloat()){
149 AliWarning("GRP/GRP/Data entry: missing value for the beam energy ! Using 0.");
150 AliError("\t UNKNOWN beam type from GRP obj -> PMT gains not set in ZDC digitizer!!!\n");
153 // LHC: "multiply by 120 to get the energy in MeV"
156 if((beamType.CompareTo("P-P")) == 0){
157 //PTM gains rescaled to beam energy for p-p
159 for(Int_t j = 0; j < 5; j++){
160 fPMGain[0][j] = (661.444/beamEnergy+0.000740671)*10000000;
161 fPMGain[1][j] = (864.350/beamEnergy+0.002344)*10000000;
162 fPMGain[2][j] = (1.32312-0.000101515*beamEnergy)*10000000;
163 fPMGain[3][j] = fPMGain[0][j];
164 fPMGain[4][j] = fPMGain[1][j] ;
166 AliInfo(Form(" PMT gains for p-p @ %1.0f+%1.0f: ZN(%1.0f), ZP(%1.0f), ZEM(%1.0f)\n",
167 beamEnergy, beamEnergy, fPMGain[0][0], fPMGain[1][0], fPMGain[2][0]));
170 else if((beamType.CompareTo("A-A")) == 0){
171 // PTM gains for Pb-Pb @ 2.7_2.7 A TeV ***************
172 for(Int_t j = 0; j < 5; j++){
173 fPMGain[0][j] = 50000.;
174 fPMGain[1][j] = 100000.;
175 fPMGain[2][j] = 100000.;
176 fPMGain[3][j] = 50000.;
177 fPMGain[4][j] = 100000.;
178 fPMGain[5][j] = 100000.;
183 fADCRes[0] = 0.0000008; // ADC Resolution high gain: 200 fC/adcCh
184 fADCRes[1] = 0.0000064; // ADC Resolution low gain: 25 fC/adcCh
189 //____________________________________________________________________________
190 void AliZDCDigitizer::Exec(Option_t* /*option*/)
192 // Execute digitization
194 // ------------------------------------------------------------
195 // !!! 2nd ZDC set added
196 // *** 1st 3 arrays are digits from REAL (simulated) hits
197 // *** last 2 are copied from simulated digits
198 // --- pm[0][...] = light in ZN side C [C, Q1, Q2, Q3, Q4]
199 // --- pm[1][...] = light in ZP side C [C, Q1, Q2, Q3, Q4]
200 // --- pm[2][...] = light in ZEM [x, 1, 2, x, x]
201 // --- pm[3][...] = light in ZN side A [C, Q1, Q2, Q3, Q4] ->NEW!
202 // --- pm[4][...] = light in ZP side A [C, Q1, Q2, Q3, Q4] ->NEW!
203 // ------------------------------------------------------------
205 for(Int_t iSector1=0; iSector1<5; iSector1++)
206 for(Int_t iSector2=0; iSector2<5; iSector2++){
207 pm[iSector1][iSector2] = 0;
210 // ------------------------------------------------------------
211 // ### Out of time ADC added (22 channels)
212 // --- same codification as for signal PTMs (see above)
213 // ------------------------------------------------------------
215 for(Int_t iSector1=0; iSector1<5; iSector1++)
216 for(Int_t iSector2=0; iSector2<5; iSector2++){
217 pmoot[iSector1][iSector2] = 0;
220 // impact parameter and number of spectators
222 Int_t specNTarg = 0, specPTarg = 0;
223 Int_t specNProj = 0, specPProj = 0;
224 Float_t signalTime0 = 0.;
226 // loop over input streams
227 for(Int_t iInput = 0; iInput<fManager->GetNinputs(); iInput++){
229 // get run loader and ZDC loader
230 AliRunLoader* runLoader =
231 AliRunLoader::GetRunLoader(fManager->GetInputFolderName(iInput));
232 AliLoader* loader = runLoader->GetLoader("ZDCLoader");
233 if(!loader) continue;
236 loader->LoadSDigits();
237 TTree* treeS = loader->TreeS();
240 AliZDCSDigit* psdigit = &sdigit;
241 treeS->SetBranchAddress("ZDC", &psdigit);
244 for(Int_t iSDigit=0; iSDigit<treeS->GetEntries(); iSDigit++){
245 treeS->GetEntry(iSDigit);
247 if(!psdigit) continue;
248 if((sdigit.GetSector(1) < 0) || (sdigit.GetSector(1) > 4)){
249 AliError(Form("\nsector[0] = %d, sector[1] = %d\n",
250 sdigit.GetSector(0), sdigit.GetSector(1)));
253 // Checking if signal is inside ADC gate
254 if(iSDigit==0) signalTime0 = sdigit.GetTrackTime();
256 // Assuming a signal lenght of 20 ns, signal is in gate if
257 // signal ENDS in signalTime0+50. -> sdigit.GetTrackTime()+20<=signalTime0+50.
258 if(sdigit.GetTrackTime()<=signalTime0+30.) fIsSignalInADCGate = kTRUE;
259 if(sdigit.GetTrackTime()>signalTime0+30.){
260 fIsSignalInADCGate = kFALSE;
261 // Vedi quaderno per spiegazione approx. usata
262 // nel calcolo della fraz. di segnale perso
263 fFracLostSignal = (sdigit.GetTrackTime()-30)*(sdigit.GetTrackTime()-30)/280.;
267 Float_t sdSignal = sdigit.GetLightPM();
268 if(fIsSignalInADCGate == kFALSE){
269 AliDebug(2,Form("\t Signal time %f -> fraction %f of ZDC signal for det.(%d, %d) out of ADC gate\n",
270 sdigit.GetTrackTime(),fFracLostSignal,sdigit.GetSector(0),sdigit.GetSector(1)));
271 sdSignal = (1-fFracLostSignal)*sdSignal;
274 pm[(sdigit.GetSector(0))-1][sdigit.GetSector(1)] += sdigit.GetLightPM();
276 /*printf("\t Detector %d, Tower %d -> pm[%d][%d] = %.0f \n",
277 sdigit.GetSector(0), sdigit.GetSector(1),sdigit.GetSector(0)-1,
278 sdigit.GetSector(1), pm[sdigit.GetSector(0)-1][sdigit.GetSector(1)]);
283 loader->UnloadSDigits();
285 // get the impact parameter and the number of spectators in case of hijing
286 if(!runLoader->GetAliRun()) runLoader->LoadgAlice();
287 AliHeader* header = runLoader->GetHeader();
288 if(!header) continue;
289 AliGenEventHeader* genHeader = header->GenEventHeader();
290 if(!genHeader) continue;
291 if(!genHeader->InheritsFrom(AliGenHijingEventHeader::Class())) continue;
293 if(fSpectators2Track==kTRUE){
294 impPar = ((AliGenHijingEventHeader*) genHeader)->ImpactParameter();
295 specNProj = ((AliGenHijingEventHeader*) genHeader)->ProjSpectatorsn();
296 specPProj = ((AliGenHijingEventHeader*) genHeader)->ProjSpectatorsp();
297 specNTarg = ((AliGenHijingEventHeader*) genHeader)->TargSpectatorsn();
298 specPTarg = ((AliGenHijingEventHeader*) genHeader)->TargSpectatorsp();
299 printf("\n\t AliZDCDigitizer: b = %1.2f fm\n"
300 " \t PROJ.: #spectator n %d, #spectator p %d\n"
301 " \t TARG.: #spectator n %d, #spectator p %d\n\n",
302 impPar, specNProj, specPProj, specNTarg, specPTarg);
307 // Applying fragmentation algorithm and adding spectator signal
308 if(fSpectators2Track==kTRUE && impPar) {
309 Int_t freeSpecNProj, freeSpecPProj;
310 Fragmentation(impPar, specNProj, specPProj, freeSpecNProj, freeSpecPProj);
311 Int_t freeSpecNTarg, freeSpecPTarg;
312 Fragmentation(impPar, specNTarg, specPTarg, freeSpecNTarg, freeSpecPTarg);
313 SpectatorSignal(1, freeSpecNProj, pm);
314 //printf(" AliZDCDigitizer -> Signal for %d PROJ free spectator n",freeSpecNProj);
315 SpectatorSignal(2, freeSpecPProj, pm);
316 //printf(" and %d free spectator p added\n",freeSpecPProj);
317 SpectatorSignal(3, freeSpecNTarg, pm);
318 //printf(" AliZDCDigitizer -> Signal for %d TARG free spectator n",freeSpecNTarg);
319 SpectatorSignal(4, freeSpecPTarg, pm);
320 //printf("and %d free spectator p added\n",freeSpecPTarg);
324 // get the output run loader and loader
325 AliRunLoader* runLoader =
326 AliRunLoader::GetRunLoader(fManager->GetOutputFolderName());
327 AliLoader* loader = runLoader->GetLoader("ZDCLoader");
329 AliError("no ZDC loader found");
333 // create the output tree
334 const char* mode = "update";
335 if(runLoader->GetEventNumber() == 0) mode = "recreate";
336 loader->LoadDigits(mode);
337 loader->MakeTree("D");
338 TTree* treeD = loader->TreeD();
340 AliZDCDigit* pdigit = &digit;
341 const Int_t kBufferSize = 4000;
342 treeD->Branch("ZDC", "AliZDCDigit", &pdigit, kBufferSize);
346 Int_t digi[2], digioot[2];
347 for(sector[0]=1; sector[0]<6; sector[0]++){
348 for(sector[1]=0; sector[1]<5; sector[1]++){
349 if((sector[0]==3) && ((sector[1]<1) || (sector[1]>2))) continue;
350 for(Int_t res=0; res<2; res++){
351 digi[res] = Phe2ADCch(sector[0], sector[1], pm[sector[0]-1][sector[1]], res)
352 + Pedestal(sector[0], sector[1], res);
354 new(pdigit) AliZDCDigit(sector, digi);
358 //printf("\t DIGIT added -> det %d quad %d - digi[0,1] = [%d, %d]\n",
359 // sector[0], sector[1], digi[0], digi[1]); // Chiara debugging!
362 } // Loop over detector
363 // Adding in-time digits for 2 reference PTM signals (after signal ch.)
364 // (for the moment the ref. signal is completely invented assuming a PMgain of 5*10^4!)
368 // Reference signal are set to 100 (high gain chain) and 800 (low gain chain)
369 if(fIsCalibration==0) {sigRef[0]=100; sigRef[1]=800;}
370 else {sigRef[0]=0; sigRef[1]=0;} // calibration -> simulation of pedestal values
372 for(Int_t iref=0; iref<2; iref++){
373 sectorRef[0] = 3*iref+1;
374 for(Int_t res=0; res<2; res++){
375 sigRef[res] += Pedestal(sectorRef[0], sectorRef[1], res);
377 new(pdigit) AliZDCDigit(sectorRef, sigRef);
381 //printf("\t RefDigit added -> det = %d, quad = %d - digi[0,1] = [%d, %d]\n",
382 // sectorRef[0], sectorRef[1], sigRef[0], sigRef[1]); // Chiara debugging!
385 // --- Adding digits for out-of-time channels after signal digits
386 for(sector[0]=1; sector[0]<6; sector[0]++){
387 for(sector[1]=0; sector[1]<5; sector[1]++){
388 if((sector[0]==3) && ((sector[1]<1) || (sector[1]>2))) continue;
389 for(Int_t res=0; res<2; res++){
390 digioot[res] = Pedestal(sector[0], sector[1], res); // out-of-time ADCs
392 new(pdigit) AliZDCDigit(sector, digioot);
396 //printf("\t DIGIToot added -> det = %d, quad = %d - digi[0,1] = [%d, %d]\n",
397 // sector[0], sector[1], digioot[0], digioot[1]); // Chiara debugging!
400 // Adding out-of-time digits for 2 reference PTM signals (after out-of-time ch.)
402 for(Int_t iref=0; iref<2; iref++){
403 sectorRef[0] = 3*iref+1;
404 for(Int_t res=0; res<2; res++){
405 sigRefoot[res] = Pedestal(sectorRef[0], sectorRef[1], res);
407 new(pdigit) AliZDCDigit(sectorRef, sigRefoot);
410 //printf("\t RefDigitoot added -> det = %d, quad = %d - digi[0,1] = [%d, %d]\n",
411 // sectorRef[0], sectorRef[1], sigRefoot[0], sigRefoot[1]); // Chiara debugging!
414 //printf("\t AliZDCDigitizer -> TreeD has %d entries\n",(Int_t) treeD->GetEntries());
416 // write the output tree
417 loader->WriteDigits("OVERWRITE");
418 loader->UnloadDigits();
422 //_____________________________________________________________________________
423 void AliZDCDigitizer::Fragmentation(Float_t impPar, Int_t specN, Int_t specP,
424 Int_t &freeSpecN, Int_t &freeSpecP) const
426 // simulate fragmentation of spectators
428 AliZDCFragment frag(impPar);
430 // Fragments generation
432 Int_t nAlpha = frag.GetNalpha();
435 Int_t ztot = frag.GetZtot();
436 Int_t ntot = frag.GetNtot();
437 frag.AttachNeutrons();
438 freeSpecN = specN-ntot-2*nAlpha;
439 freeSpecP = specP-ztot-2*nAlpha;
440 // Removing deuterons
441 Int_t ndeu = (Int_t) (freeSpecN*frag.DeuteronNumber());
444 if(freeSpecN<0) freeSpecN=0;
445 if(freeSpecP<0) freeSpecP=0;
446 AliDebug(2, Form("FreeSpn = %d, FreeSpp = %d", freeSpecN, freeSpecP));
449 //_____________________________________________________________________________
450 void AliZDCDigitizer::SpectatorSignal(Int_t SpecType, Int_t numEvents,
451 Float_t pm[5][5]) const
453 // add signal of the spectators
456 if(SpecType == 1) { // --- Signal for projectile spectator neutrons
457 hfn = "$ALICE_ROOT/ZDC/ZNCSignal.root";
459 else if(SpecType == 2) { // --- Signal for projectile spectator protons
460 hfn = "$ALICE_ROOT/ZDC/ZPCSignal.root";
462 else if(SpecType == 3) { // --- Signal for target spectator neutrons
463 hfn = "$ALICE_ROOT/ZDC/ZNASignal.root";
465 else if(SpecType == 4) { // --- Signal for target spectator protons
466 hfn = "$ALICE_ROOT/ZDC/ZPASignal.root";
469 TFile* file = TFile::Open(hfn);
470 if(!file || !file->IsOpen()) {
471 AliError((Form(" Opening file %s failed\n",hfn.Data())));
475 TNtuple* zdcSignal = (TNtuple*) file->Get("ZDCSignal");
476 Int_t nentries = (Int_t) zdcSignal->GetEntries();
479 Int_t pl, i, k, iev=0, rnd[125], volume[2];
480 for(pl=0;pl<125;pl++) rnd[pl] = 0;
481 if(numEvents > 125) {
482 AliDebug(2,Form("numEvents (%d) is larger than 125", numEvents));
485 for(pl=0;pl<numEvents;pl++){
486 rnd[pl] = (Int_t) (9999*gRandom->Rndm());
487 if(rnd[pl] >= 9999) rnd[pl] = 9998;
488 //printf(" rnd[%d] = %d\n",pl,rnd[pl]);
490 // Sorting vector in ascending order with C function QSORT
491 qsort((void*)rnd,numEvents,sizeof(Int_t),comp);
493 for(i=0; i<nentries; i++){
494 zdcSignal->GetEvent(i);
495 entry = zdcSignal->GetArgs();
496 if(entry[0] == rnd[iev]){
497 for(k=0; k<2; k++) volume[k] = (Int_t) entry[k+1];
499 Float_t lightQ = entry[7];
500 Float_t lightC = entry[8];
502 if(volume[0] != 3) { // ZN or ZP
503 pm[volume[0]-1][0] += lightC;
504 pm[volume[0]-1][volume[1]] += lightQ;
505 //printf("\n pm[%d][0] = %.0f, pm[%d][%d] = %.0f\n",(volume[0]-1),pm[volume[0]-1][0],
506 // (volume[0]-1),volume[1],pm[volume[0]-1][volume[1]]);
509 if(volume[1] == 1) pm[2][1] += lightC; // ZEM 1
510 else pm[2][2] += lightQ; // ZEM 2
511 //printf("\n pm[2][1] = %.0f, pm[2][2] = %.0f\n",pm[2][1],pm[2][2]);
514 else if(entry[0] > rnd[iev]){
519 }while(iev<numEvents);
526 //_____________________________________________________________________________
527 Int_t AliZDCDigitizer::Phe2ADCch(Int_t Det, Int_t Quad, Float_t Light,
530 // Evaluation of the ADC channel corresponding to the light yield Light
531 Int_t vADCch = (Int_t) (Light * fPMGain[Det-1][Quad] * fADCRes[Res]);
533 //printf("\t Phe2ADCch -> det %d quad %d - PMGain[%d][%d] %1.0f phe %1.2f ADC %d\n",
534 // Det,Quad,Det-1,Quad,fPMGain[Det-1][Quad],Light,vADCch);
539 //_____________________________________________________________________________
540 Int_t AliZDCDigitizer::Pedestal(Int_t Det, Int_t Quad, Int_t Res) const
542 // Returns a pedestal for detector det, PM quad, channel with res.
546 if(fIsCalibration == 0){
547 Int_t index=0, kNch=24;
549 if(Det==1) index = Quad+kNch*Res; // ZNC
550 else if(Det==2) index = (Quad+5)+kNch*Res; // ZPC
551 else if(Det==3) index = (Quad+9)+kNch*Res; // ZEM
552 else if(Det==4) index = (Quad+12)+kNch*Res; // ZNA
553 else if(Det==5) index = (Quad+17)+kNch*Res; // ZPA
555 else index = (Det-1)/3+22+kNch*Res; // Reference PMs
557 Float_t meanPed = fPedData->GetMeanPed(index);
558 Float_t pedWidth = fPedData->GetMeanPedWidth(index);
559 pedValue = gRandom->Gaus(meanPed,pedWidth);
561 /*printf("\t AliZDCDigitizer::Pedestal -> det %d quad %d res %d - Ped[%d] = %d\n",
562 Det, Quad, Res, index,(Int_t) pedValue); // Chiara debugging!
565 // To create calibration object
567 if(Res == 0) pedValue = gRandom->Gaus((35.+10.*gRandom->Rndm()),(0.5+0.2*gRandom->Rndm())); //High gain
568 else pedValue = gRandom->Gaus((250.+100.*gRandom->Rndm()),(3.5+2.*gRandom->Rndm())); //Low gain
571 return (Int_t) pedValue;
574 //_____________________________________________________________________________
575 AliCDBStorage* AliZDCDigitizer::SetStorage(const char *uri)
578 Bool_t deleteManager = kFALSE;
580 AliCDBManager *manager = AliCDBManager::Instance();
581 AliCDBStorage *defstorage = manager->GetDefaultStorage();
583 if(!defstorage || !(defstorage->Contains("ZDC"))){
584 AliWarning("No default storage set or default storage doesn't contain ZDC!");
585 manager->SetDefaultStorage(uri);
586 deleteManager = kTRUE;
589 AliCDBStorage *storage = manager->GetDefaultStorage();
592 AliCDBManager::Instance()->UnsetDefaultStorage();
593 defstorage = 0; // the storage is killed by AliCDBManager::Instance()->Destroy()
599 //_____________________________________________________________________________
600 AliZDCPedestals* AliZDCDigitizer::GetPedData() const
603 // Getting pedestal calibration object for ZDC set
605 AliCDBEntry *entry = AliCDBManager::Instance()->Get("ZDC/Calib/Pedestals");
606 if(!entry) AliFatal("No calibration data loaded!");
608 AliZDCPedestals *calibdata = dynamic_cast<AliZDCPedestals*> (entry->GetObject());
609 if(!calibdata) AliFatal("Wrong calibration object in calibration file!");
614 //_____________________________________________________________________________
615 AliZDCEnCalib* AliZDCDigitizer::GetEnCalibData() const
618 // Getting calibration object for ZDC set
620 AliCDBEntry *entry = AliCDBManager::Instance()->Get("ZDC/Calib/EnergyCalib");
621 if(!entry) AliFatal("No calibration data loaded!");
623 AliZDCEnCalib *calibdata = dynamic_cast<AliZDCEnCalib*> (entry->GetObject());
624 if(!calibdata) AliFatal("Wrong calibration object in calibration file!");
629 //_____________________________________________________________________________
630 AliZDCTowerCalib* AliZDCDigitizer::GetTowCalibData() const
633 // Getting calibration object for ZDC set
635 AliCDBEntry *entry = AliCDBManager::Instance()->Get("ZDC/Calib/TowerCalib");
636 if(!entry) AliFatal("No calibration data loaded!");
638 AliZDCTowerCalib *calibdata = dynamic_cast<AliZDCTowerCalib*> (entry->GetObject());
639 if(!calibdata) AliFatal("Wrong calibration object in calibration file!");