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()
57 // Default constructor
61 //____________________________________________________________________________
62 AliZDCDigitizer::AliZDCDigitizer(AliRunDigitizer* manager):
65 fIsCalibration=0; //By default the simulation doesn't create calib. data
66 // fIsCalibration=1; //To create pedestal calib. data
67 // Get calibration data
68 fPedData = GetPedData();
69 fCalibData = GetCalibData();
70 fRecParam = GetRecParam();
71 if(fIsCalibration!=0) printf("\n\t AliZDCDigitizer -> Creating calibration data (pedestals)\n");
75 //____________________________________________________________________________
76 AliZDCDigitizer::~AliZDCDigitizer()
83 //____________________________________________________________________________
84 AliZDCDigitizer::AliZDCDigitizer(const AliZDCDigitizer &digitizer):
89 for(Int_t i=0; i<6; i++){
90 for(Int_t j=0; j<5; j++){
91 fPMGain[i][j] = digitizer.fPMGain[i][j];
94 for(Int_t i=0; i<2; i++) fADCRes[i] = digitizer.fADCRes[i];
95 fIsCalibration = digitizer.fIsCalibration;
96 fPedData = digitizer.fPedData;
97 fCalibData = digitizer.fCalibData;
98 fRecParam = digitizer.fRecParam;
102 //____________________________________________________________________________
103 Bool_t AliZDCDigitizer::Init()
105 // Initialize the digitizer
106 // NB -> PM gain vs. HV & ADC resolutions will move to DCDB ***************
107 for(Int_t j = 0; j < 5; j++){
108 fPMGain[0][j] = 50000.;
109 fPMGain[1][j] = 100000.;
110 fPMGain[2][j] = 100000.;
111 fPMGain[3][j] = 50000.;
112 fPMGain[4][j] = 100000.;
113 fPMGain[5][j] = 100000.;
116 fADCRes[0] = 0.0000008; // ADC Resolution high gain: 200 fC/adcCh
117 fADCRes[1] = 0.0000064; // ADC Resolution low gain: 25 fC/adcCh
122 //____________________________________________________________________________
123 void AliZDCDigitizer::Exec(Option_t* /*option*/)
125 // Execute digitization
127 // ------------------------------------------------------------
128 // !!! 2nd ZDC set added
129 // *** 1st 3 arrays are digits from REAL (simulated) hits
130 // *** last 2 are copied from simulated digits
131 // --- pm[0][...] = light in ZN right [C, Q1, Q2, Q3, Q4]
132 // --- pm[1][...] = light in ZP right [C, Q1, Q2, Q3, Q4]
133 // --- pm[2][...] = light in ZEM [x, 1, 2, x, x]
134 // --- pm[3][...] = light in ZN left [C, Q1, Q2, Q3, Q4] ->NEW!
135 // --- pm[4][...] = light in ZP left [C, Q1, Q2, Q3, Q4] ->NEW!
136 // ------------------------------------------------------------
138 for(Int_t iSector1=0; iSector1<5; iSector1++)
139 for(Int_t iSector2=0; iSector2<5; iSector2++){
140 pm[iSector1][iSector2] = 0;
143 // ------------------------------------------------------------
144 // ### Out of time ADC added (22 channels)
145 // --- same codification as for signal PTMs (see above)
146 // ------------------------------------------------------------
148 for(Int_t iSector1=0; iSector1<5; iSector1++)
149 for(Int_t iSector2=0; iSector2<5; iSector2++){
150 pmoot[iSector1][iSector2] = 0;
153 // impact parameter and number of spectators
158 // loop over input streams
159 for(Int_t iInput = 0; iInput<fManager->GetNinputs(); iInput++){
161 // get run loader and ZDC loader
162 AliRunLoader* runLoader =
163 AliRunLoader::GetRunLoader(fManager->GetInputFolderName(iInput));
164 AliLoader* loader = runLoader->GetLoader("ZDCLoader");
165 if(!loader) continue;
168 loader->LoadSDigits();
169 TTree* treeS = loader->TreeS();
172 AliZDCSDigit* psdigit = &sdigit;
173 treeS->SetBranchAddress("ZDC", &psdigit);
176 for(Int_t iSDigit=0; iSDigit<treeS->GetEntries(); iSDigit++){
177 treeS->GetEntry(iSDigit);
179 if(!psdigit) continue;
180 if((sdigit.GetSector(1) < 0) || (sdigit.GetSector(1) > 4)){
181 AliError(Form("\nsector[0] = %d, sector[1] = %d\n",
182 sdigit.GetSector(0), sdigit.GetSector(1)));
186 pm[(sdigit.GetSector(0))-1][sdigit.GetSector(1)] += sdigit.GetLightPM();
187 /*printf("\n\t Detector %d, Tower %d -> pm[%d][%d] = %.0f \n",
188 sdigit.GetSector(0), sdigit.GetSector(1),sdigit.GetSector(0)-1,
189 sdigit.GetSector(1), pm[sdigit.GetSector(0)-1][sdigit.GetSector(1)]); // Chiara debugging!
193 loader->UnloadSDigits();
195 // get the impact parameter and the number of spectators in case of hijing
196 if(!runLoader->GetAliRun()) runLoader->LoadgAlice();
197 AliHeader* header = runLoader->GetAliRun()->GetHeader();
198 if(!header) continue;
199 AliGenEventHeader* genHeader = header->GenEventHeader();
200 if(!genHeader) continue;
201 if(!genHeader->InheritsFrom(AliGenHijingEventHeader::Class())) continue;
202 impPar = ((AliGenHijingEventHeader*) genHeader)->ImpactParameter();
204 specN = ((AliGenHijingEventHeader*) genHeader)->ProjSpectatorsn();
205 specP = ((AliGenHijingEventHeader*) genHeader)->ProjSpectatorsp();
206 AliDebug(2, Form("\n AliZDCDigitizer -> b = %f fm, Nspecn = %d, Nspecp = %d\n",
207 impPar, specN, specP));
208 printf("\n\t AliZDCDigitizer -> b = %f fm, # generated spectator n = %d,"
209 " # generated spectator p = %d\n", impPar, specN, specP);
214 Int_t freeSpecN, freeSpecP;
215 Fragmentation(impPar, specN, specP, freeSpecN, freeSpecP);
216 printf("\n\t AliZDCDigitizer -> Adding signal for %d free spectator n\n",freeSpecN);
217 SpectatorSignal(1, freeSpecN, pm);
218 printf("\t AliZDCDigitizer -> Adding signal for %d free spectator p\n\n",freeSpecP);
219 SpectatorSignal(2, freeSpecP, pm);
223 // get the output run loader and loader
224 AliRunLoader* runLoader =
225 AliRunLoader::GetRunLoader(fManager->GetOutputFolderName());
226 AliLoader* loader = runLoader->GetLoader("ZDCLoader");
228 AliError("no ZDC loader found");
232 // create the output tree
233 const char* mode = "update";
234 if(runLoader->GetEventNumber() == 0) mode = "recreate";
235 loader->LoadDigits(mode);
236 loader->MakeTree("D");
237 TTree* treeD = loader->TreeD();
239 AliZDCDigit* pdigit = &digit;
240 const Int_t kBufferSize = 4000;
241 treeD->Branch("ZDC", "AliZDCDigit", &pdigit, kBufferSize);
244 Int_t sector[2], sectorL[2];
245 Int_t digi[2], digiL[2], digioot[2];
246 for(sector[0]=1; sector[0]<=3; sector[0]++){
247 for(sector[1]=0; sector[1]<5; sector[1]++){
248 if((sector[0]==3) && ((sector[1]<1) || (sector[1]>2))) continue;
249 for(Int_t res=0; res<2; res++){
250 digi[res] = Phe2ADCch(sector[0], sector[1], pm[sector[0]-1][sector[1]], res)
251 + Pedestal(sector[0], sector[1], res);
253 /*printf("\t DIGIT added -> det %d quad %d - digi[0,1] = [%d, %d]\n",
254 sector[0], sector[1], digi[0], digi[1]); // Chiara debugging!
257 new(pdigit) AliZDCDigit(sector, digi);
260 // --- Adding digits for 2nd ZDC set (left side w.r.t. IP) ---
261 // --- they are copied from right ZDC digits
262 if(sector[0]==1 || sector[0]==2){
263 sectorL[0] = sector[0]+3;
264 sectorL[1] = sector[1];
265 for(Int_t res=0; res<2; res++){
266 digiL[res] = Phe2ADCch(sectorL[0], sectorL[1], pm[sector[0]-1][sector[1]], res)
267 + Pedestal(sectorL[0], sectorL[1], res);
269 /*printf("\t DIGIT added -> det %d quad %d - digi[0,1] = [%d, %d]\n",
270 sectorL[0], sectorL[1], digiL[0], digiL[1]); // Chiara debugging!
273 new(pdigit) AliZDCDigit(sectorL, digiL);
277 } // Loop over detector
278 // Adding in-time digits for 2 reference PTM signals (after signal ch.)
279 // (for the moment the ref. signal is completely invented assuming a PMgain of 5*10^4!)
283 if(fIsCalibration==0) {sigRef[0]=100; sigRef[1]=800;}
284 else {sigRef[0]=0; sigRef[1]=0;}
286 for(Int_t iref=0; iref<2; iref++){
287 sectorRef[0] = 3*iref+1;
288 for(Int_t res=0; res<2; res++){
289 sigRef[res] += Pedestal(sectorRef[0], sectorRef[1], res);
291 /*printf("\t RefDigit added -> det = %d, quad = %d - digi[0,1] = [%d, %d]\n",
292 sectorRef[0], sectorRef[1], sigRef[0], sigRef[1]); // Chiara debugging!
294 new(pdigit) AliZDCDigit(sectorRef, sigRef);
298 // --- Adding digits for out-of-time channels after signal digits
299 for(sector[0]=1; sector[0]<=3; sector[0]++){
300 for(sector[1]=0; sector[1]<5; sector[1]++){
301 if((sector[0]==3) && ((sector[1]<1) || (sector[1]>2))) continue;
302 for(Int_t res=0; res<2; res++){
303 digioot[res] = Pedestal(sector[0], sector[1], res); // out-of-time ADCs
305 /*printf("\t DIGIToot added -> det = %d, quad = %d - digi[0,1] = [%d, %d]\n",
306 sector[0], sector[1], digioot[0], digioot[1]); // Chiara debugging!
309 new(pdigit) AliZDCDigit(sector, digioot);
312 if(sector[0]==1 || sector[0]==2){
313 sectorL[0] = sector[0]+3;
314 sectorL[1] = sector[1];
315 for(Int_t res=0; res<2; res++){
316 digioot[res] = Pedestal(sectorL[0], sectorL[1], res); // out-of-time ADCs
318 /*printf("\t DIGIToot added -> det = %d, quad = %d - digi[0,1] = [%d, %d]\n",
319 sectorL[0], sectorL[1], digioot[0], digioot[1]); // Chiara debugging!
322 new(pdigit) AliZDCDigit(sectorL, digioot);
328 // Adding out-of-time digits for 2 reference PTM signals (after out-of-time ch.)
330 for(Int_t iref=0; iref<2; iref++){
331 sectorRef[0] = 3*iref+1;
332 for(Int_t res=0; res<2; res++){
333 sigRefoot[res] = Pedestal(sectorRef[0], sectorRef[1], res);
335 /*printf("\t RefDigitoot added -> det = %d, quad = %d - digi[0,1] = [%d, %d]\n",
336 sectorRef[0], sectorRef[1], sigRefoot[0], sigRefoot[1]); // Chiara debugging!
338 new(pdigit) AliZDCDigit(sectorRef, sigRefoot);
342 //printf("\t AliZDCDigitizer -> TreeD has %d entries\n",(Int_t) treeD->GetEntries());
344 // write the output tree
345 loader->WriteDigits("OVERWRITE");
346 loader->UnloadDigits();
350 //_____________________________________________________________________________
351 void AliZDCDigitizer::Fragmentation(Float_t impPar, Int_t specN, Int_t specP,
352 Int_t &freeSpecN, Int_t &freeSpecP) const
354 // simulate fragmentation of spectators
356 Int_t zz[100], nn[100];
357 AliZDCFragment frag(impPar);
358 for(Int_t j=0; j<=99; j++){
363 // Fragments generation
365 frag.GenerateIMF(zz, nAlpha);
370 frag.AttachNeutrons(zz, nn, ztot, ntot);
371 freeSpecN = specN-ntot-2*nAlpha;
372 freeSpecP = specP-ztot-2*nAlpha;
373 // Removing deuterons
374 Int_t ndeu = (Int_t) (frag.DeuteronNumber());
377 if(freeSpecN<0) freeSpecN=0;
378 if(freeSpecP<0) freeSpecP=0;
379 AliDebug(2, Form("FreeSpn = %d, FreeSpp = %d", freeSpecN, freeSpecP));
382 //_____________________________________________________________________________
383 void AliZDCDigitizer::SpectatorSignal(Int_t SpecType, Int_t numEvents,
384 Float_t pm[3][5]) const
386 // add signal of the spectators
389 if(SpecType == 1) { // --- Signal for spectator neutrons
390 file = TFile::Open("$ALICE_ROOT/ZDC/ZNsignalntu.root");
391 } else if(SpecType == 2) { // --- Signal for spectator protons
392 file = TFile::Open("$ALICE_ROOT/ZDC/ZPsignalntu.root");
394 if(!file || !file->IsOpen()) {
395 AliError("Opening of file failed");
399 TNtuple* zdcSignal = (TNtuple*) file->Get("ZDCSignal");
400 Int_t nentries = (Int_t) zdcSignal->GetEntries();
402 Float_t *entry, hitsSpec[7];
403 Int_t pl, i, j, k, iev=0, rnd[125], volume[2];
404 for(pl=0;pl<125;pl++) rnd[pl] = 0;
405 if(numEvents > 125) {
406 AliWarning(Form("numEvents (%d) is larger than 125", numEvents));
409 for(pl=0;pl<numEvents;pl++){
410 rnd[pl] = (Int_t) (9999*gRandom->Rndm());
411 if(rnd[pl] >= 9999) rnd[pl] = 9998;
412 //printf(" rnd[%d] = %d\n",pl,rnd[pl]);
414 // Sorting vector in ascending order with C function QSORT
415 qsort((void*)rnd,numEvents,sizeof(Int_t),comp);
417 for(i=0; i<nentries; i++){
418 zdcSignal->GetEvent(i);
419 entry = zdcSignal->GetArgs();
420 if(entry[0] == rnd[iev]){
421 for(k=0; k<2; k++) volume[k] = (Int_t) entry[k+1];
422 for(j=0; j<7; j++) hitsSpec[j] = entry[j+3];
424 Float_t lightQ = hitsSpec[4];
425 Float_t lightC = hitsSpec[5];
426 AliDebug(3, Form("SpectatorSignal -> vol = (%d, %d), lightQ = %.0f, lightC = %.0f",
427 volume[0], volume[1], lightQ, lightC));
428 //printf("\n Volume = (%d, %d), lightQ = %.0f, lightC = %.0f",
429 // volume[0], volume[1], lightQ, lightC);
430 if(volume[0] != 3) { // ZN or ZP
431 pm[volume[0]-1][0] += lightC;
432 pm[volume[0]-1][volume[1]] += lightQ;
433 //printf("\n pm[%d][0] = %.0f, pm[%d][%d] = %.0f\n",(volume[0]-1),pm[volume[0]-1][0],
434 // (volume[0]-1),volume[1],pm[volume[0]-1][volume[1]]);
437 if(volume[1] == 1) pm[2][1] += lightC; // ZEM 1
438 else pm[2][2] += lightQ; // ZEM 2
439 //printf("\n pm[2][1] = %.0f, pm[2][2] = %.0f\n",pm[2][1],pm[2][2]);
442 else if(entry[0] > rnd[iev]){
447 }while(iev<numEvents);
454 //_____________________________________________________________________________
455 Int_t AliZDCDigitizer::Phe2ADCch(Int_t Det, Int_t Quad, Float_t Light,
458 // Evaluation of the ADC channel corresponding to the light yield Light
459 Int_t vADCch = (Int_t) (Light * fPMGain[Det-1][Quad] * fADCRes[Res]);
460 //printf("\t Phe2ADCch -> det %d quad %d - phe %.0f ADC %d\n", Det,Quad,Light,vADCch);
465 //_____________________________________________________________________________
466 Int_t AliZDCDigitizer::Pedestal(Int_t Det, Int_t Quad, Int_t Res) const
468 // Returns a pedestal for detector det, PM quad, channel with res.
473 if(fIsCalibration == 0){
474 Float_t meanPed, pedWidth;
477 if(Det==1) index = Quad+24*Res; // ZN1
478 else if(Det==2) index = (Quad+5)+24*Res; // ZP1
479 else if(Det==3) index = (Quad+9)+24*Res; // ZEM
480 else if(Det==4) index = (Quad+12)+24*Res; // ZN2
481 else if(Det==5) index = (Quad+17)+24*Res; // ZP2
483 else index = (Det-1)/3+22+24*Res; // Reference PMs
485 meanPed = fPedData->GetMeanPed(index);
486 pedWidth = fPedData->GetMeanPedWidth(index);
487 pedValue = gRandom->Gaus(meanPed,pedWidth);
489 /*printf("\t AliZDCDigitizer::Pedestal -> det %d quad %d res %d - Ped[%d] = %d\n",
490 Det, Quad, Res, index,(Int_t) pedValue); // Chiara debugging!
493 // To create calibration object
495 if(Res == 0) pedValue = gRandom->Gaus((35.+10.*gRandom->Rndm()),(0.5+0.2*gRandom->Rndm())); //High gain
496 else pedValue = gRandom->Gaus((250.+100.*gRandom->Rndm()),(3.5+2.*gRandom->Rndm())); //Low gain
499 return (Int_t) pedValue;
502 //_____________________________________________________________________________
503 AliCDBStorage* AliZDCDigitizer::SetStorage(const char *uri)
506 Bool_t deleteManager = kFALSE;
508 AliCDBManager *manager = AliCDBManager::Instance();
509 AliCDBStorage *defstorage = manager->GetDefaultStorage();
511 if(!defstorage || !(defstorage->Contains("ZDC"))){
512 AliWarning("No default storage set or default storage doesn't contain ZDC!");
513 manager->SetDefaultStorage(uri);
514 deleteManager = kTRUE;
517 AliCDBStorage *storage = manager->GetDefaultStorage();
520 AliCDBManager::Instance()->UnsetDefaultStorage();
521 defstorage = 0; // the storage is killed by AliCDBManager::Instance()->Destroy()
527 //_____________________________________________________________________________
528 AliZDCPedestals* AliZDCDigitizer::GetPedData() const
531 // Getting pedestal calibration object for ZDC set
533 AliCDBEntry *entry = AliCDBManager::Instance()->Get("ZDC/Calib/Pedestals");
534 if(!entry) AliFatal("No calibration data loaded!");
536 AliZDCPedestals *calibdata = dynamic_cast<AliZDCPedestals*> (entry->GetObject());
537 if(!calibdata) AliFatal("Wrong calibration object in calibration file!");
542 //_____________________________________________________________________________
543 AliZDCCalib* AliZDCDigitizer::GetCalibData() const
546 // Getting calibration object for ZDC set
548 AliCDBEntry *entry = AliCDBManager::Instance()->Get("ZDC/Calib/Calib");
549 if(!entry) AliFatal("No calibration data loaded!");
551 AliZDCCalib *calibdata = dynamic_cast<AliZDCCalib*> (entry->GetObject());
552 if(!calibdata) AliFatal("Wrong calibration object in calibration file!");
557 //_____________________________________________________________________________
558 AliZDCRecParam* AliZDCDigitizer::GetRecParam() const
561 // Getting calibration object for ZDC set
563 AliCDBEntry *entry = AliCDBManager::Instance()->Get("ZDC/Calib/RecParam");
564 if(!entry) AliFatal("No calibration data loaded!");
566 AliZDCRecParam *calibdata = dynamic_cast<AliZDCRecParam*> (entry->GetObject());
567 if(!calibdata) AliFatal("Wrong calibration object in calibration file!");