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