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