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>
37 // --- AliRoot header files ---
40 #include "AliVZEROhit.h"
41 #include "AliRunLoader.h"
42 #include "AliLoader.h"
43 #include "AliGRPObject.h"
44 #include "AliDigitizationInput.h"
45 #include "AliCDBManager.h"
46 #include "AliCDBStorage.h"
47 #include "AliCDBEntry.h"
48 #include "AliVZEROCalibData.h"
49 #include "AliCTPTimeParams.h"
50 #include "AliLHCClockPhase.h"
51 #include "AliVZEROdigit.h"
52 #include "AliVZERODigitizer.h"
53 #include "AliVZEROSDigit.h"
55 ClassImp(AliVZERODigitizer)
57 AliVZERODigitizer::AliVZERODigitizer()
59 fCalibData(GetCalibData()),
60 fPhotoCathodeEfficiency(0.18),
65 fSinglePhESpectrum(NULL),
70 // default constructor
71 // Initialize OCDB and containers used in the digitization
76 //____________________________________________________________________________
77 AliVZERODigitizer::AliVZERODigitizer(AliVZERO *vzero, DigiTask_t task)
79 fCalibData(GetCalibData()),
80 fPhotoCathodeEfficiency(0.18),
85 fSinglePhESpectrum(NULL),
91 // Initialize OCDB and containers used in the digitization
96 //____________________________________________________________________________
97 AliVZERODigitizer::AliVZERODigitizer(AliDigitizationInput* digInput)
98 :AliDigitizer(digInput),
99 fCalibData(GetCalibData()),
100 fPhotoCathodeEfficiency(0.18),
105 fSinglePhESpectrum(NULL),
111 // Initialize OCDB and containers used in the digitization
116 //____________________________________________________________________________
117 AliVZERODigitizer::~AliVZERODigitizer()
135 if (fSinglePhESpectrum) {
136 delete fSinglePhESpectrum;
137 fSinglePhESpectrum = NULL;
140 for(Int_t i = 0 ; i < 64; ++i) {
141 if (fTime[i]) delete [] fTime[i];
145 //_____________________________________________________________________________
146 Bool_t AliVZERODigitizer::Init()
148 // Initialises the digitizer
149 // Initialize OCDB and containers used in the digitization
151 // check if the digitizer was already initialized
152 if (fSignalShape) return kTRUE;
154 fSignalShape = new TF1("VZEROSignalShape",this,&AliVZERODigitizer::SignalShape,0,200,6,"AliVZERODigitizer","SignalShape");
155 // fSignalShape->SetParameters(0,1.57345e1,-4.25603e-1,2.9,6.40982,3.69339e-01);
156 // fSignalShape->SetParameters(1.34330e+00,1.13007e+02,-4.95705e-01,
157 // 3.68911e+00,1.01040e+00, 3.94675e-01);
158 fSignalShape->SetParameters(-1.07335e+00,2.16002e+01,-1.26133e-01,
159 1.41619e+00,5.50334e-01,3.86111e-01);
160 fPMResponse = new TF1("VZEROPMResponse",this,&AliVZERODigitizer::PMResponse,-kPMRespTime,2.*kPMRespTime,0,"AliVZERODigitizer","PMResponse");
161 fSinglePhESpectrum = new TF1("VZEROSinglePhESpectrum",this,&AliVZERODigitizer::SinglePhESpectrum,0,20,0,"AliVZERODigitizer","SinglePhESpectrum");
163 // Now get the CTP L0->L1 delay
164 AliCDBEntry *entry = AliCDBManager::Instance()->Get("GRP/CTP/CTPtiming");
165 if (!entry) AliFatal("CTP timing parameters are not found in OCDB !");
166 AliCTPTimeParams *ctpParams = (AliCTPTimeParams*)entry->GetObject();
167 Float_t l1Delay = (Float_t)ctpParams->GetDelayL1L0()*25.0;
169 AliCDBEntry *entry1 = AliCDBManager::Instance()->Get("GRP/CTP/TimeAlign");
170 if (!entry1) AliFatal("CTP time-alignment is not found in OCDB !");
171 AliCTPTimeParams *ctpTimeAlign = (AliCTPTimeParams*)entry1->GetObject();
172 l1Delay += ((Float_t)ctpTimeAlign->GetDelayL1L0()*25.0);
174 AliCDBEntry *entry2 = AliCDBManager::Instance()->Get("VZERO/Calib/TimeDelays");
175 if (!entry2) AliFatal("VZERO time delays are not found in OCDB !");
176 TH1F *delays = (TH1F*)entry2->GetObject();
178 AliCDBEntry *entry3 = AliCDBManager::Instance()->Get("GRP/Calib/LHCClockPhase");
179 if (!entry3) AliFatal("LHC clock-phase shift is not found in OCDB !");
180 AliLHCClockPhase *phase = (AliLHCClockPhase*)entry3->GetObject();
182 for(Int_t i = 0 ; i < 64; ++i) {
184 for(Int_t j = 0; j < kNClocks; ++j) fAdc[i][j] = 0;
185 fLeadingTime[i] = fTimeWidth[i] = 0;
187 fPmGain[i] = fCalibData->GetGain(i);
189 fAdcPedestal[i][0] = fCalibData->GetPedestal(i);
190 fAdcSigma[i][0] = fCalibData->GetSigma(i);
191 fAdcPedestal[i][1] = fCalibData->GetPedestal(i+64);
192 fAdcSigma[i][1] = fCalibData->GetSigma(i+64);
194 Int_t board = AliVZEROCalibData::GetBoardNumber(i);
195 fNBins[i] = TMath::Nint(((Float_t)(fCalibData->GetMatchWindow(board)+1)*25.0+
196 (Float_t)kMaxTDCWidth*fCalibData->GetWidthResolution(board))/
197 fCalibData->GetTimeResolution(board));
198 fNBinsLT[i] = TMath::Nint(((Float_t)(fCalibData->GetMatchWindow(board)+1)*25.0)/
199 fCalibData->GetTimeResolution(board));
200 fBinSize[i] = fCalibData->GetTimeResolution(board);
201 fHptdcOffset[i] = (((Float_t)fCalibData->GetRollOver(board)-
202 (Float_t)fCalibData->GetTriggerCountOffset(board))*25.0+
203 fCalibData->GetTimeOffset(i)-
205 phase->GetMeanPhase()+
206 delays->GetBinContent(i+1)+
208 fClockOffset[i] = (((Float_t)fCalibData->GetRollOver(board)-
209 (Float_t)fCalibData->GetTriggerCountOffset(board))*25.0+
210 fCalibData->GetTimeOffset(i)-
214 fTime[i] = new Float_t[fNBins[i]];
215 memset(fTime[i],0,fNBins[i]*sizeof(Float_t));
221 //____________________________________________________________________________
222 void AliVZERODigitizer::Digitize(Option_t* /*option*/)
224 // Creates digits from hits
227 if (fVZERO && !fDigInput) {
228 AliLoader *loader = fVZERO->GetLoader();
230 AliError("Can not get VZERO Loader via AliVZERO object!");
233 AliRunLoader* runLoader = AliRunLoader::Instance();
234 for (Int_t iEvent = 0; iEvent < runLoader->GetNumberOfEvents(); ++iEvent) {
235 runLoader->GetEvent(iEvent);
236 if (fTask == kHits2Digits) {
243 WriteSDigits(loader);
247 else if (fDigInput) {
250 AliRunLoader *currentLoader = AliRunLoader::GetRunLoader(fDigInput->GetOutputFolderName());
251 AliLoader *loader = currentLoader->GetLoader("VZEROLoader");
253 AliError("Cannot get VZERO Loader via RunDigitizer!");
259 AliFatal("Invalid digitization task! Exiting!");
263 //____________________________________________________________________________
264 void AliVZERODigitizer::AddDigit(Int_t pmnumber, Float_t time, Float_t width, Bool_t integrator, Short_t *chargeADC, Int_t *labels)
269 TClonesArray &ldigits = *fDigits;
271 new(ldigits[fNdigits++]) AliVZEROdigit(pmnumber,time,width,integrator,chargeADC,labels);
274 //____________________________________________________________________________
275 void AliVZERODigitizer::AddSDigit(Int_t pmnumber, Int_t nbins, Float_t *charges, Int_t *labels)
280 TClonesArray &ldigits = *fDigits;
282 new(ldigits[fNdigits++]) AliVZEROSDigit(pmnumber,nbins,charges,labels);
285 //____________________________________________________________________________
286 void AliVZERODigitizer::ResetDigits()
292 if (fDigits) fDigits->Clear();
295 //____________________________________________________________________________
296 AliVZEROCalibData* AliVZERODigitizer::GetCalibData() const
299 AliCDBManager *man = AliCDBManager::Instance();
301 AliCDBEntry *entry=0;
303 entry = man->Get("VZERO/Calib/Data");
306 // AliWarning("Load of calibration data from default storage failed!");
307 // AliWarning("Calibration data will be loaded from local storage ($ALICE_ROOT)");
308 // Int_t runNumber = man->GetRun();
309 // entry = man->GetStorage("local://$ALICE_ROOT/OCDB")
310 // ->Get("VZERO/Calib/Data",runNumber);
314 // Retrieval of data in directory VZERO/Calib/Data:
317 AliVZEROCalibData *calibdata = 0;
319 if (entry) calibdata = (AliVZEROCalibData*) entry->GetObject();
320 if (!calibdata) AliFatal("No calibration data from calibration database !");
326 double AliVZERODigitizer::SignalShape(double *x, double *par)
328 // this function simulates the time
329 // of arrival of the photons at the
332 if (xx <= par[0]) return 0;
333 Double_t a = 1./TMath::Power((xx-par[0])/par[1],1./par[2]);
334 if (xx <= par[3]) return a;
335 Double_t b = 1./TMath::Power((xx-par[3])/par[4],1./par[5]);
336 Double_t f = a*b/(a+b);
337 AliDebug(100,Form("x=%f func=%f",xx,f));
341 double AliVZERODigitizer::PMResponse(double *x, double * /* par */)
343 // this function describes the
344 // PM time response to a single
346 Double_t xx = x[0]+kPMRespTime;
347 return xx*xx*TMath::Exp(-xx*xx/(kPMRespTime*kPMRespTime));
350 double AliVZERODigitizer::SinglePhESpectrum(double *x, double * /* par */)
352 // this function describes the
353 // PM amplitude response to a single
356 if (xx < 0) return 0;
357 return (TMath::Poisson(xx,kPMNbOfSecElec)+kPMTransparency*TMath::Poisson(xx,1.0));
360 Int_t AliVZERODigitizer::Cell2Pmt(Int_t cell) const
362 // The method maps the scintillator
363 // indexes to the PM ones
364 if (cell < 0 || cell >= 80) {
365 AliError(Form("Wrong VZERO cell index %d",cell));
368 if (cell < 16) return cell;
369 if (cell < 48) return 8 + cell/2;
373 void AliVZERODigitizer::DigitizeHits()
375 // Digitize the hits to the level of
376 // SDigits (fTime arrays)
378 for(Int_t i = 0 ; i < 64; ++i) {
379 memset(fTime[i],0,fNBins[i]*sizeof(Float_t));
380 fLabels[i][0] = fLabels[i][1] = fLabels[i][2] = -1;
382 Float_t integral = fPMResponse->Integral(-kPMRespTime,2.*kPMRespTime);
383 Float_t meansPhE = fSinglePhESpectrum->Mean(0,20);
385 AliLoader* loader = fVZERO->GetLoader();
387 AliError("Can not get VZERO Loader!");
391 TTree* treeH = loader->TreeH();
393 AliError("Cannot get TreeH!");
396 TClonesArray* hits = fVZERO->Hits();
398 // Now makes Digits from hits
399 Int_t nTracks = (Int_t) treeH->GetEntries();
400 for (Int_t iTrack = 0; iTrack < nTracks; iTrack++) {
402 treeH->GetEvent(iTrack);
403 Int_t nHits = hits->GetEntriesFast();
404 for (Int_t iHit = 0; iHit < nHits; iHit++) {
405 AliVZEROhit* hit = (AliVZEROhit *)hits->UncheckedAt(iHit);
406 Int_t nPhot = hit->Nphot();
407 Int_t cell = hit->Cell();
408 Int_t pmt = Cell2Pmt(cell);
409 if (pmt < 0) continue;
410 Int_t trackLabel = hit->GetTrack();
411 for(Int_t l = 0; l < 3; ++l) {
412 if (fLabels[pmt][l] < 0) {
413 fLabels[pmt][l] = trackLabel;
417 Float_t dt_scintillator = gRandom->Gaus(0,kIntTimeRes);
418 Float_t t = dt_scintillator + 1e9*hit->Tof();
419 if (pmt < 32) t += kV0CDelayCables;
420 t += fHptdcOffset[pmt];
422 Float_t prob = fCalibData->GetLightYields(pmt)*fPhotoCathodeEfficiency; // Optical losses included!
424 nPhE = (Int_t)gRandom->Gaus(prob*Float_t(nPhot)+0.5,
425 sqrt(Float_t(nPhot)*prob*(1.-prob)));
427 nPhE = gRandom->Binomial(nPhot,prob);
428 Float_t charge = TMath::Qe()*fPmGain[pmt]*fBinSize[pmt]/integral;
429 for (Int_t iPhE = 0; iPhE < nPhE; ++iPhE) {
430 Float_t tPhE = t + fSignalShape->GetRandom(0,fBinSize[pmt]*Float_t(fNBins[pmt]));
431 Float_t gainVar = fSinglePhESpectrum->GetRandom(0,20)/meansPhE;
432 Int_t firstBin = TMath::Max(0,(Int_t)((tPhE-kPMRespTime)/fBinSize[pmt]));
433 Int_t lastBin = TMath::Min(fNBins[pmt]-1,(Int_t)((tPhE+2.*kPMRespTime)/fBinSize[pmt]));
434 for(Int_t iBin = firstBin; iBin <= lastBin; ++iBin) {
435 Float_t tempT = fBinSize[pmt]*(0.5+iBin)-tPhE;
436 fTime[pmt][iBin] += gainVar*charge*fPMResponse->Eval(tempT);
441 loader->UnloadHits();
445 void AliVZERODigitizer::DigitizeSDigits()
447 // Digitize the fTime arrays (SDigits) to the level of
448 // Digits (fAdc arrays)
449 for(Int_t i = 0 ; i < 64; ++i) {
450 for(Int_t j = 0; j < kNClocks; ++j) fAdc[i][j] = 0;
451 fLeadingTime[i] = fTimeWidth[i] = 0;
454 Float_t maximum = 0.9*fSignalShape->GetMaximum(0,200); // Not exact, one needs to do this on the convoluted
455 Float_t integral2 = fSignalShape->Integral(0,200); // function. Anyway the effect is small <10% on the 2.5 ADC thr
456 for (Int_t ipmt = 0; ipmt < 64; ++ipmt) {
457 Float_t thr = fCalibData->GetCalibDiscriThr(ipmt,kFALSE)*kChargePerADC*maximum*fBinSize[ipmt]/integral2;
458 Bool_t ltFound = kFALSE, ttFound = kFALSE;
459 for (Int_t iBin = 0; iBin < fNBins[ipmt]; ++iBin) {
460 Float_t t = fBinSize[ipmt]*Float_t(iBin);
461 if (fTime[ipmt][iBin] > thr) {
462 if (!ltFound && (iBin < fNBinsLT[ipmt])) {
464 fLeadingTime[ipmt] = t;
471 fTimeWidth[ipmt] = t - fLeadingTime[ipmt];
475 Float_t tadc = t - fClockOffset[ipmt];
476 Int_t clock = kNClocks/2 - Int_t(tadc/25.0);
477 if (clock >= 0 && clock < kNClocks)
478 fAdc[ipmt][clock] += fTime[ipmt][iBin]/kChargePerADC;
480 AliDebug(1,Form("Channel %d Offset %f Time %f",ipmt,fClockOffset[ipmt],fLeadingTime[ipmt]));
481 Int_t board = AliVZEROCalibData::GetBoardNumber(ipmt);
482 if (ltFound && ttFound) {
483 fTimeWidth[ipmt] = fCalibData->GetWidthResolution(board)*
484 Float_t(Int_t(fTimeWidth[ipmt]/fCalibData->GetWidthResolution(board)));
485 if (fTimeWidth[ipmt] < Float_t(kMinTDCWidth)*fCalibData->GetWidthResolution(board))
486 fTimeWidth[ipmt] = Float_t(kMinTDCWidth)*fCalibData->GetWidthResolution(board);
487 if (fTimeWidth[ipmt] > Float_t(kMaxTDCWidth)*fCalibData->GetWidthResolution(board))
488 fTimeWidth[ipmt] = Float_t(kMaxTDCWidth)*fCalibData->GetWidthResolution(board);
492 fEvenOrOdd = gRandom->Integer(2);
493 for (Int_t j=0; j<64; ++j){
494 for (Int_t iClock = 0; iClock < kNClocks; ++iClock) {
495 Int_t integrator = (iClock + fEvenOrOdd) % 2;
496 AliDebug(1,Form("ADC %d %d %f",j,iClock,fAdc[j][iClock]));
497 fAdc[j][iClock] += gRandom->Gaus(fAdcPedestal[j][integrator], fAdcSigma[j][integrator]);
503 void AliVZERODigitizer::WriteDigits(AliLoader *loader)
505 // Take fAdc arrays filled by the previous
506 // method and produce and add digits to the digit Tree
508 loader->LoadDigits("UPDATE");
510 if (!loader->TreeD()) loader->MakeTree("D");
511 loader->MakeDigitsContainer();
512 TTree* treeD = loader->TreeD();
514 treeD->Branch("VZERODigit", &fDigits);
516 Short_t *chargeADC = new Short_t[kNClocks];
517 for (Int_t i=0; i<64; i++) {
518 for (Int_t j = 0; j < kNClocks; ++j) {
519 Int_t tempadc = Int_t(fAdc[i][j]);
520 if (tempadc > 1023) tempadc = 1023;
521 chargeADC[j] = tempadc;
523 AddDigit(i, fLeadingTime[i], fTimeWidth[i], Bool_t((10+fEvenOrOdd)%2), chargeADC, fLabels[i]);
528 loader->WriteDigits("OVERWRITE");
529 loader->UnloadDigits();
533 void AliVZERODigitizer::WriteSDigits(AliLoader *loader)
535 // Take fTime arrays filled by the previous
536 // method and produce and add sdigits to the sdigit Tree
538 loader->LoadSDigits("UPDATE");
540 if (!loader->TreeS()) loader->MakeTree("S");
541 loader->MakeSDigitsContainer();
542 TTree* treeS = loader->TreeS();
544 treeS->Branch("VZEROSDigit", &fDigits);
546 for (Int_t ipmt = 0; ipmt < 64; ++ipmt) {
547 AddSDigit(ipmt,fNBins[ipmt],fTime[ipmt],fLabels[ipmt]);
551 loader->WriteSDigits("OVERWRITE");
552 loader->UnloadSDigits();
556 void AliVZERODigitizer::ReadSDigits()
558 // Read SDigits which are then to precessed
559 // in the following method
560 for(Int_t i = 0 ; i < 64; ++i) {
561 memset(fTime[i],0,fNBins[i]*sizeof(Float_t));
562 fLabels[i][0] = fLabels[i][1] = fLabels[i][2] = -1;
565 // Loop over input files
566 Int_t nFiles= fDigInput->GetNinputs();
567 for (Int_t inputFile = 0; inputFile < nFiles; inputFile++) {
568 // Get the current loader
569 AliRunLoader* currentLoader =
570 AliRunLoader::GetRunLoader(fDigInput->GetInputFolderName(inputFile));
572 AliLoader *loader = currentLoader->GetLoader("VZEROLoader");
573 loader->LoadSDigits("READ");
575 // Get the tree of summable digits
576 TTree* sdigitsTree = loader->TreeS();
578 AliError("No sdigit tree from digInput");
583 TBranch* sdigitsBranch = sdigitsTree->GetBranch("VZEROSDigit");
584 if (!sdigitsBranch) {
585 AliError("Failed to get sdigit branch");
589 // Set the branch address
590 TClonesArray *sdigitsArray = NULL;
591 sdigitsBranch->SetAddress(&sdigitsArray);
593 // Sum contributions from the sdigits
594 // Get number of entries in the tree
595 Int_t nentries = Int_t(sdigitsBranch->GetEntries());
596 for (Int_t entry = 0; entry < nentries; ++entry) {
597 sdigitsBranch->GetEntry(entry);
598 // Get the number of sdigits
599 Int_t nsdigits = sdigitsArray->GetEntries();
600 for (Int_t sdigit = 0; sdigit < nsdigits; sdigit++) {
601 AliVZEROSDigit* sDigit = static_cast<AliVZEROSDigit*>(sdigitsArray->UncheckedAt(sdigit));
602 Int_t pmNumber = sDigit->PMNumber();
603 Int_t nbins = sDigit->GetNBins();
604 if (nbins != fNBins[pmNumber]) {
605 AliError(Form("Incompatible number of bins between digitizer (%d) and sdigit (%d) for PM %d! Skipping sdigit!",
606 fNBins[pmNumber],nbins,pmNumber));
610 Float_t *charges = sDigit->GetCharges();
611 for(Int_t iBin = 0; iBin < nbins; ++iBin) fTime[pmNumber][iBin] += charges[iBin];
613 Int_t *labels = sDigit->GetTracks();
615 for(Int_t i = 0; i < 3; ++i) {
616 if (fLabels[pmNumber][i] < 0) {
617 if (labels[j] < 0) break;
618 fLabels[pmNumber][i] = labels[j];
624 loader->UnloadSDigits();
628 //____________________________________________________________________
630 AliVZERODigitizer::DigitsArray()
632 // Initialize digit array if not already and
633 // return pointer to it.
635 fDigits = new TClonesArray("AliVZEROdigit", 64);
641 //____________________________________________________________________
643 AliVZERODigitizer::SDigitsArray()
645 // Initialize sdigit array if not already and
646 // return pointer to it.
648 fDigits = new TClonesArray("AliVZEROSDigit", 64);