]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ZDC/AliZDCDigitizer.cxx
Upsilon 8 TeV
[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         for(Int_t iH = 0; iH < listOfHeaders->GetEntries(); ++iH) {
305           AliGenEventHeader *currHeader = dynamic_cast <AliGenEventHeader *> (listOfHeaders->At(iH));
306           if (currHeader && currHeader->InheritsFrom(AliGenHijingEventHeader::Class())) {
307             hijingHeader = dynamic_cast <AliGenHijingEventHeader*> (currHeader);
308             break;
309           }
310         }
311       }
312       else{
313         printf(" No list of headers from generator \n");
314       }
315     }
316     if(!hijingHeader){ 
317         printf(" No HIJING header found in list of headers from generator\n");
318     }
319     
320     if(hijingHeader && fSpectators2Track==kTRUE){
321       impPar = hijingHeader->ImpactParameter();
322       specNProj = hijingHeader->ProjSpectatorsn();
323       specPProj = hijingHeader->ProjSpectatorsp();
324       specNTarg = hijingHeader->TargSpectatorsn();
325       specPTarg = hijingHeader->TargSpectatorsp();
326       AliInfo(Form("\t AliZDCDigitizer: b = %1.2f fm\n"
327       " \t    PROJECTILE:  #spectator n %d, #spectator p %d\n"
328       " \t    TARGET:  #spectator n %d, #spectator p %d\n", 
329       impPar, specNProj, specPProj, specNTarg, specPTarg));
330     
331       // Applying fragmentation algorithm and adding spectator signal
332       Int_t freeSpecNProj=0, freeSpecPProj=0;
333       if(specNProj!=0 || specPProj!=0) Fragmentation(impPar, specNProj, specPProj, freeSpecNProj, freeSpecPProj);
334       Int_t freeSpecNTarg=0, freeSpecPTarg=0;
335       if(specNTarg!=0 || specPTarg!=0) Fragmentation(impPar, specNTarg, specPTarg, freeSpecNTarg, freeSpecPTarg);
336       if(freeSpecNProj!=0) SpectatorSignal(1, freeSpecNProj, pm);
337       if(freeSpecPProj!=0) SpectatorSignal(2, freeSpecPProj, pm);
338       AliInfo(Form("\t AliZDCDigitizer -> Adding spectator signal for PROJECTILE: %d free  n and %d free p\n",freeSpecNProj,freeSpecPProj));
339       if(freeSpecNTarg!=0) SpectatorSignal(3, freeSpecNTarg, pm);
340       if(freeSpecPTarg!=0) SpectatorSignal(4, freeSpecPTarg, pm);
341       AliInfo(Form("\t AliZDCDigitizer -> Adding spectator signal for TARGET: %d free  n and %d free p\n",freeSpecNTarg,freeSpecPTarg));
342     }
343   }
344
345
346   // get the output run loader and loader
347   AliRunLoader* runLoader = 
348     AliRunLoader::GetRunLoader(fDigInput->GetOutputFolderName());
349   AliLoader* loader = runLoader->GetLoader("ZDCLoader");
350   if(!loader) {
351     AliError("no ZDC loader found");
352     return;
353   }
354
355   // create the output tree
356   const char* mode = "update";
357   if(runLoader->GetEventNumber() == 0) mode = "recreate";
358   loader->LoadDigits(mode);
359   loader->MakeTree("D");
360   TTree* treeD = loader->TreeD();
361   AliZDCDigit digit;
362   AliZDCDigit* pdigit = &digit;
363   const Int_t kBufferSize = 4000;
364   treeD->Branch("ZDC", "AliZDCDigit", &pdigit, kBufferSize);
365
366   // Create digits
367   Int_t sector[2];
368   Int_t digi[2], digioot[2];
369   for(sector[0]=1; sector[0]<6; sector[0]++){
370     for(sector[1]=0; sector[1]<5; sector[1]++){
371         if((sector[0]==3) && ((sector[1]<1) || (sector[1]>2))) continue;
372         for(Int_t res=0; res<2; res++){
373            digi[res] = Phe2ADCch(sector[0], sector[1], pm[sector[0]-1][sector[1]], res) 
374                     + Pedestal(sector[0], sector[1], res);
375         }
376         new(pdigit) AliZDCDigit(sector, digi);
377         treeD->Fill();
378
379         //Ch. debug
380         //printf("\t DIGIT added -> det %d quad %d - digi[0,1] = [%d, %d]\n",
381         //     sector[0], sector[1], digi[0], digi[1]); // Chiara debugging!
382         
383     }
384   } // Loop over detector
385   // Adding in-time digits for 2 reference PTM signals (after signal ch.)
386   // (for the moment the ref. signal is completely invented assuming a PMgain of 5*10^4!)
387   Int_t sectorRef[2];
388   sectorRef[1] = 5;
389   Int_t sigRef[2];
390   // Reference signal are set to 100 (high gain chain) and 800 (low gain chain)
391   if(fIsCalibration==0) {sigRef[0]=100;  sigRef[1]=800;}
392   else {sigRef[0]=0;  sigRef[1]=0;} // calibration -> simulation of pedestal values
393   //
394   for(Int_t iref=0; iref<2; iref++){
395      sectorRef[0] = 3*iref+1;
396      for(Int_t res=0; res<2; res++){
397        sigRef[res] += Pedestal(sectorRef[0], sectorRef[1], res);
398      }
399      new(pdigit) AliZDCDigit(sectorRef, sigRef);
400      treeD->Fill();     
401      
402      //Ch. debug
403      //printf("\t RefDigit added -> det = %d, quad = %d - digi[0,1] = [%d, %d]\n",
404      //    sectorRef[0], sectorRef[1], sigRef[0], sigRef[1]); // Chiara debugging!     
405   }
406   //
407   // --- Adding digits for out-of-time channels after signal digits
408   for(sector[0]=1; sector[0]<6; sector[0]++){
409     for(sector[1]=0; sector[1]<5; sector[1]++){
410         if((sector[0]==3) && ((sector[1]<1) || (sector[1]>2))) continue;
411         for(Int_t res=0; res<2; res++){
412            digioot[res] = Pedestal(sector[0], sector[1], res); // out-of-time ADCs
413         }
414         new(pdigit) AliZDCDigit(sector, digioot);
415         treeD->Fill();
416
417         //Ch. debug
418         //printf("\t DIGIToot added -> det = %d, quad = %d - digi[0,1] = [%d, %d]\n",
419         //     sector[0], sector[1], digioot[0], digioot[1]); // Chiara debugging!      
420     }
421   }
422   // Adding out-of-time digits for 2 reference PTM signals (after out-of-time ch.)
423   Int_t sigRefoot[2];
424   for(Int_t iref=0; iref<2; iref++){
425      sectorRef[0] = 3*iref+1;
426      for(Int_t res=0; res<2; res++){
427        sigRefoot[res] = Pedestal(sectorRef[0], sectorRef[1], res);
428      }
429      new(pdigit) AliZDCDigit(sectorRef, sigRefoot);
430      treeD->Fill();
431      //Ch. debug
432      //printf("\t RefDigitoot added -> det = %d, quad = %d - digi[0,1] = [%d, %d]\n",
433      //    sectorRef[0], sectorRef[1], sigRefoot[0], sigRefoot[1]); // Chiara debugging!
434      
435   }
436   //printf("\t AliZDCDigitizer -> TreeD has %d entries\n",(Int_t) treeD->GetEntries());
437
438   // write the output tree
439   loader->WriteDigits("OVERWRITE");
440   loader->UnloadDigits();
441 }
442
443
444 //_____________________________________________________________________________
445 void AliZDCDigitizer::ReadPMTGains()
446 {
447 // Read PMT gain from an external file
448
449   char *fname = gSystem->ExpandPathName("$ALICE_ROOT/ZDC/PMTGainsdata.txt");
450   FILE *fdata = fopen(fname,"r");
451   if(fdata==NULL){
452      AliWarning(" Can't open file $ALICE_ROOT/ZDC/PMTGainsdata.txt to read ZDC PMT Gains\n");
453      AliWarning("  -> ZDC signal will be pedestal!!!!!!!!!!!!\n\n");
454      return;
455   }
456   int read=1;
457   Float_t data[5];
458   Int_t beam[12], det[12];
459   Float_t gain[12], aEne[12], bEne[12];
460   for(int ir=0; ir<12; ir++){
461     for(int ic=0; ic<5; ic++){
462        read = fscanf(fdata,"%f ",&data[ic]);
463        if(read==0) AliDebug(3, " Error in reading PMT gains from external file ");
464     }
465     beam[ir] = data[0];
466     det[ir] = data[1];
467     gain[ir] = data[2];
468     aEne[ir] = data[3];
469     bEne[ir] = data[4];
470   }
471   fclose(fdata);
472   
473   if(((fBeamType.CompareTo("P-P")) == 0)){
474     for(int i=0; i<12; i++){
475       if(beam[i]==0 && fBeamEnergy!=0.){
476         if(det[i]!=31 && det[i]!=32){
477           for(Int_t j=0; j<5; j++) fPMGain[det[i]-1][j] = gain[i]*(aEne[i]/fBeamEnergy+bEne[i]);
478         }
479         else if(det[i] == 31) fPMGain[2][1] = gain[i]*(aEne[i]-fBeamEnergy*bEne[i]);
480         else if(det[i] == 32) fPMGain[2][2] = gain[i]*(aEne[i]-fBeamEnergy*bEne[i]);
481       }
482     }
483     //
484     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",
485         fBeamEnergy, fBeamEnergy, fPMGain[0][0], fPMGain[1][0], fPMGain[2][1], fPMGain[3][0], fPMGain[4][0]);     
486   }
487   else if(((fBeamType.CompareTo("A-A")) == 0)){
488     for(int i=0; i<12; i++){
489       if(beam[i]==1){
490         Float_t scalGainFactor = fBeamEnergy/2760.;
491         if(det[i]!=31 && det[i]!=32){
492           for(Int_t j=0; j<5; j++) fPMGain[det[i]-1][j] = gain[i]/(aEne[i]*scalGainFactor);
493         }
494         else{
495           for(int iq=1; iq<3; iq++) fPMGain[2][iq] = gain[i]/(aEne[i]*scalGainFactor);
496         }
497       }
498      }  
499      //
500      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",
501         fBeamEnergy, fBeamEnergy, fPMGain[0][0], fPMGain[1][0], fPMGain[2][1]);
502   }
503   else if(((fBeamType.CompareTo("p-A")) == 0) || ((fBeamType.CompareTo("P-A")) == 0) 
504        || ((fBeamType.CompareTo("A-p")) == 0) || ((fBeamType.CompareTo("A-P")) == 0)){
505     for(int i=0; i<12; i++){
506       if(beam[i]==0 && fBeamEnergy!=0.){
507         if(det[i]==1 || det[i]==2){
508           for(Int_t j=0; j<5; j++) fPMGain[det[i]-1][j] = gain[i]*(aEne[i]/fBeamEnergy+bEne[i]);
509         }
510       }
511       if(beam[i]==1){
512         Float_t scalGainFactor = fBeamEnergy/2760.;
513         Float_t npartScalingFactor = 208./15.;
514         if(det[i]==4 || det[i]==5){
515           for(Int_t j=0; j<5; j++) fPMGain[det[i]-1][j] = npartScalingFactor*gain[i]/(aEne[i]*scalGainFactor);
516         }
517         else if(det[i]==31 || det[i]==32){
518           for(int iq=1; iq<3; iq++) fPMGain[2][iq] = npartScalingFactor*gain[i]/(aEne[i]*scalGainFactor);
519         }
520       }
521     }
522     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",
523         fPMGain[0][0], fPMGain[1][0], fPMGain[2][1], fPMGain[3][0], fPMGain[4][0]);
524   }
525 }
526
527 //_____________________________________________________________________________
528 void AliZDCDigitizer::CalculatePMTGains()
529 {
530 // Calculate PMT gain according to beam type and beam energy
531   if( ((fBeamType.CompareTo("P-P")) == 0) ||  ((fBeamType.CompareTo("p-p"))) ){
532     // PTM gains rescaled to beam energy for p-p
533     // New correction coefficients for PMT gains needed
534     // to reproduce experimental spectra (from Grazia Jul 2010)
535     if(fBeamEnergy != 0){
536       for(Int_t j = 0; j < 5; j++){
537           fPMGain[0][j] = 1.515831*(661.444/fBeamEnergy+0.000740671)*10000000;
538           fPMGain[1][j] = 0.674234*(864.350/fBeamEnergy+0.00234375)*10000000;
539           fPMGain[3][j] = 1.350938*(661.444/fBeamEnergy+0.000740671)*10000000; 
540           fPMGain[4][j] = 0.678597*(864.350/fBeamEnergy+0.00234375)*10000000;
541       }
542       fPMGain[2][1] = 0.869654*(1.32312-0.000101515*fBeamEnergy)*10000000;
543       fPMGain[2][2] = 1.030883*(1.32312-0.000101515*fBeamEnergy)*10000000;
544       //
545       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",
546         fBeamEnergy, fBeamEnergy, fPMGain[0][0], fPMGain[1][0], fPMGain[2][1], fPMGain[3][0], fPMGain[4][0]);
547      
548     }
549   }
550   else if(((fBeamType.CompareTo("A-A")) == 0)){
551     // PTM gains for Pb-Pb @ 2.7+2.7 A TeV ***************
552     // rescaled for Pb-Pb @ 1.38+1.38 A TeV ***************
553     // Values corrected after 2010 Pb-Pb data taking (7/2/2011 - Ch.)
554     // Experimental data compared to EMD simulation for single nucleon peaks:
555     // ZN gains must be divided by 4, ZP gains by 10!
556     Float_t scalGainFactor = fBeamEnergy/2760.;
557     for(Int_t j = 0; j < 5; j++){
558        fPMGain[0][j] = 50000./(4*scalGainFactor);  // ZNC                
559        fPMGain[1][j] = 100000./(5*scalGainFactor); // ZPC       
560        fPMGain[2][j] = 100000./scalGainFactor;     // ZEM
561        fPMGain[3][j] = 50000./(4*scalGainFactor);  // ZNA                
562        fPMGain[4][j] = 100000./(5*scalGainFactor); // ZPA    
563     }
564     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",
565         fBeamEnergy, fBeamEnergy, fPMGain[0][0], fPMGain[1][0], fPMGain[2][1]);
566   }
567   else if(((fBeamType.CompareTo("p-A")) == 0) || ((fBeamType.CompareTo("P-A"))) ){
568     // PTM gains for Pb-Pb @ 1.38+1.38 A TeV on side A
569     // PTM gains rescaled to beam energy for p-p on side C
570     // WARNING! Energies are set by hand for 2011 pA RUN!!!
571     Float_t scalGainFactor = fBeamEnergy/2760.;
572     Float_t npartScalingFactor = 208./15.;
573     
574     for(Int_t j = 0; j < 5; j++){
575        fPMGain[0][j] = 1.515831*(661.444/fBeamEnergy+0.000740671)*10000000; //ZNC (p)
576        fPMGain[1][j] = 0.674234*(864.350/fBeamEnergy+0.00234375)*10000000;  //ZPC (p)
577        if(j<2) fPMGain[2][j] = npartScalingFactor*100000./scalGainFactor;          // ZEM (Pb)
578        // Npart max scales from 400 in Pb-Pb to ~8 in pPb -> *40.
579        fPMGain[3][j] = npartScalingFactor*50000/(4*scalGainFactor);  // ZNA (Pb)             
580        fPMGain[4][j] = npartScalingFactor*100000/(5*scalGainFactor); // ZPA (Pb)  
581     }
582     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",
583         fPMGain[0][0], fPMGain[1][0], fPMGain[2][1], fPMGain[3][0], fPMGain[4][0]);
584   }
585   else if(((fBeamType.CompareTo("A-p")) == 0) || ((fBeamType.CompareTo("A-P"))) ){
586     // PTM gains for Pb-Pb @ 1.38+1.38 A TeV on side 
587     // PTM gains rescaled to beam energy for p-p on side C
588     // WARNING! Energies are set by hand for 2011 pA RUN!!!
589     Float_t scalGainFactor = fBeamEnergy/2760.;
590     Float_t npartScalingFactor = 208./15.;
591     
592     for(Int_t j = 0; j < 5; j++){
593        fPMGain[3][j] = 1.350938*(661.444/fBeamEnergy+0.000740671)*10000000;  //ZNA (p)
594        fPMGain[4][j] = 0.678597*(864.350/fBeamEnergy+0.00234375)*10000000;   //ZPA (p)
595        // Npart max scales from 400 in Pb-Pb to ~8 in pPb -> *40.
596        fPMGain[1][j] = npartScalingFactor*50000/(4*scalGainFactor);  // ZNC (Pb)             
597        fPMGain[2][j] = npartScalingFactor*100000/(5*scalGainFactor); // ZPC (Pb)  
598     }
599     fPMGain[2][1] = 0.869654*(1.32312-0.000101515*fBeamEnergy)*10000000; // ZEM (pp)
600     fPMGain[2][2] = 1.030883*(1.32312-0.000101515*fBeamEnergy)*10000000; // ZEM (pp)
601     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",
602         fPMGain[0][0], fPMGain[1][0], fPMGain[2][1], fPMGain[3][0], fPMGain[4][0]);
603   }
604 }
605
606 //_____________________________________________________________________________
607 void AliZDCDigitizer::Fragmentation(Float_t impPar, Int_t specN, Int_t specP,
608                                     Int_t &freeSpecN, Int_t &freeSpecP) const
609 {
610 // simulate fragmentation of spectators
611
612   AliZDCFragment frag(impPar);
613
614   // Fragments generation
615   frag.GenerateIMF();
616   Int_t nAlpha = frag.GetNalpha();
617
618   // Attach neutrons
619   frag.AttachNeutrons();
620   Int_t ztot = frag.GetZtot();
621   Int_t ntot = frag.GetNtot();
622   
623   // Removing fragments and alpha pcs
624   freeSpecN = specN-ntot-2*nAlpha;
625   freeSpecP = specP-ztot-2*nAlpha;
626   
627   // Removing deuterons
628   Int_t ndeu = (Int_t) (freeSpecN*frag.DeuteronNumber());
629   freeSpecN -= ndeu;
630   freeSpecP -= ndeu;
631   //
632   if(freeSpecN<0) freeSpecN=0;
633   if(freeSpecP<0) freeSpecP=0;
634   AliDebug(2, Form("FreeSpn = %d, FreeSpp = %d", freeSpecN, freeSpecP));
635 }
636
637 //_____________________________________________________________________________
638 void AliZDCDigitizer::SpectatorSignal(Int_t SpecType, Int_t numEvents, 
639                                       Float_t pm[5][5]) const
640 {
641 // add signal of the spectators
642   
643   TFile *specSignalFile = TFile::Open("$ALICE_ROOT/ZDC/SpectatorSignal.root");
644   if(!specSignalFile || !specSignalFile->IsOpen()) {
645     AliError((" Opening file $ALICE_ROOT/ZDC/SpectatorSignal.root failed\n"));
646     return;
647   }
648
649   TNtuple* zdcSignal=0x0;
650   
651   Float_t sqrtS = 2*fBeamEnergy;
652   //
653   if(TMath::Abs(sqrtS-5500) < 100.){
654     AliInfo(" Extracting signal from SpectatorSignal/energy5500 ");
655     specSignalFile->cd("energy5500");
656     //
657     if(SpecType == 1) {    // --- Signal for projectile spectator neutrons
658       specSignalFile->GetObject("energy5500/ZNCSignal;1",zdcSignal);
659       if(!zdcSignal) AliError("  PROBLEM!!! Can't retrieve ZNCSignal from SpectatorSignal.root file");
660     } 
661     else if(SpecType == 2) { // --- Signal for projectile spectator protons
662       specSignalFile->GetObject("energy5500/ZPCSignal;1",zdcSignal);
663       if(!zdcSignal) AliError("  PROBLEM!!! Can't retrieve ZPCSignal from SpectatorSignal.root file");
664     }
665     else if(SpecType == 3) { // --- Signal for target spectator neutrons
666       specSignalFile->GetObject("energy5500/ZNASignal;1",zdcSignal);
667       if(!zdcSignal) AliError("  PROBLEM!!! Can't retrieve ZNASignal from SpectatorSignal.root file");
668     }
669     else if(SpecType == 4) { // --- Signal for target spectator protons
670       specSignalFile->GetObject("energy5500/ZPASignal;1",zdcSignal);
671       if(!zdcSignal) AliError("  PROBLEM!!! Can't retrieve ZPASignal from SpectatorSignal.root file");
672     }
673   }
674   else if(TMath::Abs(sqrtS-2760) < 100.){
675     AliInfo(" Extracting signal from SpectatorSignal/energy2760 ");
676     specSignalFile->cd("energy2760");
677     //
678     if(SpecType == 1) {    // --- Signal for projectile spectator neutrons
679       specSignalFile->GetObject("energy2760/ZNCSignal;1",zdcSignal);
680       if(!zdcSignal) AliError("  PROBLEM!!! Can't retrieve ZNCSignal from SpectatorSignal.root file");
681     } 
682     else if(SpecType == 2) { // --- Signal for projectile spectator protons
683       specSignalFile->GetObject("energy2760/ZPCSignal;1",zdcSignal);
684       if(!zdcSignal) AliError("  PROBLEM!!! Can't retrieve ZPCSignal from SpectatorSignal.root file");
685     }
686     else if(SpecType == 3) { // --- Signal for target spectator neutrons
687       specSignalFile->GetObject("energy2760/ZNASignal;1",zdcSignal);
688       if(!zdcSignal) AliError("  PROBLEM!!! Can't retrieve ZNASignal from SpectatorSignal.root file");
689     }
690     else if(SpecType == 4) { // --- Signal for target spectator protons
691       specSignalFile->GetObject("energy2760/ZPASignal;1",zdcSignal);
692       if(!zdcSignal) AliError("  PROBLEM!!! Can't retrieve ZPASignal from SpectatorSignal.root file");
693     }
694   }
695   
696   if(!zdcSignal){
697     printf("\n No spectator signal available for ZDC digitization\n");
698     return;
699   } 
700   
701   Int_t nentries = (Int_t) zdcSignal->GetEntries();
702   
703   Float_t *entry;
704   Int_t pl, i, k, iev=0, rnd[125], volume[2];
705   for(pl=0;pl<125;pl++) rnd[pl] = 0;
706   if(numEvents > 125) {
707     AliDebug(2,Form("numEvents (%d) is larger than 125", numEvents));
708     numEvents = 125;
709   }
710   for(pl=0;pl<numEvents;pl++){
711      rnd[pl] = (Int_t) (9999*gRandom->Rndm());
712      if(rnd[pl] >= 9999) rnd[pl] = 9998;
713      //printf(" rnd[%d] = %d\n",pl,rnd[pl]);     
714   }
715   // Sorting vector in ascending order with C function QSORT 
716   qsort((void*)rnd,numEvents,sizeof(Int_t),comp);
717   do{
718      for(i=0; i<nentries; i++){  
719         zdcSignal->GetEvent(i);
720         entry = zdcSignal->GetArgs();
721         if(entry[0] == rnd[iev]){
722           for(k=0; k<2; k++) volume[k] = (Int_t) entry[k+1];
723           //
724           Float_t lightQ = entry[7];
725           Float_t lightC = entry[8];
726           //
727           if(volume[0] != 3) {  // ZN or ZP
728             pm[volume[0]-1][0] += lightC;
729             pm[volume[0]-1][volume[1]] += lightQ;
730             //printf("\n   pm[%d][0] = %.0f, pm[%d][%d] = %.0f\n",(volume[0]-1),pm[volume[0]-1][0],
731             //  (volume[0]-1),volume[1],pm[volume[0]-1][volume[1]]);
732           } 
733           else { 
734             if(volume[1] == 1) pm[2][1] += lightC; // ZEM 1
735             else               pm[2][2] += lightQ; // ZEM 2
736             //printf("\n   pm[2][1] = %.0f, pm[2][2] = %.0f\n",pm[2][1],pm[2][2]);
737           }
738         }
739         else if(entry[0] > rnd[iev]){
740           iev++;
741           continue;
742         }
743      }
744   }while(iev<numEvents);
745   
746   specSignalFile->Close();
747   delete specSignalFile;
748 }
749
750
751 //_____________________________________________________________________________
752 Int_t AliZDCDigitizer::Phe2ADCch(Int_t Det, Int_t Quad, Float_t Light, 
753                                  Int_t Res) const
754 {
755   // Evaluation of the ADC channel corresponding to the light yield Light
756   Int_t vADCch = (Int_t) (Light * fPMGain[Det-1][Quad] * fADCRes[Res]);
757   // Ch. debug
758   //printf("\t Phe2ADCch -> det %d quad %d - PMGain[%d][%d] %1.0f phe %1.2f  ADC %d\n", 
759   //    Det,Quad,Det-1,Quad,fPMGain[Det-1][Quad],Light,vADCch);
760
761   return vADCch;
762 }
763
764 //_____________________________________________________________________________
765 Int_t AliZDCDigitizer::Pedestal(Int_t Det, Int_t Quad, Int_t Res) const
766 {
767   // Returns a pedestal for detector det, PM quad, channel with res.
768   //
769   Float_t pedValue;
770   // Normal run
771   if(fIsCalibration == 0){
772     Int_t index=0, kNch=24;
773     if(Quad!=5){
774       if(Det==1)        index = Quad+kNch*Res;    // ZNC
775       else if(Det==2)   index = (Quad+5)+kNch*Res;  // ZPC
776       else if(Det==3)   index = (Quad+9)+kNch*Res;  // ZEM
777       else if(Det==4)   index = (Quad+12)+kNch*Res; // ZNA
778       else if(Det==5)   index = (Quad+17)+kNch*Res; // ZPA
779     }
780     else index = (Det-1)/3+22+kNch*Res; // Reference PMs
781     //
782     Float_t meanPed = fPedData->GetMeanPed(index);
783     Float_t pedWidth = fPedData->GetMeanPedWidth(index);
784     pedValue = gRandom->Gaus(meanPed,pedWidth);
785     //
786     /*printf("\t  AliZDCDigitizer::Pedestal -> det %d quad %d res %d - Ped[%d] = %d\n",
787         Det, Quad, Res, index,(Int_t) pedValue); // Chiara debugging!
788     */
789   }
790   // To create calibration object
791   else{
792     if(Res == 0) pedValue = gRandom->Gaus((35.+10.*gRandom->Rndm()),(0.5+0.2*gRandom->Rndm())); //High gain
793     else  pedValue = gRandom->Gaus((250.+100.*gRandom->Rndm()),(3.5+2.*gRandom->Rndm())); //Low gain
794   }
795
796   return (Int_t) pedValue;
797 }
798
799 //_____________________________________________________________________________
800 AliCDBStorage* AliZDCDigitizer::SetStorage(const char *uri) 
801 {
802
803   Bool_t deleteManager = kFALSE;
804   
805   AliCDBManager *manager = AliCDBManager::Instance();
806   AliCDBStorage *defstorage = manager->GetDefaultStorage();
807   
808   if(!defstorage || !(defstorage->Contains("ZDC"))){ 
809      AliWarning("No default storage set or default storage doesn't contain ZDC!");
810      manager->SetDefaultStorage(uri);
811      deleteManager = kTRUE;
812   }
813  
814   AliCDBStorage *storage = manager->GetDefaultStorage();
815
816   if(deleteManager){
817     AliCDBManager::Instance()->UnsetDefaultStorage();
818     defstorage = 0;   // the storage is killed by AliCDBManager::Instance()->Destroy()
819   }
820
821   return storage; 
822 }
823
824 //_____________________________________________________________________________
825 AliZDCPedestals* AliZDCDigitizer::GetPedData() const
826 {
827
828   // Getting pedestal calibration object for ZDC set
829
830   AliCDBEntry  *entry = AliCDBManager::Instance()->Get("ZDC/Calib/Pedestals");
831   if(!entry) AliFatal("No calibration data loaded!");  
832
833   AliZDCPedestals *calibdata = dynamic_cast<AliZDCPedestals*>  (entry->GetObject());
834   if(!calibdata)  AliFatal("Wrong calibration object in calibration  file!");
835
836   return calibdata;
837 }
838