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