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