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