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 if(fBeamEnergy != 0){
156 for(Int_t j = 0; j < 5; j++){
157 fPMGain[0][j] = 1.515831*(661.444/fBeamEnergy+0.000740671)*10000000;
158 fPMGain[1][j] = 0.674234*(864.350/fBeamEnergy+0.00234375)*10000000;
159 fPMGain[3][j] = 1.350938*(661.444/fBeamEnergy+0.000740671)*10000000;
160 fPMGain[4][j] = 0.678597*(864.350/fBeamEnergy+0.00234375)*10000000;
162 fPMGain[2][1] = 0.869654*(1.32312-0.000101515*fBeamEnergy)*10000000;
163 fPMGain[2][2] = 1.030883*(1.32312-0.000101515*fBeamEnergy)*10000000;
164 AliInfo(Form(" PMT gains for p-p @ %1.0f+%1.0f GeV: ZN(%1.0f), ZP(%1.0f), ZEM(%1.0f)\n",
165 fBeamEnergy, fBeamEnergy, fPMGain[0][0], fPMGain[1][0], fPMGain[2][1]));
168 else if((beamType.CompareTo("A-A")) == 0){
169 // PTM gains for Pb-Pb @ 2.7+2.7 A TeV ***************
170 // rescaled for Pb-Pb @ 1.38+1.38 A TeV ***************
171 // New correction coefficients for PMT gains needed
172 // to reproduce experimental spectra (from Grazia Jul 2010)
173 Float_t scalGainFactor = fBeamEnergy/2760.;
174 for(Int_t j = 0; j < 5; j++){
175 fPMGain[0][j] = 50000./scalGainFactor;
176 fPMGain[1][j] = 100000./scalGainFactor;
177 fPMGain[2][j] = 100000./scalGainFactor;
178 fPMGain[3][j] = 50000./scalGainFactor;
179 fPMGain[4][j] = 100000./scalGainFactor;
181 AliInfo(Form(" PMT gains for Pb-Pb @ %1.0f+%1.0f A GeV: ZN(%1.0f), ZP(%1.0f), ZEM(%1.0f)\n",
182 fBeamEnergy, fBeamEnergy, fPMGain[0][0], fPMGain[1][0], fPMGain[2][1]));
186 fADCRes[0] = 0.0000008; // ADC Resolution high gain: 200 fC/adcCh
187 fADCRes[1] = 0.0000064; // ADC Resolution low gain: 25 fC/adcCh
192 //____________________________________________________________________________
193 void AliZDCDigitizer::Exec(Option_t* /*option*/)
195 // Execute digitization
197 // ------------------------------------------------------------
198 // !!! 2nd ZDC set added
199 // *** 1st 3 arrays are digits from REAL (simulated) hits
200 // *** last 2 are copied from simulated digits
201 // --- pm[0][...] = light in ZN side C [C, Q1, Q2, Q3, Q4]
202 // --- pm[1][...] = light in ZP side C [C, Q1, Q2, Q3, Q4]
203 // --- pm[2][...] = light in ZEM [x, 1, 2, x, x]
204 // --- pm[3][...] = light in ZN side A [C, Q1, Q2, Q3, Q4] ->NEW!
205 // --- pm[4][...] = light in ZP side A [C, Q1, Q2, Q3, Q4] ->NEW!
206 // ------------------------------------------------------------
208 for(Int_t iSector1=0; iSector1<5; iSector1++)
209 for(Int_t iSector2=0; iSector2<5; iSector2++){
210 pm[iSector1][iSector2] = 0;
213 // ------------------------------------------------------------
214 // ### Out of time ADC added (22 channels)
215 // --- same codification as for signal PTMs (see above)
216 // ------------------------------------------------------------
218 for(Int_t iSector1=0; iSector1<5; iSector1++)
219 for(Int_t iSector2=0; iSector2<5; iSector2++){
220 pmoot[iSector1][iSector2] = 0;
223 // impact parameter and number of spectators
225 Int_t specNTarg = 0, specPTarg = 0;
226 Int_t specNProj = 0, specPProj = 0;
227 Float_t signalTime0 = 0.;
229 // loop over input streams
230 for(Int_t iInput = 0; iInput<fManager->GetNinputs(); iInput++){
232 // get run loader and ZDC loader
233 AliRunLoader* runLoader =
234 AliRunLoader::GetRunLoader(fManager->GetInputFolderName(iInput));
235 AliLoader* loader = runLoader->GetLoader("ZDCLoader");
236 if(!loader) continue;
239 loader->LoadSDigits();
240 TTree* treeS = loader->TreeS();
243 AliZDCSDigit* psdigit = &sdigit;
244 treeS->SetBranchAddress("ZDC", &psdigit);
247 for(Int_t iSDigit=0; iSDigit<treeS->GetEntries(); iSDigit++){
248 treeS->GetEntry(iSDigit);
250 if(!psdigit) continue;
251 if((sdigit.GetSector(1) < 0) || (sdigit.GetSector(1) > 4)){
252 AliError(Form("\nsector[0] = %d, sector[1] = %d\n",
253 sdigit.GetSector(0), sdigit.GetSector(1)));
256 // Checking if signal is inside ADC gate
257 if(iSDigit==0) signalTime0 = sdigit.GetTrackTime();
259 // Assuming a signal lenght of 20 ns, signal is in gate if
260 // signal ENDS in signalTime0+50. -> sdigit.GetTrackTime()+20<=signalTime0+50.
261 if(sdigit.GetTrackTime()<=signalTime0+30.) fIsSignalInADCGate = kTRUE;
262 if(sdigit.GetTrackTime()>signalTime0+30.){
263 fIsSignalInADCGate = kFALSE;
264 // Vedi quaderno per spiegazione approx. usata
265 // nel calcolo della fraz. di segnale perso
266 fFracLostSignal = (sdigit.GetTrackTime()-30)*(sdigit.GetTrackTime()-30)/280.;
270 Float_t sdSignal = sdigit.GetLightPM();
271 if(fIsSignalInADCGate == kFALSE){
272 AliDebug(2,Form("\t Signal time %f -> fraction %f of ZDC signal for det.(%d, %d) out of ADC gate\n",
273 sdigit.GetTrackTime(),fFracLostSignal,sdigit.GetSector(0),sdigit.GetSector(1)));
274 sdSignal = (1-fFracLostSignal)*sdSignal;
277 pm[(sdigit.GetSector(0))-1][sdigit.GetSector(1)] += sdigit.GetLightPM();
279 /*printf("\t Detector %d, Tower %d -> pm[%d][%d] = %.0f \n",
280 sdigit.GetSector(0), sdigit.GetSector(1),sdigit.GetSector(0)-1,
281 sdigit.GetSector(1), pm[sdigit.GetSector(0)-1][sdigit.GetSector(1)]);
286 loader->UnloadSDigits();
288 // get the impact parameter and the number of spectators in case of hijing
289 if(!runLoader->GetAliRun()) runLoader->LoadgAlice();
290 AliHeader* header = runLoader->GetHeader();
291 if(!header) continue;
292 AliGenEventHeader* genHeader = header->GenEventHeader();
293 if(!genHeader) continue;
294 if(!genHeader->InheritsFrom(AliGenHijingEventHeader::Class())) continue;
296 if(fSpectators2Track==kTRUE){
297 impPar = ((AliGenHijingEventHeader*) genHeader)->ImpactParameter();
298 specNProj = ((AliGenHijingEventHeader*) genHeader)->ProjSpectatorsn();
299 specPProj = ((AliGenHijingEventHeader*) genHeader)->ProjSpectatorsp();
300 specNTarg = ((AliGenHijingEventHeader*) genHeader)->TargSpectatorsn();
301 specPTarg = ((AliGenHijingEventHeader*) genHeader)->TargSpectatorsp();
302 printf("\n\t AliZDCDigitizer: b = %1.2f fm\n"
303 " \t PROJ.: #spectator n %d, #spectator p %d\n"
304 " \t TARG.: #spectator n %d, #spectator p %d\n",
305 impPar, specNProj, specPProj, specNTarg, specPTarg);
310 // Applying fragmentation algorithm and adding spectator signal
311 if(fSpectators2Track==kTRUE && impPar) {
312 Int_t freeSpecNProj, freeSpecPProj;
313 Fragmentation(impPar, specNProj, specPProj, freeSpecNProj, freeSpecPProj);
314 Int_t freeSpecNTarg, freeSpecPTarg;
315 Fragmentation(impPar, specNTarg, specPTarg, freeSpecNTarg, freeSpecPTarg);
316 SpectatorSignal(1, freeSpecNProj, pm);
317 printf("\t AliZDCDigitizer -> Adding signal for %d PROJ free spectator n",freeSpecNProj);
318 SpectatorSignal(2, freeSpecPProj, pm);
319 printf(" and %d free spectator p\n",freeSpecPProj);
320 SpectatorSignal(3, freeSpecNTarg, pm);
321 printf("\t AliZDCDigitizer -> Adding signal for %d TARG free spectator n",freeSpecNTarg);
322 SpectatorSignal(4, freeSpecPTarg, pm);
323 printf(" and %d free spectator p\n\n",freeSpecPTarg);
327 // get the output run loader and loader
328 AliRunLoader* runLoader =
329 AliRunLoader::GetRunLoader(fManager->GetOutputFolderName());
330 AliLoader* loader = runLoader->GetLoader("ZDCLoader");
332 AliError("no ZDC loader found");
336 // create the output tree
337 const char* mode = "update";
338 if(runLoader->GetEventNumber() == 0) mode = "recreate";
339 loader->LoadDigits(mode);
340 loader->MakeTree("D");
341 TTree* treeD = loader->TreeD();
343 AliZDCDigit* pdigit = &digit;
344 const Int_t kBufferSize = 4000;
345 treeD->Branch("ZDC", "AliZDCDigit", &pdigit, kBufferSize);
349 Int_t digi[2], digioot[2];
350 for(sector[0]=1; sector[0]<6; sector[0]++){
351 for(sector[1]=0; sector[1]<5; sector[1]++){
352 if((sector[0]==3) && ((sector[1]<1) || (sector[1]>2))) continue;
353 for(Int_t res=0; res<2; res++){
354 digi[res] = Phe2ADCch(sector[0], sector[1], pm[sector[0]-1][sector[1]], res)
355 + Pedestal(sector[0], sector[1], res);
357 new(pdigit) AliZDCDigit(sector, digi);
361 //printf("\t DIGIT added -> det %d quad %d - digi[0,1] = [%d, %d]\n",
362 // sector[0], sector[1], digi[0], digi[1]); // Chiara debugging!
365 } // Loop over detector
366 // Adding in-time digits for 2 reference PTM signals (after signal ch.)
367 // (for the moment the ref. signal is completely invented assuming a PMgain of 5*10^4!)
371 // Reference signal are set to 100 (high gain chain) and 800 (low gain chain)
372 if(fIsCalibration==0) {sigRef[0]=100; sigRef[1]=800;}
373 else {sigRef[0]=0; sigRef[1]=0;} // calibration -> simulation of pedestal values
375 for(Int_t iref=0; iref<2; iref++){
376 sectorRef[0] = 3*iref+1;
377 for(Int_t res=0; res<2; res++){
378 sigRef[res] += Pedestal(sectorRef[0], sectorRef[1], res);
380 new(pdigit) AliZDCDigit(sectorRef, sigRef);
384 //printf("\t RefDigit added -> det = %d, quad = %d - digi[0,1] = [%d, %d]\n",
385 // sectorRef[0], sectorRef[1], sigRef[0], sigRef[1]); // Chiara debugging!
388 // --- Adding digits for out-of-time channels after signal digits
389 for(sector[0]=1; sector[0]<6; sector[0]++){
390 for(sector[1]=0; sector[1]<5; sector[1]++){
391 if((sector[0]==3) && ((sector[1]<1) || (sector[1]>2))) continue;
392 for(Int_t res=0; res<2; res++){
393 digioot[res] = Pedestal(sector[0], sector[1], res); // out-of-time ADCs
395 new(pdigit) AliZDCDigit(sector, digioot);
399 //printf("\t DIGIToot added -> det = %d, quad = %d - digi[0,1] = [%d, %d]\n",
400 // sector[0], sector[1], digioot[0], digioot[1]); // Chiara debugging!
403 // Adding out-of-time digits for 2 reference PTM signals (after out-of-time ch.)
405 for(Int_t iref=0; iref<2; iref++){
406 sectorRef[0] = 3*iref+1;
407 for(Int_t res=0; res<2; res++){
408 sigRefoot[res] = Pedestal(sectorRef[0], sectorRef[1], res);
410 new(pdigit) AliZDCDigit(sectorRef, sigRefoot);
413 //printf("\t RefDigitoot added -> det = %d, quad = %d - digi[0,1] = [%d, %d]\n",
414 // sectorRef[0], sectorRef[1], sigRefoot[0], sigRefoot[1]); // Chiara debugging!
417 //printf("\t AliZDCDigitizer -> TreeD has %d entries\n",(Int_t) treeD->GetEntries());
419 // write the output tree
420 loader->WriteDigits("OVERWRITE");
421 loader->UnloadDigits();
425 //_____________________________________________________________________________
426 void AliZDCDigitizer::Fragmentation(Float_t impPar, Int_t specN, Int_t specP,
427 Int_t &freeSpecN, Int_t &freeSpecP) const
429 // simulate fragmentation of spectators
431 AliZDCFragment frag(impPar);
433 // Fragments generation
435 Int_t nAlpha = frag.GetNalpha();
438 frag.AttachNeutrons();
439 Int_t ztot = frag.GetZtot();
440 Int_t ntot = frag.GetNtot();
442 // Removing fragments and alpha pcs
443 freeSpecN = specN-ntot-2*nAlpha;
444 freeSpecP = specP-ztot-2*nAlpha;
446 // Removing deuterons
447 Int_t ndeu = (Int_t) (freeSpecN*frag.DeuteronNumber());
451 if(freeSpecN<0) freeSpecN=0;
452 if(freeSpecP<0) freeSpecP=0;
453 AliDebug(2, Form("FreeSpn = %d, FreeSpp = %d", freeSpecN, freeSpecP));
456 //_____________________________________________________________________________
457 void AliZDCDigitizer::SpectatorSignal(Int_t SpecType, Int_t numEvents,
458 Float_t pm[5][5]) const
460 // add signal of the spectators
462 TFile *specSignalFile = TFile::Open("$ALICE_ROOT/ZDC/SpectatorSignal.root");
463 if(!specSignalFile || !specSignalFile->IsOpen()) {
464 AliError((" Opening file $ALICE_ROOT/ZDC/SpectatorSignal.root failed\n"));
468 TNtuple* zdcSignal=0x0;
470 Float_t sqrtS = 2*fBeamEnergy;
472 if(TMath::Abs(sqrtS-5500) < 100.){
473 specSignalFile->cd("energy5500");
475 if(SpecType == 1) { // --- Signal for projectile spectator neutrons
476 specSignalFile->GetObject("energy5500/ZNCSignal;1",zdcSignal);
477 if(!zdcSignal) AliError(" PROBLEM!!! Can't retrieve ZNCSignal from SpectatorSignal.root file");
479 else if(SpecType == 2) { // --- Signal for projectile spectator protons
480 specSignalFile->GetObject("energy5500/ZPCSignal;1",zdcSignal);
481 if(!zdcSignal) AliError(" PROBLEM!!! Can't retrieve ZPCSignal from SpectatorSignal.root file");
483 else if(SpecType == 3) { // --- Signal for target spectator neutrons
484 specSignalFile->GetObject("energy5500/ZNASignal;1",zdcSignal);
485 if(!zdcSignal) AliError(" PROBLEM!!! Can't retrieve ZNASignal from SpectatorSignal.root file");
487 else if(SpecType == 4) { // --- Signal for target spectator protons
488 specSignalFile->GetObject("energy5500/ZPASignal;1",zdcSignal);
489 if(!zdcSignal) AliError(" PROBLEM!!! Can't retrieve ZPASignal from SpectatorSignal.root file");
492 else if(TMath::Abs(sqrtS-2760) < 100.){
493 specSignalFile->cd("energy2760");
495 if(SpecType == 1) { // --- Signal for projectile spectator neutrons
496 specSignalFile->GetObject("energy2760/ZNCSignal;1",zdcSignal);
497 if(!zdcSignal) AliError(" PROBLEM!!! Can't retrieve ZNCSignal from SpectatorSignal.root file");
499 else if(SpecType == 2) { // --- Signal for projectile spectator protons
500 specSignalFile->GetObject("energy2760/ZPCSignal;1",zdcSignal);
501 if(!zdcSignal) AliError(" PROBLEM!!! Can't retrieve ZPCSignal from SpectatorSignal.root file");
503 else if(SpecType == 3) { // --- Signal for target spectator neutrons
504 specSignalFile->GetObject("energy2760/ZNASignal;1",zdcSignal);
505 if(!zdcSignal) AliError(" PROBLEM!!! Can't retrieve ZNASignal from SpectatorSignal.root file");
507 else if(SpecType == 4) { // --- Signal for target spectator protons
508 specSignalFile->GetObject("energy2760/ZPASignal;1",zdcSignal);
509 if(!zdcSignal) AliError(" PROBLEM!!! Can't retrieve ZPASignal from SpectatorSignal.root file");
514 printf("\n No spectator signal available for ZDC digitization\n");
518 Int_t nentries = (Int_t) zdcSignal->GetEntries();
521 Int_t pl, i, k, iev=0, rnd[125], volume[2];
522 for(pl=0;pl<125;pl++) rnd[pl] = 0;
523 if(numEvents > 125) {
524 AliDebug(2,Form("numEvents (%d) is larger than 125", numEvents));
527 for(pl=0;pl<numEvents;pl++){
528 rnd[pl] = (Int_t) (9999*gRandom->Rndm());
529 if(rnd[pl] >= 9999) rnd[pl] = 9998;
530 //printf(" rnd[%d] = %d\n",pl,rnd[pl]);
532 // Sorting vector in ascending order with C function QSORT
533 qsort((void*)rnd,numEvents,sizeof(Int_t),comp);
535 for(i=0; i<nentries; i++){
536 zdcSignal->GetEvent(i);
537 entry = zdcSignal->GetArgs();
538 if(entry[0] == rnd[iev]){
539 for(k=0; k<2; k++) volume[k] = (Int_t) entry[k+1];
541 Float_t lightQ = entry[7];
542 Float_t lightC = entry[8];
544 if(volume[0] != 3) { // ZN or ZP
545 pm[volume[0]-1][0] += lightC;
546 pm[volume[0]-1][volume[1]] += lightQ;
547 //printf("\n pm[%d][0] = %.0f, pm[%d][%d] = %.0f\n",(volume[0]-1),pm[volume[0]-1][0],
548 // (volume[0]-1),volume[1],pm[volume[0]-1][volume[1]]);
551 if(volume[1] == 1) pm[2][1] += lightC; // ZEM 1
552 else pm[2][2] += lightQ; // ZEM 2
553 //printf("\n pm[2][1] = %.0f, pm[2][2] = %.0f\n",pm[2][1],pm[2][2]);
556 else if(entry[0] > rnd[iev]){
561 }while(iev<numEvents);
563 specSignalFile->Close();
564 delete specSignalFile;
568 //_____________________________________________________________________________
569 Int_t AliZDCDigitizer::Phe2ADCch(Int_t Det, Int_t Quad, Float_t Light,
572 // Evaluation of the ADC channel corresponding to the light yield Light
573 Int_t vADCch = (Int_t) (Light * fPMGain[Det-1][Quad] * fADCRes[Res]);
575 //printf("\t Phe2ADCch -> det %d quad %d - PMGain[%d][%d] %1.0f phe %1.2f ADC %d\n",
576 // Det,Quad,Det-1,Quad,fPMGain[Det-1][Quad],Light,vADCch);
581 //_____________________________________________________________________________
582 Int_t AliZDCDigitizer::Pedestal(Int_t Det, Int_t Quad, Int_t Res) const
584 // Returns a pedestal for detector det, PM quad, channel with res.
588 if(fIsCalibration == 0){
589 Int_t index=0, kNch=24;
591 if(Det==1) index = Quad+kNch*Res; // ZNC
592 else if(Det==2) index = (Quad+5)+kNch*Res; // ZPC
593 else if(Det==3) index = (Quad+9)+kNch*Res; // ZEM
594 else if(Det==4) index = (Quad+12)+kNch*Res; // ZNA
595 else if(Det==5) index = (Quad+17)+kNch*Res; // ZPA
597 else index = (Det-1)/3+22+kNch*Res; // Reference PMs
599 Float_t meanPed = fPedData->GetMeanPed(index);
600 Float_t pedWidth = fPedData->GetMeanPedWidth(index);
601 pedValue = gRandom->Gaus(meanPed,pedWidth);
603 /*printf("\t AliZDCDigitizer::Pedestal -> det %d quad %d res %d - Ped[%d] = %d\n",
604 Det, Quad, Res, index,(Int_t) pedValue); // Chiara debugging!
607 // To create calibration object
609 if(Res == 0) pedValue = gRandom->Gaus((35.+10.*gRandom->Rndm()),(0.5+0.2*gRandom->Rndm())); //High gain
610 else pedValue = gRandom->Gaus((250.+100.*gRandom->Rndm()),(3.5+2.*gRandom->Rndm())); //Low gain
613 return (Int_t) pedValue;
616 //_____________________________________________________________________________
617 AliCDBStorage* AliZDCDigitizer::SetStorage(const char *uri)
620 Bool_t deleteManager = kFALSE;
622 AliCDBManager *manager = AliCDBManager::Instance();
623 AliCDBStorage *defstorage = manager->GetDefaultStorage();
625 if(!defstorage || !(defstorage->Contains("ZDC"))){
626 AliWarning("No default storage set or default storage doesn't contain ZDC!");
627 manager->SetDefaultStorage(uri);
628 deleteManager = kTRUE;
631 AliCDBStorage *storage = manager->GetDefaultStorage();
634 AliCDBManager::Instance()->UnsetDefaultStorage();
635 defstorage = 0; // the storage is killed by AliCDBManager::Instance()->Destroy()
641 //_____________________________________________________________________________
642 AliZDCPedestals* AliZDCDigitizer::GetPedData() const
645 // Getting pedestal calibration object for ZDC set
647 AliCDBEntry *entry = AliCDBManager::Instance()->Get("ZDC/Calib/Pedestals");
648 if(!entry) AliFatal("No calibration data loaded!");
650 AliZDCPedestals *calibdata = dynamic_cast<AliZDCPedestals*> (entry->GetObject());
651 if(!calibdata) AliFatal("Wrong calibration object in calibration file!");