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