1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
16 ///_________________________________________________________________________
18 /// This class constructs Digits out of Hits
22 // --- Standard library ---
24 // --- ROOT system ---
27 // --- AliRoot header files ---
28 #include "AliVZEROConst.h"
31 #include "AliVZEROhit.h"
32 #include "AliRunLoader.h"
33 #include "AliLoader.h"
34 #include "AliRunDigitizer.h"
35 #include "AliCDBManager.h"
36 #include "AliCDBStorage.h"
37 #include "AliCDBEntry.h"
38 #include "AliVZEROCalibData.h"
40 #include "AliVZEROdigit.h"
41 #include "AliVZERODigitizer.h"
43 ClassImp(AliVZERODigitizer)
45 AliVZERODigitizer::AliVZERODigitizer()
47 // default constructor
52 fPhotoCathodeEfficiency = 0.18;
54 fPMGain = TMath::Power((fPMVoltage / 112.5) ,7.04277);
56 fCalibData = GetCalibData();
59 //____________________________________________________________________________
60 AliVZERODigitizer::AliVZERODigitizer(AliRunDigitizer* manager)
61 :AliDigitizer(manager)
69 fPhotoCathodeEfficiency = 0.18;
71 fPMGain = TMath::Power( (fPMVoltage / 112.5) ,7.04277 );
73 fCalibData = GetCalibData();
77 //____________________________________________________________________________
78 AliVZERODigitizer::~AliVZERODigitizer()
89 //_____________________________________________________________________________
90 Bool_t AliVZERODigitizer::Init()
92 // Initialises the digitizer
94 // Initialises the Digit array
95 fDigits = new TClonesArray ("AliVZEROdigit", 1000);
100 //____________________________________________________________________________
101 void AliVZERODigitizer::Exec(Option_t* /*option*/)
103 // Creates digits from hits
105 Int_t map[80]; // 48 values on V0C + 32 on V0A
106 Int_t adc[64]; // 32 PMs on V0C + 32 PMs on V0A
107 Float_t time[80], time2[64], adc_gain[80];
109 Float_t cPM = fPhotoCathodeEfficiency * fPMGain;
111 // Retrieval of ADC gain values from local CDB
112 // I use only the first 64th values of the calibration array in CDB
113 // as I have no beam burst structure - odd or even beam burst number
115 // Reminder : We have 16 scintillating cells mounted on 8 PMs
116 // on Ring 3 and Ring 4 in V0C - added to produce ADC outputs
119 for(Int_t i=0; i<16; i++) { adc_gain[i] = fCalibData->GetGain(i) ; };
121 for(Int_t j=16; j<48; j=j+2) {
123 adc_gain[j] = fCalibData->GetGain(i) ;
124 adc_gain[j+1] = fCalibData->GetGain(i) ; }
125 for(Int_t i=48; i<80; i++) { adc_gain[i] = fCalibData->GetGain(i-16) ; }
127 // for(Int_t i=0; i<80; i++) { printf(" i = %d gain = %f\n\n", i, adc_gain[i] );}
129 AliRunLoader* outRunLoader =
130 AliRunLoader::GetRunLoader(fManager->GetOutputFolderName());
132 Error("Exec", "Can not get output Run Loader");
135 AliLoader* outLoader = outRunLoader->GetLoader("VZEROLoader");
137 Error("Exec", "Can not get output VZERO Loader");
140 outLoader->LoadDigits("update");
141 if (!outLoader->TreeD()) outLoader->MakeTree("D");
142 outLoader->MakeDigitsContainer();
143 TTree* treeD = outLoader->TreeD();
144 Int_t bufsize = 16000;
145 treeD->Branch("VZERODigit", &fDigits, bufsize);
147 for (Int_t iInput = 0; iInput < fManager->GetNinputs(); iInput++) {
148 AliRunLoader* runLoader =
149 AliRunLoader::GetRunLoader(fManager->GetInputFolderName(iInput));
150 AliLoader* loader = runLoader->GetLoader("VZEROLoader");
152 Error("Exec", "Can not get VZERO Loader for input %d", iInput);
155 if (!runLoader->GetAliRun()) runLoader->LoadgAlice();
157 AliVZERO* vzero = (AliVZERO*) runLoader->GetAliRun()->GetDetector("VZERO");
159 Error("Exec", "No VZERO detector for input %d", iInput);
163 TTree* treeH = loader->TreeH();
165 Error("Exec", "Cannot get TreeH for input %d", iInput);
168 Float_t timeV0 = 1e12;
169 for(Int_t i=0; i<80; i++) { map[i] = 0; time[i] = 0.0; }
171 TClonesArray* hits = vzero->Hits();
173 // Now makes Digits from hits
175 Int_t nTracks = (Int_t) treeH->GetEntries();
176 for (Int_t iTrack = 0; iTrack < nTracks; iTrack++) {
178 treeH->GetEvent(iTrack);
179 Int_t nHits = hits->GetEntriesFast();
180 for (Int_t iHit = 0; iHit < nHits; iHit++) {
181 AliVZEROhit* hit = (AliVZEROhit *)hits->UncheckedAt(iHit);
182 Int_t nPhot = hit->Nphot();
183 Int_t cell = hit->Cell();
185 Float_t dt_scintillator = gRandom->Gaus(0,0.7);
186 time[cell] = dt_scintillator + 1e9*hit->Tof();
187 if(time[cell] < timeV0) timeV0 = time[cell];
191 loader->UnloadHits();
195 // Now builds the scintillator cell response (80 cells i.e. 80 responses)
197 for (Int_t i=0; i<80; i++) {
198 Float_t q1 = Float_t ( map[i] )* cPM * kQe;
199 Float_t noise = gRandom->Gaus(10.5,3.22);
200 Float_t pmResponse = q1/kC*TMath::Power(ktheta/kthau,1/(1-ktheta/kthau))
202 map[i] = Int_t( pmResponse * adc_gain[i]);
206 // Now transforms 80 cell responses into 64 photomultiplier responses
208 for (Int_t j=0; j<32; j++){
212 for (Int_t j=48; j<80; j++){
214 time2[j-16]= time[j];}
216 for (Int_t j=0; j<16; j++){
217 adc[16+j] = map [16 + 2*j]+ map [16 + 2*j + 1];
218 time2[16+j] = TMath::Min(time [16 + 2*j], time [16 + 2*j + 1]);}
220 // Now add digits to the digit Tree
222 for (Int_t i=0; i<64; i++) {
224 // printf(" Event, cell, adc, tof = %d %d %d %f\n",
225 // outRunLoader->GetEventNumber(),i, map[i], time[i]*100.0);
226 // multiply by 10 to have 100 ps per channel :
227 AddDigit(i, adc[i], Int_t(time2[i]*10.0) );
233 outLoader->WriteDigits("OVERWRITE");
234 outLoader->UnloadDigits();
238 //____________________________________________________________________________
239 void AliVZERODigitizer::AddDigit(Int_t PMnumber, Int_t adc, Int_t time)
244 TClonesArray &ldigits = *fDigits;
245 new(ldigits[fNdigits++]) AliVZEROdigit(PMnumber,adc,time);
247 //____________________________________________________________________________
248 void AliVZERODigitizer::ResetDigit()
254 if (fDigits) fDigits->Delete();
257 //____________________________________________________________________________
258 AliVZEROCalibData* AliVZERODigitizer::GetCalibData() const
261 AliCDBManager *man = AliCDBManager::Instance();
263 AliCDBEntry *entry=0;
265 entry = man->Get("VZERO/Calib/Data");
268 AliWarning("Load of calibration data from default storage failed!");
269 AliWarning("Calibration data will be loaded from local storage ($ALICE_ROOT)");
270 Int_t runNumber = man->GetRun();
271 entry = man->GetStorage("local://$ALICE_ROOT")
272 ->Get("VZERO/Calib/Data",runNumber);
276 // Retrieval of data in directory VZERO/Calib/Data:
279 AliVZEROCalibData *calibdata = 0;
281 if (entry) calibdata = (AliVZEROCalibData*) entry->GetObject();
282 if (!calibdata) AliError("No calibration data from calibration database !");