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