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