#include <TH2F.h>
#include <TParticle.h>
#include "AliESDEvent.h"
#include <AliESDtrack.h>
#include <AliMCEvent.h>
#include <AliLog.h>
#include <AliKFParticle.h>
#include <AliKFVertex.h>
#include "AliHFEpriVtx.h"
ClassImp(AliHFEpriVtx)
AliHFEpriVtx::AliHFEpriVtx():
fESD1(0x0)
,fMCEvent(0x0)
,fNtrackswoPid(0)
,fHNtrackswoPid(0x0)
,fNESDprimVtxContributor(0x0)
,fNESDprimVtxIndices(0x0)
,fDiffDCAvsPt(0x0)
,fDiffDCAvsNt(0x0)
,fNsectrk2prim(0)
,fPVxRe(-999.)
,fPVyRe(-999.)
,fPVzRe(-999.)
{
Init();
}
AliHFEpriVtx::AliHFEpriVtx(const AliHFEpriVtx &p):
TObject(p)
,fESD1(0x0)
,fMCEvent(0x0)
,fNtrackswoPid(p.fNtrackswoPid)
,fHNtrackswoPid(0x0)
,fNESDprimVtxContributor(0x0)
,fNESDprimVtxIndices(0x0)
,fDiffDCAvsPt(0x0)
,fDiffDCAvsNt(0x0)
,fNsectrk2prim(p.fNsectrk2prim)
,fPVxRe(p.fPVxRe)
,fPVyRe(p.fPVyRe)
,fPVzRe(p.fPVzRe)
{
}
AliHFEpriVtx&
AliHFEpriVtx::operator=(const AliHFEpriVtx &)
{
AliInfo("Not yet implemented.");
return *this;
}
AliHFEpriVtx::~AliHFEpriVtx()
{
AliInfo("Analysis Done.");
}
void AliHFEpriVtx::Init()
{
fNtrackswoPid = 0;
for (int i=0; i<10; i++){
fPrimVtx[i].fNtrackCount = 0 ;
fPrimVtx[i].fNprimVtxContributorCount = 0 ;
}
}
void AliHFEpriVtx::CreateHistograms(TString hnopt)
{
fkSourceLabel[kAll]="all";
fkSourceLabel[kDirectCharm]="directCharm";
fkSourceLabel[kDirectBeauty]="directBeauty";
fkSourceLabel[kBeautyCharm]="beauty2charm";
fkSourceLabel[kGamma]="gamma";
fkSourceLabel[kPi0]="pi0";
fkSourceLabel[kElse]="others";
fkSourceLabel[kBeautyGamma]="beauty22gamma";
fkSourceLabel[kBeautyPi0]="beauty22pi0";
fkSourceLabel[kBeautyElse]="beauty22others";
TString hname;
for (Int_t isource = 0; isource < 10; isource++ ){
hname=hnopt+"ntracks_"+fkSourceLabel[isource];
fPrimVtx[isource].fNtracks = new TH1F(hname,hname,50,0,50);
hname=hnopt+"nPrimVtxContributor_"+fkSourceLabel[isource];
fPrimVtx[isource].fNprimVtxContributor = new TH1F(hname,hname,100,0,100);
hname=hnopt+"PtElec_"+fkSourceLabel[isource];
fPrimVtx[isource].fPtElec = new TH1F(hname,hname,250,0,50);
hname=hnopt+"PtElecContributor_"+fkSourceLabel[isource];
fPrimVtx[isource].fPtElecContributor = new TH1F(hname,hname,250,0,50);
}
hname=hnopt+"ntrackswopid";
fHNtrackswoPid = new TH1F(hname,hname,50,0,50);
hname=hnopt+"nESDprimVtxContributor";
fNESDprimVtxContributor = new TH1I(hname,hname,100,0,100);
hname=hnopt+"nESDprimVtxIndices";
fNESDprimVtxIndices= new TH1I(hname,hname,100,0,100);
hname=hnopt+"diffDCAvsPt";
fDiffDCAvsPt = new TH2F(hname,hname,250,0,50,500,0,1);
hname=hnopt+"diffDCAvsNt";
fDiffDCAvsNt = new TH2F(hname,hname,100,0,100,500,0,1);
}
void AliHFEpriVtx::CountNtracks(Int_t sourcePart, Int_t recpid, Double_t recprob)
{
fNtrackswoPid++;
if (!recpid && recprob>0.5)
fPrimVtx[kAll].fNtrackCount++;
if(sourcePart<0) return;
fPrimVtx[sourcePart].fNtrackCount++;
}
void AliHFEpriVtx::FillNtracks()
{
fHNtrackswoPid->Fill(fNtrackswoPid);
for (int i=0; i<10; i++){
fPrimVtx[i].fNtracks->Fill(fPrimVtx[i].fNtrackCount);
}
}
Int_t AliHFEpriVtx::GetMCPID(AliESDtrack const *track)
{
AliMCParticle *mctrack = NULL;
if(!(mctrack = dynamic_cast<AliMCParticle *>(fMCEvent->GetTrack(TMath::Abs(track->GetLabel()))))) return 0;
TParticle *mcpart = mctrack->Particle();
if ( !mcpart ) return 0;
Int_t pdgCode = mcpart->GetPdgCode();
return pdgCode;
}
void AliHFEpriVtx::GetNPriVxtContributor()
{
const AliESDVertex *primvtx = fESD1->GetPrimaryVertex();
fNESDprimVtxContributor->Fill(primvtx->GetNContributors());
fNESDprimVtxIndices->Fill(primvtx->GetNIndices());
}
void AliHFEpriVtx::CountPriVxtElecContributor(AliESDtrack *ESDelectron, Int_t sourcePart, Int_t recpid, Double_t recprob)
{
if (recpid || recprob<0.5) return;
Int_t elecTrkID = ESDelectron->GetID();
AliMCParticle *mctrack = NULL;
if(!(mctrack = dynamic_cast<AliMCParticle *>(fMCEvent->GetTrack(TMath::Abs(ESDelectron->GetLabel()))))) return;
TParticle *mcpart = mctrack->Particle();
if(!mcpart){
AliDebug(1, "no mc particle, return\n");
return;
}
AliKFParticle::SetField(fESD1->GetMagneticField());
AliKFParticle kfElectron(*ESDelectron,11);
AliKFVertex kfESDprimary;
const AliESDVertex *primvtx = fESD1->GetPrimaryVertex();
Int_t n=primvtx->GetNIndices();
if (n>0 && primvtx->GetStatus()){
kfESDprimary = AliKFVertex(*primvtx);
Double_t dcaBFrmElec = kfElectron.GetDistanceFromVertexXY(kfESDprimary);
UShort_t *priIndex = primvtx->GetIndices();
fPrimVtx[kAll].fPtElec->Fill(mcpart->Pt());
if(sourcePart>=0) fPrimVtx[sourcePart].fPtElec->Fill(mcpart->Pt());
for (Int_t i=0;i<n;i++){
Int_t idx = Int_t(priIndex[i]);
if (idx == elecTrkID){
fPrimVtx[kAll].fNprimVtxContributorCount++;
fPrimVtx[kAll].fPtElecContributor->Fill(mcpart->Pt());
if(sourcePart<0) continue;
fPrimVtx[sourcePart].fNprimVtxContributorCount++;
fPrimVtx[sourcePart].fPtElecContributor->Fill(mcpart->Pt());
kfESDprimary -= kfElectron;
Double_t dcaAFrmElec = kfElectron.GetDistanceFromVertexXY(kfESDprimary);
fDiffDCAvsPt->Fill(mcpart->Pt(),dcaBFrmElec-dcaAFrmElec);
fDiffDCAvsNt->Fill(n,dcaBFrmElec-dcaAFrmElec);
}
}
}
}
void AliHFEpriVtx::FillNprimVtxContributor() const
{
for (int i=0; i<10; i++){
fPrimVtx[i].fNprimVtxContributor->Fill(fPrimVtx[i].fNprimVtxContributorCount);
}
}
Double_t AliHFEpriVtx::GetDistanceFromRecalVertexXY(const AliESDtrack * const ESDelectron)
{
Int_t elecTrkID = ESDelectron->GetID();
AliKFParticle::SetField(fESD1->GetMagneticField());
AliKFParticle kfElectron(*ESDelectron,11);
AliKFVertex kfESDprimary;
const AliESDVertex *primvtx = fESD1->GetPrimaryVertex();
Int_t n=primvtx->GetNIndices();
if (n>0 && primvtx->GetStatus()){
kfESDprimary = AliKFVertex(*primvtx);
UShort_t *priIndex = primvtx->GetIndices();
for (Int_t i=0;i<n;i++){
Int_t idx = Int_t(priIndex[i]);
if (idx == elecTrkID){
kfESDprimary -= kfElectron;
Double_t dcaAFrmElec = kfElectron.GetDistanceFromVertexXY(kfESDprimary);
return dcaAFrmElec;
}
}
}
return -1;
}
void AliHFEpriVtx::RecalcPrimvtx(Int_t nkftrk, const Int_t * const trkid, const AliKFParticle * const kftrk)
{
const AliESDVertex *primvtx = fESD1->GetPrimaryVertex();
AliKFVertex kfESDprimary;
Int_t n = primvtx->GetNIndices();
fNsectrk2prim = 0;
fPVxRe = -999.;
fPVyRe = -999.;
fPVyRe = -999.;
if (n>0 && primvtx->GetStatus()){
kfESDprimary = AliKFVertex(*primvtx);
UShort_t *priIndex = primvtx->GetIndices();
for (Int_t j=0; j<nkftrk; j++){
for (Int_t i=0;i<n;i++){
Int_t idx = Int_t(priIndex[i]);
if (idx == trkid[j]){
kfESDprimary -= kftrk[j];
fNsectrk2prim++;
}
}
}
}
fPVxRe = kfESDprimary.GetX();
fPVyRe = kfESDprimary.GetY();
fPVzRe = kfESDprimary.GetZ();
}
void AliHFEpriVtx::RecalcPrimvtx(const AliESDtrack * const ESDelectron)
{
Int_t elecTrkID = ESDelectron->GetID();
AliKFParticle::SetField(fESD1->GetMagneticField());
AliKFParticle kfElectron(*ESDelectron,11);
const AliESDVertex *primvtx = fESD1->GetPrimaryVertex();
AliKFVertex kfESDprimary;
Int_t n = primvtx->GetNIndices();
fPVxRe = -999.;
fPVyRe = -999.;
fPVyRe = -999.;
if (n>0 && primvtx->GetStatus()){
kfESDprimary = AliKFVertex(*primvtx);
UShort_t *priIndex = primvtx->GetIndices();
for (Int_t i=0;i<n;i++){
Int_t idx = Int_t(priIndex[i]);
if (idx == elecTrkID){
kfESDprimary -= kfElectron;
}
}
}
fPVxRe = kfESDprimary.GetX();
fPVyRe = kfESDprimary.GetY();
fPVzRe = kfESDprimary.GetZ();
}