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 **************************************************************************/
18 ///_________________________________________________________________________
20 /// This class constructs Digits out of Hits
24 // --- Standard library ---
26 // --- ROOT system ---
30 #include <TGeoManager.h>
31 #include <TGeoPhysicalNode.h>
32 #include <AliGeomManager.h>
35 // --- AliRoot header files ---
36 #include "AliVZEROConst.h"
39 #include "AliVZEROhit.h"
40 #include "AliRunLoader.h"
41 #include "AliLoader.h"
42 #include "AliGRPObject.h"
43 #include "AliRunDigitizer.h"
44 #include "AliCDBManager.h"
45 #include "AliCDBStorage.h"
46 #include "AliCDBEntry.h"
47 #include "AliVZEROCalibData.h"
49 #include "AliVZEROdigit.h"
50 #include "AliVZERODigitizer.h"
52 ClassImp(AliVZERODigitizer)
54 AliVZERODigitizer::AliVZERODigitizer()
56 fCalibData(GetCalibData()),
57 fPhotoCathodeEfficiency(0.18),
59 fPMGain(TMath::Power((fPMVoltage / 112.5) ,7.04277)),
66 // default constructor
71 // fPhotoCathodeEfficiency = 0.18;
72 // fPMVoltage = 768.0;
73 // fPMGain = TMath::Power((fPMVoltage / 112.5) ,7.04277);
75 // fCalibData = GetCalibData();
78 //____________________________________________________________________________
79 AliVZERODigitizer::AliVZERODigitizer(AliRunDigitizer* manager)
80 :AliDigitizer(manager),
81 fCalibData(GetCalibData()),
82 fPhotoCathodeEfficiency(0.18),
84 fPMGain(TMath::Power((fPMVoltage / 112.5) ,7.04277)),
96 // fPhotoCathodeEfficiency = 0.18;
97 // fPMVoltage = 768.0;
98 // fPMGain = TMath::Power( (fPMVoltage / 112.5) ,7.04277 );
100 // fCalibData = GetCalibData();
104 //____________________________________________________________________________
105 AliVZERODigitizer::~AliVZERODigitizer()
116 //_____________________________________________________________________________
117 Bool_t AliVZERODigitizer::Init()
119 // Initialises the digitizer
121 // Initialises the Digit array
122 fDigits = new TClonesArray ("AliVZEROdigit", 1000);
124 // TGeoHMatrix *im = AliGeomManager::GetMatrix("VZERO/V0C");
131 //____________________________________________________________________________
132 void AliVZERODigitizer::Exec(Option_t* /*option*/)
134 // Creates digits from hits
136 Int_t map[80]; // 48 values on V0C + 32 on V0A
137 // Int_t pmNumber[80];
138 Int_t adc[64]; // 32 PMs on V0C + 32 PMs on V0A
139 Float_t time[80], time_ref[80], time2[64];
140 Float_t adc_gain[80];
141 Float_t adc_pedestal[64],adc_sigma[64];
143 Float_t pmGain_smeared[64];
146 // Smearing of the PM tubes intrinsic characteristics
148 for(Int_t ii=0; ii<64; ii++) {
149 pmGain_smeared[ii] = gRandom->Gaus(fPMGain, fPMGain/5); }
151 // Retrieval of ADC gain values and pedestal information from local CDB
152 // I use only the first 64th values of the calibration array in CDB
153 // as I have no beam burst structure - odd or even beam burst number
155 // Reminder : We have 16 scintillating cells mounted on 8 PMs
156 // on Ring 3 and Ring 4 in V0C - added to produce ADC outputs
159 for(Int_t i=0; i<16; i++) {
160 adc_gain[i] = fCalibData->GetGain(i);
161 cPM[i] = fPhotoCathodeEfficiency * pmGain_smeared[i];}
163 for(Int_t j=16; j<48; j=j+2) {
165 adc_gain[j] = fCalibData->GetGain(i);
166 adc_gain[j+1] = fCalibData->GetGain(i);
167 cPM[j] = fPhotoCathodeEfficiency * pmGain_smeared[i];
168 cPM[j+1] = fPhotoCathodeEfficiency * pmGain_smeared[i]; }
170 for(Int_t i=48; i<80; i++){
171 adc_gain[i] = fCalibData->GetGain(i-16);
172 cPM[i] = fPhotoCathodeEfficiency * pmGain_smeared[i-16];};
174 for(Int_t i=0; i<64; i++){ adc_pedestal[i] = fCalibData->GetPedestal(i);
175 adc_sigma[i] = fCalibData->GetSigma(i); };
177 // for(Int_t i=0; i<64; i++) { printf(" i = %d pedestal = %f sigma = %f \n\n",
178 // i, adc_pedestal[i], adc_sigma[i] );}
180 AliRunLoader* outRunLoader =
181 AliRunLoader::GetRunLoader(fManager->GetOutputFolderName());
183 Error("Exec", "Can not get output Run Loader");
186 AliLoader* outLoader = outRunLoader->GetLoader("VZEROLoader");
188 Error("Exec", "Can not get output VZERO Loader");
191 outLoader->LoadDigits("update");
192 if (!outLoader->TreeD()) outLoader->MakeTree("D");
193 outLoader->MakeDigitsContainer();
194 TTree* treeD = outLoader->TreeD();
195 Int_t bufsize = 16000;
196 treeD->Branch("VZERODigit", &fDigits, bufsize);
198 for (Int_t iInput = 0; iInput < fManager->GetNinputs(); iInput++) {
199 AliRunLoader* runLoader =
200 AliRunLoader::GetRunLoader(fManager->GetInputFolderName(iInput));
201 AliLoader* loader = runLoader->GetLoader("VZEROLoader");
203 Error("Exec", "Can not get VZERO Loader for input %d", iInput);
206 if (!runLoader->GetAliRun()) runLoader->LoadgAlice();
208 AliVZERO* vzero = (AliVZERO*) runLoader->GetAliRun()->GetDetector("VZERO");
210 Error("Exec", "No VZERO detector for input %d", iInput);
214 TTree* treeH = loader->TreeH();
216 Error("Exec", "Cannot get TreeH for input %d", iInput);
219 for(Int_t i=0; i<80; i++) {map[i] = 0; time[i] = 0.0;}
221 TClonesArray* hits = vzero->Hits();
223 // Now makes Digits from hits
225 Int_t nTracks = (Int_t) treeH->GetEntries();
226 for (Int_t iTrack = 0; iTrack < nTracks; iTrack++) {
227 for (Int_t i=0; i<80; i++) {time_ref[i] = 999999.0;}
229 treeH->GetEvent(iTrack);
230 Int_t nHits = hits->GetEntriesFast();
231 for (Int_t iHit = 0; iHit < nHits; iHit++) {
232 AliVZEROhit* hit = (AliVZEROhit *)hits->UncheckedAt(iHit);
233 Int_t nPhot = hit->Nphot();
234 Int_t cell = hit->Cell();
236 Float_t dt_scintillator = gRandom->Gaus(0,0.7);
237 Float_t t = dt_scintillator + 1e9*hit->Tof();
239 if(t < time_ref[cell]) time_ref[cell] = t;
240 time[cell] = TMath::Min(t,time_ref[cell]); }
244 loader->UnloadHits();
248 // Now builds the scintillator cell response (80 cells i.e. 80 responses)
250 for (Int_t i=0; i<80; i++) {
251 Float_t q1 = Float_t ( map[i] )* cPM[i] * kQe;
252 Float_t noise = gRandom->Gaus(10.5,3.22);
253 Float_t pmResponse = q1/kC*TMath::Power(ktheta/kthau,1/(1-ktheta/kthau))
255 if(fCollisionMode >0) adc_gain[i] = adc_gain[i]/70.0; // reduce dynamics in Ion Collision Mode
256 map[i] = Int_t( pmResponse * adc_gain[i]);
257 Float_t MIP = 1.0/fCalibData->GetMIPperADC(GetPMNumber(i));
258 if(fCollisionMode >0) MIP=2.0;
259 // printf("cell = %d, ADC = %d, TDC = %f \n",i,map[i], time[i]*10.0 );
260 if(map[i] > (int(( MIP/2 ) + 0.5)) )
261 {map[i] = Int_t(gRandom->Gaus(map[i], (int(( MIP/6 ) + 0.5)) ));}
264 // Now transforms 80 cell responses into 64 photomultiplier responses
265 // Also adds the ADC pedestals taken out of the calibration data base
267 for (Int_t j=0; j<16; j++){
268 adc[j] = static_cast<Int_t>(map [j] + gRandom->Gaus(adc_pedestal[j], adc_sigma[j]));
271 for (Int_t j=48; j<80; j++){
272 adc[j-16] = static_cast<Int_t>(map [j]
273 + gRandom->Gaus(adc_pedestal[j-16],adc_sigma[j-16]));
274 time2[j-16]= time[j]; }
276 for (Int_t j=0; j<16; j++){
277 adc[16+j] = static_cast<Int_t>(map [16+2*j]+ map [16+2*j+1]
278 + gRandom->Gaus(adc_pedestal[16+j], adc_sigma[16+j]));
279 Float_t min_time = TMath::Min(time [16+2*j],time [16+2*j+1]);
280 time2[16+j] = min_time;
281 if(min_time==0.0){time2[16+j]=TMath::Max(time[16+2*j],time[16+2*j+1]);}
285 // Now add digits to the digit Tree
287 for (Int_t i=0; i<64; i++) {
289 // printf(" Event, cell, adc, tof = %d %d %d %f\n",
290 // outRunLoader->GetEventNumber(),i, map[i], time2[i]*10.0);
291 // multiply by 10 to have 100 ps per channel :
292 AddDigit(i, adc[i], Int_t((time2[i]*10.0) +0.5)) ;}
296 outLoader->WriteDigits("OVERWRITE");
297 outLoader->UnloadDigits();
301 //____________________________________________________________________________
302 void AliVZERODigitizer::AddDigit(Int_t PMnumber, Int_t adc, Int_t time)
307 TClonesArray &ldigits = *fDigits;
308 new(ldigits[fNdigits++]) AliVZEROdigit(PMnumber,adc,time);
310 //____________________________________________________________________________
311 void AliVZERODigitizer::ResetDigit()
317 if (fDigits) fDigits->Delete();
320 //____________________________________________________________________________
321 void AliVZERODigitizer::GetCollisionMode()
323 // Retrieves the collision mode from GRP data
325 // Initialization of the GRP entry
327 Int_t run = AliCDBManager::Instance()->GetRun();
329 // printf("\n ++++++ Run Number retrieved as %d \n",run);
331 AliCDBEntry* entry = AliCDBManager::Instance()->Get("GRP/GRP/Data",run);
332 AliGRPObject* grpData = 0x0;
335 TMap* m = dynamic_cast<TMap*>(entry->GetObject()); // old GRP entry
338 grpData = new AliGRPObject();
339 grpData->ReadValuesFromMap(m);
342 grpData = dynamic_cast<AliGRPObject*>(entry->GetObject()); // new GRP entry
345 AliCDBManager::Instance()->UnloadFromCache("GRP/GRP/Data");
348 if(!grpData) AliError("No GRP entry found in OCDB!");
350 // Retrieval of collision mode
352 TString beamType = grpData->GetBeamType();
353 if(beamType==AliGRPObject::GetInvalidString()){
354 AliError("GRP/GRP/Data entry: missing value for the beam type !");
355 AliError("\t VZERO cannot retrieve beam type\n");
359 if( (beamType.CompareTo("P-P") ==0) || (beamType.CompareTo("p-p") ==0) ){
362 else if( (beamType.CompareTo("Pb-Pb") ==0) || (beamType.CompareTo("A-A") ==0) ){
366 fBeamEnergy = grpData->GetBeamEnergy();
367 if(fBeamEnergy==AliGRPObject::GetInvalidFloat()) {
368 AliError("GRP/GRP/Data entry: missing value for the beam energy ! Using 0");
372 // printf("\n ++++++ Beam type and collision mode retrieved as %s %d @ %1.3f GeV ++++++\n\n",beamType.Data(), fCollisionMode, fBeamEnergy);
376 //____________________________________________________________________________
377 AliVZEROCalibData* AliVZERODigitizer::GetCalibData() const
380 AliCDBManager *man = AliCDBManager::Instance();
382 AliCDBEntry *entry=0;
384 entry = man->Get("VZERO/Calib/Data");
387 // AliWarning("Load of calibration data from default storage failed!");
388 // AliWarning("Calibration data will be loaded from local storage ($ALICE_ROOT)");
389 // Int_t runNumber = man->GetRun();
390 // entry = man->GetStorage("local://$ALICE_ROOT/OCDB")
391 // ->Get("VZERO/Calib/Data",runNumber);
395 // Retrieval of data in directory VZERO/Calib/Data:
398 AliVZEROCalibData *calibdata = 0;
400 if (entry) calibdata = (AliVZEROCalibData*) entry->GetObject();
401 if (!calibdata) AliFatal("No calibration data from calibration database !");
407 //____________________________________________________________________________
408 Int_t AliVZERODigitizer::GetPMNumber(Int_t cell) const
412 Int_t pmNumber[80] = { 0, 1, 2, 3, 4, 5, 6, 7,
413 8, 9, 10, 11, 12, 13, 14, 15,
414 16, 16, 17, 17, 18, 18, 19, 19, 20, 20, 21, 21, 22, 22, 23, 23,
415 24, 24, 25, 25, 26, 26, 27, 27, 28, 28, 29, 29, 30, 30, 31, 31,
416 32, 33, 34, 35, 36, 37, 38, 39,
417 40, 41, 42, 43, 44, 45, 46, 47,
418 48, 49, 50, 51, 52, 53, 54, 55,
419 56, 57, 58, 59, 60, 61, 62, 63 };
421 return pmNumber[cell];