#include <TMath.h>
#include <TTree.h>
#include <TMap.h>
#include <TGeoManager.h>
#include <TGeoPhysicalNode.h>
#include <AliGeomManager.h>
#include <TRandom.h>
#include <TF1.h>
#include <TH1F.h>
#include "AliRun.h"
#include "AliDetector.h"
#include "AliAD.h"
#include "AliADhit.h"
#include "AliADConst.h"
#include "AliRunLoader.h"
#include "AliLoader.h"
#include "AliGRPObject.h"
#include "AliDigitizationInput.h"
#include "AliCDBManager.h"
#include "AliCDBStorage.h"
#include "AliCDBEntry.h"
#include "AliADCalibData.h"
#include "AliCTPTimeParams.h"
#include "AliLHCClockPhase.h"
#include "AliADdigit.h"
#include "AliADDigitizer.h"
#include "AliADSDigit.h"
ClassImp(AliADDigitizer)
AliADDigitizer::AliADDigitizer()
:AliDigitizer(),
fCalibData(GetCalibData()),
fNdigits(0),
fDigits(0),
fSignalShape(NULL),
fPMResponse(NULL),
fSinglePhESpectrum(NULL),
fEvenOrOdd(kFALSE),
fTask(kHits2Digits),
fAD(NULL)
{
Init();
}
AliADDigitizer::AliADDigitizer(AliAD *AD, DigiTask_t task)
:AliDigitizer(),
fCalibData(GetCalibData()),
fNdigits(0),
fDigits(0),
fSignalShape(NULL),
fPMResponse(NULL),
fSinglePhESpectrum(NULL),
fEvenOrOdd(kFALSE),
fTask(task),
fAD(AD)
{
Init();
}
AliADDigitizer::AliADDigitizer(AliDigitizationInput* digInput)
:AliDigitizer(digInput),
fCalibData(GetCalibData()),
fNdigits(0),
fDigits(0),
fSignalShape(NULL),
fPMResponse(NULL),
fSinglePhESpectrum(NULL),
fEvenOrOdd(kFALSE),
fTask(kHits2Digits),
fAD(NULL)
{
Init();
}
AliADDigitizer::~AliADDigitizer()
{
if (fDigits) {
fDigits->Delete();
delete fDigits;
fDigits=0;
}
if (fSignalShape) {
delete fSignalShape;
fSignalShape = NULL;
}
if (fPMResponse) {
delete fPMResponse;
fPMResponse = NULL;
}
if (fSinglePhESpectrum) {
delete fSinglePhESpectrum;
fSinglePhESpectrum = NULL;
}
for(Int_t i = 0 ; i < 16; ++i) {
if (fTime[i]) delete [] fTime[i];
}
}
Bool_t AliADDigitizer::Init()
{
if (fSignalShape) return kTRUE;
fSignalShape = new TF1("ADSignalShape",this,&AliADDigitizer::SignalShape,0,200,6,"AliADDigitizer","SignalShape");
fSignalShape->SetParameters(-1.07335e+00,2.16002e+01,-1.26133e-01,
1.41619e+00,5.50334e-01,3.86111e-01);
fPMResponse = new TF1("ADPMResponse",this,&AliADDigitizer::PMResponse,-kPMRespTime,2.*kPMRespTime,0,"AliADDigitizer","PMResponse");
fSinglePhESpectrum = new TF1("ADSinglePhESpectrum",this,&AliADDigitizer::SinglePhESpectrum,0,20,0,"AliADDigitizer","SinglePhESpectrum");
AliCDBEntry *entry = AliCDBManager::Instance()->Get("GRP/CTP/CTPtiming");
if (!entry) AliFatal("CTP timing parameters are not found in OCDB !");
AliCTPTimeParams *ctpParams = (AliCTPTimeParams*)entry->GetObject();
Float_t l1Delay = (Float_t)ctpParams->GetDelayL1L0()*25.0;
AliCDBEntry *entry1 = AliCDBManager::Instance()->Get("GRP/CTP/TimeAlign");
if (!entry1) AliFatal("CTP time-alignment is not found in OCDB !");
AliCTPTimeParams *ctpTimeAlign = (AliCTPTimeParams*)entry1->GetObject();
l1Delay += ((Float_t)ctpTimeAlign->GetDelayL1L0()*25.0);
AliCDBEntry *entry2 = AliCDBManager::Instance()->Get("VZERO/Calib/TimeDelays");
if (!entry2) AliFatal("AD time delays are not found in OCDB !");
TH1F *delays = (TH1F*)entry2->GetObject();
AliCDBEntry *entry3 = AliCDBManager::Instance()->Get("GRP/Calib/LHCClockPhase");
if (!entry3) AliFatal("LHC clock-phase shift is not found in OCDB !");
AliLHCClockPhase *phase = (AliLHCClockPhase*)entry3->GetObject();
for(Int_t i = 0 ; i < 16; ++i) {
for(Int_t j = 0; j < kNClocks; ++j) fAdc[i][j] = 0;
fLeadingTime[i] = fTimeWidth[i] = 0;
fPmGain[i] = fCalibData->GetGain(i);
fAdcPedestal[i][0] = fCalibData->GetPedestal(i);
fAdcSigma[i][0] = fCalibData->GetSigma(i);
fAdcPedestal[i][1] = fCalibData->GetPedestal(i+16);
fAdcSigma[i][1] = fCalibData->GetSigma(i+16);
Int_t board = AliADCalibData::GetBoardNumber(i);
fNBins[i] = TMath::Nint(((Float_t)(fCalibData->GetMatchWindow(board)+1)*25.0+
(Float_t)kMaxTDCWidth*fCalibData->GetWidthResolution(board))/
fCalibData->GetTimeResolution(board));
fNBinsLT[i] = TMath::Nint(((Float_t)(fCalibData->GetMatchWindow(board)+1)*25.0)/
fCalibData->GetTimeResolution(board));
fBinSize[i] = fCalibData->GetTimeResolution(board);
fHptdcOffset[i] = (((Float_t)fCalibData->GetRollOver(board)-
(Float_t)fCalibData->GetTriggerCountOffset(board))*25.0+
fCalibData->GetTimeOffset(i)-
l1Delay-
phase->GetMeanPhase()+
delays->GetBinContent(i+1)+
kV0Offset);
fClockOffset[i] = (((Float_t)fCalibData->GetRollOver(board)-
(Float_t)fCalibData->GetTriggerCountOffset(board))*25.0+
fCalibData->GetTimeOffset(i)-
l1Delay+
kV0Offset);
fTime[i] = new Float_t[fNBins[i]];
memset(fTime[i],0,fNBins[i]*sizeof(Float_t));
}
return kTRUE;
}
void AliADDigitizer::Digitize(Option_t* )
{
fNdigits = 0;
if (fAD && !fDigInput) {
AliLoader *loader = fAD->GetLoader();
if (!loader) {
AliError("Can not get AD Loader via AliAD object!");
return;
}
AliRunLoader* runLoader = AliRunLoader::Instance();
for (Int_t iEvent = 0; iEvent < runLoader->GetNumberOfEvents(); ++iEvent) {
runLoader->GetEvent(iEvent);
if (fTask == kHits2Digits) {
DigitizeHits();
DigitizeSDigits();
WriteDigits(loader);
}
else {
DigitizeHits();
WriteSDigits(loader);
}
}
}
else if (fDigInput) {
ReadSDigits();
DigitizeSDigits();
AliRunLoader *currentLoader = AliRunLoader::GetRunLoader(fDigInput->GetOutputFolderName());
AliLoader *loader = currentLoader->GetLoader("ADLoader");
if (!loader) {
AliError("Cannot get AD Loader via RunDigitizer!");
return;
}
WriteDigits(loader);
}
else {
AliFatal("Invalid digitization task! Exiting!");
}
}
void AliADDigitizer::DigitizeHits()
{
for(Int_t i = 0 ; i < 16; ++i) {
memset(fTime[i],0,fNBins[i]*sizeof(Float_t));
fLabels[i][0] = fLabels[i][1] = fLabels[i][2] = -1;
}
Float_t integral = fPMResponse->Integral(-kPMRespTime,2.*kPMRespTime);
Float_t meansPhE = fSinglePhESpectrum->Mean(0,20);
AliLoader* loader = fAD->GetLoader();
if (!loader) {
AliError("Can not get AD Loader!");
return;
}
loader->LoadHits();
TTree* treeH = loader->TreeH();
if (!treeH) {
AliError("Cannot get TreeH!");
return;
}
TClonesArray* hits = fAD->Hits();
Int_t nTracks = (Int_t) treeH->GetEntries();
for (Int_t iTrack = 0; iTrack < nTracks; iTrack++) {
fAD->ResetHits();
treeH->GetEvent(iTrack);
Int_t nHits = hits->GetEntriesFast();
for (Int_t iHit = 0; iHit < nHits; iHit++) {
AliADhit* hit = (AliADhit *)hits->UncheckedAt(iHit);
Int_t nPhot = hit->GetNphot();
Int_t pmt = hit->GetCell();
if (pmt < 0) continue;
Int_t trackLabel = hit->GetTrack();
for(Int_t l = 0; l < 3; ++l) {
if (fLabels[pmt][l] < 0) {
fLabels[pmt][l] = trackLabel;
break;
}
}
Float_t dt_scintillator = gRandom->Gaus(0,kIntTimeRes);
Float_t t = dt_scintillator + hit->GetTof();
t += fHptdcOffset[pmt];
Int_t nPhE;
Float_t prob = fCalibData->GetLightYields(pmt)*kPhotoCathodeEfficiency;
if (nPhot > 100)
nPhE = (Int_t)gRandom->Gaus(prob*Float_t(nPhot)+0.5,
sqrt(Float_t(nPhot)*prob*(1.-prob)));
else
nPhE = gRandom->Binomial(nPhot,prob);
Float_t charge = TMath::Qe()*fPmGain[pmt]*fBinSize[pmt]/integral;
for (Int_t iPhE = 0; iPhE < nPhE; ++iPhE) {
Float_t tPhE = t + fSignalShape->GetRandom(0,fBinSize[pmt]*Float_t(fNBins[pmt]));
Float_t gainVar = fSinglePhESpectrum->GetRandom(0,20)/meansPhE;
Int_t firstBin = TMath::Max(0,(Int_t)((tPhE-kPMRespTime)/fBinSize[pmt]));
Int_t lastBin = TMath::Min(fNBins[pmt]-1,(Int_t)((tPhE+2.*kPMRespTime)/fBinSize[pmt]));
for(Int_t iBin = firstBin; iBin <= lastBin; ++iBin) {
Float_t tempT = fBinSize[pmt]*(0.5+iBin)-tPhE;
fTime[pmt][iBin] += gainVar*charge*fPMResponse->Eval(tempT);
}
}
}
}
loader->UnloadHits();
}
void AliADDigitizer::DigitizeSDigits()
{
for(Int_t i = 0 ; i < 16; ++i) {
for(Int_t j = 0; j < kNClocks; ++j) fAdc[i][j] = 0;
fLeadingTime[i] = fTimeWidth[i] = 0;
}
Float_t maximum = 0.9*fSignalShape->GetMaximum(0,200);
Float_t integral2 = fSignalShape->Integral(0,200);
for (Int_t ipmt = 0; ipmt < 16; ++ipmt) {
Float_t thr = fCalibData->GetCalibDiscriThr(ipmt,kFALSE)*kChargePerADC*maximum*fBinSize[ipmt]/integral2;
Bool_t ltFound = kFALSE, ttFound = kFALSE;
for (Int_t iBin = 0; iBin < fNBins[ipmt]; ++iBin) {
Float_t t = fBinSize[ipmt]*Float_t(iBin);
if (fTime[ipmt][iBin] > thr) {
if (!ltFound && (iBin < fNBinsLT[ipmt])) {
ltFound = kTRUE;
fLeadingTime[ipmt] = t;
}
}
else {
if (ltFound) {
if (!ttFound) {
ttFound = kTRUE;
fTimeWidth[ipmt] = t - fLeadingTime[ipmt];
}
}
}
Float_t tadc = t - fClockOffset[ipmt];
Int_t clock = kNClocks/2 - Int_t(tadc/25.0);
if (clock >= 0 && clock < kNClocks)
fAdc[ipmt][clock] += fTime[ipmt][iBin]/kChargePerADC;
}
AliDebug(1,Form("Channel %d Offset %f Time %f",ipmt,fClockOffset[ipmt],fLeadingTime[ipmt]));
Int_t board = AliADCalibData::GetBoardNumber(ipmt);
if (ltFound && ttFound) {
fTimeWidth[ipmt] = fCalibData->GetWidthResolution(board)*
Float_t(Int_t(fTimeWidth[ipmt]/fCalibData->GetWidthResolution(board)));
if (fTimeWidth[ipmt] < Float_t(kMinTDCWidth)*fCalibData->GetWidthResolution(board))
fTimeWidth[ipmt] = Float_t(kMinTDCWidth)*fCalibData->GetWidthResolution(board);
if (fTimeWidth[ipmt] > Float_t(kMaxTDCWidth)*fCalibData->GetWidthResolution(board))
fTimeWidth[ipmt] = Float_t(kMaxTDCWidth)*fCalibData->GetWidthResolution(board);
}
}
fEvenOrOdd = gRandom->Integer(2);
for (Int_t j=0; j<16; ++j){
for (Int_t iClock = 0; iClock < kNClocks; ++iClock) {
Int_t integrator = (iClock + fEvenOrOdd) % 2;
AliDebug(1,Form("ADC %d %d %f",j,iClock,fAdc[j][iClock]));
fAdc[j][iClock] += gRandom->Gaus(fAdcPedestal[j][integrator], fAdcSigma[j][integrator]);
}
}
}
void AliADDigitizer::ReadSDigits()
{
for(Int_t i = 0 ; i < 16; ++i) {
memset(fTime[i],0,fNBins[i]*sizeof(Float_t));
fLabels[i][0] = fLabels[i][1] = fLabels[i][2] = -1;
}
Int_t nFiles= fDigInput->GetNinputs();
for (Int_t inputFile = 0; inputFile < nFiles; inputFile++) {
AliRunLoader* currentLoader =
AliRunLoader::GetRunLoader(fDigInput->GetInputFolderName(inputFile));
AliLoader *loader = currentLoader->GetLoader("ADLoader");
loader->LoadSDigits("READ");
TTree* sdigitsTree = loader->TreeS();
if (!sdigitsTree) {
AliError("No sdigit tree from digInput");
continue;
}
TBranch* sdigitsBranch = sdigitsTree->GetBranch("ADSDigit");
if (!sdigitsBranch) {
AliError("Failed to get sdigit branch");
return;
}
TClonesArray *sdigitsArray = NULL;
sdigitsBranch->SetAddress(&sdigitsArray);
Int_t nentries = Int_t(sdigitsBranch->GetEntries());
for (Int_t entry = 0; entry < nentries; ++entry) {
sdigitsBranch->GetEntry(entry);
Int_t nsdigits = sdigitsArray->GetEntries();
for (Int_t sdigit = 0; sdigit < nsdigits; sdigit++) {
AliADSDigit* sDigit = static_cast<AliADSDigit*>(sdigitsArray->UncheckedAt(sdigit));
Int_t pmNumber = sDigit->PMNumber();
Int_t nbins = sDigit->GetNBins();
if (nbins != fNBins[pmNumber]) {
AliError(Form("Incompatible number of bins between digitizer (%d) and sdigit (%d) for PM %d! Skipping sdigit!",
fNBins[pmNumber],nbins,pmNumber));
continue;
}
Float_t *charges = sDigit->GetCharges();
for(Int_t iBin = 0; iBin < nbins; ++iBin) fTime[pmNumber][iBin] += charges[iBin];
Int_t *labels = sDigit->GetTracks();
Int_t j = 0;
for(Int_t i = 0; i < 3; ++i) {
if (fLabels[pmNumber][i] < 0) {
if (labels[j] < 0) break;
fLabels[pmNumber][i] = labels[j];
j++;
}
}
}
}
loader->UnloadSDigits();
}
}
void AliADDigitizer::WriteDigits(AliLoader *loader)
{
loader->LoadDigits("UPDATE");
if (!loader->TreeD()) loader->MakeTree("D");
loader->MakeDigitsContainer();
TTree* treeD = loader->TreeD();
DigitsArray();
treeD->Branch("ADDigit", &fDigits);
Short_t *chargeADC = new Short_t[kNClocks];
for (Int_t i=0; i<16; i++) {
for (Int_t j = 0; j < kNClocks; ++j) {
Int_t tempadc = Int_t(fAdc[i][j]);
if (tempadc > 1023) tempadc = 1023;
chargeADC[j] = tempadc;
}
AddDigit(i, fLeadingTime[i], fTimeWidth[i], Bool_t((10+fEvenOrOdd)%2), chargeADC, fLabels[i]);
}
delete [] chargeADC;
treeD->Fill();
loader->WriteDigits("OVERWRITE");
loader->UnloadDigits();
ResetDigits();
}
void AliADDigitizer::WriteSDigits(AliLoader *loader)
{
loader->LoadSDigits("UPDATE");
if (!loader->TreeS()) loader->MakeTree("S");
loader->MakeSDigitsContainer();
TTree* treeS = loader->TreeS();
SDigitsArray();
treeS->Branch("ADSDigit", &fDigits);
for (Int_t ipmt = 0; ipmt < 16; ++ipmt) {
AddSDigit(ipmt,fNBins[ipmt],fTime[ipmt],fLabels[ipmt]);
}
treeS->Fill();
loader->WriteSDigits("OVERWRITE");
loader->UnloadSDigits();
ResetDigits();
}
void AliADDigitizer::AddDigit(Int_t pmnumber, Float_t time, Float_t width, Bool_t integrator, Short_t *chargeADC, Int_t *labels)
{
TClonesArray &ldigits = *fDigits;
new(ldigits[fNdigits++]) AliADdigit(pmnumber,time,width,integrator,chargeADC,labels);
}
void AliADDigitizer::AddSDigit(Int_t pmnumber, Int_t nbins, Float_t *charges, Int_t *labels)
{
TClonesArray &ldigits = *fDigits;
new(ldigits[fNdigits++]) AliADSDigit(pmnumber,nbins,charges,labels);
}
void AliADDigitizer::ResetDigits()
{
fNdigits = 0;
if (fDigits) fDigits->Clear();
}
AliADCalibData* AliADDigitizer::GetCalibData() const
{
AliADCalibData *calibdata = new AliADCalibData();
return calibdata;
}
double AliADDigitizer::SignalShape(double *x, double *par)
{
Double_t xx = x[0];
if (xx <= par[0]) return 0;
Double_t a = 1./TMath::Power((xx-par[0])/par[1],1./par[2]);
if (xx <= par[3]) return a;
Double_t b = 1./TMath::Power((xx-par[3])/par[4],1./par[5]);
Double_t f = a*b/(a+b);
AliDebug(100,Form("x=%f func=%f",xx,f));
return f;
}
double AliADDigitizer::PMResponse(double *x, double * )
{
Double_t xx = x[0]+kPMRespTime;
return xx*xx*TMath::Exp(-xx*xx/(kPMRespTime*kPMRespTime));
}
double AliADDigitizer::SinglePhESpectrum(double *x, double * )
{
Double_t xx = x[0];
if (xx < 0) return 0;
return (TMath::Poisson(xx,kPMNbOfSecElec)+kPMTransparency*TMath::Poisson(xx,1.0));
}
TClonesArray* AliADDigitizer::DigitsArray()
{
if (!fDigits) {
fDigits = new TClonesArray("AliADdigit", 16);
fNdigits = 0;
}
return fDigits;
}
TClonesArray* AliADDigitizer::SDigitsArray()
{
if (!fDigits) {
fDigits = new TClonesArray("AliADSDigit", 16);
fNdigits = 0;
}
return fDigits;
}