]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ZDC/AliZDCDigitizer.cxx
Minor fixes in the event tag to take into account the new way of storing the trigger...
[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   // Get calibration data
65   fCalibData = GetCalibData(); 
66
67 }
68
69 //____________________________________________________________________________
70 AliZDCDigitizer::~AliZDCDigitizer()
71 {
72 // Destructor
73
74 }
75
76
77 //____________________________________________________________________________
78 Bool_t AliZDCDigitizer::Init()
79 {
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.;
89    }
90    // ADC Caen V965
91   fADCRes[0] = 0.0000008; // ADC Resolution high gain: 200 fC/adcCh
92   fADCRes[1] = 0.0000064; // ADC Resolution low gain:  25  fC/adcCh
93
94   return kTRUE;
95 }
96
97 //____________________________________________________________________________
98 void AliZDCDigitizer::Exec(Option_t* /*option*/)
99 {
100   // Execute digitization
101
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!
110   
111   for (Int_t iSector1=0; iSector1<5; iSector1++) 
112     for (Int_t iSector2=0; iSector2<5; iSector2++) {
113       pm[iSector1][iSector2] = 0;
114     }
115
116   // impact parameter and number of spectators
117   Float_t impPar = -1;
118   Int_t specN = 0;
119   Int_t specP = 0;
120
121   // loop over input streams
122   for (Int_t iInput = 0; iInput<fManager->GetNinputs(); iInput++) {
123
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;
129
130     // load sdigits
131     loader->LoadSDigits();
132     TTree* treeS = loader->TreeS();
133     if (!treeS) continue;
134     AliZDCSDigit sdigit;
135     AliZDCSDigit* psdigit = &sdigit;
136     treeS->SetBranchAddress("ZDC", &psdigit);
137
138     // loop over sdigits
139     for (Int_t iSDigit=0; iSDigit<treeS->GetEntries(); iSDigit++) {
140       treeS->GetEntry(iSDigit);
141       //
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)));
146         continue;
147       }
148       //
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!
153           */
154     }
155
156     loader->UnloadSDigits();
157
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();
166     // 
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);
172   }
173
174   // add spectators
175   if (impPar >= 0) {
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);
182   }
183
184
185   // get the output run loader and loader
186   AliRunLoader* runLoader = 
187     AliRunLoader::GetRunLoader(fManager->GetOutputFolderName());
188   AliLoader* loader = runLoader->GetLoader("ZDCLoader");
189   if (!loader) {
190     AliError("no ZDC loader found");
191     return;
192   }
193
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();
200   AliZDCDigit digit;
201   AliZDCDigit* pdigit = &digit;
202   const Int_t kBufferSize = 4000;
203   treeD->Branch("ZDC", "AliZDCDigit", &pdigit, kBufferSize);
204
205   // Create digits
206   Int_t sector[2], sectorL[2];
207   Int_t digi[2], digiL[2];
208   for (sector[0]=1; sector[0]<=3; sector[0]++)
209     for (sector[1]=0; sector[1]<5; sector[1]++)  {
210         if ((sector[0]==3) && ((sector[1]<1) || (sector[1]>2))) continue;
211         for (Int_t res=0; res<2; res++){
212           digi[res] = Phe2ADCch(sector[0], sector[1], pm[sector[0]-1][sector[1]], res) 
213                     + Pedestal(sector[0], sector[1], res);
214         }
215         //printf("\t DIGIT added -> det = %d, quad = %d - digi[0,1] = [%d, %d]\n",
216         //      sector[0], sector[1], digi[0], digi[1]); // Chiara debugging!
217         new(pdigit) AliZDCDigit(sector, digi);
218         treeD->Fill();
219         //
220         // --- Adding digits for 2nd ZDC set (left side w.r.t. IP) ---
221         // --- they are copied from right ZDC digits
222         if(sector[0]==1 || sector[0]==2){
223            sectorL[0] = sector[0]+3;
224            sectorL[1] = sector[1];
225            for (Int_t res=0; res<2; res++){
226              digiL[res] = Phe2ADCch(sectorL[0], sectorL[1], pm[sector[0]-1][sector[1]], res) 
227                     + Pedestal(sectorL[0], sectorL[1], res);
228            }
229            //printf("\t DIGIT added -> det = %d, quad = %d - digi[0,1] = [%d, %d]\n",
230            //      sectorL[0], sectorL[1], digiL[0], digiL[1]); // Chiara debugging!
231           new(pdigit) AliZDCDigit(sectorL, digiL);
232         }
233         //
234         treeD->Fill();
235     }
236
237   // write the output tree
238   loader->WriteDigits("OVERWRITE");
239   loader->UnloadDigits();
240 }
241
242
243 //_____________________________________________________________________________
244 void AliZDCDigitizer::Fragmentation(Float_t impPar, Int_t specN, Int_t specP,
245                                     Int_t &freeSpecN, Int_t &freeSpecP) const
246 {
247 // simulate fragmentation of spectators
248
249   Int_t zz[100], nn[100];
250   AliZDCFragment frag(impPar);
251   for (Int_t j=0; j<=99; j++){
252      zz[j] =0;
253      nn[j] =0;
254   }
255
256   // Fragments generation
257   Int_t nAlpha;
258   frag.GenerateIMF(zz, nAlpha);
259
260   // Attach neutrons
261   Int_t ztot=0;
262   Int_t ntot=0;
263   frag.AttachNeutrons(zz, nn, ztot, ntot);
264   freeSpecN = specN-ntot-2*nAlpha;
265   freeSpecP = specP-ztot-2*nAlpha;
266   if(freeSpecN<0) freeSpecN=0;
267   if(freeSpecP<0) freeSpecP=0;
268   AliDebug(2, Form("FreeSpn = %d, FreeSpp = %d", freeSpecN, freeSpecP));
269 }
270
271 //_____________________________________________________________________________
272 void AliZDCDigitizer::SpectatorSignal(Int_t SpecType, Int_t numEvents, 
273                                       Float_t pm[3][5]) const
274 {
275 // add signal of the spectators
276  
277   TFile* file = NULL;
278   if (SpecType == 1) {          // --- Signal for spectator neutrons
279     file = TFile::Open("$ALICE/$ALICE_LEVEL/ZDC/ZNsignalntu.root");
280   } else if (SpecType == 2) {   // --- Signal for spectator protons
281     file = TFile::Open("$ALICE/$ALICE_LEVEL/ZDC/ZPsignalntu.root");
282   }
283   if (!file || !file->IsOpen()) {
284     AliError("Opening of file failed");
285     return;
286   }
287
288   TNtuple* zdcSignal = (TNtuple*) file->Get("ZDCSignal");
289   Int_t nentries = (Int_t) zdcSignal->GetEntries();
290   
291   Float_t *entry, hitsSpec[7];
292   Int_t pl, i, j, k, iev=0, rnd[125], volume[2];
293   for(pl=0;pl<125;pl++) rnd[pl] = 0;
294   if (numEvents > 125) {
295     AliWarning(Form("numEvents (%d) is larger than 125", numEvents));
296     numEvents = 125;
297   }
298   for(pl=0;pl<numEvents;pl++){
299      rnd[pl] = (Int_t) (9999*gRandom->Rndm());
300      if(rnd[pl] >= 9999) rnd[pl] = 9998;
301      //printf(" rnd[%d] = %d\n",pl,rnd[pl]);     
302   }
303   // Sorting vector in ascending order with C function QSORT 
304   qsort((void*)rnd,numEvents,sizeof(Int_t),comp);
305   do{
306      for(i=0; i<nentries; i++){  
307         zdcSignal->GetEvent(i);
308         entry = zdcSignal->GetArgs();
309         if(entry[0] == rnd[iev]){
310           for(k=0; k<2; k++) volume[k] = (Int_t) entry[k+1];
311           for(j=0; j<7; j++) hitsSpec[j] = entry[j+3];
312           //
313           Float_t lightQ = hitsSpec[4];
314           Float_t lightC = hitsSpec[5];
315           AliDebug(3, Form("SpectatorSignal -> vol = (%d, %d), lightQ = %.0f, lightC = %.0f",
316                            volume[0], volume[1], lightQ, lightC));
317           //printf("\n   Volume = (%d, %d), lightQ = %.0f, lightC = %.0f",
318           //                 volume[0], volume[1], lightQ, lightC);
319           if (volume[0] < 3) {  // ZN or ZP
320             pm[volume[0]-1][0] += lightC;
321             pm[volume[0]-1][volume[1]] += lightQ;
322             //printf("\n   pm[%d][0] = %.0f, pm[%d][%d] = %.0f\n",(volume[0]-1),pm[volume[0]-1][0],
323             //  (volume[0]-1),volume[1],pm[volume[0]-1][volume[1]]);
324           } 
325           else { 
326             if (volume[1] == 1) pm[2][1] += lightC; // ZEM 1
327             else                pm[2][2] += lightQ; // ZEM 2
328             //printf("\n   pm[2][1] = %.0f, pm[2][2] = %.0f\n",pm[2][1],pm[2][2]);
329           }
330         }
331         else if(entry[0] > rnd[iev]){
332           iev++;
333           continue;
334         }
335      }
336   }while(iev<numEvents);
337   
338   file->Close();
339   delete file;
340 }
341
342
343 //_____________________________________________________________________________
344 Int_t AliZDCDigitizer::Phe2ADCch(Int_t Det, Int_t Quad, Float_t Light, 
345                                  Int_t Res) const
346 {
347   // Evaluation of the ADC channel corresponding to the light yield Light
348   Int_t vADCch = (Int_t) (Light * fPMGain[Det-1][Quad] * fADCRes[Res]);
349   //printf("\t Phe2ADCch -> det %d quad %d - phe %.0f  ADC %d\n", Det,Quad,Light,ADCch);
350
351   return vADCch;
352 }
353
354 //_____________________________________________________________________________
355 Int_t AliZDCDigitizer::Pedestal(Int_t Det, Int_t Quad, Int_t Res) const
356 {
357   // Returns a pedestal for detector det, PM quad, channel with res.
358   //
359   Float_t meanPed;
360   if(Det != 3) meanPed = fCalibData->GetMeanPed(10*(Det-1)+Quad+5*Res);
361   else         meanPed = fCalibData->GetMeanPed(10*(Det-1)+Quad+1*Res);
362   
363   //printf("\t Pedestal -> det = %d, quad = %d, res = %d - Ped[%d] = %d\n",
364   //    Det, Quad, Res,10*(Det-1)+Quad+5*Res,(Int_t) meanPed); // Chiara debugging!
365   
366   // To create calibration object
367   /*Float_t meanPed;
368   meanPed = gRandom->Gaus((40.+10.*gRandom->Rndm()),5.);
369   */
370
371   return (Int_t) meanPed;
372 }
373
374 //_____________________________________________________________________________
375 AliCDBStorage* AliZDCDigitizer::SetStorage(const char *uri) 
376 {
377
378   Bool_t deleteManager = kFALSE;
379   
380   AliCDBManager *manager = AliCDBManager::Instance();
381   AliCDBStorage *defstorage = manager->GetDefaultStorage();
382   
383   if(!defstorage || !(defstorage->Contains("ZDC"))){ 
384      AliWarning("No default storage set or default storage doesn't contain ZDC!");
385      manager->SetDefaultStorage(uri);
386      deleteManager = kTRUE;
387   }
388  
389   AliCDBStorage *storage = manager->GetDefaultStorage();
390
391   if(deleteManager){
392     AliCDBManager::Instance()->UnsetDefaultStorage();
393     defstorage = 0;   // the storage is killed by AliCDBManager::Instance()->Destroy()
394   }
395
396   return storage; 
397 }
398
399 //_____________________________________________________________________________
400 AliZDCCalibData* AliZDCDigitizer::GetCalibData() const
401 {
402
403   // Getting calibration object for ZDC set
404
405   AliCDBEntry  *entry = AliCDBManager::Instance()->Get("ZDC/Calib/Data");
406   AliZDCCalibData *calibdata = (AliZDCCalibData*) entry->GetObject();
407
408   if (!calibdata)  AliWarning("No calibration data from calibration database !");
409
410   return calibdata;
411 }
412