]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ZDC/AliZDCDigitizer.cxx
Updated DA for STANDALONE PEDESTAL runs
[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("\n\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 
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, # generated spectator n = %d," 
202     " # generated spectator p = %d\n", impPar, specN, specP);
203   }
204
205   // add spectators
206   if(impPar >= 0) {
207     Int_t freeSpecN, freeSpecP;
208     Fragmentation(impPar, specN, specP, freeSpecN, freeSpecP);
209     printf("\n\t AliZDCDigitizer -> Adding signal for %d free spectator n\n",freeSpecN);
210     SpectatorSignal(1, freeSpecN, pm);
211     printf("\t AliZDCDigitizer -> Adding signal for %d free spectator p\n\n",freeSpecP);
212     SpectatorSignal(2, freeSpecP, pm);
213   }
214
215
216   // get the output run loader and loader
217   AliRunLoader* runLoader = 
218     AliRunLoader::GetRunLoader(fManager->GetOutputFolderName());
219   AliLoader* loader = runLoader->GetLoader("ZDCLoader");
220   if(!loader) {
221     AliError("no ZDC loader found");
222     return;
223   }
224
225   // create the output tree
226   const char* mode = "update";
227   if(runLoader->GetEventNumber() == 0) mode = "recreate";
228   loader->LoadDigits(mode);
229   loader->MakeTree("D");
230   TTree* treeD = loader->TreeD();
231   AliZDCDigit digit;
232   AliZDCDigit* pdigit = &digit;
233   const Int_t kBufferSize = 4000;
234   treeD->Branch("ZDC", "AliZDCDigit", &pdigit, kBufferSize);
235
236   // Create digits
237   Int_t sector[2], sectorL[2];
238   Int_t digi[2], digiL[2], digioot[2];
239   for(sector[0]=1; sector[0]<=3; sector[0]++){
240     for(sector[1]=0; sector[1]<5; sector[1]++){
241         if((sector[0]==3) && ((sector[1]<1) || (sector[1]>2))) continue;
242         for(Int_t res=0; res<2; res++){
243            digi[res] = Phe2ADCch(sector[0], sector[1], pm[sector[0]-1][sector[1]], res) 
244                     + Pedestal(sector[0], sector[1], res);
245         }
246         /*printf("\t DIGIT added -> det = %d, quad = %d - digi[0,1] = [%d, %d]\n",
247              sector[0], sector[1], digi[0], digi[1]); // Chiara debugging!
248         */
249         //
250         new(pdigit) AliZDCDigit(sector, digi);
251         treeD->Fill();
252         //
253         // --- Adding digits for 2nd ZDC set (left side w.r.t. IP) ---
254         // --- they are copied from right ZDC digits
255         if(sector[0]==1 || sector[0]==2){
256            sectorL[0] = sector[0]+3;
257            sectorL[1] = sector[1];
258            for(Int_t res=0; res<2; res++){
259              digiL[res] = Phe2ADCch(sectorL[0], sectorL[1], pm[sector[0]-1][sector[1]], res) 
260                     + Pedestal(sectorL[0], sectorL[1], res);
261            }
262            /*printf("\t DIGIT added -> det = %d, quad = %d - digi[0,1] = [%d, %d]\n",
263                  sectorL[0], sectorL[1], digiL[0], digiL[1]); // Chiara debugging!
264            */
265            //
266            new(pdigit) AliZDCDigit(sectorL, digiL);
267            treeD->Fill();
268         }
269     }
270   } // Loop over detector
271   // Adding in-time digits for 2 reference PTM signals (after signal ch.)
272   // (for the moment the ref. signal is completely invented assuming a PMgain of 5*10^4!)
273   Int_t sectorRef[2];
274   sectorRef[1] = 5;
275   Int_t sigRef[2] = {100, 800};
276   for(Int_t iref=0; iref<2; iref++){
277      sectorRef[0] = 3*iref+1;
278      for(Int_t res=0; res<2; res++){
279        sigRef[res] += Pedestal(sectorRef[0], sectorRef[1], res);
280      }
281      /*printf("\t RefDigit added -> det = %d, quad = %d - digi[0,1] = [%d, %d]\n",
282          sectorRef[0], sectorRef[1], sigRef[0], sigRef[1]); // Chiara debugging!
283      */
284      new(pdigit) AliZDCDigit(sectorRef, sigRef);
285      treeD->Fill();     
286   }
287   //
288   // --- Adding digits for out-of-time channels after signal digits
289   for(sector[0]=1; sector[0]<=3; sector[0]++){
290     for(sector[1]=0; sector[1]<5; sector[1]++){
291         if((sector[0]==3) && ((sector[1]<1) || (sector[1]>2))) continue;
292         for(Int_t res=0; res<2; res++){
293            digioot[res] = Pedestal(sector[0], sector[1], res); // out-of-time ADCs
294         }
295         /*printf("\t DIGIToot added -> det = %d, quad = %d - digi[0,1] = [%d, %d]\n",
296              sector[0], sector[1], digioot[0], digioot[1]); // Chiara debugging!
297         */
298         //
299         new(pdigit) AliZDCDigit(sector, digioot);
300         treeD->Fill();
301         //
302         if(sector[0]==1 || sector[0]==2){
303            sectorL[0] = sector[0]+3;
304            sectorL[1] = sector[1];
305            for(Int_t res=0; res<2; res++){
306              digioot[res] = Pedestal(sectorL[0], sectorL[1], res); // out-of-time ADCs
307            }
308            /*printf("\t DIGIToot added -> det = %d, quad = %d - digi[0,1] = [%d, %d]\n",
309                  sectorL[0], sectorL[1], digioot[0], digioot[1]); // Chiara debugging!
310            */
311            //
312            new(pdigit) AliZDCDigit(sectorL, digioot);
313            treeD->Fill();
314         }
315         //
316     }
317   }
318   // Adding out-of-time digits for 2 reference PTM signals (after out-of-time ch.)
319   Int_t sigRefoot[2];
320   for(Int_t iref=0; iref<2; iref++){
321      sectorRef[0] = 3*iref+1;
322      for(Int_t res=0; res<2; res++){
323        sigRefoot[res] = Pedestal(sectorRef[0], sectorRef[1], res);
324      }
325      /*printf("\t RefDigitoot added -> det = %d, quad = %d - digi[0,1] = [%d, %d]\n",
326          sectorRef[0], sectorRef[1], sigRefoot[0], sigRefoot[1]); // Chiara debugging!
327      */
328      new(pdigit) AliZDCDigit(sectorRef, sigRefoot);
329      treeD->Fill();
330      
331   }
332   //printf("\t AliZDCDigitizer -> TreeD has %d entries\n",(Int_t) treeD->GetEntries());
333
334   // write the output tree
335   loader->WriteDigits("OVERWRITE");
336   loader->UnloadDigits();
337 }
338
339
340 //_____________________________________________________________________________
341 void AliZDCDigitizer::Fragmentation(Float_t impPar, Int_t specN, Int_t specP,
342                                     Int_t &freeSpecN, Int_t &freeSpecP) const
343 {
344 // simulate fragmentation of spectators
345
346   Int_t zz[100], nn[100];
347   AliZDCFragment frag(impPar);
348   for(Int_t j=0; j<=99; j++){
349      zz[j] =0;
350      nn[j] =0;
351   }
352
353   // Fragments generation
354   Int_t nAlpha;
355   frag.GenerateIMF(zz, nAlpha);
356
357   // Attach neutrons
358   Int_t ztot=0;
359   Int_t ntot=0;
360   frag.AttachNeutrons(zz, nn, ztot, ntot);
361   freeSpecN = specN-ntot-2*nAlpha;
362   freeSpecP = specP-ztot-2*nAlpha;
363   // Removing deuterons
364   Int_t ndeu = (Int_t) (frag.DeuteronNumber());
365   freeSpecN -= ndeu;
366   //
367   if(freeSpecN<0) freeSpecN=0;
368   if(freeSpecP<0) freeSpecP=0;
369   AliDebug(2, Form("FreeSpn = %d, FreeSpp = %d", freeSpecN, freeSpecP));
370 }
371
372 //_____________________________________________________________________________
373 void AliZDCDigitizer::SpectatorSignal(Int_t SpecType, Int_t numEvents, 
374                                       Float_t pm[3][5]) const
375 {
376 // add signal of the spectators
377  
378   TFile* file = NULL;
379   if(SpecType == 1) {           // --- Signal for spectator neutrons
380     file = TFile::Open("$ALICE_ROOT/ZDC/ZNsignalntu.root");
381   } else if(SpecType == 2) {    // --- Signal for spectator protons
382     file = TFile::Open("$ALICE_ROOT/ZDC/ZPsignalntu.root");
383   }
384   if(!file || !file->IsOpen()) {
385     AliError("Opening of file failed");
386     return;
387   }
388
389   TNtuple* zdcSignal = (TNtuple*) file->Get("ZDCSignal");
390   Int_t nentries = (Int_t) zdcSignal->GetEntries();
391   
392   Float_t *entry, hitsSpec[7];
393   Int_t pl, i, j, k, iev=0, rnd[125], volume[2];
394   for(pl=0;pl<125;pl++) rnd[pl] = 0;
395   if(numEvents > 125) {
396     AliWarning(Form("numEvents (%d) is larger than 125", numEvents));
397     numEvents = 125;
398   }
399   for(pl=0;pl<numEvents;pl++){
400      rnd[pl] = (Int_t) (9999*gRandom->Rndm());
401      if(rnd[pl] >= 9999) rnd[pl] = 9998;
402      //printf(" rnd[%d] = %d\n",pl,rnd[pl]);     
403   }
404   // Sorting vector in ascending order with C function QSORT 
405   qsort((void*)rnd,numEvents,sizeof(Int_t),comp);
406   do{
407      for(i=0; i<nentries; i++){  
408         zdcSignal->GetEvent(i);
409         entry = zdcSignal->GetArgs();
410         if(entry[0] == rnd[iev]){
411           for(k=0; k<2; k++) volume[k] = (Int_t) entry[k+1];
412           for(j=0; j<7; j++) hitsSpec[j] = entry[j+3];
413           //
414           Float_t lightQ = hitsSpec[4];
415           Float_t lightC = hitsSpec[5];
416           AliDebug(3, Form("SpectatorSignal -> vol = (%d, %d), lightQ = %.0f, lightC = %.0f",
417                            volume[0], volume[1], lightQ, lightC));
418           //printf("\n   Volume = (%d, %d), lightQ = %.0f, lightC = %.0f",
419           //                 volume[0], volume[1], lightQ, lightC);
420           if(volume[0] < 3) {  // ZN or ZP
421             pm[volume[0]-1][0] += lightC;
422             pm[volume[0]-1][volume[1]] += lightQ;
423             //printf("\n   pm[%d][0] = %.0f, pm[%d][%d] = %.0f\n",(volume[0]-1),pm[volume[0]-1][0],
424             //  (volume[0]-1),volume[1],pm[volume[0]-1][volume[1]]);
425           } 
426           else { 
427             if(volume[1] == 1) pm[2][1] += lightC; // ZEM 1
428             else                pm[2][2] += lightQ; // ZEM 2
429             //printf("\n   pm[2][1] = %.0f, pm[2][2] = %.0f\n",pm[2][1],pm[2][2]);
430           }
431         }
432         else if(entry[0] > rnd[iev]){
433           iev++;
434           continue;
435         }
436      }
437   }while(iev<numEvents);
438   
439   file->Close();
440   delete file;
441 }
442
443
444 //_____________________________________________________________________________
445 Int_t AliZDCDigitizer::Phe2ADCch(Int_t Det, Int_t Quad, Float_t Light, 
446                                  Int_t Res) const
447 {
448   // Evaluation of the ADC channel corresponding to the light yield Light
449   Int_t vADCch = (Int_t) (Light * fPMGain[Det-1][Quad] * fADCRes[Res]);
450   //printf("\t Phe2ADCch -> det %d quad %d - phe %.0f  ADC %d\n", Det,Quad,Light,ADCch);
451
452   return vADCch;
453 }
454
455 //_____________________________________________________________________________
456 Int_t AliZDCDigitizer::Pedestal(Int_t Det, Int_t Quad, Int_t Res) const
457 {
458   // Returns a pedestal for detector det, PM quad, channel with res.
459   //
460   Float_t PedValue;
461   
462   // Normal run
463   if(fIsCalibration == 0){
464     Float_t meanPed, Pedwidth;
465     Int_t index=0;
466     if(Quad!=5){
467       if(Det==1 || Det==2)      index = 10*(Det-1)+Quad+5*Res;   // ZN1, ZP1
468       else if(Det==3)           index = 10*(Det-1)+(Quad-1)+Res; // ZEM
469       else if(Det==4 || Det==5) index = 10*(Det-2)+Quad+5*Res+4; // ZN2, ZP2
470     }
471     else index = 10*(Quad-1)+(Det-1)*1/3+2*Res+4; // Reference PMs
472     //
473     meanPed = fCalibData->GetMeanPed(index);
474     Pedwidth = fCalibData->GetMeanPedWidth(index);
475     PedValue = gRandom->Gaus(meanPed,Pedwidth);
476     //
477     /*printf("\t Pedestal -> det = %d, quad = %d, res = %d - Ped[%d] = %d\n",
478         Det, Quad, Res, index,(Int_t) PedValue); // Chiara debugging!
479     */
480   }
481   // To create calibration object
482   else{
483     if(Res == 0) PedValue = gRandom->Gaus((35.+10.*gRandom->Rndm()),(0.5+0.2*gRandom->Rndm())); //High gain
484     else  PedValue = gRandom->Gaus((250.+100.*gRandom->Rndm()),(3.5+2.*gRandom->Rndm())); //Low gain
485   }
486
487   return (Int_t) PedValue;
488 }
489
490 //_____________________________________________________________________________
491 AliCDBStorage* AliZDCDigitizer::SetStorage(const char *uri) 
492 {
493
494   Bool_t deleteManager = kFALSE;
495   
496   AliCDBManager *manager = AliCDBManager::Instance();
497   AliCDBStorage *defstorage = manager->GetDefaultStorage();
498   
499   if(!defstorage || !(defstorage->Contains("ZDC"))){ 
500      AliWarning("No default storage set or default storage doesn't contain ZDC!");
501      manager->SetDefaultStorage(uri);
502      deleteManager = kTRUE;
503   }
504  
505   AliCDBStorage *storage = manager->GetDefaultStorage();
506
507   if(deleteManager){
508     AliCDBManager::Instance()->UnsetDefaultStorage();
509     defstorage = 0;   // the storage is killed by AliCDBManager::Instance()->Destroy()
510   }
511
512   return storage; 
513 }
514
515 //_____________________________________________________________________________
516 AliZDCCalibData* AliZDCDigitizer::GetCalibData() const
517 {
518
519   // Getting calibration object for ZDC set
520
521   AliCDBEntry  *entry = AliCDBManager::Instance()->Get("ZDC/Calib/Data");
522   if(!entry) AliFatal("No calibration data loaded!");  
523
524   AliZDCCalibData *calibdata = dynamic_cast<AliZDCCalibData*>  (entry->GetObject());
525   if(!calibdata)  AliFatal("Wrong calibration object in calibration  file!");
526
527   return calibdata;
528 }
529