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