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 ClassImp(AliZDCDigitizer)
54 //____________________________________________________________________________
55 AliZDCDigitizer::AliZDCDigitizer() :
57 fIsSignalInADCGate(kFALSE),
60 fSpectators2Track(kFALSE),
63 // Default constructor
64 for(Int_t i=0; i<2; i++) fADCRes[i]=0.;
67 //____________________________________________________________________________
68 AliZDCDigitizer::AliZDCDigitizer(AliRunDigitizer* manager):
69 AliDigitizer(manager),
70 fIsCalibration(0), //By default the simulation doesn't create calib. data
71 fIsSignalInADCGate(kFALSE),
73 fPedData(GetPedData()),
74 fSpectators2Track(kFALSE),
77 // Get calibration data
78 if(fIsCalibration!=0) printf("\n\t AliZDCDigitizer -> Creating calibration data (pedestals)\n");
79 for(Int_t i=0; i<5; i++){
80 for(Int_t j=0; j<5; j++)
85 //____________________________________________________________________________
86 AliZDCDigitizer::~AliZDCDigitizer()
93 //____________________________________________________________________________
94 AliZDCDigitizer::AliZDCDigitizer(const AliZDCDigitizer &digitizer):
96 fIsCalibration(digitizer.fIsCalibration),
97 fIsSignalInADCGate(digitizer.fIsSignalInADCGate),
98 fFracLostSignal(digitizer.fFracLostSignal),
99 fPedData(digitizer.fPedData),
100 fSpectators2Track(digitizer.fSpectators2Track),
101 fBeamEnergy(digitizer.fBeamEnergy)
105 for(Int_t i=0; i<5; i++){
106 for(Int_t j=0; j<5; j++){
107 fPMGain[i][j] = digitizer.fPMGain[i][j];
110 for(Int_t i=0; i<2; i++) fADCRes[i] = digitizer.fADCRes[i];
115 //____________________________________________________________________________
116 Bool_t AliZDCDigitizer::Init()
118 // Initialize the digitizer
120 AliCDBEntry* entry = AliCDBManager::Instance()->Get("GRP/GRP/Data");
121 if(!entry) AliFatal("No calibration data loaded!");
122 AliGRPObject* grpData = 0x0;
124 TMap* m = dynamic_cast<TMap*>(entry->GetObject()); // old GRP entry
127 grpData = new AliGRPObject();
128 grpData->ReadValuesFromMap(m);
131 grpData = dynamic_cast<AliGRPObject*>(entry->GetObject()); // new GRP entry
134 AliCDBManager::Instance()->UnloadFromCache("GRP/GRP/Data");
137 AliError("No GRP entry found in OCDB! \n ");
141 TString beamType = grpData->GetBeamType();
142 if(beamType==AliGRPObject::GetInvalidString()){
143 AliError("\t UNKNOWN beam type from GRP obj -> PMT gains not set in ZDC digitizer!!!\n");
146 fBeamEnergy = grpData->GetBeamEnergy();
147 if(fBeamEnergy==AliGRPObject::GetInvalidFloat()){
148 AliWarning("GRP/GRP/Data entry: missing value for the beam energy ! Using 0.");
149 AliError("\t UNKNOWN beam type from GRP obj -> PMT gains not set in ZDC digitizer!!!\n");
153 if((beamType.CompareTo("P-P")) == 0 || (beamType.CompareTo("p-p")) == 0){
154 //PTM gains rescaled to beam energy for p-p
155 // New correction coefficients for PMT gains needed
156 // to reproduce experimental spectra (from Grazia Jul 2010)
157 if(fBeamEnergy != 0){
158 for(Int_t j = 0; j < 5; j++){
159 fPMGain[0][j] = 1.515831*(661.444/fBeamEnergy+0.000740671)*10000000;
160 fPMGain[1][j] = 0.674234*(864.350/fBeamEnergy+0.00234375)*10000000;
161 fPMGain[3][j] = 1.350938*(661.444/fBeamEnergy+0.000740671)*10000000;
162 fPMGain[4][j] = 0.678597*(864.350/fBeamEnergy+0.00234375)*10000000;
164 fPMGain[2][1] = 0.869654*(1.32312-0.000101515*fBeamEnergy)*10000000;
165 fPMGain[2][2] = 1.030883*(1.32312-0.000101515*fBeamEnergy)*10000000;
166 AliInfo(Form(" PMT gains for p-p @ %1.0f+%1.0f GeV: ZN(%1.0f), ZP(%1.0f), ZEM(%1.0f)\n",
167 fBeamEnergy, fBeamEnergy, fPMGain[0][0], fPMGain[1][0], fPMGain[2][1]));
170 else if((beamType.CompareTo("A-A")) == 0){
171 // PTM gains for Pb-Pb @ 2.7+2.7 A TeV ***************
172 // rescaled for Pb-Pb @ 1.38+1.38 A TeV ***************
173 // Values corrected after 2010 Pb-Pb data taking (7/2/2011 - Ch.)
174 // Experimental data compared to EMD simulation for single nucleon peaks:
175 // ZN gains must be divided by 4, ZP gains by 10!
176 Float_t scalGainFactor = fBeamEnergy/2760.;
177 for(Int_t j = 0; j < 5; j++){
178 fPMGain[0][j] = 50000./(4*scalGainFactor); // ZNC
179 fPMGain[1][j] = 100000./(10*scalGainFactor); // ZPC
180 fPMGain[2][j] = 100000./scalGainFactor; // ZEM
181 fPMGain[3][j] = 50000./(4*scalGainFactor); // ZNA
182 fPMGain[4][j] = 100000./(10*scalGainFactor); // ZPA
184 AliInfo(Form(" PMT gains for Pb-Pb @ %1.0f+%1.0f A GeV: ZN(%1.0f), ZP(%1.0f), ZEM(%1.0f)\n",
185 fBeamEnergy, fBeamEnergy, fPMGain[0][0], fPMGain[1][0], fPMGain[2][1]));
189 fADCRes[0] = 0.0000008; // ADC Resolution high gain: 200 fC/adcCh
190 fADCRes[1] = 0.0000064; // ADC Resolution low gain: 25 fC/adcCh
195 //____________________________________________________________________________
196 void AliZDCDigitizer::Exec(Option_t* /*option*/)
198 // Execute digitization
200 // ------------------------------------------------------------
201 // !!! 2nd ZDC set added
202 // *** 1st 3 arrays are digits from REAL (simulated) hits
203 // *** last 2 are copied from simulated digits
204 // --- pm[0][...] = light in ZN side C [C, Q1, Q2, Q3, Q4]
205 // --- pm[1][...] = light in ZP side C [C, Q1, Q2, Q3, Q4]
206 // --- pm[2][...] = light in ZEM [x, 1, 2, x, x]
207 // --- pm[3][...] = light in ZN side A [C, Q1, Q2, Q3, Q4] ->NEW!
208 // --- pm[4][...] = light in ZP side A [C, Q1, Q2, Q3, Q4] ->NEW!
209 // ------------------------------------------------------------
211 for(Int_t iSector1=0; iSector1<5; iSector1++)
212 for(Int_t iSector2=0; iSector2<5; iSector2++){
213 pm[iSector1][iSector2] = 0;
216 // ------------------------------------------------------------
217 // ### Out of time ADC added (22 channels)
218 // --- same codification as for signal PTMs (see above)
219 // ------------------------------------------------------------
221 for(Int_t iSector1=0; iSector1<5; iSector1++)
222 for(Int_t iSector2=0; iSector2<5; iSector2++){
223 pmoot[iSector1][iSector2] = 0;
226 // impact parameter and number of spectators
228 Int_t specNTarg = 0, specPTarg = 0;
229 Int_t specNProj = 0, specPProj = 0;
230 Float_t signalTime0 = 0.;
232 // loop over input streams
233 for(Int_t iInput = 0; iInput<fManager->GetNinputs(); iInput++){
235 // get run loader and ZDC loader
236 AliRunLoader* runLoader =
237 AliRunLoader::GetRunLoader(fManager->GetInputFolderName(iInput));
238 AliLoader* loader = runLoader->GetLoader("ZDCLoader");
239 if(!loader) continue;
242 loader->LoadSDigits();
243 TTree* treeS = loader->TreeS();
246 AliZDCSDigit* psdigit = &sdigit;
247 treeS->SetBranchAddress("ZDC", &psdigit);
250 for(Int_t iSDigit=0; iSDigit<treeS->GetEntries(); iSDigit++){
251 treeS->GetEntry(iSDigit);
253 if(!psdigit) continue;
254 if((sdigit.GetSector(1) < 0) || (sdigit.GetSector(1) > 4)){
255 AliError(Form("\nsector[0] = %d, sector[1] = %d\n",
256 sdigit.GetSector(0), sdigit.GetSector(1)));
259 // Checking if signal is inside ADC gate
260 if(iSDigit==0) signalTime0 = sdigit.GetTrackTime();
262 // Assuming a signal lenght of 20 ns, signal is in gate if
263 // signal ENDS in signalTime0+50. -> sdigit.GetTrackTime()+20<=signalTime0+50.
264 if(sdigit.GetTrackTime()<=signalTime0+30.) fIsSignalInADCGate = kTRUE;
265 if(sdigit.GetTrackTime()>signalTime0+30.){
266 fIsSignalInADCGate = kFALSE;
267 // Vedi quaderno per spiegazione approx. usata
268 // nel calcolo della fraz. di segnale perso
269 fFracLostSignal = (sdigit.GetTrackTime()-30)*(sdigit.GetTrackTime()-30)/280.;
273 Float_t sdSignal = sdigit.GetLightPM();
274 if(fIsSignalInADCGate == kFALSE){
275 AliDebug(2,Form("\t Signal time %f -> fraction %f of ZDC signal for det.(%d, %d) out of ADC gate\n",
276 sdigit.GetTrackTime(),fFracLostSignal,sdigit.GetSector(0),sdigit.GetSector(1)));
277 sdSignal = (1-fFracLostSignal)*sdSignal;
280 pm[(sdigit.GetSector(0))-1][sdigit.GetSector(1)] += sdigit.GetLightPM();
282 /*printf("\t Detector %d, Tower %d -> pm[%d][%d] = %.0f \n",
283 sdigit.GetSector(0), sdigit.GetSector(1),sdigit.GetSector(0)-1,
284 sdigit.GetSector(1), pm[sdigit.GetSector(0)-1][sdigit.GetSector(1)]);
289 loader->UnloadSDigits();
291 // get the impact parameter and the number of spectators in case of hijing
292 if(!runLoader->GetAliRun()) runLoader->LoadgAlice();
293 AliHeader* header = runLoader->GetHeader();
294 if(!header) continue;
295 AliGenEventHeader* genHeader = header->GenEventHeader();
296 if(!genHeader) continue;
297 if(!genHeader->InheritsFrom(AliGenHijingEventHeader::Class())) continue;
299 if(fSpectators2Track==kTRUE){
300 impPar = ((AliGenHijingEventHeader*) genHeader)->ImpactParameter();
301 specNProj = ((AliGenHijingEventHeader*) genHeader)->ProjSpectatorsn();
302 specPProj = ((AliGenHijingEventHeader*) genHeader)->ProjSpectatorsp();
303 specNTarg = ((AliGenHijingEventHeader*) genHeader)->TargSpectatorsn();
304 specPTarg = ((AliGenHijingEventHeader*) genHeader)->TargSpectatorsp();
305 printf("\n\t AliZDCDigitizer: b = %1.2f fm\n"
306 " \t PROJ.: #spectator n %d, #spectator p %d\n"
307 " \t TARG.: #spectator n %d, #spectator p %d\n",
308 impPar, specNProj, specPProj, specNTarg, specPTarg);
313 // Applying fragmentation algorithm and adding spectator signal
314 if(fSpectators2Track==kTRUE && impPar) {
315 Int_t freeSpecNProj, freeSpecPProj;
316 Fragmentation(impPar, specNProj, specPProj, freeSpecNProj, freeSpecPProj);
317 Int_t freeSpecNTarg, freeSpecPTarg;
318 Fragmentation(impPar, specNTarg, specPTarg, freeSpecNTarg, freeSpecPTarg);
319 SpectatorSignal(1, freeSpecNProj, pm);
320 printf("\t AliZDCDigitizer -> Adding signal for %d PROJ free spectator n",freeSpecNProj);
321 SpectatorSignal(2, freeSpecPProj, pm);
322 printf(" and %d free spectator p\n",freeSpecPProj);
323 SpectatorSignal(3, freeSpecNTarg, pm);
324 printf("\t AliZDCDigitizer -> Adding signal for %d TARG free spectator n",freeSpecNTarg);
325 SpectatorSignal(4, freeSpecPTarg, pm);
326 printf(" and %d free spectator p\n\n",freeSpecPTarg);
330 // get the output run loader and loader
331 AliRunLoader* runLoader =
332 AliRunLoader::GetRunLoader(fManager->GetOutputFolderName());
333 AliLoader* loader = runLoader->GetLoader("ZDCLoader");
335 AliError("no ZDC loader found");
339 // create the output tree
340 const char* mode = "update";
341 if(runLoader->GetEventNumber() == 0) mode = "recreate";
342 loader->LoadDigits(mode);
343 loader->MakeTree("D");
344 TTree* treeD = loader->TreeD();
346 AliZDCDigit* pdigit = &digit;
347 const Int_t kBufferSize = 4000;
348 treeD->Branch("ZDC", "AliZDCDigit", &pdigit, kBufferSize);
352 Int_t digi[2], digioot[2];
353 for(sector[0]=1; sector[0]<6; sector[0]++){
354 for(sector[1]=0; sector[1]<5; sector[1]++){
355 if((sector[0]==3) && ((sector[1]<1) || (sector[1]>2))) continue;
356 for(Int_t res=0; res<2; res++){
357 digi[res] = Phe2ADCch(sector[0], sector[1], pm[sector[0]-1][sector[1]], res)
358 + Pedestal(sector[0], sector[1], res);
360 new(pdigit) AliZDCDigit(sector, digi);
364 //printf("\t DIGIT added -> det %d quad %d - digi[0,1] = [%d, %d]\n",
365 // sector[0], sector[1], digi[0], digi[1]); // Chiara debugging!
368 } // Loop over detector
369 // Adding in-time digits for 2 reference PTM signals (after signal ch.)
370 // (for the moment the ref. signal is completely invented assuming a PMgain of 5*10^4!)
374 // Reference signal are set to 100 (high gain chain) and 800 (low gain chain)
375 if(fIsCalibration==0) {sigRef[0]=100; sigRef[1]=800;}
376 else {sigRef[0]=0; sigRef[1]=0;} // calibration -> simulation of pedestal values
378 for(Int_t iref=0; iref<2; iref++){
379 sectorRef[0] = 3*iref+1;
380 for(Int_t res=0; res<2; res++){
381 sigRef[res] += Pedestal(sectorRef[0], sectorRef[1], res);
383 new(pdigit) AliZDCDigit(sectorRef, sigRef);
387 //printf("\t RefDigit added -> det = %d, quad = %d - digi[0,1] = [%d, %d]\n",
388 // sectorRef[0], sectorRef[1], sigRef[0], sigRef[1]); // Chiara debugging!
391 // --- Adding digits for out-of-time channels after signal digits
392 for(sector[0]=1; sector[0]<6; sector[0]++){
393 for(sector[1]=0; sector[1]<5; sector[1]++){
394 if((sector[0]==3) && ((sector[1]<1) || (sector[1]>2))) continue;
395 for(Int_t res=0; res<2; res++){
396 digioot[res] = Pedestal(sector[0], sector[1], res); // out-of-time ADCs
398 new(pdigit) AliZDCDigit(sector, digioot);
402 //printf("\t DIGIToot added -> det = %d, quad = %d - digi[0,1] = [%d, %d]\n",
403 // sector[0], sector[1], digioot[0], digioot[1]); // Chiara debugging!
406 // Adding out-of-time digits for 2 reference PTM signals (after out-of-time ch.)
408 for(Int_t iref=0; iref<2; iref++){
409 sectorRef[0] = 3*iref+1;
410 for(Int_t res=0; res<2; res++){
411 sigRefoot[res] = Pedestal(sectorRef[0], sectorRef[1], res);
413 new(pdigit) AliZDCDigit(sectorRef, sigRefoot);
416 //printf("\t RefDigitoot added -> det = %d, quad = %d - digi[0,1] = [%d, %d]\n",
417 // sectorRef[0], sectorRef[1], sigRefoot[0], sigRefoot[1]); // Chiara debugging!
420 //printf("\t AliZDCDigitizer -> TreeD has %d entries\n",(Int_t) treeD->GetEntries());
422 // write the output tree
423 loader->WriteDigits("OVERWRITE");
424 loader->UnloadDigits();
428 //_____________________________________________________________________________
429 void AliZDCDigitizer::Fragmentation(Float_t impPar, Int_t specN, Int_t specP,
430 Int_t &freeSpecN, Int_t &freeSpecP) const
432 // simulate fragmentation of spectators
434 AliZDCFragment frag(impPar);
436 // Fragments generation
438 Int_t nAlpha = frag.GetNalpha();
441 frag.AttachNeutrons();
442 Int_t ztot = frag.GetZtot();
443 Int_t ntot = frag.GetNtot();
445 // Removing fragments and alpha pcs
446 freeSpecN = specN-ntot-2*nAlpha;
447 freeSpecP = specP-ztot-2*nAlpha;
449 // Removing deuterons
450 Int_t ndeu = (Int_t) (freeSpecN*frag.DeuteronNumber());
454 if(freeSpecN<0) freeSpecN=0;
455 if(freeSpecP<0) freeSpecP=0;
456 AliDebug(2, Form("FreeSpn = %d, FreeSpp = %d", freeSpecN, freeSpecP));
459 //_____________________________________________________________________________
460 void AliZDCDigitizer::SpectatorSignal(Int_t SpecType, Int_t numEvents,
461 Float_t pm[5][5]) const
463 // add signal of the spectators
465 TFile *specSignalFile = TFile::Open("$ALICE_ROOT/ZDC/SpectatorSignal.root");
466 if(!specSignalFile || !specSignalFile->IsOpen()) {
467 AliError((" Opening file $ALICE_ROOT/ZDC/SpectatorSignal.root failed\n"));
471 TNtuple* zdcSignal=0x0;
473 Float_t sqrtS = 2*fBeamEnergy;
475 if(TMath::Abs(sqrtS-5500) < 100.){
476 specSignalFile->cd("energy5500");
478 if(SpecType == 1) { // --- Signal for projectile spectator neutrons
479 specSignalFile->GetObject("energy5500/ZNCSignal;1",zdcSignal);
480 if(!zdcSignal) AliError(" PROBLEM!!! Can't retrieve ZNCSignal from SpectatorSignal.root file");
482 else if(SpecType == 2) { // --- Signal for projectile spectator protons
483 specSignalFile->GetObject("energy5500/ZPCSignal;1",zdcSignal);
484 if(!zdcSignal) AliError(" PROBLEM!!! Can't retrieve ZPCSignal from SpectatorSignal.root file");
486 else if(SpecType == 3) { // --- Signal for target spectator neutrons
487 specSignalFile->GetObject("energy5500/ZNASignal;1",zdcSignal);
488 if(!zdcSignal) AliError(" PROBLEM!!! Can't retrieve ZNASignal from SpectatorSignal.root file");
490 else if(SpecType == 4) { // --- Signal for target spectator protons
491 specSignalFile->GetObject("energy5500/ZPASignal;1",zdcSignal);
492 if(!zdcSignal) AliError(" PROBLEM!!! Can't retrieve ZPASignal from SpectatorSignal.root file");
495 else if(TMath::Abs(sqrtS-2760) < 100.){
496 specSignalFile->cd("energy2760");
498 if(SpecType == 1) { // --- Signal for projectile spectator neutrons
499 specSignalFile->GetObject("energy2760/ZNCSignal;1",zdcSignal);
500 if(!zdcSignal) AliError(" PROBLEM!!! Can't retrieve ZNCSignal from SpectatorSignal.root file");
502 else if(SpecType == 2) { // --- Signal for projectile spectator protons
503 specSignalFile->GetObject("energy2760/ZPCSignal;1",zdcSignal);
504 if(!zdcSignal) AliError(" PROBLEM!!! Can't retrieve ZPCSignal from SpectatorSignal.root file");
506 else if(SpecType == 3) { // --- Signal for target spectator neutrons
507 specSignalFile->GetObject("energy2760/ZNASignal;1",zdcSignal);
508 if(!zdcSignal) AliError(" PROBLEM!!! Can't retrieve ZNASignal from SpectatorSignal.root file");
510 else if(SpecType == 4) { // --- Signal for target spectator protons
511 specSignalFile->GetObject("energy2760/ZPASignal;1",zdcSignal);
512 if(!zdcSignal) AliError(" PROBLEM!!! Can't retrieve ZPASignal from SpectatorSignal.root file");
517 printf("\n No spectator signal available for ZDC digitization\n");
521 Int_t nentries = (Int_t) zdcSignal->GetEntries();
524 Int_t pl, i, k, iev=0, rnd[125], volume[2];
525 for(pl=0;pl<125;pl++) rnd[pl] = 0;
526 if(numEvents > 125) {
527 AliDebug(2,Form("numEvents (%d) is larger than 125", numEvents));
530 for(pl=0;pl<numEvents;pl++){
531 rnd[pl] = (Int_t) (9999*gRandom->Rndm());
532 if(rnd[pl] >= 9999) rnd[pl] = 9998;
533 //printf(" rnd[%d] = %d\n",pl,rnd[pl]);
535 // Sorting vector in ascending order with C function QSORT
536 qsort((void*)rnd,numEvents,sizeof(Int_t),comp);
538 for(i=0; i<nentries; i++){
539 zdcSignal->GetEvent(i);
540 entry = zdcSignal->GetArgs();
541 if(entry[0] == rnd[iev]){
542 for(k=0; k<2; k++) volume[k] = (Int_t) entry[k+1];
544 Float_t lightQ = entry[7];
545 Float_t lightC = entry[8];
547 if(volume[0] != 3) { // ZN or ZP
548 pm[volume[0]-1][0] += lightC;
549 pm[volume[0]-1][volume[1]] += lightQ;
550 //printf("\n pm[%d][0] = %.0f, pm[%d][%d] = %.0f\n",(volume[0]-1),pm[volume[0]-1][0],
551 // (volume[0]-1),volume[1],pm[volume[0]-1][volume[1]]);
554 if(volume[1] == 1) pm[2][1] += lightC; // ZEM 1
555 else pm[2][2] += lightQ; // ZEM 2
556 //printf("\n pm[2][1] = %.0f, pm[2][2] = %.0f\n",pm[2][1],pm[2][2]);
559 else if(entry[0] > rnd[iev]){
564 }while(iev<numEvents);
566 specSignalFile->Close();
567 delete specSignalFile;
571 //_____________________________________________________________________________
572 Int_t AliZDCDigitizer::Phe2ADCch(Int_t Det, Int_t Quad, Float_t Light,
575 // Evaluation of the ADC channel corresponding to the light yield Light
576 Int_t vADCch = (Int_t) (Light * fPMGain[Det-1][Quad] * fADCRes[Res]);
578 //printf("\t Phe2ADCch -> det %d quad %d - PMGain[%d][%d] %1.0f phe %1.2f ADC %d\n",
579 // Det,Quad,Det-1,Quad,fPMGain[Det-1][Quad],Light,vADCch);
584 //_____________________________________________________________________________
585 Int_t AliZDCDigitizer::Pedestal(Int_t Det, Int_t Quad, Int_t Res) const
587 // Returns a pedestal for detector det, PM quad, channel with res.
591 if(fIsCalibration == 0){
592 Int_t index=0, kNch=24;
594 if(Det==1) index = Quad+kNch*Res; // ZNC
595 else if(Det==2) index = (Quad+5)+kNch*Res; // ZPC
596 else if(Det==3) index = (Quad+9)+kNch*Res; // ZEM
597 else if(Det==4) index = (Quad+12)+kNch*Res; // ZNA
598 else if(Det==5) index = (Quad+17)+kNch*Res; // ZPA
600 else index = (Det-1)/3+22+kNch*Res; // Reference PMs
602 Float_t meanPed = fPedData->GetMeanPed(index);
603 Float_t pedWidth = fPedData->GetMeanPedWidth(index);
604 pedValue = gRandom->Gaus(meanPed,pedWidth);
606 /*printf("\t AliZDCDigitizer::Pedestal -> det %d quad %d res %d - Ped[%d] = %d\n",
607 Det, Quad, Res, index,(Int_t) pedValue); // Chiara debugging!
610 // To create calibration object
612 if(Res == 0) pedValue = gRandom->Gaus((35.+10.*gRandom->Rndm()),(0.5+0.2*gRandom->Rndm())); //High gain
613 else pedValue = gRandom->Gaus((250.+100.*gRandom->Rndm()),(3.5+2.*gRandom->Rndm())); //Low gain
616 return (Int_t) pedValue;
619 //_____________________________________________________________________________
620 AliCDBStorage* AliZDCDigitizer::SetStorage(const char *uri)
623 Bool_t deleteManager = kFALSE;
625 AliCDBManager *manager = AliCDBManager::Instance();
626 AliCDBStorage *defstorage = manager->GetDefaultStorage();
628 if(!defstorage || !(defstorage->Contains("ZDC"))){
629 AliWarning("No default storage set or default storage doesn't contain ZDC!");
630 manager->SetDefaultStorage(uri);
631 deleteManager = kTRUE;
634 AliCDBStorage *storage = manager->GetDefaultStorage();
637 AliCDBManager::Instance()->UnsetDefaultStorage();
638 defstorage = 0; // the storage is killed by AliCDBManager::Instance()->Destroy()
644 //_____________________________________________________________________________
645 AliZDCPedestals* AliZDCDigitizer::GetPedData() const
648 // Getting pedestal calibration object for ZDC set
650 AliCDBEntry *entry = AliCDBManager::Instance()->Get("ZDC/Calib/Pedestals");
651 if(!entry) AliFatal("No calibration data loaded!");
653 AliZDCPedestals *calibdata = dynamic_cast<AliZDCPedestals*> (entry->GetObject());
654 if(!calibdata) AliFatal("Wrong calibration object in calibration file!");