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 "AliZDCSDigit.h"
40 #include "AliZDCDigit.h"
41 #include "AliZDCFragment.h"
42 #include "AliZDCDigitizer.h"
44 ClassImp(AliZDCDigitizer)
47 //____________________________________________________________________________
48 AliZDCDigitizer::AliZDCDigitizer()
50 // Default constructor
54 //____________________________________________________________________________
55 AliZDCDigitizer::AliZDCDigitizer(AliRunDigitizer* manager):
62 //____________________________________________________________________________
63 AliZDCDigitizer::~AliZDCDigitizer()
70 //____________________________________________________________________________
71 Bool_t AliZDCDigitizer::Init()
73 // Initialize the digitizer
75 for (Int_t i = 0; i < 3; i++) {
76 for(Int_t j = 0; j < 5; j++) {
77 fPMGain[i][j] = 100000.;
82 fADCRes[0] = 0.0000008; // ADC Resolution high gain: 200 fC/adcCh
83 fADCRes[1] = 0.0000064; // ADC Resolution low gain: 25 fC/adcCh
88 //____________________________________________________________________________
89 void AliZDCDigitizer::Exec(Option_t* /*option*/)
91 // Execute digitization
94 for (Int_t iSector1 = 0; iSector1 < 3; iSector1++)
95 for (Int_t iSector2 = 0; iSector2 < 5; iSector2++) {
96 pm[iSector1][iSector2] = 0;
98 // pm[0][...] = light in ZN [C, Q1, Q2, Q3, Q4]
99 // pm[1][...] = light in ZP [C, Q1, Q2, Q3, Q4]
100 // pm[2][...] = light in ZEM [x, 1, 2, x, x]
102 // impact parameter and number of spectators
107 // loop over input streams
108 for (Int_t iInput = 0; iInput < fManager->GetNinputs(); iInput++) {
110 // get run loader and ZDC loader
111 AliRunLoader* runLoader =
112 AliRunLoader::GetRunLoader(fManager->GetInputFolderName(iInput));
113 AliLoader* loader = runLoader->GetLoader("ZDCLoader");
114 if (!loader) continue;
117 loader->LoadSDigits();
118 TTree* treeS = loader->TreeS();
119 if (!treeS) continue;
121 AliZDCSDigit* psdigit = &sdigit;
122 treeS->SetBranchAddress("ZDC", &psdigit);
125 for (Int_t iSDigit = 0; iSDigit < treeS->GetEntries(); iSDigit++) {
126 treeS->GetEntry(iSDigit);
127 if (!psdigit) continue;
129 if ((sdigit.GetSector(1) < 0) || (sdigit.GetSector(1) > 4)) {
130 AliError(Form("\nsector[0] = %d, sector[1] = %d\n",
131 sdigit.GetSector(0), sdigit.GetSector(1)));
134 pm[sdigit.GetSector(0)-1][sdigit.GetSector(1)] += sdigit.GetLightPM();
135 //printf("\n\t sector[0] = %d, sector[1] = %d, pm[%d][%d] = %f \n",
136 // sdigit.GetSector(0), sdigit.GetSector(1),sdigit.GetSector(0)-1,
137 // sdigit.GetSector(1), pm[sdigit.GetSector(0)-1][sdigit.GetSector(1)]); // Chiara debugging!
141 loader->UnloadSDigits();
143 // get the impact parameter and the number of spectators in case of hijing
144 if (!runLoader->GetAliRun()) runLoader->LoadgAlice();
145 AliHeader* header = runLoader->GetAliRun()->GetHeader();
146 if (!header) continue;
147 AliGenEventHeader* genHeader = header->GenEventHeader();
148 if (!genHeader) continue;
149 if (!genHeader->InheritsFrom(AliGenHijingEventHeader::Class())) continue;
150 impPar = ((AliGenHijingEventHeader*) genHeader)->ImpactParameter();
151 // Until there is only 1 ZDC set the # of spectators must be divided by 2!!!
152 specN = ((AliGenHijingEventHeader*) genHeader)->Spectatorsn() / 2;
153 specP = ((AliGenHijingEventHeader*) genHeader)->Spectatorsp() / 2;
154 AliDebug(2, Form("\n b = %f fm, Nspecn = %d, Nspecp = %d\n",
155 impPar, specN, specP));
160 Int_t freeSpecN, freeSpecP;
161 Fragmentation(impPar, specN, specP, freeSpecN, freeSpecP);
162 SpectatorSignal(1, freeSpecN, pm);
163 SpectatorSignal(2, freeSpecP, pm);
166 // get the output run loader and loader
167 AliRunLoader* runLoader =
168 AliRunLoader::GetRunLoader(fManager->GetOutputFolderName());
169 AliLoader* loader = runLoader->GetLoader("ZDCLoader");
171 AliError("no ZDC loader found");
175 // create the output tree
176 const char* mode = "update";
177 if (runLoader->GetEventNumber() == 0) mode = "recreate";
178 loader->LoadDigits(mode);
179 loader->MakeTree("D");
180 TTree* treeD = loader->TreeD();
182 AliZDCDigit* pdigit = &digit;
183 const Int_t kBufferSize = 4000;
184 treeD->Branch("ZDC", "AliZDCDigit", &pdigit, kBufferSize);
189 for (sector[0] = 1; sector[0] <= 3; sector[0]++)
190 for (sector[1] = 0; sector[1] < 5; sector[1]++) {
191 if ((sector[0] == 3) && ((sector[1] < 1) || (sector[1] > 2))) continue;
192 for (Int_t res = 0; res < 2; res++){
193 digi[res] = Phe2ADCch(1, 0, pm[sector[0]-1][sector[1]], res)
195 //printf("\n\t sector[0] = %d, sector[1] = %d, digi[%d] = %d\n",
196 // sector[0], sector[1], res, digi[res]); // Chiara debugging!
198 new(pdigit) AliZDCDigit(sector, digi);
202 // write the output tree
203 loader->WriteDigits("OVERWRITE");
204 loader->UnloadDigits();
208 //_____________________________________________________________________________
209 void AliZDCDigitizer::Fragmentation(Float_t impPar, Int_t specN, Int_t specP,
210 Int_t &freeSpecN, Int_t &freeSpecP) const
212 // simulate fragmentation of spectators
214 Int_t zz[100], nn[100];
215 AliZDCFragment frag(impPar);
216 for (Int_t j=0; j<=99; j++){
221 // Fragments generation
223 frag.GenerateIMF(zz, nAlpha);
228 frag.AttachNeutrons(zz, nn, ztot, ntot);
229 freeSpecN = specN-ntot-2*nAlpha;
230 freeSpecP = specP-ztot-2*nAlpha;
231 if(freeSpecN<0) freeSpecN=0;
232 if(freeSpecP<0) freeSpecP=0;
233 AliDebug(2, Form("FreeSpn = %d, FreeSpp = %d", freeSpecN, freeSpecP));
236 //_____________________________________________________________________________
237 void AliZDCDigitizer::SpectatorSignal(Int_t SpecType, Int_t numEvents,
238 Float_t pm[3][5]) const
240 // add signal of the spectators
243 if (SpecType == 1) { // --- Signal for spectator neutrons
244 file = TFile::Open("$ALICE/$ALICE_LEVEL/ZDC/ZNsignalntu.root");
245 } else if (SpecType == 2) { // --- Signal for spectator protons
246 file = TFile::Open("$ALICE/$ALICE_LEVEL/ZDC/ZPsignalntu.root");
248 if (!file || !file->IsOpen()) {
249 AliError("Opening of file failed");
253 TNtuple* zdcSignal = (TNtuple*) file->Get("ZDCSignal");
254 Int_t nentries = (Int_t) zdcSignal->GetEntries();
256 Float_t *entry, hitsSpec[7];
257 Int_t pl, i, j, k, iev=0, rnd[125], volume[2];
258 for(pl=0;pl<125;pl++){
261 if (numEvents > 125) {
262 AliWarning(Form("numEvents (%d) is larger than 125", numEvents));
265 for(pl=0;pl<numEvents;pl++){
266 rnd[pl] = (Int_t) (9999*gRandom->Rndm());
267 if(rnd[pl] >= 9998) rnd[pl] = 9997;
268 //printf(" rnd[%d] = %d\n",pl,rnd[pl]);
270 // Sorting vector in ascending order with C function QSORT
271 qsort((void*)rnd,numEvents,sizeof(Int_t),comp);
273 for(i=0; i<nentries; i++){
274 zdcSignal->GetEvent(i);
275 entry = zdcSignal->GetArgs();
276 if(entry[0] == rnd[iev]){
277 for(k=0; k<2; k++) volume[k] = (Int_t) entry[k+1];
278 for(j=0; j<7; j++) hitsSpec[j] = entry[j+3];
280 Float_t lightQ = hitsSpec[4];
281 Float_t lightC = hitsSpec[5];
282 AliDebug(3, Form("volume = (%d, %d), lightQ = %.3f, lightC = %.3f",
283 volume[0], volume[1], lightQ, lightC));
284 if (volume[0] < 3) { // ZN or ZP
285 pm[volume[0]-1][0] += lightC;
286 pm[volume[0]-1][volume[1]] += lightQ;
288 if (volume[1] == 1) pm[2][1] += lightC;
289 else pm[2][2] += lightQ;
294 else if(entry[0] > rnd[iev]){
299 }while(iev<numEvents);
306 //_____________________________________________________________________________
307 Int_t AliZDCDigitizer::Phe2ADCch(Int_t Det, Int_t Quad, Float_t Light,
310 // Evaluation of the ADC channel corresponding to the light yield Light
312 return (Int_t) (Light * fPMGain[Det-1][Quad] * fADCRes[Res]);
315 //_____________________________________________________________________________
316 Int_t AliZDCDigitizer::Pedestal() const
318 // --- Pedestal value -> extracted from a gaussian distribution
319 // obtained from the beam test on the ZEM prototype (Aug. 2000)
321 Float_t pedMean = 50.;
322 Float_t pedWidth = 5.;
324 return (Int_t) gRandom->Gaus(pedMean, pedWidth);