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