]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ZDC/AliZDCDigitizer.cxx
Coding conventions (C.Oppedisano)
[u/mrichter/AliRoot.git] / ZDC / AliZDCDigitizer.cxx
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"
39 #include "AliCDBManager.h"
40 #include "AliCDBEntry.h"
41 #include "AliZDCSDigit.h"
42 #include "AliZDCDigit.h"
43 #include "AliZDCFragment.h"
44 #include "AliZDCDigitizer.h"
45
46 class AliCDBStorage;
47 class AliZDCCalibData;
48
49 ClassImp(AliZDCDigitizer)
50
51
52 //____________________________________________________________________________
53 AliZDCDigitizer::AliZDCDigitizer()
54 {
55 // Default constructor    
56
57 }
58
59 //____________________________________________________________________________
60 AliZDCDigitizer::AliZDCDigitizer(AliRunDigitizer* manager):
61   AliDigitizer(manager)
62 {
63
64   // Constructor    
65   fStorage =  SetStorage("local://$ALICE_ROOT");
66
67   // Get calibration data
68   fCalibData = GetCalibData(); 
69
70 }
71
72 //____________________________________________________________________________
73 AliZDCDigitizer::~AliZDCDigitizer()
74 {
75 // Destructor
76
77 }
78
79
80 //____________________________________________________________________________
81 Bool_t AliZDCDigitizer::Init()
82 {
83   // Initialize the digitizer
84   // NB -> PM gain vs. HV & ADC resolutions will move to DCDB ***************
85    for(Int_t j = 0; j < 5; j++){
86      fPMGain[0][j] = 50000.;
87      fPMGain[1][j] = 100000.;
88      fPMGain[2][j] = 100000.;
89      fPMGain[3][j] = 50000.;
90      fPMGain[4][j] = 100000.;
91      fPMGain[5][j] = 100000.;
92    }
93    // ADC Caen V965
94   fADCRes[0] = 0.0000008; // ADC Resolution high gain: 200 fC/adcCh
95   fADCRes[1] = 0.0000064; // ADC Resolution low gain:  25  fC/adcCh
96
97   return kTRUE;
98 }
99
100 //____________________________________________________________________________
101 void AliZDCDigitizer::Exec(Option_t* /*option*/)
102 {
103   // Execute digitization
104
105   Float_t pm[5][5]; // !!! 2nd ZDC set added (needed for trigger purposes!)
106   // *** 1st 3 arrays are digits from REAL (simulated) hits
107   // *** last 2 are copied from simulated digits
108   // --- pm[0][...] = light in ZN right  [C, Q1, Q2, Q3, Q4]
109   // --- pm[1][...] = light in ZP right [C, Q1, Q2, Q3, Q4]
110   // --- pm[2][...] = light in ZEM [x, 1, 2, x, x]
111   // --- pm[3][...] = light in ZN left [C, Q1, Q2, Q3, Q4] ->NEW!
112   // --- pm[4][...] = light in ZP left [C, Q1, Q2, Q3, Q4] ->NEW!
113   
114   for (Int_t iSector1=0; iSector1<5; iSector1++) 
115     for (Int_t iSector2=0; iSector2<5; iSector2++) {
116       pm[iSector1][iSector2] = 0;
117     }
118
119   // impact parameter and number of spectators
120   Float_t impPar = -1;
121   Int_t specN = 0;
122   Int_t specP = 0;
123
124   // loop over input streams
125   for (Int_t iInput = 0; iInput<fManager->GetNinputs(); iInput++) {
126
127     // get run loader and ZDC loader
128     AliRunLoader* runLoader = 
129       AliRunLoader::GetRunLoader(fManager->GetInputFolderName(iInput));
130     AliLoader* loader = runLoader->GetLoader("ZDCLoader");
131     if (!loader) continue;
132
133     // load sdigits
134     loader->LoadSDigits();
135     TTree* treeS = loader->TreeS();
136     if (!treeS) continue;
137     AliZDCSDigit sdigit;
138     AliZDCSDigit* psdigit = &sdigit;
139     treeS->SetBranchAddress("ZDC", &psdigit);
140
141     // loop over sdigits
142     for (Int_t iSDigit=0; iSDigit<treeS->GetEntries(); iSDigit++) {
143       treeS->GetEntry(iSDigit);
144       //
145       if (!psdigit) continue;
146       if ((sdigit.GetSector(1) < 0) || (sdigit.GetSector(1) > 4)) {
147         AliError(Form("\nsector[0] = %d, sector[1] = %d\n", 
148                       sdigit.GetSector(0), sdigit.GetSector(1)));
149         continue;
150       }
151       //
152       pm[(sdigit.GetSector(0))-1][sdigit.GetSector(1)] += sdigit.GetLightPM();
153       /*printf("\n\t Detector %d, Tower %d -> pm[%d][%d] = %.0f \n",
154           sdigit.GetSector(0), sdigit.GetSector(1),sdigit.GetSector(0)-1,
155           sdigit.GetSector(1), pm[sdigit.GetSector(0)-1][sdigit.GetSector(1)]); // Chiara debugging!
156           */
157     }
158
159     loader->UnloadSDigits();
160
161     // get the impact parameter and the number of spectators in case of hijing
162     if (!runLoader->GetAliRun()) runLoader->LoadgAlice();
163     AliHeader* header = runLoader->GetAliRun()->GetHeader();
164     if (!header) continue;
165     AliGenEventHeader* genHeader = header->GenEventHeader();
166     if (!genHeader) continue;
167     if (!genHeader->InheritsFrom(AliGenHijingEventHeader::Class())) continue;
168     impPar = ((AliGenHijingEventHeader*) genHeader)->ImpactParameter();
169     // 
170     specN = ((AliGenHijingEventHeader*) genHeader)->ProjSpectatorsn();
171     specP = ((AliGenHijingEventHeader*) genHeader)->ProjSpectatorsp();
172     AliDebug(2, Form("\n AliZDCDigitizer -> b = %f fm, Nspecn = %d, Nspecp = %d\n",
173                      impPar, specN, specP));
174     printf("\n\t AliZDCDigitizer -> b = %f fm, Nspecn = %d, Nspecp = %d\n", impPar, specN, specP);
175   }
176
177   // add spectators
178   if (impPar >= 0) {
179     Int_t freeSpecN, freeSpecP;
180     Fragmentation(impPar, specN, specP, freeSpecN, freeSpecP);
181     printf("\n\t AliZDCDigitizer ---- Adding signal for %d free spectator n\n",freeSpecN);
182     SpectatorSignal(1, freeSpecN, pm);
183     printf("\t AliZDCDigitizer ---- Adding signal for %d free spectator p\n\n",freeSpecP);
184     SpectatorSignal(2, freeSpecP, pm);
185   }
186
187
188   // get the output run loader and loader
189   AliRunLoader* runLoader = 
190     AliRunLoader::GetRunLoader(fManager->GetOutputFolderName());
191   AliLoader* loader = runLoader->GetLoader("ZDCLoader");
192   if (!loader) {
193     AliError("no ZDC loader found");
194     return;
195   }
196
197   // create the output tree
198   const char* mode = "update";
199   if (runLoader->GetEventNumber() == 0) mode = "recreate";
200   loader->LoadDigits(mode);
201   loader->MakeTree("D");
202   TTree* treeD = loader->TreeD();
203   AliZDCDigit digit;
204   AliZDCDigit* pdigit = &digit;
205   const Int_t kBufferSize = 4000;
206   treeD->Branch("ZDC", "AliZDCDigit", &pdigit, kBufferSize);
207
208   // Create digits
209   Int_t sector[2], sectorL[2];
210   Int_t digi[2], digiL[2];
211   for (sector[0]=1; sector[0]<=3; sector[0]++)
212     for (sector[1]=0; sector[1]<5; sector[1]++)  {
213         if ((sector[0]==3) && ((sector[1]<1) || (sector[1]>2))) continue;
214         for (Int_t res=0; res<2; res++){
215           digi[res] = Phe2ADCch(sector[0], sector[1], pm[sector[0]-1][sector[1]], res) 
216                     + Pedestal(sector[0], sector[1], res);
217         }
218         //printf("\t DIGIT added -> det = %d, quad = %d - digi[0,1] = [%d, %d]\n",
219         //      sector[0], sector[1], digi[0], digi[1]); // Chiara debugging!
220         new(pdigit) AliZDCDigit(sector, digi);
221         treeD->Fill();
222         //
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);
231            }
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!
234           new(pdigit) AliZDCDigit(sectorL, digiL);
235         }
236         //
237         treeD->Fill();
238     }
239
240   // write the output tree
241   loader->WriteDigits("OVERWRITE");
242   loader->UnloadDigits();
243 }
244
245
246 //_____________________________________________________________________________
247 void AliZDCDigitizer::Fragmentation(Float_t impPar, Int_t specN, Int_t specP,
248                                     Int_t &freeSpecN, Int_t &freeSpecP) const
249 {
250 // simulate fragmentation of spectators
251
252   Int_t zz[100], nn[100];
253   AliZDCFragment frag(impPar);
254   for (Int_t j=0; j<=99; j++){
255      zz[j] =0;
256      nn[j] =0;
257   }
258
259   // Fragments generation
260   Int_t nAlpha;
261   frag.GenerateIMF(zz, nAlpha);
262
263   // Attach neutrons
264   Int_t ztot=0;
265   Int_t ntot=0;
266   frag.AttachNeutrons(zz, nn, ztot, ntot);
267   freeSpecN = specN-ntot-2*nAlpha;
268   freeSpecP = specP-ztot-2*nAlpha;
269   if(freeSpecN<0) freeSpecN=0;
270   if(freeSpecP<0) freeSpecP=0;
271   AliDebug(2, Form("FreeSpn = %d, FreeSpp = %d", freeSpecN, freeSpecP));
272 }
273
274 //_____________________________________________________________________________
275 void AliZDCDigitizer::SpectatorSignal(Int_t SpecType, Int_t numEvents, 
276                                       Float_t pm[3][5]) const
277 {
278 // add signal of the spectators
279  
280   TFile* file = NULL;
281   if (SpecType == 1) {          // --- Signal for spectator neutrons
282     file = TFile::Open("$ALICE/$ALICE_LEVEL/ZDC/ZNsignalntu.root");
283   } else if (SpecType == 2) {   // --- Signal for spectator protons
284     file = TFile::Open("$ALICE/$ALICE_LEVEL/ZDC/ZPsignalntu.root");
285   }
286   if (!file || !file->IsOpen()) {
287     AliError("Opening of file failed");
288     return;
289   }
290
291   TNtuple* zdcSignal = (TNtuple*) file->Get("ZDCSignal");
292   Int_t nentries = (Int_t) zdcSignal->GetEntries();
293   
294   Float_t *entry, hitsSpec[7];
295   Int_t pl, i, j, k, iev=0, rnd[125], volume[2];
296   for(pl=0;pl<125;pl++) rnd[pl] = 0;
297   if (numEvents > 125) {
298     AliWarning(Form("numEvents (%d) is larger than 125", numEvents));
299     numEvents = 125;
300   }
301   for(pl=0;pl<numEvents;pl++){
302      rnd[pl] = (Int_t) (9999*gRandom->Rndm());
303      if(rnd[pl] >= 9999) rnd[pl] = 9998;
304      //printf(" rnd[%d] = %d\n",pl,rnd[pl]);     
305   }
306   // Sorting vector in ascending order with C function QSORT 
307   qsort((void*)rnd,numEvents,sizeof(Int_t),comp);
308   do{
309      for(i=0; i<nentries; i++){  
310         zdcSignal->GetEvent(i);
311         entry = zdcSignal->GetArgs();
312         if(entry[0] == rnd[iev]){
313           for(k=0; k<2; k++) volume[k] = (Int_t) entry[k+1];
314           for(j=0; j<7; j++) hitsSpec[j] = entry[j+3];
315           //
316           Float_t lightQ = hitsSpec[4];
317           Float_t lightC = hitsSpec[5];
318           AliDebug(3, Form("SpectatorSignal -> vol = (%d, %d), lightQ = %.0f, lightC = %.0f",
319                            volume[0], volume[1], lightQ, lightC));
320           //printf("\n   Volume = (%d, %d), lightQ = %.0f, lightC = %.0f",
321           //                 volume[0], volume[1], lightQ, lightC);
322           if (volume[0] < 3) {  // ZN or ZP
323             pm[volume[0]-1][0] += lightC;
324             pm[volume[0]-1][volume[1]] += lightQ;
325             //printf("\n   pm[%d][0] = %.0f, pm[%d][%d] = %.0f\n",(volume[0]-1),pm[volume[0]-1][0],
326             //  (volume[0]-1),volume[1],pm[volume[0]-1][volume[1]]);
327           } 
328           else { 
329             if (volume[1] == 1) pm[2][1] += lightC; // ZEM 1
330             else                pm[2][2] += lightQ; // ZEM 2
331             //printf("\n   pm[2][1] = %.0f, pm[2][2] = %.0f\n",pm[2][1],pm[2][2]);
332           }
333         }
334         else if(entry[0] > rnd[iev]){
335           iev++;
336           continue;
337         }
338      }
339   }while(iev<numEvents);
340   
341   file->Close();
342   delete file;
343 }
344
345
346 //_____________________________________________________________________________
347 Int_t AliZDCDigitizer::Phe2ADCch(Int_t Det, Int_t Quad, Float_t Light, 
348                                  Int_t Res) const
349 {
350   // Evaluation of the ADC channel corresponding to the light yield Light
351   Int_t vADCch = (Int_t) (Light * fPMGain[Det-1][Quad] * fADCRes[Res]);
352   //printf("\t Phe2ADCch -> det %d quad %d - phe %.0f  ADC %d\n", Det,Quad,Light,ADCch);
353
354   return vADCch;
355 }
356
357 //_____________________________________________________________________________
358 Int_t AliZDCDigitizer::Pedestal(Int_t Det, Int_t Quad, Int_t Res) const
359 {
360   // Returns a pedestal for detector det, PM quad, channel with res.
361   //
362   Float_t meanPed;
363   if(Det != 3) meanPed = fCalibData->GetMeanPed(10*(Det-1)+Quad+5*Res);
364   else         meanPed = fCalibData->GetMeanPed(10*(Det-1)+Quad+1*Res);
365   
366   //printf("\t Pedestal -> det = %d, quad = %d, res = %d - Ped[%d] = %d\n",
367   //    Det, Quad, Res,10*(Det-1)+Quad+5*Res,(Int_t) meanPed); // Chiara debugging!
368   
369   // To create calibration object
370   /*Float_t meanPed;
371   meanPed = gRandom->Gaus((40.+10.*gRandom->Rndm()),5.);
372   */
373
374   return (Int_t) meanPed;
375 }
376
377 //_____________________________________________________________________________
378 AliCDBStorage* AliZDCDigitizer::SetStorage(const char *uri) 
379 {
380
381   Bool_t deleteManager = kFALSE;
382   
383   AliCDBManager *manager = AliCDBManager::Instance();
384   AliCDBStorage *defstorage = manager->GetDefaultStorage();
385   
386   if(!defstorage || !(defstorage->Contains("ZDC"))){ 
387      AliWarning("No default storage set or default storage doesn't contain ZDC!");
388      manager->SetDefaultStorage(uri);
389      deleteManager = kTRUE;
390   }
391  
392   AliCDBStorage *storage = manager->GetDefaultStorage();
393
394   if(deleteManager){
395     AliCDBManager::Instance()->UnsetDefaultStorage();
396     defstorage = 0;   // the storage is killed by AliCDBManager::Instance()->Destroy()
397   }
398
399   return storage; 
400 }
401
402 //_____________________________________________________________________________
403 AliZDCCalibData* AliZDCDigitizer::GetCalibData() const
404 {
405   // returns pointer to AliZDCCalibData object
406   //
407   AliCDBEntry  *entry = fStorage->Get("ZDC/Calib/Data",0);  
408   AliZDCCalibData *calibdata = (AliZDCCalibData*) entry->GetObject();
409     
410   if (!calibdata)  AliWarning("No ZDC calibration data from calibration database!");
411
412   return calibdata;
413 }