]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ZDC/AliZDCDigitizer.cxx
Update to simulate pA interactions with slow nucleon tracking
[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 #include <stdlib.h>
25
26 // --- ROOT system
27 #include <TTree.h>
28 #include <TFile.h>
29 #include <TNtuple.h>
30 #include <TRandom.h>
31
32 // --- AliRoot header files
33 #include "AliLog.h"
34 #include "AliRun.h"
35 #include "AliHeader.h"
36 #include "AliGenHijingEventHeader.h"
37 #include "AliDigitizationInput.h"
38 #include "AliRunLoader.h"
39 #include "AliLoader.h"
40 #include "AliGRPObject.h"
41 #include "AliCDBManager.h"
42 #include "AliCDBEntry.h"
43 #include "AliZDCSDigit.h"
44 #include "AliZDCDigit.h"
45 #include "AliZDCFragment.h"
46 #include "AliZDCv3.h"
47 #include "AliZDCDigitizer.h"
48
49 class AliCDBStorage;
50 class AliZDCPedestals;
51
52 ClassImp(AliZDCDigitizer)
53
54
55 //____________________________________________________________________________
56 AliZDCDigitizer::AliZDCDigitizer() :
57   fIsCalibration(0), 
58   fIsSignalInADCGate(kFALSE),
59   fFracLostSignal(0.),
60   fPedData(0), 
61   fSpectators2Track(kFALSE),
62   fBeamEnergy(0.),
63   fIspASystem(kFALSE)
64 {
65   // Default constructor    
66   for(Int_t i=0; i<2; i++) fADCRes[i]=0.;
67 }
68
69 //____________________________________________________________________________
70 AliZDCDigitizer::AliZDCDigitizer(AliDigitizationInput* digInput):
71   AliDigitizer(digInput),
72   fIsCalibration(0), //By default the simulation doesn't create calib. data
73   fIsSignalInADCGate(kFALSE),
74   fFracLostSignal(0.),
75   fPedData(GetPedData()), 
76   fSpectators2Track(kFALSE),
77   fBeamEnergy(0.),
78   fIspASystem(kFALSE)
79 {
80   // Get calibration data
81   if(fIsCalibration!=0) printf("\n\t AliZDCDigitizer -> Creating calibration data (pedestals)\n");
82   for(Int_t i=0; i<5; i++){
83     for(Int_t j=0; j<5; j++)
84        fPMGain[i][j] = 0.;
85   }
86 }
87
88 //____________________________________________________________________________
89 AliZDCDigitizer::~AliZDCDigitizer()
90 {
91 // Destructor
92 // Not implemented
93 }
94
95
96 //____________________________________________________________________________
97 AliZDCDigitizer::AliZDCDigitizer(const AliZDCDigitizer &digitizer):
98   AliDigitizer(),
99   fIsCalibration(digitizer.fIsCalibration),
100   fIsSignalInADCGate(digitizer.fIsSignalInADCGate),
101   fFracLostSignal(digitizer.fFracLostSignal),
102   fPedData(digitizer.fPedData),
103   fSpectators2Track(digitizer.fSpectators2Track),
104   fBeamEnergy(digitizer.fBeamEnergy),
105   fIspASystem(digitizer.fIspASystem)
106 {
107   // Copy constructor
108
109   for(Int_t i=0; i<5; i++){
110      for(Int_t j=0; j<5; j++){
111         fPMGain[i][j]   = digitizer.fPMGain[i][j];           
112      }
113   }
114   for(Int_t i=0; i<2; i++) fADCRes[i] = digitizer.fADCRes[i];
115
116
117 }
118
119 //____________________________________________________________________________
120 Bool_t AliZDCDigitizer::Init()
121 {
122   // Initialize the digitizer
123   
124   AliCDBEntry*  entry = AliCDBManager::Instance()->Get("GRP/GRP/Data");
125   if(!entry) AliFatal("No calibration data loaded!");  
126   AliGRPObject* grpData = 0x0;
127   if(entry){
128     TMap* m = dynamic_cast<TMap*>(entry->GetObject());  // old GRP entry
129     if(m){
130       //m->Print();
131       grpData = new AliGRPObject();
132       grpData->ReadValuesFromMap(m);
133     }
134     else{
135       grpData = dynamic_cast<AliGRPObject*>(entry->GetObject());  // new GRP entry
136     }
137     entry->SetOwner(0);
138     AliCDBManager::Instance()->UnloadFromCache("GRP/GRP/Data");
139   }
140   if(!grpData){
141     AliError("No GRP entry found in OCDB! \n ");
142     return kFALSE;
143   }
144   
145   TString beamType = grpData->GetBeamType();
146   if(beamType==AliGRPObject::GetInvalidString()){
147     AliError("\t UNKNOWN beam type from GRP obj -> PMT gains not set in ZDC digitizer!!!\n");
148   }
149   
150   fBeamEnergy = grpData->GetBeamEnergy();
151   printf("   beam energy = %f GeV\n", fBeamEnergy);
152   if(fBeamEnergy==AliGRPObject::GetInvalidFloat()){
153     AliWarning("GRP/GRP/Data entry:  missing value for the beam energy ! Using 0.");
154     AliError("\t UNKNOWN beam type from GRP obj -> PMT gains not set in ZDC digitizer!!!\n");
155     fBeamEnergy = 0.;
156   }
157
158   if((beamType.CompareTo("P-P")) == 0 || (beamType.CompareTo("p-p")) == 0){
159     // PTM gains rescaled to beam energy for p-p
160     // New correction coefficients for PMT gains needed
161     // to reproduce experimental spectra (from Grazia Jul 2010)
162     if(fBeamEnergy != 0){
163       for(Int_t j = 0; j < 5; j++){
164           fPMGain[0][j] = 1.515831*(661.444/fBeamEnergy+0.000740671)*10000000;
165           fPMGain[1][j] = 0.674234*(864.350/fBeamEnergy+0.00234375)*10000000;
166           fPMGain[3][j] = 1.350938*(661.444/fBeamEnergy+0.000740671)*10000000; 
167           fPMGain[4][j] = 0.678597*(864.350/fBeamEnergy+0.00234375)*10000000;
168       }
169       fPMGain[2][1] = 0.869654*(1.32312-0.000101515*fBeamEnergy)*10000000;
170       fPMGain[2][2] = 1.030883*(1.32312-0.000101515*fBeamEnergy)*10000000;
171       //
172       AliInfo(Form("    PMT gains for p-p @ %1.0f+%1.0f GeV: ZN(%1.0f), ZP(%1.0f), ZEM(%1.0f)\n",
173         fBeamEnergy, fBeamEnergy, fPMGain[0][0], fPMGain[1][0], fPMGain[2][1]));
174     }
175     else{ // for RELDIS simulation @ sqrt(s_{NN}) = 2.76 TeV
176        Float_t scalGainFactor = 0.5;
177        for(Int_t j = 0; j < 5; j++){
178          fPMGain[0][j] = 50000./scalGainFactor;  // ZNC          
179          fPMGain[1][j] = 100000./scalGainFactor; // ZPC    
180          fPMGain[2][j] = 100000./scalGainFactor; // ZEM
181          fPMGain[3][j] = 50000./scalGainFactor;  // ZNA          
182          fPMGain[4][j] = 100000./scalGainFactor; // ZPA 
183        }
184        //
185        AliInfo(Form("   PMT gains for RELDIS simulation: ZN(%1.0f), ZP(%1.0f), ZEM(%1.0f)\n",
186           fPMGain[0][0], fPMGain[1][0], fPMGain[2][1]));
187     }
188   }
189   else if((beamType.CompareTo("A-A")) == 0){
190     // PTM gains for Pb-Pb @ 2.7+2.7 A TeV ***************
191     // rescaled for Pb-Pb @ 1.38+1.38 A TeV ***************
192     // Values corrected after 2010 Pb-Pb data taking (7/2/2011 - Ch.)
193     // Experimental data compared to EMD simulation for single nucleon peaks:
194     // ZN gains must be divided by 4, ZP gains by 10!
195     Float_t scalGainFactor = fBeamEnergy/2760.;
196     for(Int_t j = 0; j < 5; j++){
197        fPMGain[0][j] = 50000./(4*scalGainFactor);  // ZNC                
198        fPMGain[1][j] = 100000./(5*scalGainFactor); // ZPC       
199        fPMGain[2][j] = 100000./scalGainFactor;     // ZEM
200        fPMGain[3][j] = 50000./(4*scalGainFactor);  // ZNA                
201        fPMGain[4][j] = 100000./(5*scalGainFactor); // ZPA    
202     }
203     AliInfo(Form("    PMT gains for Pb-Pb @ %1.0f+%1.0f A GeV: ZN(%1.0f), ZP(%1.0f), ZEM(%1.0f)\n",
204         fBeamEnergy, fBeamEnergy, fPMGain[0][0], fPMGain[1][0], fPMGain[2][1]));
205   }
206   
207   if(fIspASystem){
208     // PTM gains for Pb-Pb @ 1.38+1.38 A TeV on side C
209     // PTM gains rescaled to beam energy for p-p on side A
210     Float_t scalGainFactor = fBeamEnergy/2760.;
211     if(TMath::Abs(scalGainFactor)<0.001) scalGainFactor=0.5;
212     Float_t fpBeamEnergy = 3500.;
213     
214     for(Int_t j = 0; j < 5; j++){
215        fPMGain[0][j] = 50000./(4*scalGainFactor);  // ZNC (Pb)           
216        fPMGain[1][j] = 100000./(5*scalGainFactor); // ZPC (Pb)  
217        fPMGain[3][j] = 1.350938*(661.444/fpBeamEnergy+0.000740671)*10000000; //ZNA (p)
218        fPMGain[4][j] = 0.678597*(864.350/fpBeamEnergy+0.00234375)*10000000;  //ZPA (p)
219     }
220     // ZEM (p)
221     fPMGain[2][1] = 0.869654*(1.32312-0.000101515*fBeamEnergy)*10000000;
222     fPMGain[2][2] = 1.030883*(1.32312-0.000101515*fBeamEnergy)*10000000;
223     AliInfo(Form("    PMT gains for p-Pb: ZNC(%1.0f), ZPC(%1.0f), ZEM(%1.0f), ZNA(%1.0f) ZPA(%1.0f)\n",
224         fPMGain[0][0], fPMGain[1][0], fPMGain[2][1], fPMGain[3][0], fPMGain[4][0]));
225   }
226     
227   // ADC Caen V965
228   fADCRes[0] = 0.0000008; // ADC Resolution high gain: 200 fC/adcCh
229   fADCRes[1] = 0.0000064; // ADC Resolution low gain:  25  fC/adcCh
230
231   return kTRUE;
232 }
233
234 //____________________________________________________________________________
235 void AliZDCDigitizer::Digitize(Option_t* /*option*/)
236 {
237   // Execute digitization
238
239   // ------------------------------------------------------------
240   // !!! 2nd ZDC set added 
241   // *** 1st 3 arrays are digits from REAL (simulated) hits
242   // *** last 2 are copied from simulated digits
243   // --- pm[0][...] = light in ZN side C  [C, Q1, Q2, Q3, Q4]
244   // --- pm[1][...] = light in ZP side C [C, Q1, Q2, Q3, Q4]
245   // --- pm[2][...] = light in ZEM [x, 1, 2, x, x]
246   // --- pm[3][...] = light in ZN side A [C, Q1, Q2, Q3, Q4] ->NEW!
247   // --- pm[4][...] = light in ZP side A [C, Q1, Q2, Q3, Q4] ->NEW!
248   // ------------------------------------------------------------
249   Float_t pm[5][5]; 
250   for(Int_t iSector1=0; iSector1<5; iSector1++) 
251     for(Int_t iSector2=0; iSector2<5; iSector2++){
252       pm[iSector1][iSector2] = 0;
253     }
254     
255   // ------------------------------------------------------------
256   // ### Out of time ADC added (22 channels)
257   // --- same codification as for signal PTMs (see above)
258   // ------------------------------------------------------------
259   Float_t pmoot[5][5];
260   for(Int_t iSector1=0; iSector1<5; iSector1++) 
261     for(Int_t iSector2=0; iSector2<5; iSector2++){
262       pmoot[iSector1][iSector2] = 0;
263     }
264
265   // impact parameter and number of spectators
266   Float_t impPar = 0;
267   Int_t specNTarg = 0, specPTarg = 0;
268   Int_t specNProj = 0, specPProj = 0;
269   Float_t signalTime0 = 0.;
270
271   // loop over input streams
272   for(Int_t iInput = 0; iInput<fDigInput->GetNinputs(); iInput++){
273
274     // get run loader and ZDC loader
275     AliRunLoader* runLoader = 
276       AliRunLoader::GetRunLoader(fDigInput->GetInputFolderName(iInput));
277     AliLoader* loader = runLoader->GetLoader("ZDCLoader");
278     if(!loader) continue;
279
280     // load sdigits
281     loader->LoadSDigits();
282     TTree* treeS = loader->TreeS();
283     if(!treeS) continue;
284     AliZDCSDigit sdigit;
285     AliZDCSDigit* psdigit = &sdigit;
286     treeS->SetBranchAddress("ZDC", &psdigit);
287
288     // loop over sdigits
289     for(Int_t iSDigit=0; iSDigit<treeS->GetEntries(); iSDigit++){
290       treeS->GetEntry(iSDigit);
291       //
292       if(!psdigit) continue;
293       if((sdigit.GetSector(1) < 0) || (sdigit.GetSector(1) > 4)){
294         AliError(Form("\nsector[0] = %d, sector[1] = %d\n", 
295                       sdigit.GetSector(0), sdigit.GetSector(1)));
296         continue;
297       }
298       // Checking if signal is inside ADC gate
299       if(iSDigit==0) signalTime0 = sdigit.GetTrackTime();
300       else{
301         // Assuming a signal lenght of 20 ns, signal is in gate if 
302         // signal ENDS in signalTime0+50. -> sdigit.GetTrackTime()+20<=signalTime0+50.
303         if(sdigit.GetTrackTime()<=signalTime0+30.) fIsSignalInADCGate = kTRUE;
304         if(sdigit.GetTrackTime()>signalTime0+30.){
305           fIsSignalInADCGate = kFALSE;
306           // Vedi quaderno per spiegazione approx. usata 
307           // nel calcolo della fraz. di segnale perso
308           fFracLostSignal = (sdigit.GetTrackTime()-30)*(sdigit.GetTrackTime()-30)/280.;
309         }
310       }
311       
312       Float_t sdSignal = sdigit.GetLightPM();
313       if(fIsSignalInADCGate == kFALSE){
314         AliDebug(2,Form("\t Signal time %f -> fraction %f of ZDC signal for det.(%d, %d) out of ADC gate\n",
315                 sdigit.GetTrackTime(),fFracLostSignal,sdigit.GetSector(0),sdigit.GetSector(1)));
316         sdSignal = (1-fFracLostSignal)*sdSignal;
317       }
318       
319       pm[(sdigit.GetSector(0))-1][sdigit.GetSector(1)] += sdigit.GetLightPM();
320       //Ch. debug
321       /*printf("\t Detector %d, Tower %d -> pm[%d][%d] = %.0f \n",
322           sdigit.GetSector(0), sdigit.GetSector(1),sdigit.GetSector(0)-1,
323           sdigit.GetSector(1), pm[sdigit.GetSector(0)-1][sdigit.GetSector(1)]); 
324       */
325       
326     }
327
328     loader->UnloadSDigits();
329
330     // get the impact parameter and the number of spectators in case of hijing
331     if(!runLoader->GetAliRun()) runLoader->LoadgAlice();
332     AliHeader* header = runLoader->GetHeader();
333     if(!header) continue;
334     AliGenEventHeader* genHeader = header->GenEventHeader();
335     if(!genHeader) continue;
336     if(!genHeader->InheritsFrom(AliGenHijingEventHeader::Class())) continue;
337     
338     if(fSpectators2Track==kTRUE && fIspASystem==kFALSE){
339       impPar = ((AliGenHijingEventHeader*) genHeader)->ImpactParameter(); 
340       specNProj = ((AliGenHijingEventHeader*) genHeader)->ProjSpectatorsn();
341       specPProj = ((AliGenHijingEventHeader*) genHeader)->ProjSpectatorsp();
342       specNTarg = ((AliGenHijingEventHeader*) genHeader)->TargSpectatorsn();
343       specPTarg = ((AliGenHijingEventHeader*) genHeader)->TargSpectatorsp();
344       /*printf("\n\t AliZDCDigitizer: b = %1.2f fm\n"
345       " \t    PROJ.:  #spectator n %d, #spectator p %d\n"
346       " \t    TARG.:  #spectator n %d, #spectator p %d\n", 
347       impPar, specNProj, specPProj, specNTarg, specPTarg);*/
348     }
349     
350   }
351
352   // Applying fragmentation algorithm and adding spectator signal
353   if(fSpectators2Track==kTRUE && impPar && fIspASystem==kFALSE) {
354     Int_t freeSpecNProj, freeSpecPProj;
355     Fragmentation(impPar, specNProj, specPProj, freeSpecNProj, freeSpecPProj);
356     Int_t freeSpecNTarg, freeSpecPTarg;
357     Fragmentation(impPar, specNTarg, specPTarg, freeSpecNTarg, freeSpecPTarg);
358     SpectatorSignal(1, freeSpecNProj, pm);
359 //    printf("\t AliZDCDigitizer -> Adding signal for %d PROJ free spectator n",freeSpecNProj);
360     SpectatorSignal(2, freeSpecPProj, pm);
361 //    printf(" and %d free spectator p\n",freeSpecPProj);
362     SpectatorSignal(3, freeSpecNTarg, pm);
363 //    printf("\t AliZDCDigitizer -> Adding signal for %d TARG free spectator n",freeSpecNTarg);
364     SpectatorSignal(4, freeSpecPTarg, pm);
365 //    printf(" and %d free spectator p\n\n",freeSpecPTarg);
366   }
367
368
369   // get the output run loader and loader
370   AliRunLoader* runLoader = 
371     AliRunLoader::GetRunLoader(fDigInput->GetOutputFolderName());
372   AliLoader* loader = runLoader->GetLoader("ZDCLoader");
373   if(!loader) {
374     AliError("no ZDC loader found");
375     return;
376   }
377
378   // create the output tree
379   const char* mode = "update";
380   if(runLoader->GetEventNumber() == 0) mode = "recreate";
381   loader->LoadDigits(mode);
382   loader->MakeTree("D");
383   TTree* treeD = loader->TreeD();
384   AliZDCDigit digit;
385   AliZDCDigit* pdigit = &digit;
386   const Int_t kBufferSize = 4000;
387   treeD->Branch("ZDC", "AliZDCDigit", &pdigit, kBufferSize);
388
389   // Create digits
390   Int_t sector[2];
391   Int_t digi[2], digioot[2];
392   for(sector[0]=1; sector[0]<6; sector[0]++){
393     for(sector[1]=0; sector[1]<5; sector[1]++){
394         if((sector[0]==3) && ((sector[1]<1) || (sector[1]>2))) continue;
395         for(Int_t res=0; res<2; res++){
396            digi[res] = Phe2ADCch(sector[0], sector[1], pm[sector[0]-1][sector[1]], res) 
397                     + Pedestal(sector[0], sector[1], res);
398         }
399         new(pdigit) AliZDCDigit(sector, digi);
400         treeD->Fill();
401
402         //Ch. debug
403         //printf("\t DIGIT added -> det %d quad %d - digi[0,1] = [%d, %d]\n",
404         //     sector[0], sector[1], digi[0], digi[1]); // Chiara debugging!
405         
406     }
407   } // Loop over detector
408   // Adding in-time digits for 2 reference PTM signals (after signal ch.)
409   // (for the moment the ref. signal is completely invented assuming a PMgain of 5*10^4!)
410   Int_t sectorRef[2];
411   sectorRef[1] = 5;
412   Int_t sigRef[2];
413   // Reference signal are set to 100 (high gain chain) and 800 (low gain chain)
414   if(fIsCalibration==0) {sigRef[0]=100;  sigRef[1]=800;}
415   else {sigRef[0]=0;  sigRef[1]=0;} // calibration -> simulation of pedestal values
416   //
417   for(Int_t iref=0; iref<2; iref++){
418      sectorRef[0] = 3*iref+1;
419      for(Int_t res=0; res<2; res++){
420        sigRef[res] += Pedestal(sectorRef[0], sectorRef[1], res);
421      }
422      new(pdigit) AliZDCDigit(sectorRef, sigRef);
423      treeD->Fill();     
424      
425      //Ch. debug
426      //printf("\t RefDigit added -> det = %d, quad = %d - digi[0,1] = [%d, %d]\n",
427      //    sectorRef[0], sectorRef[1], sigRef[0], sigRef[1]); // Chiara debugging!     
428   }
429   //
430   // --- Adding digits for out-of-time channels after signal digits
431   for(sector[0]=1; sector[0]<6; sector[0]++){
432     for(sector[1]=0; sector[1]<5; sector[1]++){
433         if((sector[0]==3) && ((sector[1]<1) || (sector[1]>2))) continue;
434         for(Int_t res=0; res<2; res++){
435            digioot[res] = Pedestal(sector[0], sector[1], res); // out-of-time ADCs
436         }
437         new(pdigit) AliZDCDigit(sector, digioot);
438         treeD->Fill();
439
440         //Ch. debug
441         //printf("\t DIGIToot added -> det = %d, quad = %d - digi[0,1] = [%d, %d]\n",
442         //     sector[0], sector[1], digioot[0], digioot[1]); // Chiara debugging!      
443     }
444   }
445   // Adding out-of-time digits for 2 reference PTM signals (after out-of-time ch.)
446   Int_t sigRefoot[2];
447   for(Int_t iref=0; iref<2; iref++){
448      sectorRef[0] = 3*iref+1;
449      for(Int_t res=0; res<2; res++){
450        sigRefoot[res] = Pedestal(sectorRef[0], sectorRef[1], res);
451      }
452      new(pdigit) AliZDCDigit(sectorRef, sigRefoot);
453      treeD->Fill();
454      //Ch. debug
455      //printf("\t RefDigitoot added -> det = %d, quad = %d - digi[0,1] = [%d, %d]\n",
456      //    sectorRef[0], sectorRef[1], sigRefoot[0], sigRefoot[1]); // Chiara debugging!
457      
458   }
459   //printf("\t AliZDCDigitizer -> TreeD has %d entries\n",(Int_t) treeD->GetEntries());
460
461   // write the output tree
462   loader->WriteDigits("OVERWRITE");
463   loader->UnloadDigits();
464 }
465
466
467 //_____________________________________________________________________________
468 void AliZDCDigitizer::Fragmentation(Float_t impPar, Int_t specN, Int_t specP,
469                                     Int_t &freeSpecN, Int_t &freeSpecP) const
470 {
471 // simulate fragmentation of spectators
472
473   AliZDCFragment frag(impPar);
474
475   // Fragments generation
476   frag.GenerateIMF();
477   Int_t nAlpha = frag.GetNalpha();
478
479   // Attach neutrons
480   frag.AttachNeutrons();
481   Int_t ztot = frag.GetZtot();
482   Int_t ntot = frag.GetNtot();
483   
484   // Removing fragments and alpha pcs
485   freeSpecN = specN-ntot-2*nAlpha;
486   freeSpecP = specP-ztot-2*nAlpha;
487   
488   // Removing deuterons
489   Int_t ndeu = (Int_t) (freeSpecN*frag.DeuteronNumber());
490   freeSpecN -= ndeu;
491   freeSpecP -= ndeu;
492   //
493   if(freeSpecN<0) freeSpecN=0;
494   if(freeSpecP<0) freeSpecP=0;
495   AliDebug(2, Form("FreeSpn = %d, FreeSpp = %d", freeSpecN, freeSpecP));
496 }
497
498 //_____________________________________________________________________________
499 void AliZDCDigitizer::SpectatorSignal(Int_t SpecType, Int_t numEvents, 
500                                       Float_t pm[5][5]) const
501 {
502 // add signal of the spectators
503   
504   TFile *specSignalFile = TFile::Open("$ALICE_ROOT/ZDC/SpectatorSignal.root");
505   if(!specSignalFile || !specSignalFile->IsOpen()) {
506     AliError((" Opening file $ALICE_ROOT/ZDC/SpectatorSignal.root failed\n"));
507     return;
508   }
509
510   TNtuple* zdcSignal=0x0;
511   
512   Float_t sqrtS = 2*fBeamEnergy;
513   //
514   if(TMath::Abs(sqrtS-5500) < 100.){
515     specSignalFile->cd("energy5500");
516     //
517     if(SpecType == 1) {    // --- Signal for projectile spectator neutrons
518       specSignalFile->GetObject("energy5500/ZNCSignal;1",zdcSignal);
519       if(!zdcSignal) AliError("  PROBLEM!!! Can't retrieve ZNCSignal from SpectatorSignal.root file");
520     } 
521     else if(SpecType == 2) { // --- Signal for projectile spectator protons
522       specSignalFile->GetObject("energy5500/ZPCSignal;1",zdcSignal);
523       if(!zdcSignal) AliError("  PROBLEM!!! Can't retrieve ZPCSignal from SpectatorSignal.root file");
524     }
525     else if(SpecType == 3) { // --- Signal for target spectator neutrons
526       specSignalFile->GetObject("energy5500/ZNASignal;1",zdcSignal);
527       if(!zdcSignal) AliError("  PROBLEM!!! Can't retrieve ZNASignal from SpectatorSignal.root file");
528     }
529     else if(SpecType == 4) { // --- Signal for target spectator protons
530       specSignalFile->GetObject("energy5500/ZPASignal;1",zdcSignal);
531       if(!zdcSignal) AliError("  PROBLEM!!! Can't retrieve ZPASignal from SpectatorSignal.root file");
532     }
533   }
534   else if(TMath::Abs(sqrtS-2760) < 100.){
535     specSignalFile->cd("energy2760");
536     //
537     if(SpecType == 1) {    // --- Signal for projectile spectator neutrons
538       specSignalFile->GetObject("energy2760/ZNCSignal;1",zdcSignal);
539       if(!zdcSignal) AliError("  PROBLEM!!! Can't retrieve ZNCSignal from SpectatorSignal.root file");
540     } 
541     else if(SpecType == 2) { // --- Signal for projectile spectator protons
542       specSignalFile->GetObject("energy2760/ZPCSignal;1",zdcSignal);
543       if(!zdcSignal) AliError("  PROBLEM!!! Can't retrieve ZPCSignal from SpectatorSignal.root file");
544     }
545     else if(SpecType == 3) { // --- Signal for target spectator neutrons
546       specSignalFile->GetObject("energy2760/ZNASignal;1",zdcSignal);
547       if(!zdcSignal) AliError("  PROBLEM!!! Can't retrieve ZNASignal from SpectatorSignal.root file");
548     }
549     else if(SpecType == 4) { // --- Signal for target spectator protons
550       specSignalFile->GetObject("energy2760/ZPASignal;1",zdcSignal);
551       if(!zdcSignal) AliError("  PROBLEM!!! Can't retrieve ZPASignal from SpectatorSignal.root file");
552     }
553   }
554   
555   if(!zdcSignal){
556     printf("\n No spectator signal available for ZDC digitization\n");
557     return;
558   } 
559   
560   Int_t nentries = (Int_t) zdcSignal->GetEntries();
561   
562   Float_t *entry;
563   Int_t pl, i, k, iev=0, rnd[125], volume[2];
564   for(pl=0;pl<125;pl++) rnd[pl] = 0;
565   if(numEvents > 125) {
566     AliDebug(2,Form("numEvents (%d) is larger than 125", numEvents));
567     numEvents = 125;
568   }
569   for(pl=0;pl<numEvents;pl++){
570      rnd[pl] = (Int_t) (9999*gRandom->Rndm());
571      if(rnd[pl] >= 9999) rnd[pl] = 9998;
572      //printf(" rnd[%d] = %d\n",pl,rnd[pl]);     
573   }
574   // Sorting vector in ascending order with C function QSORT 
575   qsort((void*)rnd,numEvents,sizeof(Int_t),comp);
576   do{
577      for(i=0; i<nentries; i++){  
578         zdcSignal->GetEvent(i);
579         entry = zdcSignal->GetArgs();
580         if(entry[0] == rnd[iev]){
581           for(k=0; k<2; k++) volume[k] = (Int_t) entry[k+1];
582           //
583           Float_t lightQ = entry[7];
584           Float_t lightC = entry[8];
585           //
586           if(volume[0] != 3) {  // ZN or ZP
587             pm[volume[0]-1][0] += lightC;
588             pm[volume[0]-1][volume[1]] += lightQ;
589             //printf("\n   pm[%d][0] = %.0f, pm[%d][%d] = %.0f\n",(volume[0]-1),pm[volume[0]-1][0],
590             //  (volume[0]-1),volume[1],pm[volume[0]-1][volume[1]]);
591           } 
592           else { 
593             if(volume[1] == 1) pm[2][1] += lightC; // ZEM 1
594             else               pm[2][2] += lightQ; // ZEM 2
595             //printf("\n   pm[2][1] = %.0f, pm[2][2] = %.0f\n",pm[2][1],pm[2][2]);
596           }
597         }
598         else if(entry[0] > rnd[iev]){
599           iev++;
600           continue;
601         }
602      }
603   }while(iev<numEvents);
604   
605   specSignalFile->Close();
606   delete specSignalFile;
607 }
608
609
610 //_____________________________________________________________________________
611 Int_t AliZDCDigitizer::Phe2ADCch(Int_t Det, Int_t Quad, Float_t Light, 
612                                  Int_t Res) const
613 {
614   // Evaluation of the ADC channel corresponding to the light yield Light
615   Int_t vADCch = (Int_t) (Light * fPMGain[Det-1][Quad] * fADCRes[Res]);
616   // Ch. debug
617   //printf("\t Phe2ADCch -> det %d quad %d - PMGain[%d][%d] %1.0f phe %1.2f  ADC %d\n", 
618   //    Det,Quad,Det-1,Quad,fPMGain[Det-1][Quad],Light,vADCch);
619
620   return vADCch;
621 }
622
623 //_____________________________________________________________________________
624 Int_t AliZDCDigitizer::Pedestal(Int_t Det, Int_t Quad, Int_t Res) const
625 {
626   // Returns a pedestal for detector det, PM quad, channel with res.
627   //
628   Float_t pedValue;
629   // Normal run
630   if(fIsCalibration == 0){
631     Int_t index=0, kNch=24;
632     if(Quad!=5){
633       if(Det==1)        index = Quad+kNch*Res;    // ZNC
634       else if(Det==2)   index = (Quad+5)+kNch*Res;  // ZPC
635       else if(Det==3)   index = (Quad+9)+kNch*Res;  // ZEM
636       else if(Det==4)   index = (Quad+12)+kNch*Res; // ZNA
637       else if(Det==5)   index = (Quad+17)+kNch*Res; // ZPA
638     }
639     else index = (Det-1)/3+22+kNch*Res; // Reference PMs
640     //
641     Float_t meanPed = fPedData->GetMeanPed(index);
642     Float_t pedWidth = fPedData->GetMeanPedWidth(index);
643     pedValue = gRandom->Gaus(meanPed,pedWidth);
644     //
645     /*printf("\t  AliZDCDigitizer::Pedestal -> det %d quad %d res %d - Ped[%d] = %d\n",
646         Det, Quad, Res, index,(Int_t) pedValue); // Chiara debugging!
647     */
648   }
649   // To create calibration object
650   else{
651     if(Res == 0) pedValue = gRandom->Gaus((35.+10.*gRandom->Rndm()),(0.5+0.2*gRandom->Rndm())); //High gain
652     else  pedValue = gRandom->Gaus((250.+100.*gRandom->Rndm()),(3.5+2.*gRandom->Rndm())); //Low gain
653   }
654
655   return (Int_t) pedValue;
656 }
657
658 //_____________________________________________________________________________
659 AliCDBStorage* AliZDCDigitizer::SetStorage(const char *uri) 
660 {
661
662   Bool_t deleteManager = kFALSE;
663   
664   AliCDBManager *manager = AliCDBManager::Instance();
665   AliCDBStorage *defstorage = manager->GetDefaultStorage();
666   
667   if(!defstorage || !(defstorage->Contains("ZDC"))){ 
668      AliWarning("No default storage set or default storage doesn't contain ZDC!");
669      manager->SetDefaultStorage(uri);
670      deleteManager = kTRUE;
671   }
672  
673   AliCDBStorage *storage = manager->GetDefaultStorage();
674
675   if(deleteManager){
676     AliCDBManager::Instance()->UnsetDefaultStorage();
677     defstorage = 0;   // the storage is killed by AliCDBManager::Instance()->Destroy()
678   }
679
680   return storage; 
681 }
682
683 //_____________________________________________________________________________
684 AliZDCPedestals* AliZDCDigitizer::GetPedData() const
685 {
686
687   // Getting pedestal calibration object for ZDC set
688
689   AliCDBEntry  *entry = AliCDBManager::Instance()->Get("ZDC/Calib/Pedestals");
690   if(!entry) AliFatal("No calibration data loaded!");  
691
692   AliZDCPedestals *calibdata = dynamic_cast<AliZDCPedestals*>  (entry->GetObject());
693   if(!calibdata)  AliFatal("Wrong calibration object in calibration  file!");
694
695   return calibdata;
696 }
697