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