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