Reconstruction of ESD data implemented
[u/mrichter/AliRoot.git] / VZERO / AliVZERODigitizer.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 ///_________________________________________________________________________
17 ///
18 /// This class constructs Digits out of Hits
19 ///
20 ///
21
22 // --- Standard library ---
23
24 // --- ROOT system ---
25 #include <TTree.h>
26 #include <TRandom.h>
27
28 // --- AliRoot header files ---
29 #include "AliVZEROConst.h"
30 #include "AliRun.h"
31 #include "AliVZERO.h"
32 #include "AliVZEROhit.h"
33 #include "AliRunLoader.h"
34 #include "AliLoader.h"
35 #include "AliRunDigitizer.h"
36 #include "AliCDBManager.h"
37 #include "AliCDBStorage.h"
38 #include "AliCDBEntry.h"
39 #include "AliVZEROCalibData.h"
40
41 #include "AliVZEROdigit.h"
42 #include "AliVZERODigitizer.h"
43
44 ClassImp(AliVZERODigitizer)
45
46  AliVZERODigitizer::AliVZERODigitizer()
47    :AliDigitizer(),
48     fCalibData(GetCalibData()),
49     fPhotoCathodeEfficiency(0.18),
50     fPMVoltage(768.0),
51     fPMGain(TMath::Power((fPMVoltage / 112.5) ,7.04277)),
52     fNdigits(0),
53     fDigits(0)
54    
55 {
56   // default constructor
57
58 //    fNdigits = 0;
59 //    fDigits  = 0;
60 //   
61 //    fPhotoCathodeEfficiency =   0.18;
62 //    fPMVoltage              =  768.0;
63 //    fPMGain = TMath::Power((fPMVoltage / 112.5) ,7.04277); 
64    
65 //   fCalibData = GetCalibData();
66 }
67
68 //____________________________________________________________________________ 
69   AliVZERODigitizer::AliVZERODigitizer(AliRunDigitizer* manager)
70                     :AliDigitizer(manager),
71                      fCalibData(GetCalibData()),
72                      fPhotoCathodeEfficiency(0.18),
73                      fPMVoltage(768.0),
74                      fPMGain(TMath::Power((fPMVoltage / 112.5) ,7.04277)),
75                      fNdigits(0),
76                      fDigits(0)
77                                         
78 {
79   // constructor
80   
81 //   fNdigits = 0;
82 //   fDigits  = 0;
83 //   
84 //   fPhotoCathodeEfficiency =   0.18;
85 //   fPMVoltage              =  768.0;
86 //   fPMGain = TMath::Power( (fPMVoltage / 112.5) ,7.04277 );
87   
88 //  fCalibData = GetCalibData();
89   
90 }
91            
92 //____________________________________________________________________________ 
93   AliVZERODigitizer::~AliVZERODigitizer()
94 {
95   // destructor
96   
97   if (fDigits) {
98     fDigits->Delete();
99     delete fDigits;
100     fDigits=0; 
101   }
102 }
103
104 //_____________________________________________________________________________
105 Bool_t AliVZERODigitizer::Init()
106 {
107   // Initialises the digitizer
108
109   // Initialises the Digit array
110   fDigits = new TClonesArray ("AliVZEROdigit", 1000);
111
112   return kTRUE;
113 }
114
115 //____________________________________________________________________________
116 void AliVZERODigitizer::Exec(Option_t* /*option*/) 
117 {   
118   // Creates digits from hits
119      
120   Int_t       map[80];    // 48 values on V0C + 32 on V0A
121   Int_t       adc[64];    // 32 PMs on V0C + 32 PMs on V0A
122   Float_t    time[80], time_ref[80], time2[64];
123   Float_t    adc_gain[80];    
124   fNdigits     =    0;  
125   Float_t cPM  = fPhotoCathodeEfficiency * fPMGain;
126
127   // Retrieval of ADC gain values from local CDB 
128   // I use only the first 64th values of the calibration array in CDB 
129   // as I have no beam burst structure - odd or even beam burst number
130   
131   // Reminder : We have 16 scintillating cells mounted on 8 PMs 
132   // on Ring 3 and Ring 4 in V0C -  added to produce  ADC outputs 
133   // on these rings... 
134    
135   for(Int_t i=0; i<16; i++) { adc_gain[i]  = fCalibData->GetGain(i) ; };
136   
137   for(Int_t j=16; j<48; j=j+2) { 
138             Int_t i=(j+17)/2;
139             adc_gain[j]   = fCalibData->GetGain(i) ; 
140             adc_gain[j+1] = fCalibData->GetGain(i) ; }
141   for(Int_t i=48; i<80; i++) { adc_gain[i]  = fCalibData->GetGain(i-16) ; }
142   
143 //  for(Int_t i=0; i<80; i++) { printf(" i = %d gain = %f\n\n", i, adc_gain[i] );} 
144             
145   AliRunLoader* outRunLoader = 
146     AliRunLoader::GetRunLoader(fManager->GetOutputFolderName());    
147   if (!outRunLoader) {
148     Error("Exec", "Can not get output Run Loader");
149     return;}
150     
151   AliLoader* outLoader = outRunLoader->GetLoader("VZEROLoader");
152   if (!outLoader) {
153     Error("Exec", "Can not get output VZERO Loader");
154     return;}
155
156   outLoader->LoadDigits("update");
157   if (!outLoader->TreeD()) outLoader->MakeTree("D");
158   outLoader->MakeDigitsContainer();
159   TTree* treeD  = outLoader->TreeD();
160   Int_t bufsize = 16000;
161   treeD->Branch("VZERODigit", &fDigits, bufsize); 
162   
163   for (Int_t iInput = 0; iInput < fManager->GetNinputs(); iInput++) {
164     AliRunLoader* runLoader = 
165       AliRunLoader::GetRunLoader(fManager->GetInputFolderName(iInput));
166     AliLoader* loader = runLoader->GetLoader("VZEROLoader");
167     if (!loader) {
168       Error("Exec", "Can not get VZERO Loader for input %d", iInput);
169       continue;}
170       
171     if (!runLoader->GetAliRun()) runLoader->LoadgAlice();
172
173     AliVZERO* vzero = (AliVZERO*) runLoader->GetAliRun()->GetDetector("VZERO");
174     if (!vzero) {
175       Error("Exec", "No VZERO detector for input %d", iInput);
176       continue;}
177       
178     loader->LoadHits();
179     TTree* treeH = loader->TreeH();
180     if (!treeH) {
181       Error("Exec", "Cannot get TreeH for input %d", iInput);
182       continue; }
183        
184     for(Int_t i=0; i<80; i++) {map[i] = 0; time[i] = 0.0;}
185               
186     TClonesArray* hits = vzero->Hits();
187              
188 //  Now makes Digits from hits
189          
190     Int_t nTracks = (Int_t) treeH->GetEntries();
191     for (Int_t iTrack = 0; iTrack < nTracks; iTrack++) {
192       for (Int_t i=0; i<80; i++) {time_ref[i] = 999999.0;}   
193       vzero->ResetHits();
194       treeH->GetEvent(iTrack);
195       Int_t nHits = hits->GetEntriesFast();
196       for (Int_t iHit = 0; iHit < nHits; iHit++) {
197         AliVZEROhit* hit = (AliVZEROhit *)hits->UncheckedAt(iHit);
198         Int_t nPhot = hit->Nphot();
199         Int_t cell  = hit->Cell();                                    
200         map[cell] += nPhot;
201         Float_t dt_scintillator = gRandom->Gaus(0,0.7);
202         Float_t t = dt_scintillator + 1e9*hit->Tof();
203         if (t > 0.0) {
204            if(t < time_ref[cell]) time_ref[cell] = t;
205            time[cell] = TMath::Min(t,time_ref[cell]);
206         }
207       }           // hit   loop      
208     }             // track loop
209
210     loader->UnloadHits();
211
212   }               // input loop
213
214 // Now builds the scintillator cell response (80 cells i.e. 80 responses)
215          
216    for (Int_t i=0; i<80; i++) {    
217        Float_t q1 = Float_t ( map[i] )* cPM * kQe;
218        Float_t noise = gRandom->Gaus(10.5,3.22);
219        Float_t pmResponse  =  q1/kC*TMath::Power(ktheta/kthau,1/(1-ktheta/kthau)) 
220         + noise*1e-3;
221        map[i] = Int_t( pmResponse * adc_gain[i]);
222      }
223       
224 // Now transforms 80 cell responses into 64 photomultiplier responses 
225         
226    for (Int_t j=0; j<32; j++){
227         adc[j]  = map [j];
228         time2[j]= time[j];}
229         
230    for (Int_t j=48; j<80; j++){
231         adc[j-16]  = map [j];
232         time2[j-16]= time[j];}
233         
234    for (Int_t j=0; j<16; j++){
235         adc[16+j] = map [16+2*j]+ map [16+2*j+1];
236         Float_t min_time = TMath::Min(time [16+2*j],time [16+2*j+1]);
237         time2[16+j] = min_time;
238         if(min_time==0.0) {
239            time2[16+j] = TMath::Max(time[16+2*j],time[16+2*j+1]);
240         }
241    }
242         
243
244 // Now add digits to the digit Tree 
245         
246    for (Int_t i=0; i<64; i++) {      
247       if(adc[i] > 0) {
248 //      printf(" Event, cell, adc, tof = %d %d %d %f\n", 
249 //                 outRunLoader->GetEventNumber(),i, map[i], time2[i]*10.0);
250 //   multiply by 10 to have 100 ps per channel :
251       AddDigit(i, adc[i], Int_t(time2[i]*10.0) );
252      } 
253    }
254     
255   treeD->Fill();
256   outLoader->WriteDigits("OVERWRITE");  
257   outLoader->UnloadDigits();     
258   ResetDigit();
259 }
260
261 //____________________________________________________________________________
262 void AliVZERODigitizer::AddDigit(Int_t PMnumber, Int_t adc, Int_t time) 
263  { 
264  
265 // Adds Digit 
266  
267   TClonesArray &ldigits = *fDigits;  
268   new(ldigits[fNdigits++]) AliVZEROdigit(PMnumber,adc,time);
269 }
270 //____________________________________________________________________________
271 void AliVZERODigitizer::ResetDigit()
272 {
273 //
274 // Clears Digits
275 //
276   fNdigits = 0;
277   if (fDigits) fDigits->Delete();
278 }
279
280 //____________________________________________________________________________
281 AliVZEROCalibData* AliVZERODigitizer::GetCalibData() const
282
283 {
284   AliCDBManager *man = AliCDBManager::Instance();
285
286   AliCDBEntry *entry=0;
287
288   entry = man->Get("VZERO/Calib/Data");
289
290   if(!entry){
291     AliWarning("Load of calibration data from default storage failed!");
292     AliWarning("Calibration data will be loaded from local storage ($ALICE_ROOT)");
293     Int_t runNumber = man->GetRun();
294     entry = man->GetStorage("local://$ALICE_ROOT")
295       ->Get("VZERO/Calib/Data",runNumber);
296         
297   }
298
299   // Retrieval of data in directory VZERO/Calib/Data:
300
301
302   AliVZEROCalibData *calibdata = 0;
303
304   if (entry) calibdata = (AliVZEROCalibData*) entry->GetObject();
305   if (!calibdata)  AliError("No calibration data from calibration database !");
306
307   return calibdata;
308
309 }
310