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