#include <AliRun.h>
#include <AliRunLoader.h>
#include "AliDigitizationInput.h"
#include <AliLoader.h>
#include "AliLog.h"
#include <AliCDBEntry.h>
#include <AliCDBManager.h>
#include "AliHMPIDDigitizer.h"
#include "AliHMPIDReconstructor.h"
#include "AliHMPIDDigit.h"
#include "AliHMPID.h"
#include "AliHMPIDParam.h"
#include <TRandom.h>
#include <TMath.h>
#include <TTree.h>
#include <TObjArray.h>
ClassImp(AliHMPIDDigitizer)
Bool_t AliHMPIDDigitizer::fgDoNoise=kTRUE;
void AliHMPIDDigitizer::Digitize(Option_t*)
{
AliDebug(1,Form("Start with %i input(s) for event %i",fDigInput->GetNinputs(),fDigInput->GetOutputEventNr()));
AliRunLoader *pInRunLoader=0;
AliLoader *pInRichLoader=0;
static TClonesArray sdigs("AliHMPIDDigit");
Int_t total=0;
for(Int_t inFileN=0;inFileN<fDigInput->GetNinputs();inFileN++){
pInRunLoader = AliRunLoader::GetRunLoader(fDigInput->GetInputFolderName(inFileN));
pInRichLoader = pInRunLoader->GetLoader("HMPIDLoader"); if(pInRichLoader==0) continue;
if (!pInRunLoader->GetAliRun()) pInRunLoader->LoadgAlice();
AliHMPID* pInRich=(AliHMPID*)pInRunLoader->GetAliRun()->GetDetector("HMPID");
pInRichLoader->LoadSDigits(); pInRichLoader->TreeS()->GetEntry(0);
AliDebug(1,Form("input %i has %i sdigits",inFileN,pInRich->SdiLst()->GetEntries()));
for(Int_t i=0;i<pInRich->SdiLst()->GetEntries();i++){
AliHMPIDDigit *pSDig=(AliHMPIDDigit*)pInRich->SdiLst()->At(i);
pSDig->AddTidOffset(fDigInput->GetMask(inFileN));
new(sdigs[total++]) AliHMPIDDigit(*pSDig);
}
pInRichLoader->UnloadSDigits(); pInRich->SdiReset();
}
AliRunLoader *pOutRunLoader = AliRunLoader::GetRunLoader(fDigInput->GetOutputFolderName());
AliLoader *pOutRichLoader = pOutRunLoader->GetLoader("HMPIDLoader");
AliRun *pArun = pOutRunLoader->GetAliRun();
AliHMPID *pOutRich = (AliHMPID*)pArun->GetDetector("HMPID");
pOutRichLoader->MakeTree("D"); pOutRich->MakeBranch("D");
Sdi2Dig(&sdigs,pOutRich->DigLst());
pOutRichLoader->TreeD()->Fill();
pOutRichLoader->WriteDigits("OVERWRITE");
sdigs.Clear();
pOutRichLoader->UnloadDigits(); pOutRich->DigReset();
}
void AliHMPIDDigitizer::Sdi2Dig(TClonesArray *pSdiLst,TObjArray *pDigLst)
{
TClonesArray *pLst[7]; Int_t iCnt[7];
for(Int_t i=0;i<7;i++){
pLst[i]=(TClonesArray*)(*pDigLst)[i];
iCnt[i]=0; if(pLst[i]->GetEntries()!=0) AliErrorClass("Some of digits lists is not empty");
}
Float_t arrNoise[7][6][80][48], arrSigmaPed[7][6][80][48];
if(fgDoNoise) {
for (Int_t iCh=AliHMPIDParam::kMinCh;iCh<=AliHMPIDParam::kMaxCh;iCh++)
for (Int_t iPc=AliHMPIDParam::kMinPc;iPc<=AliHMPIDParam::kMaxPc;iPc++)
for(Int_t iPx=AliHMPIDParam::kMinPx;iPx<=AliHMPIDParam::kMaxPx;iPx++)
for(Int_t iPy=AliHMPIDParam::kMinPy;iPy<=AliHMPIDParam::kMaxPy;iPy++){
arrNoise[iCh][iPc][iPx][iPy] = gRandom->Gaus(0,1.);
arrSigmaPed[iCh][iPc][iPx][iPy] = 1.;
}
AliCDBEntry *pDaqSigEnt = AliCDBManager::Instance()->Get("HMPID/Calib/DaqSig");
if(pDaqSigEnt){
TObjArray *pDaqSig = (TObjArray*)pDaqSigEnt->GetObject();
for(Int_t iCh=AliHMPIDParam::kMinCh;iCh<=AliHMPIDParam::kMaxCh;iCh++){
TMatrixF *pM = (TMatrixF*)pDaqSig->At(iCh);
for (Int_t iPc=AliHMPIDParam::kMinPc;iPc<=AliHMPIDParam::kMaxPc;iPc++)
for(Int_t iPx=AliHMPIDParam::kMinPx;iPx<=AliHMPIDParam::kMaxPx;iPx++)
for(Int_t iPy=AliHMPIDParam::kMinPy;iPy<=AliHMPIDParam::kMaxPy;iPy++){
Int_t padX = (iPc%2)*AliHMPIDParam::kPadPcX+iPx;
Int_t padY = (iPc/2)*AliHMPIDParam::kPadPcY+iPy;
if((*pM)(padX,padY)>0.){
arrNoise[iCh][iPc][iPx][iPy] = gRandom->Gaus(0,(*pM)(padX,padY));
arrSigmaPed[iCh][iPc][iPx][iPy] = (*pM)(padX,padY);}
else{
arrNoise[iCh][iPc][iPx][iPy] = gRandom->Gaus(0,1.);
arrSigmaPed[iCh][iPc][iPx][iPy] = 1.;}
}
}
}
}
pSdiLst->Sort();
Int_t iPad=-1,iCh=-1,iNdigPad=-1,aTids[3]={-1,-1,-1}; Float_t q=-1;
for(Int_t i=0;i<pSdiLst->GetEntries();i++){
AliHMPIDDigit *pSdig=(AliHMPIDDigit*)pSdiLst->At(i);
if(pSdig->Pad()==iPad){
q+=pSdig->Q();
iNdigPad++; if(iNdigPad<=3) aTids[iNdigPad-1]=pSdig->GetTrack(0);
continue;
}
if(i!=0 && iCh>=AliHMPIDParam::kMinCh && iCh<=AliHMPIDParam::kMaxCh){
AliHMPIDParam::Instance()->SetThreshold((TMath::Nint(arrSigmaPed[iCh][pSdig->Pc()][pSdig->PadPcX()][pSdig->PadPcY()])*AliHMPIDParam::Nsig()));
if(AliHMPIDParam::IsOverTh(q)) new((*pLst[iCh])[iCnt[iCh]++]) AliHMPIDDigit(iPad,(Int_t)q,aTids);}
iPad=pSdig->Pad(); iCh=AliHMPIDParam::A2C(iPad);
iNdigPad=1;
aTids[0]=pSdig->GetTrack(0);aTids[1]=aTids[2]=-1;
q=pSdig->Q();
if(fgDoNoise) q+=arrNoise[iCh][pSdig->Pc()][pSdig->PadPcX()][pSdig->PadPcY()];
arrNoise[iCh][pSdig->Pc()][pSdig->PadPcX()][pSdig->PadPcY()]=0;
}
if(iCh>=AliHMPIDParam::kMinCh && iCh<=AliHMPIDParam::kMaxCh){
Int_t pc = AliHMPIDParam::A2P(iPad);
Int_t px = AliHMPIDParam::A2X(iPad);
Int_t py = AliHMPIDParam::A2Y(iPad);
AliHMPIDParam::Instance()->SetThreshold((TMath::Nint(arrSigmaPed[iCh][pc][px][py])*AliHMPIDParam::Nsig()));
if(AliHMPIDParam::IsOverTh(q)) new((*pLst[iCh])[iCnt[iCh]++]) AliHMPIDDigit(iPad,(Int_t)q,aTids);
}
if(!fgDoNoise) return;
aTids[0]=aTids[1]=aTids[2]=-1;
for (Int_t iChCurr=AliHMPIDParam::kMinCh;iChCurr<=AliHMPIDParam::kMaxCh;iChCurr++){
for (Int_t iPc=AliHMPIDParam::kMinPc;iPc<=AliHMPIDParam::kMaxPc;iPc++)
for(Int_t iPx=AliHMPIDParam::kMinPx;iPx<=AliHMPIDParam::kMaxPx;iPx++)
for(Int_t iPy=AliHMPIDParam::kMinPy;iPy<=AliHMPIDParam::kMaxPy;iPy++) {
Float_t qNoise = arrNoise[iChCurr][iPc][iPx][iPy];
AliHMPIDParam::Instance()->SetThreshold((TMath::Nint(arrSigmaPed[iChCurr][iPc][iPx][iPy])*AliHMPIDParam::Nsig()));
if(AliHMPIDParam::IsOverTh(qNoise)) new((*pLst[iChCurr])[iCnt[iChCurr]++]) AliHMPIDDigit(AliHMPIDParam::Abs(iChCurr,iPc,iPx,iPy),(Int_t)qNoise,aTids);
}
}
}
AliHMPIDDigitizer.cxx:100 AliHMPIDDigitizer.cxx:101 AliHMPIDDigitizer.cxx:102 AliHMPIDDigitizer.cxx:103 AliHMPIDDigitizer.cxx:104 AliHMPIDDigitizer.cxx:105 AliHMPIDDigitizer.cxx:106 AliHMPIDDigitizer.cxx:107 AliHMPIDDigitizer.cxx:108 AliHMPIDDigitizer.cxx:109 AliHMPIDDigitizer.cxx:110 AliHMPIDDigitizer.cxx:111 AliHMPIDDigitizer.cxx:112 AliHMPIDDigitizer.cxx:113 AliHMPIDDigitizer.cxx:114 AliHMPIDDigitizer.cxx:115 AliHMPIDDigitizer.cxx:116 AliHMPIDDigitizer.cxx:117 AliHMPIDDigitizer.cxx:118 AliHMPIDDigitizer.cxx:119 AliHMPIDDigitizer.cxx:120 AliHMPIDDigitizer.cxx:121 AliHMPIDDigitizer.cxx:122 AliHMPIDDigitizer.cxx:123 AliHMPIDDigitizer.cxx:124 AliHMPIDDigitizer.cxx:125 AliHMPIDDigitizer.cxx:126 AliHMPIDDigitizer.cxx:127 AliHMPIDDigitizer.cxx:128 AliHMPIDDigitizer.cxx:129 AliHMPIDDigitizer.cxx:130 AliHMPIDDigitizer.cxx:131 AliHMPIDDigitizer.cxx:132 AliHMPIDDigitizer.cxx:133 AliHMPIDDigitizer.cxx:134 AliHMPIDDigitizer.cxx:135 AliHMPIDDigitizer.cxx:136 AliHMPIDDigitizer.cxx:137 AliHMPIDDigitizer.cxx:138 AliHMPIDDigitizer.cxx:139 AliHMPIDDigitizer.cxx:140 AliHMPIDDigitizer.cxx:141 AliHMPIDDigitizer.cxx:142 AliHMPIDDigitizer.cxx:143 AliHMPIDDigitizer.cxx:144 AliHMPIDDigitizer.cxx:145 AliHMPIDDigitizer.cxx:146 AliHMPIDDigitizer.cxx:147 AliHMPIDDigitizer.cxx:148 AliHMPIDDigitizer.cxx:149 AliHMPIDDigitizer.cxx:150 AliHMPIDDigitizer.cxx:151 AliHMPIDDigitizer.cxx:152 AliHMPIDDigitizer.cxx:153 AliHMPIDDigitizer.cxx:154 AliHMPIDDigitizer.cxx:155 AliHMPIDDigitizer.cxx:156 AliHMPIDDigitizer.cxx:157 AliHMPIDDigitizer.cxx:158 AliHMPIDDigitizer.cxx:159 AliHMPIDDigitizer.cxx:160 AliHMPIDDigitizer.cxx:161 AliHMPIDDigitizer.cxx:162 AliHMPIDDigitizer.cxx:163 AliHMPIDDigitizer.cxx:164 AliHMPIDDigitizer.cxx:165 AliHMPIDDigitizer.cxx:166 AliHMPIDDigitizer.cxx:167 AliHMPIDDigitizer.cxx:168 AliHMPIDDigitizer.cxx:169 AliHMPIDDigitizer.cxx:170 AliHMPIDDigitizer.cxx:171 AliHMPIDDigitizer.cxx:172 AliHMPIDDigitizer.cxx:173 AliHMPIDDigitizer.cxx:174