]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ZDC/AliZDCDigitizer.cxx
9fe90becae79c67756042a6cd828a6961b0e13a3
[u/mrichter/AliRoot.git] / ZDC / AliZDCDigitizer.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 //                      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"
39 #include "AliGRPObject.h"
40 #include "AliCDBManager.h"
41 #include "AliCDBEntry.h"
42 #include "AliZDCSDigit.h"
43 #include "AliZDCDigit.h"
44 #include "AliZDCFragment.h"
45 #include "AliZDCv3.h"
46 #include "AliZDCDigitizer.h"
47
48 class AliCDBStorage;
49 class AliZDCPedestals;
50
51 ClassImp(AliZDCDigitizer)
52
53
54 //____________________________________________________________________________
55 AliZDCDigitizer::AliZDCDigitizer() :
56   fIsCalibration(0), 
57   fIsSignalInADCGate(kFALSE),
58   fFracLostSignal(0.),
59   fPedData(0), 
60   fSpectators2Track(kFALSE),
61   fBeamEnergy(0.)
62 {
63   // Default constructor    
64
65 }
66
67 //____________________________________________________________________________
68 AliZDCDigitizer::AliZDCDigitizer(AliRunDigitizer* manager):
69   AliDigitizer(manager),
70   fIsCalibration(0), //By default the simulation doesn't create calib. data
71   fIsSignalInADCGate(kFALSE),
72   fFracLostSignal(0.),
73   fPedData(GetPedData()), 
74   fSpectators2Track(kFALSE),
75   fBeamEnergy(0.)
76 {
77   // Get calibration data
78   if(fIsCalibration!=0) printf("\n\t AliZDCDigitizer -> Creating calibration data (pedestals)\n");
79   for(Int_t i=0; i<6; i++){
80     for(Int_t j=0; j<5; j++)
81        fPMGain[i][j] = 0.;
82   }
83 }
84
85 //____________________________________________________________________________
86 AliZDCDigitizer::~AliZDCDigitizer()
87 {
88 // Destructor
89 // Not implemented
90 }
91
92
93 //____________________________________________________________________________
94 AliZDCDigitizer::AliZDCDigitizer(const AliZDCDigitizer &digitizer):
95   AliDigitizer(),
96   fIsCalibration(digitizer.fIsCalibration),
97   fIsSignalInADCGate(digitizer.fIsSignalInADCGate),
98   fFracLostSignal(digitizer.fFracLostSignal),
99   fPedData(digitizer.fPedData),
100   fSpectators2Track(digitizer.fSpectators2Track),
101   fBeamEnergy(digitizer.fBeamEnergy)
102 {
103   // Copy constructor
104
105   for(Int_t i=0; i<6; i++){
106      for(Int_t j=0; j<5; j++){
107         fPMGain[i][j]   = digitizer.fPMGain[i][j];           
108      }
109   }
110   for(Int_t i=0; i<2; i++) fADCRes[i] = digitizer.fADCRes[i];
111
112
113 }
114
115 //____________________________________________________________________________
116 Bool_t AliZDCDigitizer::Init()
117 {
118   // Initialize the digitizer
119   
120   AliCDBEntry*  entry = AliCDBManager::Instance()->Get("GRP/GRP/Data");
121   AliGRPObject* grpData = 0x0;
122   if(entry){
123     TMap* m = dynamic_cast<TMap*>(entry->GetObject());  // old GRP entry
124     if(m){
125       //m->Print();
126       grpData = new AliGRPObject();
127       grpData->ReadValuesFromMap(m);
128     }
129     else{
130       grpData = dynamic_cast<AliGRPObject*>(entry->GetObject());  // new GRP entry
131     }
132     entry->SetOwner(0);
133     AliCDBManager::Instance()->UnloadFromCache("GRP/GRP/Data");
134   }
135   if(!grpData) AliError("No GRP entry found in OCDB!");
136   
137   TString beamType = grpData->GetBeamType();
138   if(beamType==AliGRPObject::GetInvalidString()){
139     AliError("\t UNKNOWN beam type from GRP obj -> PMT gains not set in ZDC digitizer!!!\n");
140   }
141   
142   fBeamEnergy = grpData->GetBeamEnergy();
143   if(fBeamEnergy==AliGRPObject::GetInvalidFloat()){
144     AliWarning("GRP/GRP/Data entry:  missing value for the beam energy ! Using 0.");
145     AliError("\t UNKNOWN beam type from GRP obj -> PMT gains not set in ZDC digitizer!!!\n");
146     fBeamEnergy = 0.;
147   }
148
149   if((beamType.CompareTo("P-P")) == 0 || (beamType.CompareTo("p-p")) == 0){
150     //PTM gains rescaled to beam energy for p-p
151     if(fBeamEnergy != 0){
152     for(Int_t j = 0; j < 5; j++){
153         fPMGain[0][j] = 1.515831*(661.444/fBeamEnergy+0.000740671)*10000000;
154         fPMGain[1][j] = 0.674234*(864.350/fBeamEnergy+0.00234375)*10000000;
155         fPMGain[3][j] = 1.350938*(661.444/fBeamEnergy+0.000740671)*10000000; 
156         fPMGain[4][j] = 0.678597*(864.350/fBeamEnergy+0.00234375)*10000000;
157     }
158     fPMGain[2][1] = 0.869654*(1.32312-0.000101515*fBeamEnergy)*10000000;
159     fPMGain[2][2] = 1.030883*(1.32312-0.000101515*fBeamEnergy)*10000000;
160     AliInfo(Form("    PMT gains for p-p @ %1.0f+%1.0f GeV: ZN(%1.0f), ZP(%1.0f), ZEM(%1.0f)\n",
161         fBeamEnergy, fBeamEnergy, fPMGain[0][0], fPMGain[1][0], fPMGain[2][0]));
162     }
163   }
164   else if((beamType.CompareTo("A-A")) == 0){
165     // PTM gains for Pb-Pb @ 2.7+2.7 A TeV ***************
166     // rescaled for Pb-Pb @ 1.38+1.38 A TeV ***************
167     // New correction coefficients for PMT gains needed
168     // to reproduce experimental spectra (from Grazia Jul 2010)
169     Float_t scalGainFactor = fBeamEnergy/2760.;
170     for(Int_t j = 0; j < 5; j++){
171        fPMGain[0][j] = 50000./scalGainFactor;            
172        fPMGain[1][j] = 100000./scalGainFactor;           
173        fPMGain[3][j] = 50000./scalGainFactor;            
174        fPMGain[4][j] = 100000./scalGainFactor;           
175        fPMGain[5][j] = 100000./scalGainFactor;    
176     }
177     AliInfo(Form("    PMT gains for Pb-Pb @ %1.0f+%1.0f A GeV: ZN(%1.0f), ZP(%1.0f), ZEM(%1.0f)\n",
178         fBeamEnergy, fBeamEnergy, fPMGain[0][0], fPMGain[1][0], fPMGain[2][0]));
179   }
180     
181   // ADC Caen V965
182   fADCRes[0] = 0.0000008; // ADC Resolution high gain: 200 fC/adcCh
183   fADCRes[1] = 0.0000064; // ADC Resolution low gain:  25  fC/adcCh
184
185   return kTRUE;
186 }
187
188 //____________________________________________________________________________
189 void AliZDCDigitizer::Exec(Option_t* /*option*/)
190 {
191   // Execute digitization
192
193   // ------------------------------------------------------------
194   // !!! 2nd ZDC set added 
195   // *** 1st 3 arrays are digits from REAL (simulated) hits
196   // *** last 2 are copied from simulated digits
197   // --- pm[0][...] = light in ZN side C  [C, Q1, Q2, Q3, Q4]
198   // --- pm[1][...] = light in ZP side C [C, Q1, Q2, Q3, Q4]
199   // --- pm[2][...] = light in ZEM [x, 1, 2, x, x]
200   // --- pm[3][...] = light in ZN side A [C, Q1, Q2, Q3, Q4] ->NEW!
201   // --- pm[4][...] = light in ZP side A [C, Q1, Q2, Q3, Q4] ->NEW!
202   // ------------------------------------------------------------
203   Float_t pm[5][5]; 
204   for(Int_t iSector1=0; iSector1<5; iSector1++) 
205     for(Int_t iSector2=0; iSector2<5; iSector2++){
206       pm[iSector1][iSector2] = 0;
207     }
208     
209   // ------------------------------------------------------------
210   // ### Out of time ADC added (22 channels)
211   // --- same codification as for signal PTMs (see above)
212   // ------------------------------------------------------------
213   Float_t pmoot[5][5];
214   for(Int_t iSector1=0; iSector1<5; iSector1++) 
215     for(Int_t iSector2=0; iSector2<5; iSector2++){
216       pmoot[iSector1][iSector2] = 0;
217     }
218
219   // impact parameter and number of spectators
220   Float_t impPar = -1;
221   Int_t specNTarg = 0, specPTarg = 0;
222   Int_t specNProj = 0, specPProj = 0;
223   Float_t signalTime0 = 0.;
224
225   // loop over input streams
226   for(Int_t iInput = 0; iInput<fManager->GetNinputs(); iInput++){
227
228     // get run loader and ZDC loader
229     AliRunLoader* runLoader = 
230       AliRunLoader::GetRunLoader(fManager->GetInputFolderName(iInput));
231     AliLoader* loader = runLoader->GetLoader("ZDCLoader");
232     if(!loader) continue;
233
234     // load sdigits
235     loader->LoadSDigits();
236     TTree* treeS = loader->TreeS();
237     if(!treeS) continue;
238     AliZDCSDigit sdigit;
239     AliZDCSDigit* psdigit = &sdigit;
240     treeS->SetBranchAddress("ZDC", &psdigit);
241
242     // loop over sdigits
243     for(Int_t iSDigit=0; iSDigit<treeS->GetEntries(); iSDigit++){
244       treeS->GetEntry(iSDigit);
245       //
246       if(!psdigit) continue;
247       if((sdigit.GetSector(1) < 0) || (sdigit.GetSector(1) > 4)){
248         AliError(Form("\nsector[0] = %d, sector[1] = %d\n", 
249                       sdigit.GetSector(0), sdigit.GetSector(1)));
250         continue;
251       }
252       // Checking if signal is inside ADC gate
253       if(iSDigit==0) signalTime0 = sdigit.GetTrackTime();
254       else{
255         // Assuming a signal lenght of 20 ns, signal is in gate if 
256         // signal ENDS in signalTime0+50. -> sdigit.GetTrackTime()+20<=signalTime0+50.
257         if(sdigit.GetTrackTime()<=signalTime0+30.) fIsSignalInADCGate = kTRUE;
258         if(sdigit.GetTrackTime()>signalTime0+30.){
259           fIsSignalInADCGate = kFALSE;
260           // Vedi quaderno per spiegazione approx. usata 
261           // nel calcolo della fraz. di segnale perso
262           fFracLostSignal = (sdigit.GetTrackTime()-30)*(sdigit.GetTrackTime()-30)/280.;
263         }
264       }
265       
266       Float_t sdSignal = sdigit.GetLightPM();
267       if(fIsSignalInADCGate == kFALSE){
268         AliDebug(2,Form("\t Signal time %f -> fraction %f of ZDC signal for det.(%d, %d) out of ADC gate\n",
269                 sdigit.GetTrackTime(),fFracLostSignal,sdigit.GetSector(0),sdigit.GetSector(1)));
270         sdSignal = (1-fFracLostSignal)*sdSignal;
271       }
272       
273       pm[(sdigit.GetSector(0))-1][sdigit.GetSector(1)] += sdigit.GetLightPM();
274       //Ch. debug
275       /*printf("\t Detector %d, Tower %d -> pm[%d][%d] = %.0f \n",
276           sdigit.GetSector(0), sdigit.GetSector(1),sdigit.GetSector(0)-1,
277           sdigit.GetSector(1), pm[sdigit.GetSector(0)-1][sdigit.GetSector(1)]); 
278       */
279       
280     }
281
282     loader->UnloadSDigits();
283
284     // get the impact parameter and the number of spectators in case of hijing
285     if(!runLoader->GetAliRun()) runLoader->LoadgAlice();
286     AliHeader* header = runLoader->GetHeader();
287     if(!header) continue;
288     AliGenEventHeader* genHeader = header->GenEventHeader();
289     if(!genHeader) continue;
290     if(!genHeader->InheritsFrom(AliGenHijingEventHeader::Class())) continue;
291     
292     if(fSpectators2Track==kTRUE){
293       impPar = ((AliGenHijingEventHeader*) genHeader)->ImpactParameter(); 
294       specNProj = ((AliGenHijingEventHeader*) genHeader)->ProjSpectatorsn();
295       specPProj = ((AliGenHijingEventHeader*) genHeader)->ProjSpectatorsp();
296       specNTarg = ((AliGenHijingEventHeader*) genHeader)->TargSpectatorsn();
297       specPTarg = ((AliGenHijingEventHeader*) genHeader)->TargSpectatorsp();
298       printf("\n\t AliZDCDigitizer: b = %1.2f fm\n"
299       " \t    PROJ.:  #spectator n %d, #spectator p %d\n"
300       " \t    TARG.:  #spectator n %d, #spectator p %d\n", 
301       impPar, specNProj, specPProj, specNTarg, specPTarg);
302     }
303     
304   }
305
306   // Applying fragmentation algorithm and adding spectator signal
307   if(fSpectators2Track==kTRUE && impPar) {
308     Int_t freeSpecNProj, freeSpecPProj;
309     Fragmentation(impPar, specNProj, specPProj, freeSpecNProj, freeSpecPProj);
310     Int_t freeSpecNTarg, freeSpecPTarg;
311     Fragmentation(impPar, specNTarg, specPTarg, freeSpecNTarg, freeSpecPTarg);
312     SpectatorSignal(1, freeSpecNProj, pm);
313     printf("\t AliZDCDigitizer -> Adding signal for %d PROJ free spectator n",freeSpecNProj);
314     SpectatorSignal(2, freeSpecPProj, pm);
315     printf(" and %d free spectator p\n",freeSpecPProj);
316     SpectatorSignal(3, freeSpecNTarg, pm);
317     printf("\t AliZDCDigitizer -> Adding signal for %d TARG free spectator n",freeSpecNTarg);
318     SpectatorSignal(4, freeSpecPTarg, pm);
319     printf(" and %d free spectator p\n\n",freeSpecPTarg);
320   }
321
322
323   // get the output run loader and loader
324   AliRunLoader* runLoader = 
325     AliRunLoader::GetRunLoader(fManager->GetOutputFolderName());
326   AliLoader* loader = runLoader->GetLoader("ZDCLoader");
327   if(!loader) {
328     AliError("no ZDC loader found");
329     return;
330   }
331
332   // create the output tree
333   const char* mode = "update";
334   if(runLoader->GetEventNumber() == 0) mode = "recreate";
335   loader->LoadDigits(mode);
336   loader->MakeTree("D");
337   TTree* treeD = loader->TreeD();
338   AliZDCDigit digit;
339   AliZDCDigit* pdigit = &digit;
340   const Int_t kBufferSize = 4000;
341   treeD->Branch("ZDC", "AliZDCDigit", &pdigit, kBufferSize);
342
343   // Create digits
344   Int_t sector[2];
345   Int_t digi[2], digioot[2];
346   for(sector[0]=1; sector[0]<6; sector[0]++){
347     for(sector[1]=0; sector[1]<5; sector[1]++){
348         if((sector[0]==3) && ((sector[1]<1) || (sector[1]>2))) continue;
349         for(Int_t res=0; res<2; res++){
350            digi[res] = Phe2ADCch(sector[0], sector[1], pm[sector[0]-1][sector[1]], res) 
351                     + Pedestal(sector[0], sector[1], res);
352         }
353         new(pdigit) AliZDCDigit(sector, digi);
354         treeD->Fill();
355
356         //Ch. debug
357         //printf("\t DIGIT added -> det %d quad %d - digi[0,1] = [%d, %d]\n",
358         //     sector[0], sector[1], digi[0], digi[1]); // Chiara debugging!
359         
360     }
361   } // Loop over detector
362   // Adding in-time digits for 2 reference PTM signals (after signal ch.)
363   // (for the moment the ref. signal is completely invented assuming a PMgain of 5*10^4!)
364   Int_t sectorRef[2];
365   sectorRef[1] = 5;
366   Int_t sigRef[2];
367   // Reference signal are set to 100 (high gain chain) and 800 (low gain chain)
368   if(fIsCalibration==0) {sigRef[0]=100;  sigRef[1]=800;}
369   else {sigRef[0]=0;  sigRef[1]=0;} // calibration -> simulation of pedestal values
370   //
371   for(Int_t iref=0; iref<2; iref++){
372      sectorRef[0] = 3*iref+1;
373      for(Int_t res=0; res<2; res++){
374        sigRef[res] += Pedestal(sectorRef[0], sectorRef[1], res);
375      }
376      new(pdigit) AliZDCDigit(sectorRef, sigRef);
377      treeD->Fill();     
378      
379      //Ch. debug
380      //printf("\t RefDigit added -> det = %d, quad = %d - digi[0,1] = [%d, %d]\n",
381      //    sectorRef[0], sectorRef[1], sigRef[0], sigRef[1]); // Chiara debugging!     
382   }
383   //
384   // --- Adding digits for out-of-time channels after signal digits
385   for(sector[0]=1; sector[0]<6; sector[0]++){
386     for(sector[1]=0; sector[1]<5; sector[1]++){
387         if((sector[0]==3) && ((sector[1]<1) || (sector[1]>2))) continue;
388         for(Int_t res=0; res<2; res++){
389            digioot[res] = Pedestal(sector[0], sector[1], res); // out-of-time ADCs
390         }
391         new(pdigit) AliZDCDigit(sector, digioot);
392         treeD->Fill();
393
394         //Ch. debug
395         //printf("\t DIGIToot added -> det = %d, quad = %d - digi[0,1] = [%d, %d]\n",
396         //     sector[0], sector[1], digioot[0], digioot[1]); // Chiara debugging!      
397     }
398   }
399   // Adding out-of-time digits for 2 reference PTM signals (after out-of-time ch.)
400   Int_t sigRefoot[2];
401   for(Int_t iref=0; iref<2; iref++){
402      sectorRef[0] = 3*iref+1;
403      for(Int_t res=0; res<2; res++){
404        sigRefoot[res] = Pedestal(sectorRef[0], sectorRef[1], res);
405      }
406      new(pdigit) AliZDCDigit(sectorRef, sigRefoot);
407      treeD->Fill();
408      //Ch. debug
409      //printf("\t RefDigitoot added -> det = %d, quad = %d - digi[0,1] = [%d, %d]\n",
410      //    sectorRef[0], sectorRef[1], sigRefoot[0], sigRefoot[1]); // Chiara debugging!
411      
412   }
413   //printf("\t AliZDCDigitizer -> TreeD has %d entries\n",(Int_t) treeD->GetEntries());
414
415   // write the output tree
416   loader->WriteDigits("OVERWRITE");
417   loader->UnloadDigits();
418 }
419
420
421 //_____________________________________________________________________________
422 void AliZDCDigitizer::Fragmentation(Float_t impPar, Int_t specN, Int_t specP,
423                                     Int_t &freeSpecN, Int_t &freeSpecP) const
424 {
425 // simulate fragmentation of spectators
426
427   AliZDCFragment frag(impPar);
428
429   // Fragments generation
430   frag.GenerateIMF();
431   Int_t nAlpha = frag.GetNalpha();
432
433   // Attach neutrons
434   Int_t ztot = frag.GetZtot();
435   Int_t ntot = frag.GetNtot();
436   frag.AttachNeutrons();
437   freeSpecN = specN-ntot-2*nAlpha;
438   freeSpecP = specP-ztot-2*nAlpha;
439   // Removing deuterons
440   Int_t ndeu = (Int_t) (freeSpecN*frag.DeuteronNumber());
441   freeSpecN -= ndeu;
442   //
443   if(freeSpecN<0) freeSpecN=0;
444   if(freeSpecP<0) freeSpecP=0;
445   AliDebug(2, Form("FreeSpn = %d, FreeSpp = %d", freeSpecN, freeSpecP));
446 }
447
448 //_____________________________________________________________________________
449 void AliZDCDigitizer::SpectatorSignal(Int_t SpecType, Int_t numEvents, 
450                                       Float_t pm[5][5]) const
451 {
452 // add signal of the spectators
453   
454   TFile *specSignalFile = TFile::Open("$ALICE_ROOT/ZDC/SpectatorSignal.root");
455   if(!specSignalFile || !specSignalFile->IsOpen()) {
456     AliError((" Opening file $ALICE_ROOT/ZDC/SpectatorSignal.root failed\n"));
457     return;
458   }
459
460   TNtuple* zdcSignal;
461   
462   Float_t sqrtS = 2*fBeamEnergy;
463   //
464   if(TMath::Abs(sqrtS-5500) < 100.){
465     specSignalFile->cd("energy5500");
466     //
467     if(SpecType == 1) {    // --- Signal for projectile spectator neutrons
468       specSignalFile->GetObject("energy5500/ZNCSignal;1",zdcSignal);
469       if(!zdcSignal) AliError("  PROBLEM!!! Can't retrieve ZNCSignal from SpectatorSignal.root file");
470     } 
471     else if(SpecType == 2) { // --- Signal for projectile spectator protons
472       specSignalFile->GetObject("energy5500/ZPCSignal;1",zdcSignal);
473       if(!zdcSignal) AliError("  PROBLEM!!! Can't retrieve ZPCSignal from SpectatorSignal.root file");
474     }
475     else if(SpecType == 3) { // --- Signal for target spectator neutrons
476       specSignalFile->GetObject("energy5500/ZNASignal;1",zdcSignal);
477       if(!zdcSignal) AliError("  PROBLEM!!! Can't retrieve ZNASignal from SpectatorSignal.root file");
478     }
479     else if(SpecType == 4) { // --- Signal for target spectator protons
480       specSignalFile->GetObject("energy5500/ZPASignal;1",zdcSignal);
481       if(!zdcSignal) AliError("  PROBLEM!!! Can't retrieve ZPASignal from SpectatorSignal.root file");
482     }
483   }
484   else if(TMath::Abs(sqrtS-2760) < 100.){
485     specSignalFile->cd("energy2760");
486     //
487     if(SpecType == 1) {    // --- Signal for projectile spectator neutrons
488       specSignalFile->GetObject("energy2760/ZNCSignal;1",zdcSignal);
489       if(!zdcSignal) AliError("  PROBLEM!!! Can't retrieve ZNCSignal from SpectatorSignal.root file");
490     } 
491     else if(SpecType == 2) { // --- Signal for projectile spectator protons
492       specSignalFile->GetObject("energy2760/ZPCSignal;1",zdcSignal);
493       if(!zdcSignal) AliError("  PROBLEM!!! Can't retrieve ZPCSignal from SpectatorSignal.root file");
494     }
495     else if(SpecType == 3) { // --- Signal for target spectator neutrons
496       specSignalFile->GetObject("energy2760/ZNASignal;1",zdcSignal);
497       if(!zdcSignal) AliError("  PROBLEM!!! Can't retrieve ZNASignal from SpectatorSignal.root file");
498     }
499     else if(SpecType == 4) { // --- Signal for target spectator protons
500       specSignalFile->GetObject("energy2760/ZPASignal;1",zdcSignal);
501       if(!zdcSignal) AliError("  PROBLEM!!! Can't retrieve ZPASignal from SpectatorSignal.root file");
502     }
503   }
504   
505   Int_t nentries = (Int_t) zdcSignal->GetEntries();
506   
507   Float_t *entry;
508   Int_t pl, i, k, iev=0, rnd[125], volume[2];
509   for(pl=0;pl<125;pl++) rnd[pl] = 0;
510   if(numEvents > 125) {
511     AliDebug(2,Form("numEvents (%d) is larger than 125", numEvents));
512     numEvents = 125;
513   }
514   for(pl=0;pl<numEvents;pl++){
515      rnd[pl] = (Int_t) (9999*gRandom->Rndm());
516      if(rnd[pl] >= 9999) rnd[pl] = 9998;
517      //printf(" rnd[%d] = %d\n",pl,rnd[pl]);     
518   }
519   // Sorting vector in ascending order with C function QSORT 
520   qsort((void*)rnd,numEvents,sizeof(Int_t),comp);
521   do{
522      for(i=0; i<nentries; i++){  
523         zdcSignal->GetEvent(i);
524         entry = zdcSignal->GetArgs();
525         if(entry[0] == rnd[iev]){
526           for(k=0; k<2; k++) volume[k] = (Int_t) entry[k+1];
527           //
528           Float_t lightQ = entry[7];
529           Float_t lightC = entry[8];
530           //
531           if(volume[0] != 3) {  // ZN or ZP
532             pm[volume[0]-1][0] += lightC;
533             pm[volume[0]-1][volume[1]] += lightQ;
534             //printf("\n   pm[%d][0] = %.0f, pm[%d][%d] = %.0f\n",(volume[0]-1),pm[volume[0]-1][0],
535             //  (volume[0]-1),volume[1],pm[volume[0]-1][volume[1]]);
536           } 
537           else { 
538             if(volume[1] == 1) pm[2][1] += lightC; // ZEM 1
539             else               pm[2][2] += lightQ; // ZEM 2
540             //printf("\n   pm[2][1] = %.0f, pm[2][2] = %.0f\n",pm[2][1],pm[2][2]);
541           }
542         }
543         else if(entry[0] > rnd[iev]){
544           iev++;
545           continue;
546         }
547      }
548   }while(iev<numEvents);
549   
550   specSignalFile->Close();
551   delete specSignalFile;
552 }
553
554
555 //_____________________________________________________________________________
556 Int_t AliZDCDigitizer::Phe2ADCch(Int_t Det, Int_t Quad, Float_t Light, 
557                                  Int_t Res) const
558 {
559   // Evaluation of the ADC channel corresponding to the light yield Light
560   Int_t vADCch = (Int_t) (Light * fPMGain[Det-1][Quad] * fADCRes[Res]);
561   // Ch. debug
562   //printf("\t Phe2ADCch -> det %d quad %d - PMGain[%d][%d] %1.0f phe %1.2f  ADC %d\n", 
563   //    Det,Quad,Det-1,Quad,fPMGain[Det-1][Quad],Light,vADCch);
564
565   return vADCch;
566 }
567
568 //_____________________________________________________________________________
569 Int_t AliZDCDigitizer::Pedestal(Int_t Det, Int_t Quad, Int_t Res) const
570 {
571   // Returns a pedestal for detector det, PM quad, channel with res.
572   //
573   Float_t pedValue;
574   // Normal run
575   if(fIsCalibration == 0){
576     Int_t index=0, kNch=24;
577     if(Quad!=5){
578       if(Det==1)        index = Quad+kNch*Res;    // ZNC
579       else if(Det==2)   index = (Quad+5)+kNch*Res;  // ZPC
580       else if(Det==3)   index = (Quad+9)+kNch*Res;  // ZEM
581       else if(Det==4)   index = (Quad+12)+kNch*Res; // ZNA
582       else if(Det==5)   index = (Quad+17)+kNch*Res; // ZPA
583     }
584     else index = (Det-1)/3+22+kNch*Res; // Reference PMs
585     //
586     Float_t meanPed = fPedData->GetMeanPed(index);
587     Float_t pedWidth = fPedData->GetMeanPedWidth(index);
588     pedValue = gRandom->Gaus(meanPed,pedWidth);
589     //
590     /*printf("\t  AliZDCDigitizer::Pedestal -> det %d quad %d res %d - Ped[%d] = %d\n",
591         Det, Quad, Res, index,(Int_t) pedValue); // Chiara debugging!
592     */
593   }
594   // To create calibration object
595   else{
596     if(Res == 0) pedValue = gRandom->Gaus((35.+10.*gRandom->Rndm()),(0.5+0.2*gRandom->Rndm())); //High gain
597     else  pedValue = gRandom->Gaus((250.+100.*gRandom->Rndm()),(3.5+2.*gRandom->Rndm())); //Low gain
598   }
599
600   return (Int_t) pedValue;
601 }
602
603 //_____________________________________________________________________________
604 AliCDBStorage* AliZDCDigitizer::SetStorage(const char *uri) 
605 {
606
607   Bool_t deleteManager = kFALSE;
608   
609   AliCDBManager *manager = AliCDBManager::Instance();
610   AliCDBStorage *defstorage = manager->GetDefaultStorage();
611   
612   if(!defstorage || !(defstorage->Contains("ZDC"))){ 
613      AliWarning("No default storage set or default storage doesn't contain ZDC!");
614      manager->SetDefaultStorage(uri);
615      deleteManager = kTRUE;
616   }
617  
618   AliCDBStorage *storage = manager->GetDefaultStorage();
619
620   if(deleteManager){
621     AliCDBManager::Instance()->UnsetDefaultStorage();
622     defstorage = 0;   // the storage is killed by AliCDBManager::Instance()->Destroy()
623   }
624
625   return storage; 
626 }
627
628 //_____________________________________________________________________________
629 AliZDCPedestals* AliZDCDigitizer::GetPedData() const
630 {
631
632   // Getting pedestal calibration object for ZDC set
633
634   AliCDBEntry  *entry = AliCDBManager::Instance()->Get("ZDC/Calib/Pedestals");
635   if(!entry) AliFatal("No calibration data loaded!");  
636
637   AliZDCPedestals *calibdata = dynamic_cast<AliZDCPedestals*>  (entry->GetObject());
638   if(!calibdata)  AliFatal("Wrong calibration object in calibration  file!");
639
640   return calibdata;
641 }
642