]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ZDC/AliZDC.cxx
New files, editor for QuadSet gluing together translation and palette sub-editors.
[u/mrichter/AliRoot.git] / ZDC / AliZDC.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 //                      Zero Degree Calorimeter                              //
21 //           This class contains the basic functions for the ZDCs;           //
22 //            functions specific to one particular geometry are              //
23 //                      contained in the derived classes                     //
24 //                                                                           //
25 ///////////////////////////////////////////////////////////////////////////////
26
27 // --- ROOT system
28 #include <TBRIK.h>
29 #include <TGeometry.h>
30 #include <TNode.h>
31 #include <TTree.h>
32 #include <TFile.h>
33 #include <TSystem.h>
34
35 // --- AliRoot header files
36 #include "AliDetector.h"
37 #include "AliZDC.h"
38 #include "AliZDCHit.h"
39 #include "AliZDCSDigit.h"
40 #include "AliZDCDigit.h"
41 #include "AliZDCDigitizer.h"
42 #include "AliZDCRawStream.h"
43 #include "AliZDCCalibData.h"
44
45 #include "AliRawDataHeader.h"
46 #include "AliLoader.h"
47 #include "AliRun.h"
48 #include "AliMC.h"
49 #include "AliLog.h"
50 #include "AliDAQ.h"
51  
52 ClassImp(AliZDC)
53
54 AliZDC *gAliZDC;
55  
56 //_____________________________________________________________________________
57 AliZDC::AliZDC() :
58   AliDetector(),
59   fNoShower  (0),
60   fCalibData (0)
61
62 {
63   //
64   // Default constructor for the Zero Degree Calorimeter base class
65   //
66   
67   fIshunt = 1;
68   fNhits  = 0;
69   fHits = 0;
70   fDigits = 0;
71   fNdigits = 0;
72   
73 }
74  
75 //_____________________________________________________________________________
76 AliZDC::AliZDC(const char *name, const char *title) : 
77   AliDetector(name,title),
78   fNoShower  (0),
79   fCalibData (0)
80                  
81 {
82   //
83   // Standard constructor for the Zero Degree Calorimeter base class
84   //
85     
86   fIshunt = 1;
87   fNhits  = 0;
88   fDigits = 0;
89   fNdigits = 0;
90  
91   fHits = new TClonesArray("AliZDCHit",1000);
92   gAlice->GetMCApp()->AddHitList(fHits);
93   
94   char sensname[5],senstitle[25];
95   sprintf(sensname,"ZDC");
96   sprintf(senstitle,"ZDC dummy");
97   SetName(sensname); SetTitle(senstitle);
98
99   gAliZDC = this;
100
101 }
102
103 //____________________________________________________________________________ 
104 AliZDC::~AliZDC()
105 {
106   //
107   // ZDC destructor
108   //
109
110   fIshunt = 0;
111   gAliZDC = 0;
112
113   delete fCalibData;
114
115 }
116
117 //_____________________________________________________________________________
118 AliZDC::AliZDC(const AliZDC& ZDC) :
119   AliDetector("ZDC","ZDC")
120 {
121   // copy constructor
122     fNoShower = ZDC.fNoShower;
123     fCalibData = ZDC.fCalibData;
124     fZDCCalibFName = ZDC.fZDCCalibFName;
125 }
126
127 //_____________________________________________________________________________
128 AliZDC& AliZDC::operator=(const AliZDC& ZDC)
129 {
130   // assignement operator
131   if(this!=&ZDC){
132     fNoShower = ZDC.fNoShower;
133     fCalibData = ZDC.fCalibData;
134     fZDCCalibFName = ZDC.fZDCCalibFName;
135   } return *this;
136 }
137
138 //_____________________________________________________________________________
139 void AliZDC::AddHit(Int_t track, Int_t *vol, Float_t *hits)
140 {
141   //
142   //            Add a ZDC hit to the hit list.
143   // -> We make use of 2 array of hits:
144   // [1]  fHits (the usual one) that contains hits for each PRIMARY
145   // [2]  fStHits that contains hits for each EVENT and is used to
146   //      obtain digits at the end of each event
147   //
148   
149   static Float_t primKinEn, xImpact, yImpact, sFlag;
150
151   AliZDCHit *newquad, *curprimquad;
152   newquad = new AliZDCHit(fIshunt, track, vol, hits);
153   TClonesArray &lhits = *fHits;
154   
155   if(fNhits==0){
156       // First hit -> setting flag for primary or secondary particle
157       Int_t primary = gAlice->GetMCApp()->GetPrimary(track);     
158       if(track != primary){
159         newquad->SetSFlag(1);  // SECONDARY particle entering the ZDC
160       }
161       else if(track == primary){
162         newquad->SetSFlag(0);  // PRIMARY particle entering the ZDC
163       }  
164       sFlag     = newquad->GetSFlag();
165       primKinEn = newquad->GetPrimKinEn();
166       xImpact   = newquad->GetXImpact();
167       yImpact   = newquad->GetYImpact();
168    }
169    else{       
170       newquad->SetPrimKinEn(primKinEn);
171       newquad->SetXImpact(xImpact);
172       newquad->SetYImpact(yImpact);
173       newquad->SetSFlag(sFlag);
174    }
175  
176   Int_t j;
177   for(j=0; j<fNhits; j++){
178     // If hits are equal (same track, same volume), sum them.
179      curprimquad = (AliZDCHit*) lhits[j];
180      if(*curprimquad == *newquad){
181         *curprimquad = *curprimquad+*newquad;
182         // CH. debug
183         /*if(newquad->GetEnergy() != 0. || newquad->GetLightPMC() != 0. || 
184            newquad->GetLightPMQ() != 0.){
185           printf("\n\t --- Equal hits found\n");
186           curprimquad->Print("");
187           newquad->Print("");
188           printf("\t --- Det. %d, Quad. %d: X = %f, E = %f, LightPMC = %f, LightPMQ = %f\n",
189           curprimquad->GetVolume(0),curprimquad->GetVolume(1),curprimquad->GetXImpact(),
190           curprimquad->GetEnergy(), curprimquad->GetLightPMC(), curprimquad->GetLightPMQ());
191         }*/
192         //
193         delete newquad;
194         return;
195      } 
196   }
197
198     //Otherwise create a new hit
199     new(lhits[fNhits]) AliZDCHit(*newquad);
200     fNhits++;
201     // CH. debug
202     /*printf("\n\t New ZDC hit added! fNhits = %d\n", fNhits);
203     printf("\t Det. %d, Quad.t %d: X = %f, E = %f, LightPMC = %f, LightPMQ = %f\n",
204     newquad->GetVolume(0),newquad->GetVolume(1),newquad->GetXImpact(),
205     newquad->GetEnergy(), newquad->GetLightPMC(), newquad->GetLightPMQ());
206     */
207     delete newquad;
208 }
209
210 //_____________________________________________________________________________
211 void AliZDC::BuildGeometry()
212 {
213   //
214   // Build the ROOT TNode geometry for event display 
215   // in the Zero Degree Calorimeter
216   // This routine is dummy for the moment
217   //
218
219   TNode *node, *top;
220   TBRIK *brik;
221   const int kColorZDC  = kBlue;
222   
223   //
224   top=gAlice->GetGeometry()->GetNode("alice");
225   
226   // ZDC
227     brik = new TBRIK("S_ZDC","ZDC box","void",300,300,5);
228     top->cd();
229     node = new TNode("ZDC","ZDC","S_ZDC",0,0,600,"");
230     node->SetLineColor(kColorZDC);
231     fNodes->Add(node);
232 }
233
234 //____________________________________________________________________________
235 Float_t AliZDC::ZMin(void) const
236 {
237   // Minimum dimension of the ZDC module in z
238   return -11600.;
239 }
240
241 //____________________________________________________________________________
242 Float_t AliZDC::ZMax(void) const
243 {
244   // Maximum dimension of the ZDC module in z
245   return -11750.;
246 }
247   
248
249 //_____________________________________________________________________________
250 void AliZDC::MakeBranch(Option_t *opt)
251 {
252   //
253   // Create Tree branches for the ZDC
254   //
255
256   char branchname[10];
257   sprintf(branchname,"%s",GetName());
258
259   const char *cH = strstr(opt,"H");
260   
261   if (cH && fLoader->TreeH())
262    fHits   = new TClonesArray("AliZDCHit",1000); 
263   
264   AliDetector::MakeBranch(opt);
265 }
266
267 //_____________________________________________________________________________
268 void AliZDC::Hits2SDigits()
269 {
270   // Create summable digits from hits
271   
272   AliDebug(1,"\n        Entering AliZDC::Hits2Digits() ");
273   
274   fLoader->LoadHits("read");
275   fLoader->LoadSDigits("recreate");
276   AliRunLoader* runLoader = fLoader->GetRunLoader();
277   AliZDCSDigit sdigit;
278   AliZDCSDigit* psdigit = &sdigit;
279
280   // Event loop
281   for (Int_t iEvent = 0; iEvent < runLoader->GetNumberOfEvents(); iEvent++) {
282     Float_t pmCZN = 0, pmCZP = 0, pmQZN[4], pmQZP[4], pmZEM1 = 0, pmZEM2 = 0;
283     for (Int_t i = 0; i < 4; i++) pmQZN[i] = pmQZP[i] = 0;
284
285     runLoader->GetEvent(iEvent);
286     TTree* treeH = fLoader->TreeH();
287     Int_t ntracks = (Int_t) treeH->GetEntries();
288     ResetHits();
289
290     // Tracks loop
291     Int_t sector[2];
292     for (Int_t itrack = 0; itrack < ntracks; itrack++) {
293       treeH->GetEntry(itrack);
294       for (AliZDCHit* zdcHit = (AliZDCHit*)FirstHit(-1); zdcHit;
295                       zdcHit = (AliZDCHit*)NextHit()) { 
296                       
297         sector[0] = zdcHit->GetVolume(0);
298         sector[1] = zdcHit->GetVolume(1);
299         if ((sector[1] < 1) || (sector[1] > 4)) {
300           Error("Hits2SDigits", "sector[0] = %d, sector[1] = %d", 
301                 sector[0], sector[1]);
302           continue;
303         }
304         Float_t lightQ = zdcHit->GetLightPMQ();
305         Float_t lightC = zdcHit->GetLightPMC();
306      
307         if (sector[0] == 1) {          //ZN 
308           pmCZN += lightC;
309           pmQZN[sector[1]-1] += lightQ;
310         } else if (sector[0] == 2) {   //ZP 
311           pmCZP += lightC;
312           pmQZP[sector[1]-1] += lightQ;
313         } else if (sector[0] == 3) {   //ZEM 
314           if (sector[1] == 1) pmZEM1 += lightC;
315           else                pmZEM2 += lightQ;
316         }
317       }//Hits loop
318     }
319
320     // create the output tree
321     fLoader->MakeTree("S");
322     TTree* treeS = fLoader->TreeS();
323     const Int_t kBufferSize = 4000;
324     treeS->Branch(GetName(), "AliZDCSDigit", &psdigit, kBufferSize);
325
326     // Create sdigits for ZN
327     sector[0] = 1; // Detector = ZN
328     sector[1] = 0; // Common PM ADC
329     new(psdigit) AliZDCSDigit(sector, pmCZN);
330     if (pmCZN > 0) treeS->Fill();
331     for (Int_t j = 0; j < 4; j++) {
332       sector[1] = j+1; // Towers PM ADCs
333       new(psdigit) AliZDCSDigit(sector, pmQZN[j]);
334       if (pmQZN[j] > 0) treeS->Fill();
335     }
336   
337     // Create sdigits for ZP
338     sector[0] = 2; // Detector = ZP
339     sector[1] = 0; // Common PM ADC
340     new(psdigit) AliZDCSDigit(sector, pmCZP);
341     if (pmCZP > 0) treeS->Fill();
342     for (Int_t j = 0; j < 4; j++) {
343       sector[1] = j+1; // Towers PM ADCs
344       new(psdigit) AliZDCSDigit(sector, pmQZP[j]);
345       if (pmQZP[j] > 0) treeS->Fill();
346     }
347
348     // Create sdigits for ZEM
349     sector[0] = 3; 
350     sector[1] = 1; // Detector = ZEM1
351     new(psdigit) AliZDCSDigit(sector, pmZEM1);
352     if (pmZEM1 > 0) treeS->Fill();
353     sector[1] = 2; // Detector = ZEM2
354     new(psdigit) AliZDCSDigit(sector, pmZEM2);
355     if (pmZEM2 > 0) treeS->Fill();
356
357     // write the output tree
358     fLoader->WriteSDigits("OVERWRITE");
359   }
360
361   fLoader->UnloadHits();
362   fLoader->UnloadSDigits();
363 }
364
365 //_____________________________________________________________________________
366 AliDigitizer* AliZDC::CreateDigitizer(AliRunDigitizer* manager) const
367 {
368   // Create the digitizer for ZDC
369
370   return new AliZDCDigitizer(manager);
371 }
372
373 //_____________________________________________________________________________
374 void AliZDC::Digits2Raw()
375 {
376   // Convert ZDC digits to raw data
377
378   // Format: 22 interger values -> ZN1 (C+Q1-4), ZP1 (C+Q1-4), ZEM1, 2, ZN (C+Q1-4), ZP2 (C+Q1-4))
379   // For the CAEN module V965 we have an header, the Data Words and an End Of Block
380   // 24 channels read on 1st ADC module, 20 channels read on 2nd ADC module
381   const int knADCData1=24, knADCData2=20; 
382   UInt_t lADCHeader1; 
383   UInt_t lADCData1[knADCData1];
384   //
385   UInt_t lADCHeader2; 
386   UInt_t lADCData2[knADCData2];
387   //
388   UInt_t lADCEndBlock;
389
390   // load the digits
391   fLoader->LoadDigits("read");
392   AliZDCDigit digit;
393   AliZDCDigit* pdigit = &digit;
394   TTree* treeD = fLoader->TreeD();
395   if (!treeD) return;
396   treeD->SetBranchAddress("ZDC", &pdigit);
397   //printf("\t AliZDC::Digits2raw -> TreeD has %d entries\n",(Int_t) treeD->GetEntries());
398   //digit.Print(""); // Ch. debug
399
400
401   // Fill data array
402   // ADC header
403   UInt_t lADCHeaderGEO = 0;
404   UInt_t lADCHeaderCRATE = 0;
405   UInt_t lADCHeaderCNT1 = knADCData1;
406   UInt_t lADCHeaderCNT2 = knADCData2;
407     
408   lADCHeader1 = lADCHeaderGEO << 27 | 0x1 << 25 | lADCHeaderCRATE << 16 |
409                lADCHeaderCNT1 << 8 ;
410   //           
411   lADCHeader2 = lADCHeaderGEO << 27 | 0x1 << 25 | lADCHeaderCRATE << 16 |
412                lADCHeaderCNT2 << 8 ;
413
414   //printf("\t lADCHeader1 = %x, lADCHeader2 = %x\n",lADCHeader1, lADCHeader2);
415       
416   // ADC data word
417   UInt_t lADCDataGEO = lADCHeaderGEO;
418   UInt_t lADCDataValue1[knADCData1];
419   UInt_t lADCDataValue2[knADCData2];
420   UInt_t lADCDataOvFlw1[knADCData1];
421   UInt_t lADCDataOvFlw2[knADCData2];
422   for(Int_t i = 0; i<knADCData1 ; i++){
423     lADCDataValue1[i] = 0;
424     lADCDataOvFlw1[i] = 0;
425   }
426   for(Int_t i = 0; i<knADCData2 ; i++){
427     lADCDataValue2[i] = 0;
428     lADCDataOvFlw2[i] = 0;
429   }
430   UInt_t lADCDataChannel = 0;
431   
432   // loop over digits
433   for (Int_t iDigit = 0; iDigit<treeD->GetEntries(); iDigit++) {
434     treeD->GetEntry(iDigit);
435     if (!pdigit) continue;
436     
437     //ADC data
438     Int_t index1 = 0, index2 = 0;
439     // ADC #1 (ZN1, ZP1, ZEM1,2)
440     if(digit.GetSector(0)==1 || digit.GetSector(0)==2 || digit.GetSector(0)==3){
441       if(digit.GetSector(0)==1 || digit.GetSector(0)==2){
442         index1 = (digit.GetSector(0)-1) + digit.GetSector(1)*4; // ZN1 or ZP1
443         lADCDataChannel = (digit.GetSector(0)-1)*8 + digit.GetSector(1);
444       }
445       else if(digit.GetSector(0)==3){ // ZEM 1,2
446         index1 = 20 + (digit.GetSector(1)-1);
447         lADCDataChannel = 5 + (digit.GetSector(1)-1)*8;
448       }
449       //
450       /*printf("\t AliZDC::Digits2raw -> det %d, quad %d, index = %d, ADCch = %d\n",
451                 digit.GetSector(0),digit.GetSector(1),index1,lADCDataChannel);// Ch. debug
452       */
453       //
454       lADCDataValue1[index1] = digit.GetADCValue(0);    // High gain ADC ch.    
455       if(lADCDataValue1[index1] > 2047) lADCDataOvFlw1[index1] = 1; 
456       lADCDataValue1[index1+2] = digit.GetADCValue(1); // Low gain ADC ch.
457       if(lADCDataValue1[index1+2] > 2047) lADCDataOvFlw1[index1+2] = 1; 
458     
459       lADCData1[index1] = lADCDataGEO << 27 | lADCDataChannel << 17 | 
460                         lADCDataOvFlw1[index1] << 12 | (lADCDataValue1[index1] & 0xfff); 
461       lADCData1[index1+2] = lADCDataGEO << 27 | lADCDataChannel << 17 | 0x1 << 16 |
462                         lADCDataOvFlw1[index1+2] << 12 | (lADCDataValue1[index1+2] & 0xfff);                    
463     }
464     // ADC #2 (ZN2, ZP2)
465     else if(digit.GetSector(0)==4 || digit.GetSector(0)==5){
466       index2 = (digit.GetSector(0)-4) + digit.GetSector(1)*4; // ZN2 or ZP2
467       lADCDataChannel = (digit.GetSector(0)-4)*8 + digit.GetSector(1);
468       //
469       /*printf("\t AliZDC::Digits2raw -> det %d, quad %d, index = %d, ADCch = %d\n",
470                 digit.GetSector(0),digit.GetSector(1),index1,lADCDataChannel); // Ch. debug
471       */
472       //
473       lADCDataValue2[index2] = digit.GetADCValue(0);
474       if (lADCDataValue2[index2] > 2047) lADCDataOvFlw2[index2] = 1; 
475       lADCDataValue2[index2+2] = digit.GetADCValue(1);
476       if (lADCDataValue2[index2+2] > 2047) lADCDataOvFlw2[index2+2] = 1; 
477       //
478       lADCData2[index2] =   lADCDataGEO << 27 | lADCDataChannel << 17 | 
479                         lADCDataOvFlw2[index2] << 12 | (lADCDataValue2[index2] & 0xfff); 
480       lADCData2[index2+2] = lADCDataGEO << 27 | lADCDataChannel << 17 | 0x1 << 16 |
481                         lADCDataOvFlw2[index2+2] << 12 | (lADCDataValue2[index2+2] & 0xfff);                    
482     } 
483     if((index1<0) || (index1>23)) {
484       Error("Digits2Raw", "sector[0] = %d, sector[1] = %d", 
485             digit.GetSector(0), digit.GetSector(1));
486       continue;
487     }
488     if((index2<0) || (index2>19)) {
489       Error("Digits2Raw", "sector[0] = %d, sector[1] = %d", 
490             digit.GetSector(0), digit.GetSector(1));
491       continue;
492     }
493     
494     
495   }
496   //for(Int_t i=0;i<24;i++) printf("\t ADCData1[%d] = %x\n",i,lADCData1[i]);
497   //for(Int_t i=0;i<20;i++) printf("\t ADCData2[%d] = %x\n",i,lADCData2[i]);
498  
499   // End of Block
500   UInt_t lADCEndBlockGEO = lADCHeaderGEO;
501   UInt_t lADCEndBlockEvCount = gAlice->GetEventNrInRun();
502   
503   lADCEndBlock = lADCEndBlockGEO << 27 | 0x1 << 26 | lADCEndBlockEvCount;
504   
505   //printf("\t ADCEndBlock = %d\n",lADCEndBlock);
506
507
508   // open the output file
509   char fileName[30];
510   strcpy(fileName,AliDAQ::DdlFileName("ZDC",0));
511 #ifndef __DECCXX
512   ofstream file(fileName, ios::binary);
513 #else
514   ofstream file(fileName);
515 #endif
516
517   // write the DDL data header
518   AliRawDataHeader header;
519   header.fSize = sizeof(header) + sizeof(lADCHeader1) + sizeof(lADCData1) + 
520                 sizeof(lADCEndBlock)+ sizeof(lADCHeader2) + sizeof(lADCData2) + sizeof(lADCEndBlock);
521   /*printf("sizeof header = %d, ADCHeader1 = %d, ADCData1 = %d, ADCEndBlock = %d\n",
522           sizeof(header),sizeof(lADCHeader1),sizeof(lADCData1),sizeof(lADCEndBlock));
523   printf("sizeof header = %d, ADCHeader2 = %d, ADCData2 = %d, ADCEndBlock = %d\n",
524           sizeof(header),sizeof(lADCHeader2),sizeof(lADCData2),sizeof(lADCEndBlock));*/
525   header.SetAttribute(0);  // valid data
526   file.write((char*)(&header), sizeof(header));
527
528   // write the raw data and close the file
529   file.write((char*) &lADCHeader1, sizeof (lADCHeader1));
530   file.write((char*)(lADCData1), sizeof(lADCData1));
531   file.write((char*) &lADCEndBlock, sizeof(lADCEndBlock));
532   file.write((char*) &lADCHeader2, sizeof (lADCHeader2));
533   file.write((char*)(lADCData2), sizeof(lADCData2));
534   file.write((char*) &lADCEndBlock, sizeof(lADCEndBlock));
535   file.close();
536
537   // unload the digits
538   fLoader->UnloadDigits();
539 }
540
541 //______________________________________________________________________
542 void AliZDC::SetTreeAddress(){
543   // Set branch address for the Trees.
544   // Inputs:
545   //      none.
546   // Outputs:
547   //      none.
548   // Return:
549   //      none.
550   if (fLoader->TreeH() && (fHits == 0x0))
551     fHits   = new TClonesArray("AliZDCHit",1000);
552       
553   AliDetector::SetTreeAddress();
554 }
555  
556  
557 //Calibration methods (by Alberto Colla)
558  
559  
560 //________________________________________________________________
561 void AliZDC::CreateCalibData()
562 {
563   // 
564   //if (fCalibData) delete fCalibData; // delete previous version
565   fCalibData = new AliZDCCalibData(GetName());
566 }
567 //________________________________________________________________
568 void AliZDC::WriteCalibData(Int_t option)
569 {
570   //
571   const int kCompressLevel = 9;
572   char* fnam = GetZDCCalibFName();
573   if (!fnam || fnam[0]=='\0') {
574     fnam = gSystem->ExpandPathName("$(ALICE)/$(ALICE_LEVEL)/data/AliZDCCalib.root");
575     Warning("WriteCalibData","No File Name is provided, using default %s",fnam);
576   }
577   TFile* cdfile = TFile::Open(fnam,"UPDATE","",kCompressLevel);
578
579   // Writes Calibration Data to current directory. 
580   // User MUST take care of corresponding file opening and ->cd()... !!!
581   // By default, the object is overwritten. Use 0 option for opposite.
582   if (option) option = TObject::kOverwrite;
583   if (fCalibData) fCalibData->Write(0,option);
584   else if (fCalibData) fCalibData->Write(0,option);
585
586   cdfile->Close();
587   delete cdfile;
588 }
589
590 //________________________________________________________________
591 void AliZDC::LoadCalibData()
592 {
593   //
594   char* fnam = GetZDCCalibFName();
595   if (!fnam || fnam[0]=='\0') return; 
596   if (!gAlice->IsFileAccessible(fnam)) {
597     Error("LoadCalibData","ZDC Calibration Data file is not accessible, %s",fnam);
598     exit(1);
599   }
600   TFile* cdfile = TFile::Open(fnam);
601
602   // Loads Calibration Data from current directory. 
603   // User MUST take care of corresponding file opening and ->cd()...!!!
604   //
605   if (fCalibData) delete fCalibData; // delete previous version
606   TString dtname = "Calib_";
607   dtname += GetName();
608   fCalibData = (AliZDCCalibData*) gDirectory->Get(dtname.Data());
609   if (!fCalibData) { 
610     Error("LoadCalibData","No Calibration data found for %s",GetName());
611     exit(1);
612   }
613
614   cdfile->Close();
615   delete cdfile;
616 }
617
618
619 //Calibration methods (by Alberto Colla)