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 AliZDCCalibData;
49 ClassImp(AliZDCDigitizer)
52 //____________________________________________________________________________
53 AliZDCDigitizer::AliZDCDigitizer()
55 // Default constructor
59 //____________________________________________________________________________
60 AliZDCDigitizer::AliZDCDigitizer(AliRunDigitizer* manager):
63 fIsCalibration=0; //By default the simulation doesn't create calib. data
64 // Get calibration data
65 fCalibData = GetCalibData();
69 //____________________________________________________________________________
70 AliZDCDigitizer::~AliZDCDigitizer()
77 //____________________________________________________________________________
78 Bool_t AliZDCDigitizer::Init()
80 // Initialize the digitizer
81 // NB -> PM gain vs. HV & ADC resolutions will move to DCDB ***************
82 for(Int_t j = 0; j < 5; j++){
83 fPMGain[0][j] = 50000.;
84 fPMGain[1][j] = 100000.;
85 fPMGain[2][j] = 100000.;
86 fPMGain[3][j] = 50000.;
87 fPMGain[4][j] = 100000.;
88 fPMGain[5][j] = 100000.;
91 fADCRes[0] = 0.0000008; // ADC Resolution high gain: 200 fC/adcCh
92 fADCRes[1] = 0.0000064; // ADC Resolution low gain: 25 fC/adcCh
97 //____________________________________________________________________________
98 void AliZDCDigitizer::Exec(Option_t* /*option*/)
100 // Execute digitization
102 Float_t pm[5][5]; // !!! 2nd ZDC set added (needed for trigger purposes!)
103 // *** 1st 3 arrays are digits from REAL (simulated) hits
104 // *** last 2 are copied from simulated digits
105 // --- pm[0][...] = light in ZN right [C, Q1, Q2, Q3, Q4]
106 // --- pm[1][...] = light in ZP right [C, Q1, Q2, Q3, Q4]
107 // --- pm[2][...] = light in ZEM [x, 1, 2, x, x]
108 // --- pm[3][...] = light in ZN left [C, Q1, Q2, Q3, Q4] ->NEW!
109 // --- pm[4][...] = light in ZP left [C, Q1, Q2, Q3, Q4] ->NEW!
111 for (Int_t iSector1=0; iSector1<5; iSector1++)
112 for (Int_t iSector2=0; iSector2<5; iSector2++){
113 pm[iSector1][iSector2] = 0;
116 // impact parameter and number of spectators
121 // loop over input streams
122 for (Int_t iInput = 0; iInput<fManager->GetNinputs(); iInput++){
124 // get run loader and ZDC loader
125 AliRunLoader* runLoader =
126 AliRunLoader::GetRunLoader(fManager->GetInputFolderName(iInput));
127 AliLoader* loader = runLoader->GetLoader("ZDCLoader");
128 if (!loader) continue;
131 loader->LoadSDigits();
132 TTree* treeS = loader->TreeS();
133 if (!treeS) continue;
135 AliZDCSDigit* psdigit = &sdigit;
136 treeS->SetBranchAddress("ZDC", &psdigit);
139 for (Int_t iSDigit=0; iSDigit<treeS->GetEntries(); iSDigit++){
140 treeS->GetEntry(iSDigit);
142 if (!psdigit) continue;
143 if ((sdigit.GetSector(1) < 0) || (sdigit.GetSector(1) > 4)){
144 AliError(Form("\nsector[0] = %d, sector[1] = %d\n",
145 sdigit.GetSector(0), sdigit.GetSector(1)));
149 pm[(sdigit.GetSector(0))-1][sdigit.GetSector(1)] += sdigit.GetLightPM();
150 /*printf("\n\t Detector %d, Tower %d -> pm[%d][%d] = %.0f \n",
151 sdigit.GetSector(0), sdigit.GetSector(1),sdigit.GetSector(0)-1,
152 sdigit.GetSector(1), pm[sdigit.GetSector(0)-1][sdigit.GetSector(1)]); // Chiara debugging!
156 loader->UnloadSDigits();
158 // get the impact parameter and the number of spectators in case of hijing
159 if (!runLoader->GetAliRun()) runLoader->LoadgAlice();
160 AliHeader* header = runLoader->GetAliRun()->GetHeader();
161 if (!header) continue;
162 AliGenEventHeader* genHeader = header->GenEventHeader();
163 if (!genHeader) continue;
164 if (!genHeader->InheritsFrom(AliGenHijingEventHeader::Class())) continue;
165 impPar = ((AliGenHijingEventHeader*) genHeader)->ImpactParameter();
167 specN = ((AliGenHijingEventHeader*) genHeader)->ProjSpectatorsn();
168 specP = ((AliGenHijingEventHeader*) genHeader)->ProjSpectatorsp();
169 AliDebug(2, Form("\n AliZDCDigitizer -> b = %f fm, Nspecn = %d, Nspecp = %d\n",
170 impPar, specN, specP));
171 printf("\n\t AliZDCDigitizer -> b = %f fm, Nspecn = %d, Nspecp = %d\n", impPar, specN, specP);
176 Int_t freeSpecN, freeSpecP;
177 Fragmentation(impPar, specN, specP, freeSpecN, freeSpecP);
178 printf("\n\t AliZDCDigitizer ---- Adding signal for %d free spectator n\n",freeSpecN);
179 SpectatorSignal(1, freeSpecN, pm);
180 printf("\t AliZDCDigitizer ---- Adding signal for %d free spectator p\n\n",freeSpecP);
181 SpectatorSignal(2, freeSpecP, pm);
185 // get the output run loader and loader
186 AliRunLoader* runLoader =
187 AliRunLoader::GetRunLoader(fManager->GetOutputFolderName());
188 AliLoader* loader = runLoader->GetLoader("ZDCLoader");
190 AliError("no ZDC loader found");
194 // create the output tree
195 const char* mode = "update";
196 if (runLoader->GetEventNumber() == 0) mode = "recreate";
197 loader->LoadDigits(mode);
198 loader->MakeTree("D");
199 TTree* treeD = loader->TreeD();
201 AliZDCDigit* pdigit = &digit;
202 const Int_t kBufferSize = 4000;
203 treeD->Branch("ZDC", "AliZDCDigit", &pdigit, kBufferSize);
206 if(fIsCalibration!=0) printf("\t **** AliZDCDigitizer -> Creating calibration data (pedestals)\n");
207 Int_t sector[2], sectorL[2];
208 Int_t digi[2], digiL[2];
209 for(sector[0]=1; sector[0]<=3; sector[0]++){
210 for(sector[1]=0; sector[1]<5; sector[1]++){
211 if((sector[0]==3) && ((sector[1]<1) || (sector[1]>2))) continue;
212 for (Int_t res=0; res<2; res++){
213 digi[res] = Phe2ADCch(sector[0], sector[1], pm[sector[0]-1][sector[1]], res)
214 + Pedestal(sector[0], sector[1], res);
216 /*printf("\t DIGIT added -> det = %d, quad = %d - digi[0,1] = [%d, %d]\n",
217 sector[0], sector[1], digi[0], digi[1]); // Chiara debugging!
220 new(pdigit) AliZDCDigit(sector, digi);
223 // --- Adding digits for 2nd ZDC set (left side w.r.t. IP) ---
224 // --- they are copied from right ZDC digits
225 if(sector[0]==1 || sector[0]==2){
226 sectorL[0] = sector[0]+3;
227 sectorL[1] = sector[1];
228 for (Int_t res=0; res<2; res++){
229 digiL[res] = Phe2ADCch(sectorL[0], sectorL[1], pm[sector[0]-1][sector[1]], res)
230 + Pedestal(sectorL[0], sectorL[1], res);
232 /*printf("\t DIGIT added -> det = %d, quad = %d - digi[0,1] = [%d, %d]\n",
233 sectorL[0], sectorL[1], digiL[0], digiL[1]); // Chiara debugging!
236 new(pdigit) AliZDCDigit(sectorL, digiL);
240 //printf("\t AliZDCDigitizer -> TreeD has %d entries\n",(Int_t) treeD->GetEntries());
243 // write the output tree
244 loader->WriteDigits("OVERWRITE");
245 loader->UnloadDigits();
249 //_____________________________________________________________________________
250 void AliZDCDigitizer::Fragmentation(Float_t impPar, Int_t specN, Int_t specP,
251 Int_t &freeSpecN, Int_t &freeSpecP) const
253 // simulate fragmentation of spectators
255 Int_t zz[100], nn[100];
256 AliZDCFragment frag(impPar);
257 for (Int_t j=0; j<=99; j++){
262 // Fragments generation
264 frag.GenerateIMF(zz, nAlpha);
269 frag.AttachNeutrons(zz, nn, ztot, ntot);
270 freeSpecN = specN-ntot-2*nAlpha;
271 freeSpecP = specP-ztot-2*nAlpha;
272 if(freeSpecN<0) freeSpecN=0;
273 if(freeSpecP<0) freeSpecP=0;
274 AliDebug(2, Form("FreeSpn = %d, FreeSpp = %d", freeSpecN, freeSpecP));
277 //_____________________________________________________________________________
278 void AliZDCDigitizer::SpectatorSignal(Int_t SpecType, Int_t numEvents,
279 Float_t pm[3][5]) const
281 // add signal of the spectators
284 if (SpecType == 1) { // --- Signal for spectator neutrons
285 file = TFile::Open("$ALICE/$ALICE_LEVEL/ZDC/ZNsignalntu.root");
286 } else if (SpecType == 2) { // --- Signal for spectator protons
287 file = TFile::Open("$ALICE/$ALICE_LEVEL/ZDC/ZPsignalntu.root");
289 if (!file || !file->IsOpen()) {
290 AliError("Opening of file failed");
294 TNtuple* zdcSignal = (TNtuple*) file->Get("ZDCSignal");
295 Int_t nentries = (Int_t) zdcSignal->GetEntries();
297 Float_t *entry, hitsSpec[7];
298 Int_t pl, i, j, k, iev=0, rnd[125], volume[2];
299 for(pl=0;pl<125;pl++) rnd[pl] = 0;
300 if (numEvents > 125) {
301 AliWarning(Form("numEvents (%d) is larger than 125", numEvents));
304 for(pl=0;pl<numEvents;pl++){
305 rnd[pl] = (Int_t) (9999*gRandom->Rndm());
306 if(rnd[pl] >= 9999) rnd[pl] = 9998;
307 //printf(" rnd[%d] = %d\n",pl,rnd[pl]);
309 // Sorting vector in ascending order with C function QSORT
310 qsort((void*)rnd,numEvents,sizeof(Int_t),comp);
312 for(i=0; i<nentries; i++){
313 zdcSignal->GetEvent(i);
314 entry = zdcSignal->GetArgs();
315 if(entry[0] == rnd[iev]){
316 for(k=0; k<2; k++) volume[k] = (Int_t) entry[k+1];
317 for(j=0; j<7; j++) hitsSpec[j] = entry[j+3];
319 Float_t lightQ = hitsSpec[4];
320 Float_t lightC = hitsSpec[5];
321 AliDebug(3, Form("SpectatorSignal -> vol = (%d, %d), lightQ = %.0f, lightC = %.0f",
322 volume[0], volume[1], lightQ, lightC));
323 //printf("\n Volume = (%d, %d), lightQ = %.0f, lightC = %.0f",
324 // volume[0], volume[1], lightQ, lightC);
325 if (volume[0] < 3) { // ZN or ZP
326 pm[volume[0]-1][0] += lightC;
327 pm[volume[0]-1][volume[1]] += lightQ;
328 //printf("\n pm[%d][0] = %.0f, pm[%d][%d] = %.0f\n",(volume[0]-1),pm[volume[0]-1][0],
329 // (volume[0]-1),volume[1],pm[volume[0]-1][volume[1]]);
332 if (volume[1] == 1) pm[2][1] += lightC; // ZEM 1
333 else pm[2][2] += lightQ; // ZEM 2
334 //printf("\n pm[2][1] = %.0f, pm[2][2] = %.0f\n",pm[2][1],pm[2][2]);
337 else if(entry[0] > rnd[iev]){
342 }while(iev<numEvents);
349 //_____________________________________________________________________________
350 Int_t AliZDCDigitizer::Phe2ADCch(Int_t Det, Int_t Quad, Float_t Light,
353 // Evaluation of the ADC channel corresponding to the light yield Light
354 Int_t vADCch = (Int_t) (Light * fPMGain[Det-1][Quad] * fADCRes[Res]);
355 //printf("\t Phe2ADCch -> det %d quad %d - phe %.0f ADC %d\n", Det,Quad,Light,ADCch);
360 //_____________________________________________________________________________
361 Int_t AliZDCDigitizer::Pedestal(Int_t Det, Int_t Quad, Int_t Res) const
363 // Returns a pedestal for detector det, PM quad, channel with res.
368 if(fIsCalibration == 0){
369 Float_t meanPed, Pedwidth;
371 if(Det==1|| Det==2) index = 10*(Det-1)+Quad+5*Res; // ZN1, ZP1
372 else if(Det==3) index = 10*(Det-1)+(Quad-1)+Res; // ZEM
373 else if(Det==4|| Det==5) index = 10*(Det-2)+Quad+5*Res+4; // ZN2, ZP2
374 meanPed = fCalibData->GetMeanPed(index);
375 Pedwidth = fCalibData->GetMeanPedWidth(index);
376 PedValue = gRandom->Gaus(meanPed,Pedwidth);
378 /*printf("\t Pedestal -> det = %d, quad = %d, res = %d - Ped[%d] = %d\n",
379 Det, Quad, index,(Int_t) PedValue); // Chiara debugging!
383 // To create calibration object
384 else PedValue = gRandom->Gaus((40.+10.*gRandom->Rndm()),5.);
387 return (Int_t) PedValue;
390 //_____________________________________________________________________________
391 AliCDBStorage* AliZDCDigitizer::SetStorage(const char *uri)
394 Bool_t deleteManager = kFALSE;
396 AliCDBManager *manager = AliCDBManager::Instance();
397 AliCDBStorage *defstorage = manager->GetDefaultStorage();
399 if(!defstorage || !(defstorage->Contains("ZDC"))){
400 AliWarning("No default storage set or default storage doesn't contain ZDC!");
401 manager->SetDefaultStorage(uri);
402 deleteManager = kTRUE;
405 AliCDBStorage *storage = manager->GetDefaultStorage();
408 AliCDBManager::Instance()->UnsetDefaultStorage();
409 defstorage = 0; // the storage is killed by AliCDBManager::Instance()->Destroy()
415 //_____________________________________________________________________________
416 AliZDCCalibData* AliZDCDigitizer::GetCalibData() const
419 // Getting calibration object for ZDC set
421 AliCDBEntry *entry = AliCDBManager::Instance()->Get("ZDC/Calib/Data");
422 AliZDCCalibData *calibdata = (AliZDCCalibData*) entry->GetObject();
424 if (!calibdata) AliWarning("No calibration data from calibration database !");