]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ZDC/AliZDCDigitizer.cxx
Updated fragmentation model
[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   for(Int_t i=0; i<2; i++) fADCRes[i]=0.;
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<5; 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<5; 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   if(!entry) AliFatal("No calibration data loaded!");  
122   AliGRPObject* grpData = 0x0;
123   if(entry){
124     TMap* m = dynamic_cast<TMap*>(entry->GetObject());  // old GRP entry
125     if(m){
126       //m->Print();
127       grpData = new AliGRPObject();
128       grpData->ReadValuesFromMap(m);
129     }
130     else{
131       grpData = dynamic_cast<AliGRPObject*>(entry->GetObject());  // new GRP entry
132     }
133     entry->SetOwner(0);
134     AliCDBManager::Instance()->UnloadFromCache("GRP/GRP/Data");
135   }
136   if(!grpData){
137     AliError("No GRP entry found in OCDB! \n ");
138     return kFALSE;
139   }
140   
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");
144   }
145   
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");
150     fBeamEnergy = 0.;
151   }
152
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;
163       }
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]));
168     }
169   }
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    
183     }
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]));
186   }
187     
188   // ADC Caen V965
189   fADCRes[0] = 0.0000008; // ADC Resolution high gain: 200 fC/adcCh
190   fADCRes[1] = 0.0000064; // ADC Resolution low gain:  25  fC/adcCh
191
192   return kTRUE;
193 }
194
195 //____________________________________________________________________________
196 void AliZDCDigitizer::Exec(Option_t* /*option*/)
197 {
198   // Execute digitization
199
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   // ------------------------------------------------------------
210   Float_t pm[5][5]; 
211   for(Int_t iSector1=0; iSector1<5; iSector1++) 
212     for(Int_t iSector2=0; iSector2<5; iSector2++){
213       pm[iSector1][iSector2] = 0;
214     }
215     
216   // ------------------------------------------------------------
217   // ### Out of time ADC added (22 channels)
218   // --- same codification as for signal PTMs (see above)
219   // ------------------------------------------------------------
220   Float_t pmoot[5][5];
221   for(Int_t iSector1=0; iSector1<5; iSector1++) 
222     for(Int_t iSector2=0; iSector2<5; iSector2++){
223       pmoot[iSector1][iSector2] = 0;
224     }
225
226   // impact parameter and number of spectators
227   Float_t impPar = 0;
228   Int_t specNTarg = 0, specPTarg = 0;
229   Int_t specNProj = 0, specPProj = 0;
230   Float_t signalTime0 = 0.;
231
232   // loop over input streams
233   for(Int_t iInput = 0; iInput<fManager->GetNinputs(); iInput++){
234
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;
240
241     // load sdigits
242     loader->LoadSDigits();
243     TTree* treeS = loader->TreeS();
244     if(!treeS) continue;
245     AliZDCSDigit sdigit;
246     AliZDCSDigit* psdigit = &sdigit;
247     treeS->SetBranchAddress("ZDC", &psdigit);
248
249     // loop over sdigits
250     for(Int_t iSDigit=0; iSDigit<treeS->GetEntries(); iSDigit++){
251       treeS->GetEntry(iSDigit);
252       //
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)));
257         continue;
258       }
259       // Checking if signal is inside ADC gate
260       if(iSDigit==0) signalTime0 = sdigit.GetTrackTime();
261       else{
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.;
270         }
271       }
272       
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;
278       }
279       
280       pm[(sdigit.GetSector(0))-1][sdigit.GetSector(1)] += sdigit.GetLightPM();
281       //Ch. debug
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)]); 
285       */
286       
287     }
288
289     loader->UnloadSDigits();
290
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;
298     
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);
309     }
310     
311   }
312
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);
327   }
328
329
330   // get the output run loader and loader
331   AliRunLoader* runLoader = 
332     AliRunLoader::GetRunLoader(fManager->GetOutputFolderName());
333   AliLoader* loader = runLoader->GetLoader("ZDCLoader");
334   if(!loader) {
335     AliError("no ZDC loader found");
336     return;
337   }
338
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();
345   AliZDCDigit digit;
346   AliZDCDigit* pdigit = &digit;
347   const Int_t kBufferSize = 4000;
348   treeD->Branch("ZDC", "AliZDCDigit", &pdigit, kBufferSize);
349
350   // Create digits
351   Int_t sector[2];
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);
359         }
360         new(pdigit) AliZDCDigit(sector, digi);
361         treeD->Fill();
362
363         //Ch. debug
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!
366         
367     }
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!)
371   Int_t sectorRef[2];
372   sectorRef[1] = 5;
373   Int_t sigRef[2];
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
377   //
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);
382      }
383      new(pdigit) AliZDCDigit(sectorRef, sigRef);
384      treeD->Fill();     
385      
386      //Ch. debug
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!     
389   }
390   //
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
397         }
398         new(pdigit) AliZDCDigit(sector, digioot);
399         treeD->Fill();
400
401         //Ch. debug
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!      
404     }
405   }
406   // Adding out-of-time digits for 2 reference PTM signals (after out-of-time ch.)
407   Int_t sigRefoot[2];
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);
412      }
413      new(pdigit) AliZDCDigit(sectorRef, sigRefoot);
414      treeD->Fill();
415      //Ch. debug
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!
418      
419   }
420   //printf("\t AliZDCDigitizer -> TreeD has %d entries\n",(Int_t) treeD->GetEntries());
421
422   // write the output tree
423   loader->WriteDigits("OVERWRITE");
424   loader->UnloadDigits();
425 }
426
427
428 //_____________________________________________________________________________
429 void AliZDCDigitizer::Fragmentation(Float_t impPar, Int_t specN, Int_t specP,
430                                     Int_t &freeSpecN, Int_t &freeSpecP) const
431 {
432 // simulate fragmentation of spectators
433
434   AliZDCFragment frag(impPar);
435
436   // Fragments generation
437   frag.GenerateIMF();
438   Int_t nAlpha = frag.GetNalpha();
439
440   // Attach neutrons
441   frag.AttachNeutrons();
442   Int_t ztot = frag.GetZtot();
443   Int_t ntot = frag.GetNtot();
444   
445   // Removing fragments and alpha pcs
446   freeSpecN = specN-ntot-2*nAlpha;
447   freeSpecP = specP-ztot-2*nAlpha;
448   
449   // Removing deuterons
450   Int_t ndeu = (Int_t) (freeSpecN*frag.DeuteronNumber());
451   freeSpecN -= ndeu;
452   freeSpecP -= ndeu;
453   //
454   if(freeSpecN<0) freeSpecN=0;
455   if(freeSpecP<0) freeSpecP=0;
456   AliDebug(2, Form("FreeSpn = %d, FreeSpp = %d", freeSpecN, freeSpecP));
457 }
458
459 //_____________________________________________________________________________
460 void AliZDCDigitizer::SpectatorSignal(Int_t SpecType, Int_t numEvents, 
461                                       Float_t pm[5][5]) const
462 {
463 // add signal of the spectators
464   
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"));
468     return;
469   }
470
471   TNtuple* zdcSignal=0x0;
472   
473   Float_t sqrtS = 2*fBeamEnergy;
474   //
475   if(TMath::Abs(sqrtS-5500) < 100.){
476     specSignalFile->cd("energy5500");
477     //
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");
481     } 
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");
485     }
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");
489     }
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");
493     }
494   }
495   else if(TMath::Abs(sqrtS-2760) < 100.){
496     specSignalFile->cd("energy2760");
497     //
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");
501     } 
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");
505     }
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");
509     }
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");
513     }
514   }
515   
516   if(!zdcSignal){
517     printf("\n No spectator signal available for ZDC digitization\n");
518     return;
519   } 
520   
521   Int_t nentries = (Int_t) zdcSignal->GetEntries();
522   
523   Float_t *entry;
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));
528     numEvents = 125;
529   }
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]);     
534   }
535   // Sorting vector in ascending order with C function QSORT 
536   qsort((void*)rnd,numEvents,sizeof(Int_t),comp);
537   do{
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];
543           //
544           Float_t lightQ = entry[7];
545           Float_t lightC = entry[8];
546           //
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]]);
552           } 
553           else { 
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]);
557           }
558         }
559         else if(entry[0] > rnd[iev]){
560           iev++;
561           continue;
562         }
563      }
564   }while(iev<numEvents);
565   
566   specSignalFile->Close();
567   delete specSignalFile;
568 }
569
570
571 //_____________________________________________________________________________
572 Int_t AliZDCDigitizer::Phe2ADCch(Int_t Det, Int_t Quad, Float_t Light, 
573                                  Int_t Res) const
574 {
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]);
577   // Ch. debug
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);
580
581   return vADCch;
582 }
583
584 //_____________________________________________________________________________
585 Int_t AliZDCDigitizer::Pedestal(Int_t Det, Int_t Quad, Int_t Res) const
586 {
587   // Returns a pedestal for detector det, PM quad, channel with res.
588   //
589   Float_t pedValue;
590   // Normal run
591   if(fIsCalibration == 0){
592     Int_t index=0, kNch=24;
593     if(Quad!=5){
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
599     }
600     else index = (Det-1)/3+22+kNch*Res; // Reference PMs
601     //
602     Float_t meanPed = fPedData->GetMeanPed(index);
603     Float_t pedWidth = fPedData->GetMeanPedWidth(index);
604     pedValue = gRandom->Gaus(meanPed,pedWidth);
605     //
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!
608     */
609   }
610   // To create calibration object
611   else{
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
614   }
615
616   return (Int_t) pedValue;
617 }
618
619 //_____________________________________________________________________________
620 AliCDBStorage* AliZDCDigitizer::SetStorage(const char *uri) 
621 {
622
623   Bool_t deleteManager = kFALSE;
624   
625   AliCDBManager *manager = AliCDBManager::Instance();
626   AliCDBStorage *defstorage = manager->GetDefaultStorage();
627   
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;
632   }
633  
634   AliCDBStorage *storage = manager->GetDefaultStorage();
635
636   if(deleteManager){
637     AliCDBManager::Instance()->UnsetDefaultStorage();
638     defstorage = 0;   // the storage is killed by AliCDBManager::Instance()->Destroy()
639   }
640
641   return storage; 
642 }
643
644 //_____________________________________________________________________________
645 AliZDCPedestals* AliZDCDigitizer::GetPedData() const
646 {
647
648   // Getting pedestal calibration object for ZDC set
649
650   AliCDBEntry  *entry = AliCDBManager::Instance()->Get("ZDC/Calib/Pedestals");
651   if(!entry) AliFatal("No calibration data loaded!");  
652
653   AliZDCPedestals *calibdata = dynamic_cast<AliZDCPedestals*>  (entry->GetObject());
654   if(!calibdata)  AliFatal("Wrong calibration object in calibration  file!");
655
656   return calibdata;
657 }
658