1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
18 ///////////////////////////////////////////////////////////////////////////////
20 // ZDC digitizer class //
22 ///////////////////////////////////////////////////////////////////////////////
32 // --- AliRoot header files
35 #include "AliHeader.h"
36 #include "AliGenHijingEventHeader.h"
37 #include "AliRunDigitizer.h"
38 #include "AliRunLoader.h"
39 #include "AliCDBManager.h"
40 #include "AliCDBEntry.h"
41 #include "AliZDCSDigit.h"
42 #include "AliZDCDigit.h"
43 #include "AliZDCFragment.h"
44 #include "AliZDCDigitizer.h"
47 class AliZDCPedestals;
51 ClassImp(AliZDCDigitizer)
54 //____________________________________________________________________________
55 AliZDCDigitizer::AliZDCDigitizer() :
61 // Default constructor
65 //____________________________________________________________________________
66 AliZDCDigitizer::AliZDCDigitizer(AliRunDigitizer* manager):
67 AliDigitizer(manager),
68 fIsCalibration(0), //By default the simulation doesn't create calib. data
69 fPedData(GetPedData()),
70 fCalibData(GetCalibData()),
71 fRecParam(GetRecParam())
73 // Get calibration data
74 if(fIsCalibration!=0) printf("\n\t AliZDCDigitizer -> Creating calibration data (pedestals)\n");
78 //____________________________________________________________________________
79 AliZDCDigitizer::~AliZDCDigitizer()
86 //____________________________________________________________________________
87 AliZDCDigitizer::AliZDCDigitizer(const AliZDCDigitizer &digitizer):
89 fIsCalibration(digitizer.fIsCalibration),
90 fPedData(digitizer.fPedData),
91 fCalibData(digitizer.fCalibData),
92 fRecParam(digitizer.fRecParam)
96 for(Int_t i=0; i<6; i++){
97 for(Int_t j=0; j<5; j++){
98 fPMGain[i][j] = digitizer.fPMGain[i][j];
101 for(Int_t i=0; i<2; i++) fADCRes[i] = digitizer.fADCRes[i];
105 //____________________________________________________________________________
106 Bool_t AliZDCDigitizer::Init()
108 // Initialize the digitizer
109 // NB -> PM gain vs. HV & ADC resolutions will move to DCDB ***************
110 for(Int_t j = 0; j < 5; j++){
111 fPMGain[0][j] = 50000.;
112 fPMGain[1][j] = 100000.;
113 fPMGain[2][j] = 100000.;
114 fPMGain[3][j] = 50000.;
115 fPMGain[4][j] = 100000.;
116 fPMGain[5][j] = 100000.;
119 fADCRes[0] = 0.0000008; // ADC Resolution high gain: 200 fC/adcCh
120 fADCRes[1] = 0.0000064; // ADC Resolution low gain: 25 fC/adcCh
125 //____________________________________________________________________________
126 void AliZDCDigitizer::Exec(Option_t* /*option*/)
128 // Execute digitization
130 // ------------------------------------------------------------
131 // !!! 2nd ZDC set added
132 // *** 1st 3 arrays are digits from REAL (simulated) hits
133 // *** last 2 are copied from simulated digits
134 // --- pm[0][...] = light in ZN right [C, Q1, Q2, Q3, Q4]
135 // --- pm[1][...] = light in ZP right [C, Q1, Q2, Q3, Q4]
136 // --- pm[2][...] = light in ZEM [x, 1, 2, x, x]
137 // --- pm[3][...] = light in ZN left [C, Q1, Q2, Q3, Q4] ->NEW!
138 // --- pm[4][...] = light in ZP left [C, Q1, Q2, Q3, Q4] ->NEW!
139 // ------------------------------------------------------------
141 for(Int_t iSector1=0; iSector1<5; iSector1++)
142 for(Int_t iSector2=0; iSector2<5; iSector2++){
143 pm[iSector1][iSector2] = 0;
146 // ------------------------------------------------------------
147 // ### Out of time ADC added (22 channels)
148 // --- same codification as for signal PTMs (see above)
149 // ------------------------------------------------------------
151 for(Int_t iSector1=0; iSector1<5; iSector1++)
152 for(Int_t iSector2=0; iSector2<5; iSector2++){
153 pmoot[iSector1][iSector2] = 0;
156 // impact parameter and number of spectators
161 // loop over input streams
162 for(Int_t iInput = 0; iInput<fManager->GetNinputs(); iInput++){
164 // get run loader and ZDC loader
165 AliRunLoader* runLoader =
166 AliRunLoader::GetRunLoader(fManager->GetInputFolderName(iInput));
167 AliLoader* loader = runLoader->GetLoader("ZDCLoader");
168 if(!loader) continue;
171 loader->LoadSDigits();
172 TTree* treeS = loader->TreeS();
175 AliZDCSDigit* psdigit = &sdigit;
176 treeS->SetBranchAddress("ZDC", &psdigit);
179 for(Int_t iSDigit=0; iSDigit<treeS->GetEntries(); iSDigit++){
180 treeS->GetEntry(iSDigit);
182 if(!psdigit) continue;
183 if((sdigit.GetSector(1) < 0) || (sdigit.GetSector(1) > 4)){
184 AliError(Form("\nsector[0] = %d, sector[1] = %d\n",
185 sdigit.GetSector(0), sdigit.GetSector(1)));
189 pm[(sdigit.GetSector(0))-1][sdigit.GetSector(1)] += sdigit.GetLightPM();
190 /*printf("\n\t Detector %d, Tower %d -> pm[%d][%d] = %.0f \n",
191 sdigit.GetSector(0), sdigit.GetSector(1),sdigit.GetSector(0)-1,
192 sdigit.GetSector(1), pm[sdigit.GetSector(0)-1][sdigit.GetSector(1)]); // Chiara debugging!
196 loader->UnloadSDigits();
198 // get the impact parameter and the number of spectators in case of hijing
199 if(!runLoader->GetAliRun()) runLoader->LoadgAlice();
200 AliHeader* header = runLoader->GetAliRun()->GetHeader();
201 if(!header) continue;
202 AliGenEventHeader* genHeader = header->GenEventHeader();
203 if(!genHeader) continue;
204 if(!genHeader->InheritsFrom(AliGenHijingEventHeader::Class())) continue;
205 impPar = ((AliGenHijingEventHeader*) genHeader)->ImpactParameter();
207 specN = ((AliGenHijingEventHeader*) genHeader)->ProjSpectatorsn();
208 specP = ((AliGenHijingEventHeader*) genHeader)->ProjSpectatorsp();
209 AliDebug(2, Form("\n AliZDCDigitizer -> b = %f fm, Nspecn = %d, Nspecp = %d\n",
210 impPar, specN, specP));
211 printf("\n\t AliZDCDigitizer -> b = %f fm, # generated spectator n = %d,"
212 " # generated spectator p = %d\n", impPar, specN, specP);
217 Int_t freeSpecN, freeSpecP;
218 Fragmentation(impPar, specN, specP, freeSpecN, freeSpecP);
219 printf("\n\t AliZDCDigitizer -> Adding signal for %d free spectator n\n",freeSpecN);
220 SpectatorSignal(1, freeSpecN, pm);
221 printf("\t AliZDCDigitizer -> Adding signal for %d free spectator p\n\n",freeSpecP);
222 SpectatorSignal(2, freeSpecP, pm);
226 // get the output run loader and loader
227 AliRunLoader* runLoader =
228 AliRunLoader::GetRunLoader(fManager->GetOutputFolderName());
229 AliLoader* loader = runLoader->GetLoader("ZDCLoader");
231 AliError("no ZDC loader found");
235 // create the output tree
236 const char* mode = "update";
237 if(runLoader->GetEventNumber() == 0) mode = "recreate";
238 loader->LoadDigits(mode);
239 loader->MakeTree("D");
240 TTree* treeD = loader->TreeD();
242 AliZDCDigit* pdigit = &digit;
243 const Int_t kBufferSize = 4000;
244 treeD->Branch("ZDC", "AliZDCDigit", &pdigit, kBufferSize);
248 Int_t digi[2], digioot[2];
249 for(sector[0]=1; sector[0]<6; sector[0]++){
250 for(sector[1]=0; sector[1]<5; sector[1]++){
251 if((sector[0]==3) && ((sector[1]<1) || (sector[1]>2))) continue;
252 for(Int_t res=0; res<2; res++){
253 digi[res] = Phe2ADCch(sector[0], sector[1], pm[sector[0]-1][sector[1]], res)
254 + Pedestal(sector[0], sector[1], res);
256 /*printf("\t DIGIT added -> det %d quad %d - digi[0,1] = [%d, %d]\n",
257 sector[0], sector[1], digi[0], digi[1]); // Chiara debugging!
260 new(pdigit) AliZDCDigit(sector, digi);
263 } // Loop over detector
264 // Adding in-time digits for 2 reference PTM signals (after signal ch.)
265 // (for the moment the ref. signal is completely invented assuming a PMgain of 5*10^4!)
269 // Reference signal are set to 100 (high gain chain) and 800 (low gain chain)
270 if(fIsCalibration==0) {sigRef[0]=100; sigRef[1]=800;}
271 else {sigRef[0]=0; sigRef[1]=0;} // calibration -> simulation of pedestal values
273 for(Int_t iref=0; iref<2; iref++){
274 sectorRef[0] = 3*iref+1;
275 for(Int_t res=0; res<2; res++){
276 sigRef[res] += Pedestal(sectorRef[0], sectorRef[1], res);
278 /*printf("\t RefDigit added -> det = %d, quad = %d - digi[0,1] = [%d, %d]\n",
279 sectorRef[0], sectorRef[1], sigRef[0], sigRef[1]); // Chiara debugging!
281 new(pdigit) AliZDCDigit(sectorRef, sigRef);
285 // --- Adding digits for out-of-time channels after signal digits
286 for(sector[0]=1; sector[0]<6; sector[0]++){
287 for(sector[1]=0; sector[1]<5; sector[1]++){
288 if((sector[0]==3) && ((sector[1]<1) || (sector[1]>2))) continue;
289 for(Int_t res=0; res<2; res++){
290 digioot[res] = Pedestal(sector[0], sector[1], res); // out-of-time ADCs
292 /*printf("\t DIGIToot added -> det = %d, quad = %d - digi[0,1] = [%d, %d]\n",
293 sector[0], sector[1], digioot[0], digioot[1]); // Chiara debugging!
296 new(pdigit) AliZDCDigit(sector, digioot);
300 // Adding out-of-time digits for 2 reference PTM signals (after out-of-time ch.)
302 for(Int_t iref=0; iref<2; iref++){
303 sectorRef[0] = 3*iref+1;
304 for(Int_t res=0; res<2; res++){
305 sigRefoot[res] = Pedestal(sectorRef[0], sectorRef[1], res);
307 /*printf("\t RefDigitoot added -> det = %d, quad = %d - digi[0,1] = [%d, %d]\n",
308 sectorRef[0], sectorRef[1], sigRefoot[0], sigRefoot[1]); // Chiara debugging!
310 new(pdigit) AliZDCDigit(sectorRef, sigRefoot);
313 //printf("\t AliZDCDigitizer -> TreeD has %d entries\n",(Int_t) treeD->GetEntries());
315 // write the output tree
316 loader->WriteDigits("OVERWRITE");
317 loader->UnloadDigits();
321 //_____________________________________________________________________________
322 void AliZDCDigitizer::Fragmentation(Float_t impPar, Int_t specN, Int_t specP,
323 Int_t &freeSpecN, Int_t &freeSpecP) const
325 // simulate fragmentation of spectators
327 Int_t zz[100], nn[100];
328 AliZDCFragment frag(impPar);
329 for(Int_t j=0; j<=99; j++){
334 // Fragments generation
336 frag.GenerateIMF(zz, nAlpha);
341 frag.AttachNeutrons(zz, nn, ztot, ntot);
342 freeSpecN = specN-ntot-2*nAlpha;
343 freeSpecP = specP-ztot-2*nAlpha;
344 // Removing deuterons
345 Int_t ndeu = (Int_t) (freeSpecN*frag.DeuteronNumber());
348 if(freeSpecN<0) freeSpecN=0;
349 if(freeSpecP<0) freeSpecP=0;
350 AliDebug(2, Form("FreeSpn = %d, FreeSpp = %d", freeSpecN, freeSpecP));
353 //_____________________________________________________________________________
354 void AliZDCDigitizer::SpectatorSignal(Int_t SpecType, Int_t numEvents,
355 Float_t pm[3][5]) const
357 // add signal of the spectators
360 if(SpecType == 1) { // --- Signal for spectator neutrons
361 file = TFile::Open("$ALICE_ROOT/ZDC/ZNsignalntu.root");
362 } else if(SpecType == 2) { // --- Signal for spectator protons
363 file = TFile::Open("$ALICE_ROOT/ZDC/ZPsignalntu.root");
365 if(!file || !file->IsOpen()) {
366 AliError("Opening of file failed");
370 TNtuple* zdcSignal = (TNtuple*) file->Get("ZDCSignal");
371 Int_t nentries = (Int_t) zdcSignal->GetEntries();
373 Float_t *entry, hitsSpec[7];
374 Int_t pl, i, j, k, iev=0, rnd[125], volume[2];
375 for(pl=0;pl<125;pl++) rnd[pl] = 0;
376 if(numEvents > 125) {
377 AliWarning(Form("numEvents (%d) is larger than 125", numEvents));
380 for(pl=0;pl<numEvents;pl++){
381 rnd[pl] = (Int_t) (9999*gRandom->Rndm());
382 if(rnd[pl] >= 9999) rnd[pl] = 9998;
383 //printf(" rnd[%d] = %d\n",pl,rnd[pl]);
385 // Sorting vector in ascending order with C function QSORT
386 qsort((void*)rnd,numEvents,sizeof(Int_t),comp);
388 for(i=0; i<nentries; i++){
389 zdcSignal->GetEvent(i);
390 entry = zdcSignal->GetArgs();
391 if(entry[0] == rnd[iev]){
392 for(k=0; k<2; k++) volume[k] = (Int_t) entry[k+1];
393 for(j=0; j<7; j++) hitsSpec[j] = entry[j+3];
395 Float_t lightQ = hitsSpec[4];
396 Float_t lightC = hitsSpec[5];
397 AliDebug(3, Form("SpectatorSignal -> vol = (%d, %d), lightQ = %.0f, lightC = %.0f",
398 volume[0], volume[1], lightQ, lightC));
399 //printf("\n Volume = (%d, %d), lightQ = %.0f, lightC = %.0f",
400 // volume[0], volume[1], lightQ, lightC);
401 if(volume[0] != 3) { // ZN or ZP
402 pm[volume[0]-1][0] += lightC;
403 pm[volume[0]-1][volume[1]] += lightQ;
404 //printf("\n pm[%d][0] = %.0f, pm[%d][%d] = %.0f\n",(volume[0]-1),pm[volume[0]-1][0],
405 // (volume[0]-1),volume[1],pm[volume[0]-1][volume[1]]);
408 if(volume[1] == 1) pm[2][1] += lightC; // ZEM 1
409 else pm[2][2] += lightQ; // ZEM 2
410 //printf("\n pm[2][1] = %.0f, pm[2][2] = %.0f\n",pm[2][1],pm[2][2]);
413 else if(entry[0] > rnd[iev]){
418 }while(iev<numEvents);
425 //_____________________________________________________________________________
426 Int_t AliZDCDigitizer::Phe2ADCch(Int_t Det, Int_t Quad, Float_t Light,
429 // Evaluation of the ADC channel corresponding to the light yield Light
430 Int_t vADCch = (Int_t) (Light * fPMGain[Det-1][Quad] * fADCRes[Res]);
431 //printf("\t Phe2ADCch -> det %d quad %d - phe %.0f ADC %d\n", Det,Quad,Light,vADCch);
436 //_____________________________________________________________________________
437 Int_t AliZDCDigitizer::Pedestal(Int_t Det, Int_t Quad, Int_t Res) const
439 // Returns a pedestal for detector det, PM quad, channel with res.
443 if(fIsCalibration == 0){
444 Int_t index=0, kNch=24;
446 if(Det==1) index = Quad+kNch*Res; // ZNC
447 else if(Det==2) index = (Quad+5)+kNch*Res; // ZPC
448 else if(Det==3) index = (Quad+9)+kNch*Res; // ZEM
449 else if(Det==4) index = (Quad+12)+kNch*Res; // ZNA
450 else if(Det==5) index = (Quad+17)+kNch*Res; // ZPA
452 else index = (Det-1)/3+22+kNch*Res; // Reference PMs
454 Float_t meanPed = fPedData->GetMeanPed(index);
455 Float_t pedWidth = fPedData->GetMeanPedWidth(index);
456 pedValue = gRandom->Gaus(meanPed,pedWidth);
458 /*printf("\t AliZDCDigitizer::Pedestal -> det %d quad %d res %d - Ped[%d] = %d\n",
459 Det, Quad, Res, index,(Int_t) pedValue); // Chiara debugging!
462 // To create calibration object
464 if(Res == 0) pedValue = gRandom->Gaus((35.+10.*gRandom->Rndm()),(0.5+0.2*gRandom->Rndm())); //High gain
465 else pedValue = gRandom->Gaus((250.+100.*gRandom->Rndm()),(3.5+2.*gRandom->Rndm())); //Low gain
468 return (Int_t) pedValue;
471 //_____________________________________________________________________________
472 AliCDBStorage* AliZDCDigitizer::SetStorage(const char *uri)
475 Bool_t deleteManager = kFALSE;
477 AliCDBManager *manager = AliCDBManager::Instance();
478 AliCDBStorage *defstorage = manager->GetDefaultStorage();
480 if(!defstorage || !(defstorage->Contains("ZDC"))){
481 AliWarning("No default storage set or default storage doesn't contain ZDC!");
482 manager->SetDefaultStorage(uri);
483 deleteManager = kTRUE;
486 AliCDBStorage *storage = manager->GetDefaultStorage();
489 AliCDBManager::Instance()->UnsetDefaultStorage();
490 defstorage = 0; // the storage is killed by AliCDBManager::Instance()->Destroy()
496 //_____________________________________________________________________________
497 AliZDCPedestals* AliZDCDigitizer::GetPedData() const
500 // Getting pedestal calibration object for ZDC set
502 AliCDBEntry *entry = AliCDBManager::Instance()->Get("ZDC/Calib/Pedestals");
503 if(!entry) AliFatal("No calibration data loaded!");
505 AliZDCPedestals *calibdata = dynamic_cast<AliZDCPedestals*> (entry->GetObject());
506 if(!calibdata) AliFatal("Wrong calibration object in calibration file!");
511 //_____________________________________________________________________________
512 AliZDCCalib* AliZDCDigitizer::GetCalibData() const
515 // Getting calibration object for ZDC set
517 AliCDBEntry *entry = AliCDBManager::Instance()->Get("ZDC/Calib/Calib");
518 if(!entry) AliFatal("No calibration data loaded!");
520 AliZDCCalib *calibdata = dynamic_cast<AliZDCCalib*> (entry->GetObject());
521 if(!calibdata) AliFatal("Wrong calibration object in calibration file!");
526 //_____________________________________________________________________________
527 AliZDCRecParam* AliZDCDigitizer::GetRecParam() const
530 // Getting calibration object for ZDC set
532 AliCDBEntry *entry = AliCDBManager::Instance()->Get("ZDC/Calib/RecParam");
533 if(!entry) AliFatal("No calibration data loaded!");
535 AliZDCRecParam *calibdata = dynamic_cast<AliZDCRecParam*> (entry->GetObject());
536 if(!calibdata) AliFatal("Wrong calibration object in calibration file!");