]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ZDC/AliZDCDigitizer.cxx
for non-miscalibrated digits, set an ad-hoc conversion factor fAdC->fToT to have...
[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   // ------------------------------------------------------------
121   // !!! 2nd ZDC set added (needed for trigger purposes!)
122   // *** 1st 3 arrays are digits from REAL (simulated) hits
123   // *** last 2 are copied from simulated digits
124   // --- pm[0][...] = light in ZN right  [C, Q1, Q2, Q3, Q4]
125   // --- pm[1][...] = light in ZP right [C, Q1, Q2, Q3, Q4]
126   // --- pm[2][...] = light in ZEM [x, 1, 2, x, x]
127   // --- pm[3][...] = light in ZN left [C, Q1, Q2, Q3, Q4] ->NEW!
128   // --- pm[4][...] = light in ZP left [C, Q1, Q2, Q3, Q4] ->NEW!
129   // ------------------------------------------------------------
130   Float_t pm[5][5]; 
131   for(Int_t iSector1=0; iSector1<5; iSector1++) 
132     for(Int_t iSector2=0; iSector2<5; iSector2++){
133       pm[iSector1][iSector2] = 0;
134     }
135     
136   // ------------------------------------------------------------
137   // ### Out of time ADC added (22 channels)
138   // --- same codification as for signal PTMs (see above)
139   // ------------------------------------------------------------
140   Float_t pmoot[5][5];
141   for(Int_t iSector1=0; iSector1<5; iSector1++) 
142     for(Int_t iSector2=0; iSector2<5; iSector2++){
143       pmoot[iSector1][iSector2] = 0;
144     }
145
146   // impact parameter and number of spectators
147   Float_t impPar = -1;
148   Int_t specN = 0;
149   Int_t specP = 0;
150
151   // loop over input streams
152   for(Int_t iInput = 0; iInput<fManager->GetNinputs(); iInput++){
153
154     // get run loader and ZDC loader
155     AliRunLoader* runLoader = 
156       AliRunLoader::GetRunLoader(fManager->GetInputFolderName(iInput));
157     AliLoader* loader = runLoader->GetLoader("ZDCLoader");
158     if(!loader) continue;
159
160     // load sdigits
161     loader->LoadSDigits();
162     TTree* treeS = loader->TreeS();
163     if(!treeS) continue;
164     AliZDCSDigit sdigit;
165     AliZDCSDigit* psdigit = &sdigit;
166     treeS->SetBranchAddress("ZDC", &psdigit);
167
168     // loop over sdigits
169     for(Int_t iSDigit=0; iSDigit<treeS->GetEntries(); iSDigit++){
170       treeS->GetEntry(iSDigit);
171       //
172       if(!psdigit) continue;
173       if((sdigit.GetSector(1) < 0) || (sdigit.GetSector(1) > 4)){
174         AliError(Form("\nsector[0] = %d, sector[1] = %d\n", 
175                       sdigit.GetSector(0), sdigit.GetSector(1)));
176         continue;
177       }
178       //
179       pm[(sdigit.GetSector(0))-1][sdigit.GetSector(1)] += sdigit.GetLightPM();
180       /*printf("\n\t Detector %d, Tower %d -> pm[%d][%d] = %.0f \n",
181           sdigit.GetSector(0), sdigit.GetSector(1),sdigit.GetSector(0)-1,
182           sdigit.GetSector(1), pm[sdigit.GetSector(0)-1][sdigit.GetSector(1)]); // Chiara debugging!
183       */
184     }
185
186     loader->UnloadSDigits();
187
188     // get the impact parameter and the number of spectators in case of hijing
189     if(!runLoader->GetAliRun()) runLoader->LoadgAlice();
190     AliHeader* header = runLoader->GetAliRun()->GetHeader();
191     if(!header) continue;
192     AliGenEventHeader* genHeader = header->GenEventHeader();
193     if(!genHeader) continue;
194     if(!genHeader->InheritsFrom(AliGenHijingEventHeader::Class())) continue;
195     impPar = ((AliGenHijingEventHeader*) genHeader)->ImpactParameter();
196     // 
197     specN = ((AliGenHijingEventHeader*) genHeader)->ProjSpectatorsn();
198     specP = ((AliGenHijingEventHeader*) genHeader)->ProjSpectatorsp();
199     AliDebug(2, Form("\n AliZDCDigitizer -> b = %f fm, Nspecn = %d, Nspecp = %d\n",
200                      impPar, specN, specP));
201     printf("\n\t AliZDCDigitizer -> b = %f fm, NSpecn = %d, NSpecp = %d\n", impPar, specN, specP);
202   }
203
204   // add spectators
205   if(impPar >= 0) {
206     Int_t freeSpecN, freeSpecP;
207     Fragmentation(impPar, specN, specP, freeSpecN, freeSpecP);
208     printf("\n\t AliZDCDigitizer ---- Adding signal for %d free spectator n\n",freeSpecN);
209     SpectatorSignal(1, freeSpecN, pm);
210     printf("\t AliZDCDigitizer ---- Adding signal for %d free spectator p\n\n",freeSpecP);
211     SpectatorSignal(2, freeSpecP, pm);
212   }
213
214
215   // get the output run loader and loader
216   AliRunLoader* runLoader = 
217     AliRunLoader::GetRunLoader(fManager->GetOutputFolderName());
218   AliLoader* loader = runLoader->GetLoader("ZDCLoader");
219   if(!loader) {
220     AliError("no ZDC loader found");
221     return;
222   }
223
224   // create the output tree
225   const char* mode = "update";
226   if(runLoader->GetEventNumber() == 0) mode = "recreate";
227   loader->LoadDigits(mode);
228   loader->MakeTree("D");
229   TTree* treeD = loader->TreeD();
230   AliZDCDigit digit;
231   AliZDCDigit* pdigit = &digit;
232   const Int_t kBufferSize = 4000;
233   treeD->Branch("ZDC", "AliZDCDigit", &pdigit, kBufferSize);
234
235   // Create digits
236   Int_t sector[2], sectorL[2];
237   Int_t digi[2], digiL[2], digioot[2];
238   for(sector[0]=1; sector[0]<=3; sector[0]++){
239     for(sector[1]=0; sector[1]<5; sector[1]++){
240         if((sector[0]==3) && ((sector[1]<1) || (sector[1]>2))) continue;
241         for(Int_t res=0; res<2; res++){
242            digi[res] = Phe2ADCch(sector[0], sector[1], pm[sector[0]-1][sector[1]], res) 
243                     + Pedestal(sector[0], sector[1], res);
244         }
245         /*printf("\t DIGIT added -> det = %d, quad = %d - digi[0,1] = [%d, %d]\n",
246              sector[0], sector[1], digi[0], digi[1]); // Chiara debugging!
247         */
248         //
249         new(pdigit) AliZDCDigit(sector, digi);
250         treeD->Fill();
251         //
252         // --- Adding digits for 2nd ZDC set (left side w.r.t. IP) ---
253         // --- they are copied from right ZDC digits
254         if(sector[0]==1 || sector[0]==2){
255            sectorL[0] = sector[0]+3;
256            sectorL[1] = sector[1];
257            for(Int_t res=0; res<2; res++){
258              digiL[res] = Phe2ADCch(sectorL[0], sectorL[1], pm[sector[0]-1][sector[1]], res) 
259                     + Pedestal(sectorL[0], sectorL[1], res);
260            }
261            /*printf("\t DIGIT added -> det = %d, quad = %d - digi[0,1] = [%d, %d]\n",
262                  sectorL[0], sectorL[1], digiL[0], digiL[1]); // Chiara debugging!
263            */
264            //
265            new(pdigit) AliZDCDigit(sectorL, digiL);
266            treeD->Fill();
267         }
268         //
269         //printf("\t AliZDCDigitizer -> TreeD has %d entries\n",(Int_t) treeD->GetEntries());
270     }
271   }
272   //
273   // --- Adding digits for out-of-time channels after signal digits
274   for(sector[0]=1; sector[0]<=3; sector[0]++){
275     for(sector[1]=0; sector[1]<5; sector[1]++){
276         if((sector[0]==3) && ((sector[1]<1) || (sector[1]>2))) continue;
277         for(Int_t res=0; res<2; res++){
278            digioot[res] = Pedestal(sector[0], sector[1], res); // out-of-time ADCs
279         }
280         /*printf("\t DIGIToot added -> det = %d, quad = %d - digi[0,1] = [%d, %d]\n",
281              sector[0], sector[1], digioot[0], digioot[1]); // Chiara debugging!
282         */
283         //
284         new(pdigit) AliZDCDigit(sector, digioot);
285         treeD->Fill();
286         //
287         if(sector[0]==1 || sector[0]==2){
288            sectorL[0] = sector[0]+3;
289            sectorL[1] = sector[1];
290            for(Int_t res=0; res<2; res++){
291              digioot[res] = Pedestal(sectorL[0], sectorL[1], res); // out-of-time ADCs
292            }
293            /*printf("\t DIGIToot added -> det = %d, quad = %d - digi[0,1] = [%d, %d]\n",
294                  sectorL[0], sectorL[1], digioot[0], digioot[1]); // Chiara debugging!
295            */
296            //
297            new(pdigit) AliZDCDigit(sectorL, digioot);
298            treeD->Fill();
299         }
300         //
301         //printf("\t AliZDCDigitizer -> TreeD has %d entries\n",(Int_t) treeD->GetEntries());
302     }
303   }
304   // write the output tree
305   loader->WriteDigits("OVERWRITE");
306   loader->UnloadDigits();
307 }
308
309
310 //_____________________________________________________________________________
311 void AliZDCDigitizer::Fragmentation(Float_t impPar, Int_t specN, Int_t specP,
312                                     Int_t &freeSpecN, Int_t &freeSpecP) const
313 {
314 // simulate fragmentation of spectators
315
316   Int_t zz[100], nn[100];
317   AliZDCFragment frag(impPar);
318   for(Int_t j=0; j<=99; j++){
319      zz[j] =0;
320      nn[j] =0;
321   }
322
323   // Fragments generation
324   Int_t nAlpha;
325   frag.GenerateIMF(zz, nAlpha);
326
327   // Attach neutrons
328   Int_t ztot=0;
329   Int_t ntot=0;
330   frag.AttachNeutrons(zz, nn, ztot, ntot);
331   freeSpecN = specN-ntot-2*nAlpha;
332   freeSpecP = specP-ztot-2*nAlpha;
333   if(freeSpecN<0) freeSpecN=0;
334   if(freeSpecP<0) freeSpecP=0;
335   AliDebug(2, Form("FreeSpn = %d, FreeSpp = %d", freeSpecN, freeSpecP));
336 }
337
338 //_____________________________________________________________________________
339 void AliZDCDigitizer::SpectatorSignal(Int_t SpecType, Int_t numEvents, 
340                                       Float_t pm[3][5]) const
341 {
342 // add signal of the spectators
343  
344   TFile* file = NULL;
345   if(SpecType == 1) {           // --- Signal for spectator neutrons
346     file = TFile::Open("$ALICE_ROOT/ZDC/ZNsignalntu.root");
347   } else if(SpecType == 2) {    // --- Signal for spectator protons
348     file = TFile::Open("$ALICE_ROOT/ZDC/ZPsignalntu.root");
349   }
350   if(!file || !file->IsOpen()) {
351     AliError("Opening of file failed");
352     return;
353   }
354
355   TNtuple* zdcSignal = (TNtuple*) file->Get("ZDCSignal");
356   Int_t nentries = (Int_t) zdcSignal->GetEntries();
357   
358   Float_t *entry, hitsSpec[7];
359   Int_t pl, i, j, k, iev=0, rnd[125], volume[2];
360   for(pl=0;pl<125;pl++) rnd[pl] = 0;
361   if(numEvents > 125) {
362     AliWarning(Form("numEvents (%d) is larger than 125", numEvents));
363     numEvents = 125;
364   }
365   for(pl=0;pl<numEvents;pl++){
366      rnd[pl] = (Int_t) (9999*gRandom->Rndm());
367      if(rnd[pl] >= 9999) rnd[pl] = 9998;
368      //printf(" rnd[%d] = %d\n",pl,rnd[pl]);     
369   }
370   // Sorting vector in ascending order with C function QSORT 
371   qsort((void*)rnd,numEvents,sizeof(Int_t),comp);
372   do{
373      for(i=0; i<nentries; i++){  
374         zdcSignal->GetEvent(i);
375         entry = zdcSignal->GetArgs();
376         if(entry[0] == rnd[iev]){
377           for(k=0; k<2; k++) volume[k] = (Int_t) entry[k+1];
378           for(j=0; j<7; j++) hitsSpec[j] = entry[j+3];
379           //
380           Float_t lightQ = hitsSpec[4];
381           Float_t lightC = hitsSpec[5];
382           AliDebug(3, Form("SpectatorSignal -> vol = (%d, %d), lightQ = %.0f, lightC = %.0f",
383                            volume[0], volume[1], lightQ, lightC));
384           //printf("\n   Volume = (%d, %d), lightQ = %.0f, lightC = %.0f",
385           //                 volume[0], volume[1], lightQ, lightC);
386           if(volume[0] < 3) {  // ZN or ZP
387             pm[volume[0]-1][0] += lightC;
388             pm[volume[0]-1][volume[1]] += lightQ;
389             //printf("\n   pm[%d][0] = %.0f, pm[%d][%d] = %.0f\n",(volume[0]-1),pm[volume[0]-1][0],
390             //  (volume[0]-1),volume[1],pm[volume[0]-1][volume[1]]);
391           } 
392           else { 
393             if(volume[1] == 1) pm[2][1] += lightC; // ZEM 1
394             else                pm[2][2] += lightQ; // ZEM 2
395             //printf("\n   pm[2][1] = %.0f, pm[2][2] = %.0f\n",pm[2][1],pm[2][2]);
396           }
397         }
398         else if(entry[0] > rnd[iev]){
399           iev++;
400           continue;
401         }
402      }
403   }while(iev<numEvents);
404   
405   file->Close();
406   delete file;
407 }
408
409
410 //_____________________________________________________________________________
411 Int_t AliZDCDigitizer::Phe2ADCch(Int_t Det, Int_t Quad, Float_t Light, 
412                                  Int_t Res) const
413 {
414   // Evaluation of the ADC channel corresponding to the light yield Light
415   Int_t vADCch = (Int_t) (Light * fPMGain[Det-1][Quad] * fADCRes[Res]);
416   //printf("\t Phe2ADCch -> det %d quad %d - phe %.0f  ADC %d\n", Det,Quad,Light,ADCch);
417
418   return vADCch;
419 }
420
421 //_____________________________________________________________________________
422 Int_t AliZDCDigitizer::Pedestal(Int_t Det, Int_t Quad, Int_t Res) const
423 {
424   // Returns a pedestal for detector det, PM quad, channel with res.
425   //
426   Float_t PedValue;
427   
428   // Normal run
429   if(fIsCalibration == 0){
430     Float_t meanPed, Pedwidth;
431     Int_t index=0;
432     if(Det==1|| Det==2)         index = 10*(Det-1)+Quad+5*Res;   // ZN1, ZP1
433     else if(Det==3)             index = 10*(Det-1)+(Quad-1)+Res; // ZEM
434     else if(Det==4|| Det==5)    index = 10*(Det-2)+Quad+5*Res+4; // ZN2, ZP2
435     meanPed = fCalibData->GetMeanPed(index);
436     Pedwidth = fCalibData->GetMeanPedWidth(index);
437     PedValue = gRandom->Gaus(meanPed,Pedwidth);
438     //
439     /*printf("\t Pedestal -> det = %d, quad = %d, res = %d - Ped[%d] = %d\n",
440         Det, Quad, index,(Int_t) PedValue); // Chiara debugging!
441     */
442   }
443   
444   // To create calibration object
445   else PedValue = gRandom->Gaus((40.+10.*gRandom->Rndm()),5.);
446   
447
448   return (Int_t) PedValue;
449 }
450
451 //_____________________________________________________________________________
452 AliCDBStorage* AliZDCDigitizer::SetStorage(const char *uri) 
453 {
454
455   Bool_t deleteManager = kFALSE;
456   
457   AliCDBManager *manager = AliCDBManager::Instance();
458   AliCDBStorage *defstorage = manager->GetDefaultStorage();
459   
460   if(!defstorage || !(defstorage->Contains("ZDC"))){ 
461      AliWarning("No default storage set or default storage doesn't contain ZDC!");
462      manager->SetDefaultStorage(uri);
463      deleteManager = kTRUE;
464   }
465  
466   AliCDBStorage *storage = manager->GetDefaultStorage();
467
468   if(deleteManager){
469     AliCDBManager::Instance()->UnsetDefaultStorage();
470     defstorage = 0;   // the storage is killed by AliCDBManager::Instance()->Destroy()
471   }
472
473   return storage; 
474 }
475
476 //_____________________________________________________________________________
477 AliZDCCalibData* AliZDCDigitizer::GetCalibData() const
478 {
479
480   // Getting calibration object for ZDC set
481
482   AliCDBEntry  *entry = AliCDBManager::Instance()->Get("ZDC/Calib/Data");
483   if(!entry) AliFatal("No calibration data loaded!");  
484
485   AliZDCCalibData *calibdata = dynamic_cast<AliZDCCalibData*>  (entry->GetObject());
486   if(!calibdata)  AliFatal("Wrong calibration object in calibration  file!");
487
488   return calibdata;
489 }
490