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