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 /* $Id: AliADDigitizer.cxx $ */
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 ---
39 #include "AliDetector.h"
42 #include "AliADConst.h"
43 #include "AliRunLoader.h"
44 #include "AliLoader.h"
45 #include "AliGRPObject.h"
46 #include "AliDigitizationInput.h"
47 #include "AliCDBManager.h"
48 #include "AliCDBStorage.h"
49 #include "AliCDBEntry.h"
50 #include "AliADCalibData.h"
51 #include "AliCTPTimeParams.h"
52 #include "AliLHCClockPhase.h"
53 #include "AliADdigit.h"
54 #include "AliADDigitizer.h"
55 #include "AliADSDigit.h"
57 ClassImp(AliADDigitizer)
59 //____________________________________________________________________________
60 AliADDigitizer::AliADDigitizer()
62 fCalibData(GetCalibData()),
67 fSinglePhESpectrum(NULL),
72 // default constructor
73 // Initialize OCDB and containers used in the digitization
78 //____________________________________________________________________________
79 AliADDigitizer::AliADDigitizer(AliAD *AD, DigiTask_t task)
81 fCalibData(GetCalibData()),
86 fSinglePhESpectrum(NULL),
92 // Initialize OCDB and containers used in the digitization
97 //____________________________________________________________________________
98 AliADDigitizer::AliADDigitizer(AliDigitizationInput* digInput)
99 :AliDigitizer(digInput),
100 fCalibData(GetCalibData()),
105 fSinglePhESpectrum(NULL),
111 // Initialize OCDB and containers used in the digitization
116 //____________________________________________________________________________
117 AliADDigitizer::~AliADDigitizer()
135 if (fSinglePhESpectrum) {
136 delete fSinglePhESpectrum;
137 fSinglePhESpectrum = NULL;
140 for(Int_t i = 0 ; i < 16; ++i) {
141 if (fTime[i]) delete [] fTime[i];
145 //____________________________________________________________________________
146 Bool_t AliADDigitizer::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("ADSignalShape",this,&AliADDigitizer::SignalShape,0,200,6,"AliADDigitizer","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("ADPMResponse",this,&AliADDigitizer::PMResponse,-kPMRespTime,2.*kPMRespTime,0,"AliADDigitizer","PMResponse");
161 fSinglePhESpectrum = new TF1("ADSinglePhESpectrum",this,&AliADDigitizer::SinglePhESpectrum,0,20,0,"AliADDigitizer","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("AD 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 < 16; ++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+16);
192 fAdcSigma[i][1] = fCalibData->GetSigma(i+16);
194 Int_t board = AliADCalibData::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));
217 // AliWarning(Form("PMT %d,PM gain %f, fNBins %d, TimeBinSize %f,",i, fPmGain[i], fNBins[i],fBinSize[i]));
225 //____________________________________________________________________________
226 void AliADDigitizer::Digitize(Option_t* /*option*/)
228 // Creates digits from hits
231 if (fAD && !fDigInput) {
232 AliLoader *loader = fAD->GetLoader();
234 AliError("Can not get AD Loader via AliAD object!");
237 AliRunLoader* runLoader = AliRunLoader::Instance();
238 for (Int_t iEvent = 0; iEvent < runLoader->GetNumberOfEvents(); ++iEvent) {
239 runLoader->GetEvent(iEvent);
240 if (fTask == kHits2Digits) {
247 WriteSDigits(loader);
251 else if (fDigInput) {
254 AliRunLoader *currentLoader = AliRunLoader::GetRunLoader(fDigInput->GetOutputFolderName());
255 AliLoader *loader = currentLoader->GetLoader("ADLoader");
257 AliError("Cannot get AD Loader via RunDigitizer!");
263 AliFatal("Invalid digitization task! Exiting!");
267 //____________________________________________________________________________
268 void AliADDigitizer::DigitizeHits()
270 // Digitize the hits to the level of
271 // SDigits (fTime arrays)
273 for(Int_t i = 0 ; i < 16; ++i) {
274 memset(fTime[i],0,fNBins[i]*sizeof(Float_t));
275 fLabels[i][0] = fLabels[i][1] = fLabels[i][2] = -1;
277 Float_t integral = fPMResponse->Integral(-kPMRespTime,2.*kPMRespTime);
278 Float_t meansPhE = fSinglePhESpectrum->Mean(0,20);
280 AliLoader* loader = fAD->GetLoader();
282 AliError("Can not get AD Loader!");
286 TTree* treeH = loader->TreeH();
288 AliError("Cannot get TreeH!");
291 TClonesArray* hits = fAD->Hits();
293 // Now makes Digits from hits
294 Int_t nTracks = (Int_t) treeH->GetEntries();
295 for (Int_t iTrack = 0; iTrack < nTracks; iTrack++) {
297 treeH->GetEvent(iTrack);
298 Int_t nHits = hits->GetEntriesFast();
299 for (Int_t iHit = 0; iHit < nHits; iHit++) {
300 AliADhit* hit = (AliADhit *)hits->UncheckedAt(iHit);
301 Int_t nPhot = hit->GetNphot();
302 Int_t pmt = hit->GetCell();//One PM per cell in AD
303 if (pmt < 0) continue;
304 Int_t trackLabel = hit->GetTrack();
305 for(Int_t l = 0; l < 3; ++l) {
306 if (fLabels[pmt][l] < 0) {
307 fLabels[pmt][l] = trackLabel;
311 Float_t dt_scintillator = gRandom->Gaus(0,kIntTimeRes);
312 Float_t t = dt_scintillator + hit->GetTof();
313 //if (pmt < 32) t += kV0CDelayCables;
314 t += fHptdcOffset[pmt];
316 Float_t prob = fCalibData->GetLightYields(pmt)*kPhotoCathodeEfficiency; // Optical losses included!
318 nPhE = (Int_t)gRandom->Gaus(prob*Float_t(nPhot)+0.5,
319 sqrt(Float_t(nPhot)*prob*(1.-prob)));
321 nPhE = gRandom->Binomial(nPhot,prob);
322 Float_t charge = TMath::Qe()*fPmGain[pmt]*fBinSize[pmt]/integral;
325 for (Int_t iPhE = 0; iPhE < nPhE; ++iPhE) {
326 Float_t tPhE = t + fSignalShape->GetRandom(0,fBinSize[pmt]*Float_t(fNBins[pmt]));
327 Float_t gainVar = fSinglePhESpectrum->GetRandom(0,20)/meansPhE;
328 Int_t firstBin = TMath::Max(0,(Int_t)((tPhE-kPMRespTime)/fBinSize[pmt]));
329 Int_t lastBin = TMath::Min(fNBins[pmt]-1,(Int_t)((tPhE+2.*kPMRespTime)/fBinSize[pmt]));
330 for(Int_t iBin = firstBin; iBin <= lastBin; ++iBin) {
331 Float_t tempT = fBinSize[pmt]*(0.5+iBin)-tPhE;
332 fTime[pmt][iBin] += gainVar*charge*fPMResponse->Eval(tempT);
337 loader->UnloadHits();
340 //____________________________________________________________________________
341 void AliADDigitizer::DigitizeSDigits()
343 // Digitize the fTime arrays (SDigits) to the level of
344 // Digits (fAdc arrays)
345 for(Int_t i = 0 ; i < 16; ++i) {
346 for(Int_t j = 0; j < kNClocks; ++j) fAdc[i][j] = 0;
347 fLeadingTime[i] = fTimeWidth[i] = 0;
350 Float_t maximum = 0.9*fSignalShape->GetMaximum(0,200); // Not exact, one needs to do this on the convoluted
351 Float_t integral2 = fSignalShape->Integral(0,200); // function. Anyway the effect is small <10% on the 2.5 ADC thr
352 for (Int_t ipmt = 0; ipmt < 16; ++ipmt) {
353 Float_t thr = fCalibData->GetCalibDiscriThr(ipmt,kFALSE)*kChargePerADC*maximum*fBinSize[ipmt]/integral2;
354 //Float_t thr = 1e-25;
356 Bool_t ltFound = kFALSE, ttFound = kFALSE;
357 for (Int_t iBin = 0; iBin < fNBins[ipmt]; ++iBin) {
358 Float_t t = fBinSize[ipmt]*Float_t(iBin);
360 if (fTime[ipmt][iBin] > thr) {
361 if (!ltFound && (iBin < fNBinsLT[ipmt])) {
363 fLeadingTime[ipmt] = t;
370 fTimeWidth[ipmt] = t - fLeadingTime[ipmt];
374 Float_t tadc = t - fClockOffset[ipmt];
375 Int_t clock = kNClocks/2 - Int_t(tadc/25.0);
376 if (clock >= 0 && clock < kNClocks)
377 fAdc[ipmt][clock] += fTime[ipmt][iBin]/kChargePerADC;
379 AliDebug(1,Form("Channel %d Offset %f Time %f",ipmt,fClockOffset[ipmt],fLeadingTime[ipmt]));
380 Int_t board = AliADCalibData::GetBoardNumber(ipmt);
381 if (ltFound && ttFound) {
382 fTimeWidth[ipmt] = fCalibData->GetWidthResolution(board)*
383 Float_t(Int_t(fTimeWidth[ipmt]/fCalibData->GetWidthResolution(board)));
384 if (fTimeWidth[ipmt] < Float_t(kMinTDCWidth)*fCalibData->GetWidthResolution(board))
385 fTimeWidth[ipmt] = Float_t(kMinTDCWidth)*fCalibData->GetWidthResolution(board);
386 if (fTimeWidth[ipmt] > Float_t(kMaxTDCWidth)*fCalibData->GetWidthResolution(board))
387 fTimeWidth[ipmt] = Float_t(kMaxTDCWidth)*fCalibData->GetWidthResolution(board);
391 fEvenOrOdd = gRandom->Integer(2);
392 for (Int_t j=0; j<16; ++j){
393 for (Int_t iClock = 0; iClock < kNClocks; ++iClock) {
394 Int_t integrator = (iClock + fEvenOrOdd) % 2;
395 AliDebug(1,Form("ADC %d %d %f",j,iClock,fAdc[j][iClock]));
396 fAdc[j][iClock] += gRandom->Gaus(fAdcPedestal[j][integrator], fAdcSigma[j][integrator]);
402 //____________________________________________________________________________
403 void AliADDigitizer::ReadSDigits()
405 // Read SDigits which are then to precessed
406 // in the following method
407 for(Int_t i = 0 ; i < 16; ++i) {
408 memset(fTime[i],0,fNBins[i]*sizeof(Float_t));
409 fLabels[i][0] = fLabels[i][1] = fLabels[i][2] = -1;
412 // Loop over input files
413 Int_t nFiles= fDigInput->GetNinputs();
414 for (Int_t inputFile = 0; inputFile < nFiles; inputFile++) {
415 // Get the current loader
416 AliRunLoader* currentLoader =
417 AliRunLoader::GetRunLoader(fDigInput->GetInputFolderName(inputFile));
419 AliLoader *loader = currentLoader->GetLoader("ADLoader");
420 loader->LoadSDigits("READ");
422 // Get the tree of summable digits
423 TTree* sdigitsTree = loader->TreeS();
425 AliError("No sdigit tree from digInput");
430 TBranch* sdigitsBranch = sdigitsTree->GetBranch("ADSDigit");
431 if (!sdigitsBranch) {
432 AliError("Failed to get sdigit branch");
436 // Set the branch address
437 TClonesArray *sdigitsArray = NULL;
438 sdigitsBranch->SetAddress(&sdigitsArray);
440 // Sum contributions from the sdigits
441 // Get number of entries in the tree
442 Int_t nentries = Int_t(sdigitsBranch->GetEntries());
443 for (Int_t entry = 0; entry < nentries; ++entry) {
444 sdigitsBranch->GetEntry(entry);
445 // Get the number of sdigits
446 Int_t nsdigits = sdigitsArray->GetEntries();
448 for (Int_t sdigit = 0; sdigit < nsdigits; sdigit++) {
449 AliADSDigit* sDigit = static_cast<AliADSDigit*>(sdigitsArray->UncheckedAt(sdigit));
450 Int_t pmNumber = sDigit->PMNumber();
451 Int_t nbins = sDigit->GetNBins();
452 if (nbins != fNBins[pmNumber]) {
453 AliError(Form("Incompatible number of bins between digitizer (%d) and sdigit (%d) for PM %d! Skipping sdigit!",
454 fNBins[pmNumber],nbins,pmNumber));
458 Float_t *charges = sDigit->GetCharges();
459 for(Int_t iBin = 0; iBin < nbins; ++iBin) fTime[pmNumber][iBin] += charges[iBin];
461 Int_t *labels = sDigit->GetTracks();
463 for(Int_t i = 0; i < 3; ++i) {
464 if (fLabels[pmNumber][i] < 0) {
465 if (labels[j] < 0) break;
466 fLabels[pmNumber][i] = labels[j];
472 loader->UnloadSDigits();
477 //____________________________________________________________________________
478 void AliADDigitizer::WriteDigits(AliLoader *loader)
480 // Take fAdc arrays filled by the previous
481 // method and produce and add digits to the digit Tree
483 loader->LoadDigits("UPDATE");
485 if (!loader->TreeD()) loader->MakeTree("D");
486 loader->MakeDigitsContainer();
487 TTree* treeD = loader->TreeD();
489 treeD->Branch("ADDigit", &fDigits);
490 //fAD->MakeBranchInTree(treeD,"AD",&fDigits,1000,"");
492 Short_t *chargeADC = new Short_t[kNClocks];
493 for (Int_t i=0; i<16; i++) {
494 for (Int_t j = 0; j < kNClocks; ++j) {
495 Int_t tempadc = Int_t(fAdc[i][j]);
496 if (tempadc > 1023) tempadc = 1023;
497 chargeADC[j] = tempadc;
499 AddDigit(i, fLeadingTime[i], fTimeWidth[i], Bool_t((10+fEvenOrOdd)%2), chargeADC, fLabels[i]);
504 loader->WriteDigits("OVERWRITE");
505 loader->UnloadDigits();
509 //____________________________________________________________________________
510 void AliADDigitizer::WriteSDigits(AliLoader *loader)
512 // Take fTime arrays filled by the previous
513 // method and produce and add sdigits to the sdigit Tree
515 loader->LoadSDigits("UPDATE");
517 if (!loader->TreeS()) loader->MakeTree("S");
518 loader->MakeSDigitsContainer();
519 TTree* treeS = loader->TreeS();
521 treeS->Branch("ADSDigit", &fDigits);
522 //fAD->MakeBranchInTree(treeS,"AD",&fDigits,8000,"");
524 for (Int_t ipmt = 0; ipmt < 16; ++ipmt) {
525 AddSDigit(ipmt,fNBins[ipmt],fTime[ipmt],fLabels[ipmt]);
529 loader->WriteSDigits("OVERWRITE");
530 loader->UnloadSDigits();
536 //____________________________________________________________________________
537 void AliADDigitizer::AddDigit(Int_t pmnumber, Float_t time, Float_t width, Bool_t integrator, Short_t *chargeADC, Int_t *labels)
542 TClonesArray &ldigits = *fDigits;
544 new(ldigits[fNdigits++]) AliADdigit(pmnumber,time,width,integrator,chargeADC,labels);
547 //____________________________________________________________________________
548 void AliADDigitizer::AddSDigit(Int_t pmnumber, Int_t nbins, Float_t *charges, Int_t *labels)
553 TClonesArray &ldigits = *fDigits;
555 new(ldigits[fNdigits++]) AliADSDigit(pmnumber,nbins,charges,labels);
558 //____________________________________________________________________________
559 void AliADDigitizer::ResetDigits()
565 if (fDigits) fDigits->Clear();
568 //____________________________________________________________________________
569 AliADCalibData* AliADDigitizer::GetCalibData() const
573 AliCDBManager *man = AliCDBManager::Instance();
575 AliCDBEntry *entry=0;
577 entry = man->Get("AD/Calib/Data");
579 AliADCalibData *calibdata = 0;
581 if (entry) calibdata = (AliADCalibData*) entry->GetObject();
582 if (!calibdata) AliFatal("No calibration data from calibration database !");
584 AliADCalibData *calibdata = new AliADCalibData();
590 //____________________________________________________________________________
591 double AliADDigitizer::SignalShape(double *x, double *par)
593 // this function simulates the time
594 // of arrival of the photons at the
597 if (xx <= par[0]) return 0;
598 Double_t a = 1./TMath::Power((xx-par[0])/par[1],1./par[2]);
599 if (xx <= par[3]) return a;
600 Double_t b = 1./TMath::Power((xx-par[3])/par[4],1./par[5]);
601 Double_t f = a*b/(a+b);
602 AliDebug(100,Form("x=%f func=%f",xx,f));
606 //____________________________________________________________________________
607 double AliADDigitizer::PMResponse(double *x, double * /* par */)
609 // this function describes the
610 // PM time response to a single
612 Double_t xx = x[0]+kPMRespTime;
613 return xx*xx*TMath::Exp(-xx*xx/(kPMRespTime*kPMRespTime));
616 //____________________________________________________________________________
617 double AliADDigitizer::SinglePhESpectrum(double *x, double * /* par */)
619 // this function describes the
620 // PM amplitude response to a single
623 if (xx < 0) return 0;
624 return (TMath::Poisson(xx,kPMNbOfSecElec)+kPMTransparency*TMath::Poisson(xx,1.0));
626 //____________________________________________________________________
627 TClonesArray* AliADDigitizer::DigitsArray()
629 // Initialize digit array if not already and
630 // return pointer to it.
632 fDigits = new TClonesArray("AliADdigit", 16);
638 //____________________________________________________________________
639 TClonesArray* AliADDigitizer::SDigitsArray()
641 // Initialize sdigit array if not already and
642 // return pointer to it.
644 fDigits = new TClonesArray("AliADSDigit", 16);