]> git.uio.no Git - u/mrichter/AliRoot.git/blame - ZDC/AliZDCDigitizer.cxx
Protection against missing digit in a digits array in added (P.Hristov)
[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;
49class AliZDCRecParam;
8309c1ab 50
51ClassImp(AliZDCDigitizer)
52
53
54//____________________________________________________________________________
55AliZDCDigitizer::AliZDCDigitizer()
56{
57// Default constructor
58
59}
60
61//____________________________________________________________________________
62AliZDCDigitizer::AliZDCDigitizer(AliRunDigitizer* manager):
63 AliDigitizer(manager)
64{
abf60186 65 fIsCalibration=0; //By default the simulation doesn't create calib. data
6024ec85 66// fIsCalibration=1; //To create pedestal calib. data
78d18275 67 // Get calibration data
6024ec85 68 fPedData = GetPedData();
7f73eb6b 69 fCalibData = GetCalibData();
6024ec85 70 fRecParam = GetRecParam();
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
79
80}
81
82
cc2abffd 83//____________________________________________________________________________
84AliZDCDigitizer::AliZDCDigitizer(const AliZDCDigitizer &digitizer):
85 AliDigitizer()
86{
87 // Copy constructor
88
89 for(Int_t i=0; i<6; i++){
90 for(Int_t j=0; j<5; j++){
91 fPMGain[i][j] = digitizer.fPMGain[i][j];
92 }
93 }
94 for(Int_t i=0; i<2; i++) fADCRes[i] = digitizer.fADCRes[i];
95 fIsCalibration = digitizer.fIsCalibration;
6024ec85 96 fPedData = digitizer.fPedData;
cc2abffd 97 fCalibData = digitizer.fCalibData;
6024ec85 98 fRecParam = digitizer.fRecParam;
cc2abffd 99
100}
101
8309c1ab 102//____________________________________________________________________________
103Bool_t AliZDCDigitizer::Init()
104{
48642b09 105 // Initialize the digitizer
106 // NB -> PM gain vs. HV & ADC resolutions will move to DCDB ***************
7f73eb6b 107 for(Int_t j = 0; j < 5; j++){
108 fPMGain[0][j] = 50000.;
109 fPMGain[1][j] = 100000.;
110 fPMGain[2][j] = 100000.;
a017f01b 111 fPMGain[3][j] = 50000.;
112 fPMGain[4][j] = 100000.;
113 fPMGain[5][j] = 100000.;
7f73eb6b 114 }
115 // ADC Caen V965
116 fADCRes[0] = 0.0000008; // ADC Resolution high gain: 200 fC/adcCh
117 fADCRes[1] = 0.0000064; // ADC Resolution low gain: 25 fC/adcCh
8309c1ab 118
119 return kTRUE;
120}
121
122//____________________________________________________________________________
123void AliZDCDigitizer::Exec(Option_t* /*option*/)
124{
48642b09 125 // Execute digitization
8309c1ab 126
d79f8d50 127 // ------------------------------------------------------------
83347831 128 // !!! 2nd ZDC set added
7f73eb6b 129 // *** 1st 3 arrays are digits from REAL (simulated) hits
130 // *** last 2 are copied from simulated digits
131 // --- pm[0][...] = light in ZN right [C, Q1, Q2, Q3, Q4]
132 // --- pm[1][...] = light in ZP right [C, Q1, Q2, Q3, Q4]
48642b09 133 // --- pm[2][...] = light in ZEM [x, 1, 2, x, x]
7f73eb6b 134 // --- pm[3][...] = light in ZN left [C, Q1, Q2, Q3, Q4] ->NEW!
135 // --- pm[4][...] = light in ZP left [C, Q1, Q2, Q3, Q4] ->NEW!
d79f8d50 136 // ------------------------------------------------------------
137 Float_t pm[5][5];
8574b308 138 for(Int_t iSector1=0; iSector1<5; iSector1++)
139 for(Int_t iSector2=0; iSector2<5; iSector2++){
8309c1ab 140 pm[iSector1][iSector2] = 0;
141 }
d79f8d50 142
143 // ------------------------------------------------------------
144 // ### Out of time ADC added (22 channels)
145 // --- same codification as for signal PTMs (see above)
146 // ------------------------------------------------------------
147 Float_t pmoot[5][5];
8574b308 148 for(Int_t iSector1=0; iSector1<5; iSector1++)
149 for(Int_t iSector2=0; iSector2<5; iSector2++){
d79f8d50 150 pmoot[iSector1][iSector2] = 0;
151 }
8309c1ab 152
153 // impact parameter and number of spectators
154 Float_t impPar = -1;
155 Int_t specN = 0;
156 Int_t specP = 0;
157
158 // loop over input streams
8574b308 159 for(Int_t iInput = 0; iInput<fManager->GetNinputs(); iInput++){
8309c1ab 160
161 // get run loader and ZDC loader
162 AliRunLoader* runLoader =
163 AliRunLoader::GetRunLoader(fManager->GetInputFolderName(iInput));
164 AliLoader* loader = runLoader->GetLoader("ZDCLoader");
8574b308 165 if(!loader) continue;
8309c1ab 166
167 // load sdigits
168 loader->LoadSDigits();
169 TTree* treeS = loader->TreeS();
8574b308 170 if(!treeS) continue;
8309c1ab 171 AliZDCSDigit sdigit;
172 AliZDCSDigit* psdigit = &sdigit;
173 treeS->SetBranchAddress("ZDC", &psdigit);
174
175 // loop over sdigits
8574b308 176 for(Int_t iSDigit=0; iSDigit<treeS->GetEntries(); iSDigit++){
8309c1ab 177 treeS->GetEntry(iSDigit);
7f73eb6b 178 //
8574b308 179 if(!psdigit) continue;
180 if((sdigit.GetSector(1) < 0) || (sdigit.GetSector(1) > 4)){
8309c1ab 181 AliError(Form("\nsector[0] = %d, sector[1] = %d\n",
182 sdigit.GetSector(0), sdigit.GetSector(1)));
183 continue;
184 }
7f73eb6b 185 //
48642b09 186 pm[(sdigit.GetSector(0))-1][sdigit.GetSector(1)] += sdigit.GetLightPM();
7f73eb6b 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,
189 sdigit.GetSector(1), pm[sdigit.GetSector(0)-1][sdigit.GetSector(1)]); // Chiara debugging!
abf60186 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();
8309c1ab 197 AliHeader* header = runLoader->GetAliRun()->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));
d27e20dd 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
7f73eb6b 244 Int_t sector[2], sectorL[2];
d79f8d50 245 Int_t digi[2], digiL[2], digioot[2];
abf60186 246 for(sector[0]=1; sector[0]<=3; sector[0]++){
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 }
58e12fee 253 /*printf("\t DIGIT added -> det %d quad %d - digi[0,1] = [%d, %d]\n",
abf60186 254 sector[0], sector[1], digi[0], digi[1]); // Chiara debugging!
d79f8d50 255 */
abf60186 256 //
257 new(pdigit) AliZDCDigit(sector, digi);
8309c1ab 258 treeD->Fill();
7f73eb6b 259 //
260 // --- Adding digits for 2nd ZDC set (left side w.r.t. IP) ---
261 // --- they are copied from right ZDC digits
262 if(sector[0]==1 || sector[0]==2){
263 sectorL[0] = sector[0]+3;
264 sectorL[1] = sector[1];
8574b308 265 for(Int_t res=0; res<2; res++){
7f73eb6b 266 digiL[res] = Phe2ADCch(sectorL[0], sectorL[1], pm[sector[0]-1][sector[1]], res)
267 + Pedestal(sectorL[0], sectorL[1], res);
268 }
58e12fee 269 /*printf("\t DIGIT added -> det %d quad %d - digi[0,1] = [%d, %d]\n",
abf60186 270 sectorL[0], sectorL[1], digiL[0], digiL[1]); // Chiara debugging!
271 */
272 //
273 new(pdigit) AliZDCDigit(sectorL, digiL);
274 treeD->Fill();
7f73eb6b 275 }
8309c1ab 276 }
83347831 277 } // Loop over detector
278 // Adding in-time digits for 2 reference PTM signals (after signal ch.)
279 // (for the moment the ref. signal is completely invented assuming a PMgain of 5*10^4!)
280 Int_t sectorRef[2];
281 sectorRef[1] = 5;
9d38ced8 282 Int_t sigRef[2];
283 if(fIsCalibration==0) {sigRef[0]=100; sigRef[1]=800;}
284 else {sigRef[0]=0; sigRef[1]=0;}
285 //
83347831 286 for(Int_t iref=0; iref<2; iref++){
287 sectorRef[0] = 3*iref+1;
288 for(Int_t res=0; res<2; res++){
289 sigRef[res] += Pedestal(sectorRef[0], sectorRef[1], res);
290 }
291 /*printf("\t RefDigit added -> det = %d, quad = %d - digi[0,1] = [%d, %d]\n",
292 sectorRef[0], sectorRef[1], sigRef[0], sigRef[1]); // Chiara debugging!
293 */
294 new(pdigit) AliZDCDigit(sectorRef, sigRef);
295 treeD->Fill();
abf60186 296 }
d79f8d50 297 //
298 // --- Adding digits for out-of-time channels after signal digits
299 for(sector[0]=1; sector[0]<=3; sector[0]++){
300 for(sector[1]=0; sector[1]<5; sector[1]++){
301 if((sector[0]==3) && ((sector[1]<1) || (sector[1]>2))) continue;
8574b308 302 for(Int_t res=0; res<2; res++){
d79f8d50 303 digioot[res] = Pedestal(sector[0], sector[1], res); // out-of-time ADCs
304 }
305 /*printf("\t DIGIToot added -> det = %d, quad = %d - digi[0,1] = [%d, %d]\n",
306 sector[0], sector[1], digioot[0], digioot[1]); // Chiara debugging!
307 */
308 //
309 new(pdigit) AliZDCDigit(sector, digioot);
310 treeD->Fill();
311 //
312 if(sector[0]==1 || sector[0]==2){
313 sectorL[0] = sector[0]+3;
314 sectorL[1] = sector[1];
8574b308 315 for(Int_t res=0; res<2; res++){
d79f8d50 316 digioot[res] = Pedestal(sectorL[0], sectorL[1], res); // out-of-time ADCs
317 }
318 /*printf("\t DIGIToot added -> det = %d, quad = %d - digi[0,1] = [%d, %d]\n",
319 sectorL[0], sectorL[1], digioot[0], digioot[1]); // Chiara debugging!
320 */
321 //
322 new(pdigit) AliZDCDigit(sectorL, digioot);
323 treeD->Fill();
324 }
325 //
d79f8d50 326 }
327 }
83347831 328 // Adding out-of-time digits for 2 reference PTM signals (after out-of-time ch.)
329 Int_t sigRefoot[2];
330 for(Int_t iref=0; iref<2; iref++){
331 sectorRef[0] = 3*iref+1;
332 for(Int_t res=0; res<2; res++){
333 sigRefoot[res] = Pedestal(sectorRef[0], sectorRef[1], res);
334 }
335 /*printf("\t RefDigitoot added -> det = %d, quad = %d - digi[0,1] = [%d, %d]\n",
336 sectorRef[0], sectorRef[1], sigRefoot[0], sigRefoot[1]); // Chiara debugging!
337 */
338 new(pdigit) AliZDCDigit(sectorRef, sigRefoot);
339 treeD->Fill();
340
341 }
342 //printf("\t AliZDCDigitizer -> TreeD has %d entries\n",(Int_t) treeD->GetEntries());
343
8309c1ab 344 // write the output tree
345 loader->WriteDigits("OVERWRITE");
346 loader->UnloadDigits();
347}
348
349
350//_____________________________________________________________________________
351void AliZDCDigitizer::Fragmentation(Float_t impPar, Int_t specN, Int_t specP,
352 Int_t &freeSpecN, Int_t &freeSpecP) const
353{
354// simulate fragmentation of spectators
355
356 Int_t zz[100], nn[100];
357 AliZDCFragment frag(impPar);
8574b308 358 for(Int_t j=0; j<=99; j++){
8309c1ab 359 zz[j] =0;
360 nn[j] =0;
361 }
362
363 // Fragments generation
364 Int_t nAlpha;
365 frag.GenerateIMF(zz, nAlpha);
366
367 // Attach neutrons
368 Int_t ztot=0;
369 Int_t ntot=0;
370 frag.AttachNeutrons(zz, nn, ztot, ntot);
371 freeSpecN = specN-ntot-2*nAlpha;
372 freeSpecP = specP-ztot-2*nAlpha;
3848c666 373 // Removing deuterons
a0607c1f 374 Int_t ndeu = (Int_t) (freeSpecN*frag.DeuteronNumber());
3848c666 375 freeSpecN -= ndeu;
376 //
8309c1ab 377 if(freeSpecN<0) freeSpecN=0;
378 if(freeSpecP<0) freeSpecP=0;
379 AliDebug(2, Form("FreeSpn = %d, FreeSpp = %d", freeSpecN, freeSpecP));
380}
381
382//_____________________________________________________________________________
383void AliZDCDigitizer::SpectatorSignal(Int_t SpecType, Int_t numEvents,
384 Float_t pm[3][5]) const
385{
386// add signal of the spectators
387
388 TFile* file = NULL;
8574b308 389 if(SpecType == 1) { // --- Signal for spectator neutrons
e689b8ea 390 file = TFile::Open("$ALICE_ROOT/ZDC/ZNsignalntu.root");
8574b308 391 } else if(SpecType == 2) { // --- Signal for spectator protons
e689b8ea 392 file = TFile::Open("$ALICE_ROOT/ZDC/ZPsignalntu.root");
8309c1ab 393 }
8574b308 394 if(!file || !file->IsOpen()) {
8309c1ab 395 AliError("Opening of file failed");
396 return;
397 }
398
399 TNtuple* zdcSignal = (TNtuple*) file->Get("ZDCSignal");
400 Int_t nentries = (Int_t) zdcSignal->GetEntries();
401
402 Float_t *entry, hitsSpec[7];
403 Int_t pl, i, j, k, iev=0, rnd[125], volume[2];
48642b09 404 for(pl=0;pl<125;pl++) rnd[pl] = 0;
8574b308 405 if(numEvents > 125) {
8309c1ab 406 AliWarning(Form("numEvents (%d) is larger than 125", numEvents));
407 numEvents = 125;
408 }
409 for(pl=0;pl<numEvents;pl++){
410 rnd[pl] = (Int_t) (9999*gRandom->Rndm());
85a96223 411 if(rnd[pl] >= 9999) rnd[pl] = 9998;
8309c1ab 412 //printf(" rnd[%d] = %d\n",pl,rnd[pl]);
413 }
414 // Sorting vector in ascending order with C function QSORT
415 qsort((void*)rnd,numEvents,sizeof(Int_t),comp);
416 do{
417 for(i=0; i<nentries; i++){
418 zdcSignal->GetEvent(i);
419 entry = zdcSignal->GetArgs();
420 if(entry[0] == rnd[iev]){
421 for(k=0; k<2; k++) volume[k] = (Int_t) entry[k+1];
422 for(j=0; j<7; j++) hitsSpec[j] = entry[j+3];
48642b09 423 //
8309c1ab 424 Float_t lightQ = hitsSpec[4];
425 Float_t lightC = hitsSpec[5];
78d18275 426 AliDebug(3, Form("SpectatorSignal -> vol = (%d, %d), lightQ = %.0f, lightC = %.0f",
8309c1ab 427 volume[0], volume[1], lightQ, lightC));
48642b09 428 //printf("\n Volume = (%d, %d), lightQ = %.0f, lightC = %.0f",
429 // volume[0], volume[1], lightQ, lightC);
acf9550c 430 if(volume[0] != 3) { // ZN or ZP
8309c1ab 431 pm[volume[0]-1][0] += lightC;
432 pm[volume[0]-1][volume[1]] += lightQ;
48642b09 433 //printf("\n pm[%d][0] = %.0f, pm[%d][%d] = %.0f\n",(volume[0]-1),pm[volume[0]-1][0],
434 // (volume[0]-1),volume[1],pm[volume[0]-1][volume[1]]);
435 }
436 else {
8574b308 437 if(volume[1] == 1) pm[2][1] += lightC; // ZEM 1
acf9550c 438 else pm[2][2] += lightQ; // ZEM 2
48642b09 439 //printf("\n pm[2][1] = %.0f, pm[2][2] = %.0f\n",pm[2][1],pm[2][2]);
440 }
8309c1ab 441 }
442 else if(entry[0] > rnd[iev]){
443 iev++;
444 continue;
445 }
446 }
447 }while(iev<numEvents);
448
449 file->Close();
450 delete file;
451}
452
453
454//_____________________________________________________________________________
455Int_t AliZDCDigitizer::Phe2ADCch(Int_t Det, Int_t Quad, Float_t Light,
456 Int_t Res) const
457{
48642b09 458 // Evaluation of the ADC channel corresponding to the light yield Light
8a2624cc 459 Int_t vADCch = (Int_t) (Light * fPMGain[Det-1][Quad] * fADCRes[Res]);
58e12fee 460 //printf("\t Phe2ADCch -> det %d quad %d - phe %.0f ADC %d\n", Det,Quad,Light,vADCch);
7f73eb6b 461
8a2624cc 462 return vADCch;
48642b09 463}
8309c1ab 464
48642b09 465//_____________________________________________________________________________
466Int_t AliZDCDigitizer::Pedestal(Int_t Det, Int_t Quad, Int_t Res) const
467{
8a2624cc 468 // Returns a pedestal for detector det, PM quad, channel with res.
469 //
65b16e87 470 Float_t pedValue;
abf60186 471 // Normal run
472 if(fIsCalibration == 0){
c35ed519 473 Int_t index=0, kNch=24;
83347831 474 if(Quad!=5){
c35ed519 475 if(Det==1) index = Quad+kNch*Res; // ZNC
476 else if(Det==2) index = (Quad+5)+kNch*Res; // ZPC
477 else if(Det==3) index = (Quad+9)+kNch*Res; // ZEM
478 else if(Det==4) index = (Quad+12)+kNch*Res; // ZNA
479 else if(Det==5) index = (Quad+17)+kNch*Res; // ZPA
83347831 480 }
c35ed519 481 else index = (Det-1)/3+22+kNch*Res; // Reference PMs
83347831 482 //
c35ed519 483 Float_t meanPed = fPedData->GetMeanPed(index);
484 Float_t pedWidth = fPedData->GetMeanPedWidth(index);
65b16e87 485 pedValue = gRandom->Gaus(meanPed,pedWidth);
abf60186 486 //
58e12fee 487 /*printf("\t AliZDCDigitizer::Pedestal -> det %d quad %d res %d - Ped[%d] = %d\n",
65b16e87 488 Det, Quad, Res, index,(Int_t) pedValue); // Chiara debugging!
abf60186 489 */
490 }
7f73eb6b 491 // To create calibration object
83347831 492 else{
65b16e87 493 if(Res == 0) pedValue = gRandom->Gaus((35.+10.*gRandom->Rndm()),(0.5+0.2*gRandom->Rndm())); //High gain
494 else pedValue = gRandom->Gaus((250.+100.*gRandom->Rndm()),(3.5+2.*gRandom->Rndm())); //Low gain
83347831 495 }
dbc8d7fc 496
65b16e87 497 return (Int_t) pedValue;
8309c1ab 498}
499
500//_____________________________________________________________________________
78d18275 501AliCDBStorage* AliZDCDigitizer::SetStorage(const char *uri)
8309c1ab 502{
8309c1ab 503
78d18275 504 Bool_t deleteManager = kFALSE;
48642b09 505
78d18275 506 AliCDBManager *manager = AliCDBManager::Instance();
507 AliCDBStorage *defstorage = manager->GetDefaultStorage();
48642b09 508
78d18275 509 if(!defstorage || !(defstorage->Contains("ZDC"))){
510 AliWarning("No default storage set or default storage doesn't contain ZDC!");
511 manager->SetDefaultStorage(uri);
512 deleteManager = kTRUE;
513 }
514
515 AliCDBStorage *storage = manager->GetDefaultStorage();
516
517 if(deleteManager){
518 AliCDBManager::Instance()->UnsetDefaultStorage();
519 defstorage = 0; // the storage is killed by AliCDBManager::Instance()->Destroy()
520 }
521
522 return storage;
523}
48642b09 524
78d18275 525//_____________________________________________________________________________
6024ec85 526AliZDCPedestals* AliZDCDigitizer::GetPedData() const
527{
528
529 // Getting pedestal calibration object for ZDC set
530
531 AliCDBEntry *entry = AliCDBManager::Instance()->Get("ZDC/Calib/Pedestals");
532 if(!entry) AliFatal("No calibration data loaded!");
533
534 AliZDCPedestals *calibdata = dynamic_cast<AliZDCPedestals*> (entry->GetObject());
535 if(!calibdata) AliFatal("Wrong calibration object in calibration file!");
536
537 return calibdata;
538}
539
540//_____________________________________________________________________________
541AliZDCCalib* AliZDCDigitizer::GetCalibData() const
542{
543
544 // Getting calibration object for ZDC set
545
546 AliCDBEntry *entry = AliCDBManager::Instance()->Get("ZDC/Calib/Calib");
547 if(!entry) AliFatal("No calibration data loaded!");
548
549 AliZDCCalib *calibdata = dynamic_cast<AliZDCCalib*> (entry->GetObject());
550 if(!calibdata) AliFatal("Wrong calibration object in calibration file!");
551
552 return calibdata;
553}
554
555//_____________________________________________________________________________
556AliZDCRecParam* AliZDCDigitizer::GetRecParam() const
78d18275 557{
4fda3ba1 558
559 // Getting calibration object for ZDC set
560
6024ec85 561 AliCDBEntry *entry = AliCDBManager::Instance()->Get("ZDC/Calib/RecParam");
457a440d 562 if(!entry) AliFatal("No calibration data loaded!");
563
6024ec85 564 AliZDCRecParam *calibdata = dynamic_cast<AliZDCRecParam*> (entry->GetObject());
457a440d 565 if(!calibdata) AliFatal("Wrong calibration object in calibration file!");
8309c1ab 566
78d18275 567 return calibdata;
8309c1ab 568}
4fda3ba1 569