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