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