]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ZDC/ZDCsim/AliZDCDigitizer.cxx
Removing annoying prints
[u/mrichter/AliRoot.git] / ZDC / ZDCsim / 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("A-A")) == 0)){
480     for(int i=0; i<12; i++){
481       if(beam[i]==1){
482         Float_t scalGainFactor = fBeamEnergy/2760.;
483         if(det[i]!=31 && det[i]!=32){
484           for(Int_t j=0; j<5; j++) fPMGain[det[i]-1][j] = gain[i]/(aEne[i]*scalGainFactor);
485         }
486         else{
487           for(int iq=1; iq<3; iq++) fPMGain[2][iq] = gain[i]/(aEne[i]*scalGainFactor);
488         }
489       }
490      }  
491      //
492      printf("\n    AliZDCDigitizer::ReadPMTGains -> ZDC PMT gains for %s @ %1.0f+%1.0f A GeV: ZN(%1.0f), ZP(%1.0f), ZEM(%1.0f)\n",
493         fBeamType.Data(),fBeamEnergy, fBeamEnergy, fPMGain[0][0], fPMGain[1][0], fPMGain[2][1]);
494   }
495   else if(((fBeamType.CompareTo("p-A")) == 0) || ((fBeamType.CompareTo("P-A")) == 0)){
496     for(int i=0; i<12; i++){
497       if(beam[i]==0 && fBeamEnergy!=0.){
498         if(det[i]==1 || det[i]==2){ // Setting ZNC/ZPC gains as for p-p
499           for(Int_t j=0; j<5; j++) fPMGain[det[i]-1][j] = gain[i]*(aEne[i]/fBeamEnergy+bEne[i]);
500         }
501       }
502       if(beam[i]==1){
503         Float_t scalGainFactor = fBeamEnergy/2760.;
504         Float_t npartScalingFactor = 208./15.;
505         if(det[i]==4 || det[i]==5){ // Setting ZNA/ZPA gains as for Pb-Pb
506           for(Int_t j=0; j<5; j++) fPMGain[det[i]-1][j] = npartScalingFactor*gain[i]/(aEne[i]*scalGainFactor);
507         }
508         else if(det[i]==31 || det[i]==32){ // Setting ZEM gains as for Pb-Pb
509           for(int iq=1; iq<3; iq++) fPMGain[2][iq] = npartScalingFactor*gain[i]/(aEne[i]*scalGainFactor);
510         }
511       }
512     }
513     printf("\n    AliZDCDigitizer::ReadPMTGains -> ZDC PMT gains for %s: ZNC(%1.0f), ZPC(%1.0f), ZEM(%1.0f), ZNA(%1.0f) ZPA(%1.0f)\n",
514         fBeamType.Data(),fPMGain[0][0], fPMGain[1][0], fPMGain[2][1], fPMGain[3][0], fPMGain[4][0]);
515   }
516   else if(((fBeamType.CompareTo("A-p")) == 0) || ((fBeamType.CompareTo("A-P")) == 0)){
517     for(int i=0; i<12; i++){
518       if(beam[i]==0 && fBeamEnergy!=0.){ // Setting ZNA/ZPA gains as for p-p
519         if(det[i]==4 || det[i]==5){
520           for(Int_t j=0; j<5; j++) fPMGain[det[i]-1][j] = gain[i]*(aEne[i]/fBeamEnergy+bEne[i]);
521         }
522         // Setting ZEM gains as for p-p
523         else if(det[i] == 31) fPMGain[2][1] = gain[i]*(aEne[i]-fBeamEnergy*bEne[i]);
524         else if(det[i] == 32) fPMGain[2][2] = gain[i]*(aEne[i]-fBeamEnergy*bEne[i]);
525       }
526       if(beam[i]==1){
527         Float_t scalGainFactor = fBeamEnergy/2760.;
528         Float_t npartScalingFactor = 208./15.;
529         if(det[i]==1 || det[i]==2){ // Setting ZNC/ZPC gains as for Pb-Pb
530           for(Int_t j=0; j<5; j++) fPMGain[det[i]-1][j] = npartScalingFactor*gain[i]/(aEne[i]*scalGainFactor);
531         }
532       }
533     }
534     printf("\n    AliZDCDigitizer::ReadPMTGains -> ZDC PMT gains for %s: ZNC(%1.0f), ZPC(%1.0f), ZEM(%1.0f), ZNA(%1.0f) ZPA(%1.0f)\n",
535         fBeamType.Data(),fPMGain[0][0], fPMGain[1][0], fPMGain[2][1], fPMGain[3][0], fPMGain[4][0]);
536   }
537   else{ // Setting pp gains as default value
538     for(int i=0; i<12; i++){
539       if(beam[i]==0 && fBeamEnergy!=0.){
540         if(det[i]!=31 && det[i]!=32){
541           for(Int_t j=0; j<5; j++) fPMGain[det[i]-1][j] = gain[i]*(aEne[i]/fBeamEnergy+bEne[i]);
542         }
543         else if(det[i] == 31) fPMGain[2][1] = gain[i]*(aEne[i]-fBeamEnergy*bEne[i]);
544         else if(det[i] == 32) fPMGain[2][2] = gain[i]*(aEne[i]-fBeamEnergy*bEne[i]);
545       }
546     }
547     //
548     printf("\n    AliZDCDigitizer::ReadPMTGains -> ZDC PMT gains for %s @ %1.0f+%1.0f GeV: ZNC(%1.0f), ZPC(%1.0f), ZEM(%1.0f), ZNA(%1.0f) ZPA(%1.0f)\n",
549         fBeamType.Data(),fBeamEnergy, fBeamEnergy, fPMGain[0][0], fPMGain[1][0], fPMGain[2][1], fPMGain[3][0], fPMGain[4][0]);     
550   }
551
552 }
553
554 //_____________________________________________________________________________
555 void AliZDCDigitizer::CalculatePMTGains()
556 {
557 // Calculate PMT gain according to beam type and beam energy
558   if( ((fBeamType.CompareTo("P-P")) == 0) ||  ((fBeamType.CompareTo("p-p"))) ){
559     // PTM gains rescaled to beam energy for p-p
560     // New correction coefficients for PMT gains needed
561     // to reproduce experimental spectra (from Grazia Jul 2010)
562     if(fBeamEnergy != 0){
563       for(Int_t j = 0; j < 5; j++){
564           fPMGain[0][j] = 1.515831*(661.444/fBeamEnergy+0.000740671)*10000000;
565           fPMGain[1][j] = 0.674234*(864.350/fBeamEnergy+0.00234375)*10000000;
566           fPMGain[3][j] = 1.350938*(661.444/fBeamEnergy+0.000740671)*10000000; 
567           fPMGain[4][j] = 0.678597*(864.350/fBeamEnergy+0.00234375)*10000000;
568       }
569       fPMGain[2][1] = 0.869654*(1.32312-0.000101515*fBeamEnergy)*10000000;
570       fPMGain[2][2] = 1.030883*(1.32312-0.000101515*fBeamEnergy)*10000000;
571       //
572       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",
573         fBeamEnergy, fBeamEnergy, fPMGain[0][0], fPMGain[1][0], fPMGain[2][1], fPMGain[3][0], fPMGain[4][0]);
574      
575     }
576   }
577   else if(((fBeamType.CompareTo("A-A")) == 0)){
578     // PTM gains for Pb-Pb @ 2.7+2.7 A TeV ***************
579     // rescaled for Pb-Pb @ 1.38+1.38 A TeV ***************
580     // Values corrected after 2010 Pb-Pb data taking (7/2/2011 - Ch.)
581     // Experimental data compared to EMD simulation for single nucleon peaks:
582     // ZN gains must be divided by 4, ZP gains by 10!
583     Float_t scalGainFactor = fBeamEnergy/2760.;
584     for(Int_t j = 0; j < 5; j++){
585        fPMGain[0][j] = 50000./(4*scalGainFactor);  // ZNC                
586        fPMGain[1][j] = 100000./(5*scalGainFactor); // ZPC       
587        fPMGain[2][j] = 100000./scalGainFactor;     // ZEM
588        fPMGain[3][j] = 50000./(4*scalGainFactor);  // ZNA                
589        fPMGain[4][j] = 100000./(5*scalGainFactor); // ZPA    
590     }
591     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",
592         fBeamEnergy, fBeamEnergy, fPMGain[0][0], fPMGain[1][0], fPMGain[2][1]);
593   }
594   else if(((fBeamType.CompareTo("p-A")) == 0) || ((fBeamType.CompareTo("P-A"))) ){
595     // PTM gains for Pb-Pb @ 1.38+1.38 A TeV on side A
596     // PTM gains rescaled to beam energy for p-p on side C
597     // WARNING! Energies are set by hand for 2011 pA RUN!!!
598     Float_t scalGainFactor = fBeamEnergy/2760.;
599     Float_t npartScalingFactor = 208./15.;
600     
601     for(Int_t j = 0; j < 5; j++){
602        fPMGain[0][j] = 1.515831*(661.444/fBeamEnergy+0.000740671)*10000000; //ZNC (p)
603        fPMGain[1][j] = 0.674234*(864.350/fBeamEnergy+0.00234375)*10000000;  //ZPC (p)
604        if(j<2) fPMGain[2][j] = npartScalingFactor*100000./scalGainFactor;          // ZEM (Pb)
605        // Npart max scales from 400 in Pb-Pb to ~8 in pPb -> *40.
606        fPMGain[3][j] = npartScalingFactor*50000/(4*scalGainFactor);  // ZNA (Pb)             
607        fPMGain[4][j] = npartScalingFactor*100000/(5*scalGainFactor); // ZPA (Pb)  
608     }
609     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",
610         fPMGain[0][0], fPMGain[1][0], fPMGain[2][1], fPMGain[3][0], fPMGain[4][0]);
611   }
612   else if(((fBeamType.CompareTo("A-p")) == 0) || ((fBeamType.CompareTo("A-P"))) ){
613     // PTM gains for Pb-Pb @ 1.38+1.38 A TeV on side A
614     // PTM gains rescaled to beam energy for p-p on side C
615      Float_t scalGainFactor = fBeamEnergy/2760.;
616     Float_t npartScalingFactor = 208./15.;
617     
618     for(Int_t j = 0; j < 5; j++){
619        fPMGain[3][j] = 1.350938*(661.444/fBeamEnergy+0.000740671)*10000000;  //ZNA (p)
620        fPMGain[4][j] = 0.678597*(864.350/fBeamEnergy+0.00234375)*10000000;   //ZPA (p)
621        // Npart max scales from 400 in Pb-Pb to ~8 in pPb -> *40.
622        fPMGain[1][j] = npartScalingFactor*50000/(4*scalGainFactor);  // ZNC (Pb)             
623        fPMGain[2][j] = npartScalingFactor*100000/(5*scalGainFactor); // ZPC (Pb)  
624     }
625     fPMGain[2][1] = 0.869654*(1.32312-0.000101515*fBeamEnergy)*10000000; // ZEM (pp)
626     fPMGain[2][2] = 1.030883*(1.32312-0.000101515*fBeamEnergy)*10000000; // ZEM (pp)
627     printf("\n    AliZDCDigitizer::CalculatePMTGains -> ZDC PMT gains for Pb-p: ZNC(%1.0f), ZPC(%1.0f), ZEM(%1.0f), ZNA(%1.0f) ZPA(%1.0f)\n",
628         fPMGain[0][0], fPMGain[1][0], fPMGain[2][1], fPMGain[3][0], fPMGain[4][0]);
629   }
630 }
631
632 //_____________________________________________________________________________
633 void AliZDCDigitizer::Fragmentation(Float_t impPar, Int_t specN, Int_t specP,
634                                     Int_t &freeSpecN, Int_t &freeSpecP) const
635 {
636 // simulate fragmentation of spectators
637
638   AliZDCFragment frag(impPar);
639
640   // Fragments generation
641   frag.GenerateIMF();
642   Int_t nAlpha = frag.GetNalpha();
643
644   // Attach neutrons
645   frag.AttachNeutrons();
646   Int_t ztot = frag.GetZtot();
647   Int_t ntot = frag.GetNtot();
648   
649   // Removing fragments and alpha pcs
650   freeSpecN = specN-ntot-2*nAlpha;
651   freeSpecP = specP-ztot-2*nAlpha;
652   
653   // Removing deuterons
654   Int_t ndeu = (Int_t) (freeSpecN*frag.DeuteronNumber());
655   freeSpecN -= ndeu;
656   freeSpecP -= ndeu;
657   //
658   if(freeSpecN<0) freeSpecN=0;
659   if(freeSpecP<0) freeSpecP=0;
660   AliDebug(2, Form("FreeSpn = %d, FreeSpp = %d", freeSpecN, freeSpecP));
661 }
662
663 //_____________________________________________________________________________
664 void AliZDCDigitizer::SpectatorSignal(Int_t SpecType, Int_t numEvents, Float_t pm[5][5]) 
665 {
666 // add signal of the spectators
667   
668   if(!fSpectatorData) fSpectatorData = TFile::Open("$ALICE_ROOT/ZDC/SpectatorSignal.root");
669   if(!fSpectatorData || !fSpectatorData->IsOpen()) {
670     AliError((" Opening file $ALICE_ROOT/ZDC/SpectatorSignal.root failed\n"));
671     return;
672   }
673
674   TNtuple* zdcSignal=0x0;
675   
676   Float_t sqrtS = 2*fBeamEnergy;
677   //
678   if(TMath::Abs(sqrtS-5500) < 100.){
679     AliInfo(" Extracting signal from SpectatorSignal/energy5500 ");
680     fSpectatorData->cd("energy5500");
681     //
682     if(SpecType == 1) {    // --- Signal for projectile spectator neutrons
683       fSpectatorData->GetObject("energy5500/ZNCSignal;1",zdcSignal);
684       if(!zdcSignal) AliError("  PROBLEM!!! Can't retrieve ZNCSignal from SpectatorSignal.root file");
685     } 
686     else if(SpecType == 2) { // --- Signal for projectile spectator protons
687       fSpectatorData->GetObject("energy5500/ZPCSignal;1",zdcSignal);
688       if(!zdcSignal) AliError("  PROBLEM!!! Can't retrieve ZPCSignal from SpectatorSignal.root file");
689     }
690     else if(SpecType == 3) { // --- Signal for target spectator neutrons
691       fSpectatorData->GetObject("energy5500/ZNASignal;1",zdcSignal);
692       if(!zdcSignal) AliError("  PROBLEM!!! Can't retrieve ZNASignal from SpectatorSignal.root file");
693     }
694     else if(SpecType == 4) { // --- Signal for target spectator protons
695       fSpectatorData->GetObject("energy5500/ZPASignal;1",zdcSignal);
696       if(!zdcSignal) AliError("  PROBLEM!!! Can't retrieve ZPASignal from SpectatorSignal.root file");
697     }
698   }
699   else if(TMath::Abs(sqrtS-2760) < 100.){
700     AliInfo(" Extracting signal from SpectatorSignal/energy2760 ");
701     fSpectatorData->cd("energy2760");
702     //
703     if(SpecType == 1) {    // --- Signal for projectile spectator neutrons
704       fSpectatorData->GetObject("energy2760/ZNCSignal;1",zdcSignal);
705       if(!zdcSignal) AliError("  PROBLEM!!! Can't retrieve ZNCSignal from SpectatorSignal.root file");
706     } 
707     else if(SpecType == 2) { // --- Signal for projectile spectator protons
708       fSpectatorData->GetObject("energy2760/ZPCSignal;1",zdcSignal);
709       if(!zdcSignal) AliError("  PROBLEM!!! Can't retrieve ZPCSignal from SpectatorSignal.root file");
710     }
711     else if(SpecType == 3) { // --- Signal for target spectator neutrons
712       fSpectatorData->GetObject("energy2760/ZNASignal;1",zdcSignal);
713       if(!zdcSignal) AliError("  PROBLEM!!! Can't retrieve ZNASignal from SpectatorSignal.root file");
714     }
715     else if(SpecType == 4) { // --- Signal for target spectator protons
716       fSpectatorData->GetObject("energy2760/ZPASignal;1",zdcSignal);
717       if(!zdcSignal) AliError("  PROBLEM!!! Can't retrieve ZPASignal from SpectatorSignal.root file");
718     }
719   }
720   
721   if(!zdcSignal){
722     printf("\n No spectator signal available for ZDC digitization\n");
723     return;
724   } 
725   
726   Int_t nentries = (Int_t) zdcSignal->GetEntries();
727   
728   Float_t *entry;
729   Int_t pl, i, k, iev=0, rnd[125], volume[2];
730   for(pl=0;pl<125;pl++) rnd[pl] = 0;
731   if(numEvents > 125) {
732     AliDebug(2,Form("numEvents (%d) is larger than 125", numEvents));
733     numEvents = 125;
734   }
735   for(pl=0;pl<numEvents;pl++){
736      rnd[pl] = (Int_t) (9999*gRandom->Rndm());
737      if(rnd[pl] >= 9999) rnd[pl] = 9998;
738      //printf(" rnd[%d] = %d\n",pl,rnd[pl]);     
739   }
740   // Sorting vector in ascending order with C function QSORT 
741   qsort((void*)rnd,numEvents,sizeof(Int_t),comp);
742   do{
743      for(i=0; i<nentries; i++){  
744         zdcSignal->GetEvent(i);
745         entry = zdcSignal->GetArgs();
746         if(entry[0] == rnd[iev]){
747           for(k=0; k<2; k++) volume[k] = (Int_t) entry[k+1];
748           //
749           Float_t lightQ = entry[7];
750           Float_t lightC = entry[8];
751           //
752           if(volume[0] != 3) {  // ZN or ZP
753             pm[volume[0]-1][0] += lightC;
754             pm[volume[0]-1][volume[1]] += lightQ;
755             //printf("\n   pm[%d][0] = %.0f, pm[%d][%d] = %.0f\n",(volume[0]-1),pm[volume[0]-1][0],
756             //  (volume[0]-1),volume[1],pm[volume[0]-1][volume[1]]);
757           } 
758           else { 
759             if(volume[1] == 1) pm[2][1] += lightC; // ZEM 1
760             else               pm[2][2] += lightQ; // ZEM 2
761             //printf("\n   pm[2][1] = %.0f, pm[2][2] = %.0f\n",pm[2][1],pm[2][2]);
762           }
763         }
764         else if(entry[0] > rnd[iev]){
765           iev++;
766           continue;
767         }
768      }
769   }while(iev<numEvents);
770   
771 }
772
773
774 //_____________________________________________________________________________
775 Int_t AliZDCDigitizer::Phe2ADCch(Int_t Det, Int_t Quad, Float_t Light, 
776                                  Int_t Res) const
777 {
778   // Evaluation of the ADC channel corresponding to the light yield Light
779   Int_t vADCch = (Int_t) (Light * fPMGain[Det-1][Quad] * fADCRes[Res]);
780   // Ch. debug
781   //printf("\t Phe2ADCch -> det %d quad %d - PMGain[%d][%d] %1.0f phe %1.2f  ADC %d\n", 
782   //    Det,Quad,Det-1,Quad,fPMGain[Det-1][Quad],Light,vADCch);
783
784   return vADCch;
785 }
786
787 //_____________________________________________________________________________
788 Int_t AliZDCDigitizer::Pedestal(Int_t Det, Int_t Quad, Int_t Res) const
789 {
790   // Returns a pedestal for detector det, PM quad, channel with res.
791   //
792   Float_t pedValue=0.;
793   // Normal run
794   if(fIsCalibration == 0){
795     Int_t index=0, kNch=24;
796     if(Quad!=5){
797       if(Det==1)        index = Quad+kNch*Res;    // ZNC
798       else if(Det==2)   index = (Quad+5)+kNch*Res;  // ZPC
799       else if(Det==3)   index = (Quad+9)+kNch*Res;  // ZEM
800       else if(Det==4)   index = (Quad+12)+kNch*Res; // ZNA
801       else if(Det==5)   index = (Quad+17)+kNch*Res; // ZPA
802     }
803     else index = (Det-1)/3+22+kNch*Res; // Reference PMs
804     //
805     if(fPedData){
806       Float_t meanPed = fPedData->GetMeanPed(index);
807       Float_t pedWidth = fPedData->GetMeanPedWidth(index);
808       pedValue = gRandom->Gaus(meanPed,pedWidth);
809       //
810       /*printf("\t  AliZDCDigitizer::Pedestal -> det %d quad %d res %d - Ped[%d] = %d\n",
811         Det, Quad, Res, index,(Int_t) pedValue); // Chiara debugging!
812       */
813     }
814     else{
815       printf("  AliZDCDigitizer::Pedestal -> No valid pedestal calibration object loaded!\n\n");
816     }
817   }
818   // To create calibration object
819   else{
820     if(Res == 0) pedValue = gRandom->Gaus((35.+10.*gRandom->Rndm()),(0.5+0.2*gRandom->Rndm())); //High gain
821     else  pedValue = gRandom->Gaus((250.+100.*gRandom->Rndm()),(3.5+2.*gRandom->Rndm())); //Low gain
822   }
823
824   return (Int_t) pedValue;
825 }
826
827 //_____________________________________________________________________________
828 AliCDBStorage* AliZDCDigitizer::SetStorage(const char *uri) 
829 {
830
831   Bool_t deleteManager = kFALSE;
832   
833   AliCDBManager *manager = AliCDBManager::Instance();
834   AliCDBStorage *defstorage = manager->GetDefaultStorage();
835   
836   if(!defstorage || !(defstorage->Contains("ZDC"))){ 
837      AliWarning("No default storage set or default storage doesn't contain ZDC!");
838      manager->SetDefaultStorage(uri);
839      deleteManager = kTRUE;
840   }
841  
842   AliCDBStorage *storage = manager->GetDefaultStorage();
843
844   if(deleteManager){
845     AliCDBManager::Instance()->UnsetDefaultStorage();
846     defstorage = 0;   // the storage is killed by AliCDBManager::Instance()->Destroy()
847   }
848
849   return storage; 
850 }
851
852 //_____________________________________________________________________________
853 AliZDCPedestals* AliZDCDigitizer::GetPedData() const
854 {
855
856   // Getting pedestal calibration object for ZDC set
857   AliZDCPedestals *calibdata = 0x0;
858   AliCDBEntry  *entry = AliCDBManager::Instance()->Get("ZDC/Calib/Pedestals");
859   if(!entry) AliFatal("No calibration data loaded!");  
860   else{
861     calibdata = dynamic_cast<AliZDCPedestals*>  (entry->GetObject());
862     if(!calibdata)  AliFatal("Wrong calibration object in calibration  file!");
863   }
864   return calibdata;
865 }
866