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");
154 if((beamType.CompareTo("P-P")) == 0){
155 //PTM gains rescaled to beam energy for p-p
157 for(Int_t j = 0; j < 5; j++){
158 fPMGain[0][j] = (661.444/beamEnergy+0.000740671)*10000000;
159 fPMGain[1][j] = (864.350/beamEnergy+0.002344)*10000000;
160 fPMGain[2][j] = (1.32312-0.000101515*beamEnergy)*10000000;
161 fPMGain[3][j] = fPMGain[0][j];
162 fPMGain[4][j] = fPMGain[1][j] ;
164 AliInfo(Form(" PMT gains for p-p @ %1.0f+%1.0f: ZN(%1.0f), ZP(%1.0f), ZEM(%1.0f)\n",
165 beamEnergy, beamEnergy, fPMGain[0][0], fPMGain[1][0], fPMGain[2][0]));
168 else if((beamType.CompareTo("A-A")) == 0){
169 // PTM gains for Pb-Pb @ 2.7_2.7 A TeV ***************
170 for(Int_t j = 0; j < 5; j++){
171 fPMGain[0][j] = 50000.;
172 fPMGain[1][j] = 100000.;
173 fPMGain[2][j] = 100000.;
174 fPMGain[3][j] = 50000.;
175 fPMGain[4][j] = 100000.;
176 fPMGain[5][j] = 100000.;
181 fADCRes[0] = 0.0000008; // ADC Resolution high gain: 200 fC/adcCh
182 fADCRes[1] = 0.0000064; // ADC Resolution low gain: 25 fC/adcCh
187 //____________________________________________________________________________
188 void AliZDCDigitizer::Exec(Option_t* /*option*/)
190 // Execute digitization
192 // ------------------------------------------------------------
193 // !!! 2nd ZDC set added
194 // *** 1st 3 arrays are digits from REAL (simulated) hits
195 // *** last 2 are copied from simulated digits
196 // --- pm[0][...] = light in ZN side C [C, Q1, Q2, Q3, Q4]
197 // --- pm[1][...] = light in ZP side C [C, Q1, Q2, Q3, Q4]
198 // --- pm[2][...] = light in ZEM [x, 1, 2, x, x]
199 // --- pm[3][...] = light in ZN side A [C, Q1, Q2, Q3, Q4] ->NEW!
200 // --- pm[4][...] = light in ZP side A [C, Q1, Q2, Q3, Q4] ->NEW!
201 // ------------------------------------------------------------
203 for(Int_t iSector1=0; iSector1<5; iSector1++)
204 for(Int_t iSector2=0; iSector2<5; iSector2++){
205 pm[iSector1][iSector2] = 0;
208 // ------------------------------------------------------------
209 // ### Out of time ADC added (22 channels)
210 // --- same codification as for signal PTMs (see above)
211 // ------------------------------------------------------------
213 for(Int_t iSector1=0; iSector1<5; iSector1++)
214 for(Int_t iSector2=0; iSector2<5; iSector2++){
215 pmoot[iSector1][iSector2] = 0;
218 // impact parameter and number of spectators
220 Int_t specNTarg = 0, specPTarg = 0;
221 Int_t specNProj = 0, specPProj = 0;
222 Float_t signalTime0 = 0.;
224 // loop over input streams
225 for(Int_t iInput = 0; iInput<fManager->GetNinputs(); iInput++){
227 // get run loader and ZDC loader
228 AliRunLoader* runLoader =
229 AliRunLoader::GetRunLoader(fManager->GetInputFolderName(iInput));
230 AliLoader* loader = runLoader->GetLoader("ZDCLoader");
231 if(!loader) continue;
234 loader->LoadSDigits();
235 TTree* treeS = loader->TreeS();
238 AliZDCSDigit* psdigit = &sdigit;
239 treeS->SetBranchAddress("ZDC", &psdigit);
242 for(Int_t iSDigit=0; iSDigit<treeS->GetEntries(); iSDigit++){
243 treeS->GetEntry(iSDigit);
245 if(!psdigit) continue;
246 if((sdigit.GetSector(1) < 0) || (sdigit.GetSector(1) > 4)){
247 AliError(Form("\nsector[0] = %d, sector[1] = %d\n",
248 sdigit.GetSector(0), sdigit.GetSector(1)));
251 // Checking if signal is inside ADC gate
252 if(iSDigit==0) signalTime0 = sdigit.GetTrackTime();
254 // Assuming a signal lenght of 20 ns, signal is in gate if
255 // signal ENDS in signalTime0+50. -> sdigit.GetTrackTime()+20<=signalTime0+50.
256 if(sdigit.GetTrackTime()<=signalTime0+30.) fIsSignalInADCGate = kTRUE;
257 if(sdigit.GetTrackTime()>signalTime0+30.){
258 fIsSignalInADCGate = kFALSE;
259 // Vedi quaderno per spiegazione approx. usata
260 // nel calcolo della fraz. di segnale perso
261 fFracLostSignal = (sdigit.GetTrackTime()-30)*(sdigit.GetTrackTime()-30)/280.;
265 Float_t sdSignal = sdigit.GetLightPM();
266 if(fIsSignalInADCGate == kFALSE){
267 AliDebug(2,Form("\t Signal time %f -> fraction %f of ZDC signal for det.(%d, %d) out of ADC gate\n",
268 sdigit.GetTrackTime(),fFracLostSignal,sdigit.GetSector(0),sdigit.GetSector(1)));
269 sdSignal = (1-fFracLostSignal)*sdSignal;
272 pm[(sdigit.GetSector(0))-1][sdigit.GetSector(1)] += sdigit.GetLightPM();
274 /*printf("\t Detector %d, Tower %d -> pm[%d][%d] = %.0f \n",
275 sdigit.GetSector(0), sdigit.GetSector(1),sdigit.GetSector(0)-1,
276 sdigit.GetSector(1), pm[sdigit.GetSector(0)-1][sdigit.GetSector(1)]);
281 loader->UnloadSDigits();
283 // get the impact parameter and the number of spectators in case of hijing
284 if(!runLoader->GetAliRun()) runLoader->LoadgAlice();
285 AliHeader* header = runLoader->GetHeader();
286 if(!header) continue;
287 AliGenEventHeader* genHeader = header->GenEventHeader();
288 if(!genHeader) continue;
289 if(!genHeader->InheritsFrom(AliGenHijingEventHeader::Class())) continue;
291 if(fSpectators2Track==kTRUE){
292 impPar = ((AliGenHijingEventHeader*) genHeader)->ImpactParameter();
293 specNProj = ((AliGenHijingEventHeader*) genHeader)->ProjSpectatorsn();
294 specPProj = ((AliGenHijingEventHeader*) genHeader)->ProjSpectatorsp();
295 specNTarg = ((AliGenHijingEventHeader*) genHeader)->TargSpectatorsn();
296 specPTarg = ((AliGenHijingEventHeader*) genHeader)->TargSpectatorsp();
297 printf("\n\t AliZDCDigitizer: b = %1.2f fm\n"
298 " \t PROJ.: #spectator n %d, #spectator p %d\n"
299 " \t TARG.: #spectator n %d, #spectator p %d\n\n",
300 impPar, specNProj, specPProj, specNTarg, specPTarg);
305 // Applying fragmentation algorithm and adding spectator signal
306 if(fSpectators2Track==kTRUE && impPar) {
307 Int_t freeSpecNProj, freeSpecPProj;
308 Fragmentation(impPar, specNProj, specPProj, freeSpecNProj, freeSpecPProj);
309 Int_t freeSpecNTarg, freeSpecPTarg;
310 Fragmentation(impPar, specNTarg, specPTarg, freeSpecNTarg, freeSpecPTarg);
311 SpectatorSignal(1, freeSpecNProj, pm);
312 //printf(" AliZDCDigitizer -> Signal for %d PROJ free spectator n",freeSpecNProj);
313 SpectatorSignal(2, freeSpecPProj, pm);
314 //printf(" and %d free spectator p added\n",freeSpecPProj);
315 SpectatorSignal(3, freeSpecNTarg, pm);
316 //printf(" AliZDCDigitizer -> Signal for %d TARG free spectator n",freeSpecNTarg);
317 SpectatorSignal(4, freeSpecPTarg, pm);
318 //printf("and %d free spectator p added\n",freeSpecPTarg);
322 // get the output run loader and loader
323 AliRunLoader* runLoader =
324 AliRunLoader::GetRunLoader(fManager->GetOutputFolderName());
325 AliLoader* loader = runLoader->GetLoader("ZDCLoader");
327 AliError("no ZDC loader found");
331 // create the output tree
332 const char* mode = "update";
333 if(runLoader->GetEventNumber() == 0) mode = "recreate";
334 loader->LoadDigits(mode);
335 loader->MakeTree("D");
336 TTree* treeD = loader->TreeD();
338 AliZDCDigit* pdigit = &digit;
339 const Int_t kBufferSize = 4000;
340 treeD->Branch("ZDC", "AliZDCDigit", &pdigit, kBufferSize);
344 Int_t digi[2], digioot[2];
345 for(sector[0]=1; sector[0]<6; sector[0]++){
346 for(sector[1]=0; sector[1]<5; sector[1]++){
347 if((sector[0]==3) && ((sector[1]<1) || (sector[1]>2))) continue;
348 for(Int_t res=0; res<2; res++){
349 digi[res] = Phe2ADCch(sector[0], sector[1], pm[sector[0]-1][sector[1]], res)
350 + Pedestal(sector[0], sector[1], res);
352 new(pdigit) AliZDCDigit(sector, digi);
356 //printf("\t DIGIT added -> det %d quad %d - digi[0,1] = [%d, %d]\n",
357 // sector[0], sector[1], digi[0], digi[1]); // Chiara debugging!
360 } // Loop over detector
361 // Adding in-time digits for 2 reference PTM signals (after signal ch.)
362 // (for the moment the ref. signal is completely invented assuming a PMgain of 5*10^4!)
366 // Reference signal are set to 100 (high gain chain) and 800 (low gain chain)
367 if(fIsCalibration==0) {sigRef[0]=100; sigRef[1]=800;}
368 else {sigRef[0]=0; sigRef[1]=0;} // calibration -> simulation of pedestal values
370 for(Int_t iref=0; iref<2; iref++){
371 sectorRef[0] = 3*iref+1;
372 for(Int_t res=0; res<2; res++){
373 sigRef[res] += Pedestal(sectorRef[0], sectorRef[1], res);
375 new(pdigit) AliZDCDigit(sectorRef, sigRef);
379 //printf("\t RefDigit added -> det = %d, quad = %d - digi[0,1] = [%d, %d]\n",
380 // sectorRef[0], sectorRef[1], sigRef[0], sigRef[1]); // Chiara debugging!
383 // --- Adding digits for out-of-time channels after signal digits
384 for(sector[0]=1; sector[0]<6; sector[0]++){
385 for(sector[1]=0; sector[1]<5; sector[1]++){
386 if((sector[0]==3) && ((sector[1]<1) || (sector[1]>2))) continue;
387 for(Int_t res=0; res<2; res++){
388 digioot[res] = Pedestal(sector[0], sector[1], res); // out-of-time ADCs
390 new(pdigit) AliZDCDigit(sector, digioot);
394 //printf("\t DIGIToot added -> det = %d, quad = %d - digi[0,1] = [%d, %d]\n",
395 // sector[0], sector[1], digioot[0], digioot[1]); // Chiara debugging!
398 // Adding out-of-time digits for 2 reference PTM signals (after out-of-time ch.)
400 for(Int_t iref=0; iref<2; iref++){
401 sectorRef[0] = 3*iref+1;
402 for(Int_t res=0; res<2; res++){
403 sigRefoot[res] = Pedestal(sectorRef[0], sectorRef[1], res);
405 new(pdigit) AliZDCDigit(sectorRef, sigRefoot);
408 //printf("\t RefDigitoot added -> det = %d, quad = %d - digi[0,1] = [%d, %d]\n",
409 // sectorRef[0], sectorRef[1], sigRefoot[0], sigRefoot[1]); // Chiara debugging!
412 //printf("\t AliZDCDigitizer -> TreeD has %d entries\n",(Int_t) treeD->GetEntries());
414 // write the output tree
415 loader->WriteDigits("OVERWRITE");
416 loader->UnloadDigits();
420 //_____________________________________________________________________________
421 void AliZDCDigitizer::Fragmentation(Float_t impPar, Int_t specN, Int_t specP,
422 Int_t &freeSpecN, Int_t &freeSpecP) const
424 // simulate fragmentation of spectators
426 AliZDCFragment frag(impPar);
428 // Fragments generation
430 Int_t nAlpha = frag.GetNalpha();
433 Int_t ztot = frag.GetZtot();
434 Int_t ntot = frag.GetNtot();
435 frag.AttachNeutrons();
436 freeSpecN = specN-ntot-2*nAlpha;
437 freeSpecP = specP-ztot-2*nAlpha;
438 // Removing deuterons
439 Int_t ndeu = (Int_t) (freeSpecN*frag.DeuteronNumber());
442 if(freeSpecN<0) freeSpecN=0;
443 if(freeSpecP<0) freeSpecP=0;
444 AliDebug(2, Form("FreeSpn = %d, FreeSpp = %d", freeSpecN, freeSpecP));
447 //_____________________________________________________________________________
448 void AliZDCDigitizer::SpectatorSignal(Int_t SpecType, Int_t numEvents,
449 Float_t pm[5][5]) const
451 // add signal of the spectators
454 if(SpecType == 1) { // --- Signal for projectile spectator neutrons
455 hfn = "$ALICE_ROOT/ZDC/ZNCSignal.root";
457 else if(SpecType == 2) { // --- Signal for projectile spectator protons
458 hfn = "$ALICE_ROOT/ZDC/ZPCSignal.root";
460 else if(SpecType == 3) { // --- Signal for target spectator neutrons
461 hfn = "$ALICE_ROOT/ZDC/ZNASignal.root";
463 else if(SpecType == 4) { // --- Signal for target spectator protons
464 hfn = "$ALICE_ROOT/ZDC/ZPASignal.root";
467 TFile* file = TFile::Open(hfn);
468 if(!file || !file->IsOpen()) {
469 AliError((Form(" Opening file %s failed\n",hfn.Data())));
473 TNtuple* zdcSignal = (TNtuple*) file->Get("ZDCSignal");
474 Int_t nentries = (Int_t) zdcSignal->GetEntries();
477 Int_t pl, i, k, iev=0, rnd[125], volume[2];
478 for(pl=0;pl<125;pl++) rnd[pl] = 0;
479 if(numEvents > 125) {
480 AliDebug(2,Form("numEvents (%d) is larger than 125", numEvents));
483 for(pl=0;pl<numEvents;pl++){
484 rnd[pl] = (Int_t) (9999*gRandom->Rndm());
485 if(rnd[pl] >= 9999) rnd[pl] = 9998;
486 //printf(" rnd[%d] = %d\n",pl,rnd[pl]);
488 // Sorting vector in ascending order with C function QSORT
489 qsort((void*)rnd,numEvents,sizeof(Int_t),comp);
491 for(i=0; i<nentries; i++){
492 zdcSignal->GetEvent(i);
493 entry = zdcSignal->GetArgs();
494 if(entry[0] == rnd[iev]){
495 for(k=0; k<2; k++) volume[k] = (Int_t) entry[k+1];
497 Float_t lightQ = entry[7];
498 Float_t lightC = entry[8];
500 if(volume[0] != 3) { // ZN or ZP
501 pm[volume[0]-1][0] += lightC;
502 pm[volume[0]-1][volume[1]] += lightQ;
503 //printf("\n pm[%d][0] = %.0f, pm[%d][%d] = %.0f\n",(volume[0]-1),pm[volume[0]-1][0],
504 // (volume[0]-1),volume[1],pm[volume[0]-1][volume[1]]);
507 if(volume[1] == 1) pm[2][1] += lightC; // ZEM 1
508 else pm[2][2] += lightQ; // ZEM 2
509 //printf("\n pm[2][1] = %.0f, pm[2][2] = %.0f\n",pm[2][1],pm[2][2]);
512 else if(entry[0] > rnd[iev]){
517 }while(iev<numEvents);
524 //_____________________________________________________________________________
525 Int_t AliZDCDigitizer::Phe2ADCch(Int_t Det, Int_t Quad, Float_t Light,
528 // Evaluation of the ADC channel corresponding to the light yield Light
529 Int_t vADCch = (Int_t) (Light * fPMGain[Det-1][Quad] * fADCRes[Res]);
531 //printf("\t Phe2ADCch -> det %d quad %d - PMGain[%d][%d] %1.0f phe %1.2f ADC %d\n",
532 // Det,Quad,Det-1,Quad,fPMGain[Det-1][Quad],Light,vADCch);
537 //_____________________________________________________________________________
538 Int_t AliZDCDigitizer::Pedestal(Int_t Det, Int_t Quad, Int_t Res) const
540 // Returns a pedestal for detector det, PM quad, channel with res.
544 if(fIsCalibration == 0){
545 Int_t index=0, kNch=24;
547 if(Det==1) index = Quad+kNch*Res; // ZNC
548 else if(Det==2) index = (Quad+5)+kNch*Res; // ZPC
549 else if(Det==3) index = (Quad+9)+kNch*Res; // ZEM
550 else if(Det==4) index = (Quad+12)+kNch*Res; // ZNA
551 else if(Det==5) index = (Quad+17)+kNch*Res; // ZPA
553 else index = (Det-1)/3+22+kNch*Res; // Reference PMs
555 Float_t meanPed = fPedData->GetMeanPed(index);
556 Float_t pedWidth = fPedData->GetMeanPedWidth(index);
557 pedValue = gRandom->Gaus(meanPed,pedWidth);
559 /*printf("\t AliZDCDigitizer::Pedestal -> det %d quad %d res %d - Ped[%d] = %d\n",
560 Det, Quad, Res, index,(Int_t) pedValue); // Chiara debugging!
563 // To create calibration object
565 if(Res == 0) pedValue = gRandom->Gaus((35.+10.*gRandom->Rndm()),(0.5+0.2*gRandom->Rndm())); //High gain
566 else pedValue = gRandom->Gaus((250.+100.*gRandom->Rndm()),(3.5+2.*gRandom->Rndm())); //Low gain
569 return (Int_t) pedValue;
572 //_____________________________________________________________________________
573 AliCDBStorage* AliZDCDigitizer::SetStorage(const char *uri)
576 Bool_t deleteManager = kFALSE;
578 AliCDBManager *manager = AliCDBManager::Instance();
579 AliCDBStorage *defstorage = manager->GetDefaultStorage();
581 if(!defstorage || !(defstorage->Contains("ZDC"))){
582 AliWarning("No default storage set or default storage doesn't contain ZDC!");
583 manager->SetDefaultStorage(uri);
584 deleteManager = kTRUE;
587 AliCDBStorage *storage = manager->GetDefaultStorage();
590 AliCDBManager::Instance()->UnsetDefaultStorage();
591 defstorage = 0; // the storage is killed by AliCDBManager::Instance()->Destroy()
597 //_____________________________________________________________________________
598 AliZDCPedestals* AliZDCDigitizer::GetPedData() const
601 // Getting pedestal calibration object for ZDC set
603 AliCDBEntry *entry = AliCDBManager::Instance()->Get("ZDC/Calib/Pedestals");
604 if(!entry) AliFatal("No calibration data loaded!");
606 AliZDCPedestals *calibdata = dynamic_cast<AliZDCPedestals*> (entry->GetObject());
607 if(!calibdata) AliFatal("Wrong calibration object in calibration file!");
612 //_____________________________________________________________________________
613 AliZDCEnCalib* AliZDCDigitizer::GetEnCalibData() const
616 // Getting calibration object for ZDC set
618 AliCDBEntry *entry = AliCDBManager::Instance()->Get("ZDC/Calib/EnergyCalib");
619 if(!entry) AliFatal("No calibration data loaded!");
621 AliZDCEnCalib *calibdata = dynamic_cast<AliZDCEnCalib*> (entry->GetObject());
622 if(!calibdata) AliFatal("Wrong calibration object in calibration file!");
627 //_____________________________________________________________________________
628 AliZDCTowerCalib* AliZDCDigitizer::GetTowCalibData() const
631 // Getting calibration object for ZDC set
633 AliCDBEntry *entry = AliCDBManager::Instance()->Get("ZDC/Calib/TowerCalib");
634 if(!entry) AliFatal("No calibration data loaded!");
636 AliZDCTowerCalib *calibdata = dynamic_cast<AliZDCTowerCalib*> (entry->GetObject());
637 if(!calibdata) AliFatal("Wrong calibration object in calibration file!");