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