#include "hiterator.h"
#include "hcategory.h"
#include "hpidefficiencycalculation.h"
#include "hpidtrackcand.h"
#include "hpidparticlesim.h"
#include "hpidgeanttrackset.h"
#include "hpairsim.h"
#include "hshowerhittof.h"
#include "htofhit.h"
#include "hpidcandidate.h"
#include "piddef.h"
#include "hiterator.h"
#include "hcategory.h"
#include "TNtuple.h"
#include "TVector3.h"
#include "hgeantkine.h"
#include "hkicktrack.h"
#include <iostream>
#include "piddef.h"
#include "hpidfl.h"
#include "hpairgeantdata.h"
#include "hrichhitsim.h"
#include "hmdcsegsim.h"
#include "htofhitsim.h"
#include "hshowerhittoftrack.h"
#include "richdef.h"
#include "hmdcdef.h"
#include "piddef.h"
#include "tofdef.h"
#include "hpidtrackcandsim.h"
#include "hmdctrkcand.h"
#include "showertofinodef.h"
#include "hmdctrackgdef.h"
#include "hmetamatch.h"
#include "hgeantmdc.h"
#include "hgeanttof.h"
#include "hgeantshower.h"
using namespace std;
ClassImp(HPidEfficiencyCalculation)
HPidEfficiencyCalculation::HPidEfficiencyCalculation(Int_t SELECTED_MOMENTUM_ALG,
Int_t o_sourceID,
Int_t o_targetmedium,
const Char_t* cOutFileName)
:HSUDummyRec(kFALSE){
clearMembers();
sOutFileName = cOutFileName;
selectedMomAlg = SELECTED_MOMENTUM_ALG;
sourceID=o_sourceID;
targetmedium=o_targetmedium;
}
HPidEfficiencyCalculation::~HPidEfficiencyCalculation(){
if(NULL!=effNtuple){
printf("Ntuple not null!!!\n");
}
}
void HPidEfficiencyCalculation::clearMembers(){
pItGeant = NULL;
pItGeant1 = NULL;
pItPart = NULL;
pItPairs = NULL;
pItMdcSeg = NULL;
pItRichHit = NULL;
pItTOFHit = NULL;
pItShowerHitTof= NULL;
effNtuple = NULL;
resetNtuple();
}
Bool_t HPidEfficiencyCalculation::init(void){
pItMdcGeant = getIterator(catMdcGeantRaw,kTRUE);
if( NULL==pItMdcGeant){
return kFALSE;
}
pItMdcTrkCand = getIterator(catMdcTrkCand,kTRUE);
if( NULL==pItMdcTrkCand){
return kFALSE;
}
pItMetaMatch = getIterator(catMetaMatch,kTRUE);
if( NULL==pItMetaMatch){
return kFALSE;
}
pItPidTrackCand = getIterator(catPidTrackCand,kTRUE);
if( NULL==pItPidTrackCand){
return kFALSE;
}
pItGeant = getIterator(catGeantKine,kTRUE);
if( NULL==pItGeant){
return kFALSE;
}
pItGeant1 = getIterator(catGeantKine,kTRUE);
if( NULL==pItGeant1){
return kFALSE;
}
pItPart = getIterator(catPidPart,kTRUE);
if( NULL==pItPart){
return kFALSE;
}
pItPairs = getIterator(catPairs,kTRUE);
if( NULL==pItPairs){
return kFALSE;
}
pItMdcSeg = getIterator(catMdcSeg,kTRUE);
if( NULL==pItMdcSeg){
return kFALSE;
}
pItRichHit = getIterator(catRichHit,kTRUE);
if( NULL==pItRichHit){
return kFALSE;
}
pItTOFHit = getIterator(catTofHit,kTRUE);
if( NULL==pItTOFHit){
return kFALSE;
}
pItTOFGeant = getIterator(catTofGeantRaw,kTRUE);
if( NULL==pItTOFGeant){
return kFALSE;
}
pItShowerGeant = getIterator(catShowerGeantRaw,kTRUE);
if( NULL==pItShowerGeant){
return kFALSE;
}
pItShowerHitTof = getIterator(catShowerHitTofTrack,kTRUE);
if( NULL==pItShowerHitTof){
return kFALSE;
}
if(openOutFile(sOutFileName.Data()) == NULL){
return kFALSE;
}
Char_t cBuffer[200000];
snprintf(cBuffer,sizeof(cBuffer),"%s:%s:%s:%s:%s:%s:%s:%s:%s:%s:%s",
EFF_NTUPLE_VARS_1,EFF_NTUPLE_VARS_2,EFF_NTUPLE_VARS_3,
EFF_NTUPLE_VARS_4,EFF_NTUPLE_VARS_5,EFF_NTUPLE_VARS_6,
EFF_NTUPLE_VARS_7,EFF_NTUPLE_VARS_8,EFF_NTUPLE_VARS_9,
EFF_NTUPLE_VARS_10,EFF_NTUPLE_VARS_11);
printf("Ntuple =%s\n",cBuffer);
if((effNtuple=new TNtuple("E","Ntuple with efficiency data",cBuffer))==NULL){
Error("init", "Cannot create NTuple");
return kFALSE;
}
setInitOk();
return kTRUE;
}
Int_t HPidEfficiencyCalculation::checkGeantTrackSet(const HPidGeantTrackSet* pGTS)
{
Int_t nCorrTracks=pGTS->getNCorrelatedTrackIds();
for(Int_t i=0;i<nCorrTracks-1;i++)
{
if(pGTS->getCorrelationFlag(i)<pGTS->getCorrelationFlag(i+1))
{
}
}
return 0;
}
Int_t HPidEfficiencyCalculation::execute(void) {
pItGeant->Reset();
const HPidHitData *pHitData;
const HPidTrackData *pTrackData;
HPidParticleSim *pPart;
HPidTrackCandSim *pPidTrackCandSim;
HMetaMatch *pMetaMatch;
HPairSim *pPairSim;
HGeantKine *pKine;
Int_t countPart=0;
Int_t s;
s=0;
while((pKine=(HGeantKine*)pItGeant->Next()))
{
countPart++;
if(pKine->getID()!=2 && pKine->getID()!=3)
{
continue;
}
resetNtuple();
if(!isEmbeddedTrack(pKine))
{
continue;
}
if(isEmbeddedTrack(pKine))
{
eNtupleVars[E_isEmbeddedTrack]=kTRUE;
}
if(!isInAcceptance(pKine))
{
continue;
}
if(isInAcceptance(pKine))
{
eNtupleVars[E_isAccepted]=kTRUE;
}
if(isCloseToOtherParticle(pKine))
{
eNtupleVars[E_isCloseToOtherParticle]=kTRUE;
if(isCloseToEmbeddedParticle(pKine))
{
eNtupleVars[E_isCloseToEmbeddedParticle]=kTRUE;
}
}
eNtupleVars[E_track_number]=pKine->getTrack();
eNtupleVars[E_parent_track_number]=pKine->getParentTrack();
eNtupleVars[E_geant_id] = pKine->getID();
Int_t parent, medium, mechanism;
pKine->getCreator(parent, medium, mechanism);
eNtupleVars[E_medium] = medium;
eNtupleVars[E_mechanism] = mechanism;
eNtupleVars[E_geant_mom] = pKine->getTotalMomentum();
pKine->getMomentum(eNtupleVars[E_geant_mom_x],eNtupleVars[E_geant_mom_y],eNtupleVars[E_geant_mom_z]);
TVector3 geantMomentum(eNtupleVars[E_geant_mom_x],eNtupleVars[E_geant_mom_y],eNtupleVars[E_geant_mom_z]);
eNtupleVars[E_geant_theta] = geantMomentum.Theta();
eNtupleVars[E_geant_phi] = geantMomentum.Phi();
eNtupleVars[E_metaIndex]=findCorrespondingMETAHit(pKine->getTrack());
eNtupleVars[E_ringIndex]=findCorrespondingRichRing(pKine->getTrack());
pItPidTrackCand->Reset();
while((pPidTrackCandSim = (HPidTrackCandSim*)pItPidTrackCand->Next()))
{
const HPidGeantTrackSet* pTrackSetCand = pPidTrackCandSim->getGeantTrackSet();
for(UInt_t c = 0; c<pTrackSetCand->getNCorrelatedTrackIds();c++)
{
if(pTrackSetCand->getGeantTrackID(c)==eNtupleVars[E_track_number])
{
if(pTrackSetCand->getCorrelationFlag(c)>=eNtupleVars[E_PidTrackCandIndex])
{
eNtupleVars[E_PidTrackCandIndex]=pTrackSetCand->getCorrelationFlag(c);
eNtupleVars[E_RKIndex]=pPidTrackCandSim->getTrackData()->getRKTrackInd();
eNtupleVars[E_KickIndex]=pPidTrackCandSim->getTrackData()->getKickTrackInd();
break;
}
}
}
}
eNtupleVars[E_innerSegmentIndex]=findCorrespondingMDCSegment(pKine->getTrack(),0);
eNtupleVars[E_outerSegmentIndex]=findCorrespondingMDCSegment(pKine->getTrack(),1);
pItMetaMatch->Reset();
Int_t idx=-1;
while((pMetaMatch = (HMetaMatch*)pItMetaMatch->Next()))
{
idx++;
if(pMetaMatch->getSystem()!=eNtupleVars[E_system]) continue;
if(pMetaMatch->getMetaHitInd()!=eNtupleVars[E_metaIndex]) continue;
if(pMetaMatch->getNCandForRich()<=0) continue;
if(pMetaMatch->getARichInd(0)!=eNtupleVars[E_ringIndex] &&
pMetaMatch->getARichInd(1)!=eNtupleVars[E_ringIndex] &&
pMetaMatch->getARichInd(2)!=eNtupleVars[E_ringIndex]) continue;
eNtupleVars[E_matchIndex]=idx;
break;
}
UInt_t maximumCorrelation=0;
pItPart->Reset();
while((pPart = (HPidParticleSim*)pItPart->Next()) && eNtupleVars[E_isIdentified]==kFALSE)
{
const HPidGeantTrackSet* pTrackSet = pPart->getGeantTrackSet();
for(UInt_t c = 0; c<pTrackSet->getNCorrelatedTrackIds();c++)
{
if(pTrackSet->getGeantTrackID(c)==eNtupleVars[E_track_number])
{
eNtupleVars[E_isPartOfATrack]=kTRUE;
if(isSeenInSingleDetector(pTrackSet->getCorrelationFlag(c)))
{
break;
}
if(pTrackSet->getCorrelationFlag(c)>=maximumCorrelation || isHiddenLepton(pTrackSet))
{
eNtupleVars[E_geant_flag]=pTrackSet->getCorrelationFlag(c);
maximumCorrelation=pTrackSet->getCorrelationFlag(c);
eNtupleVars[E_isFullyReconstructed]=kFALSE;
if(eNtupleVars[E_geant_flag]>=76 || isHiddenLepton(pTrackSet))
{
eNtupleVars[E_isFullyReconstructed] = kTRUE;
}
pHitData = pPart->getHitData();
pTrackData = pPart->getTrackData();
eNtupleVars[E_charge] = pPart->getCharge();
eNtupleVars[E_rec_mom] = pPart->P();
eNtupleVars[E_beta] = pPart->getBetaExp();
eNtupleVars[E_sector] = pHitData->getSector();
eNtupleVars[E_theta] = pPart->thetaDeg();
eNtupleVars[E_phi] = pPart->phiDeg();
eNtupleVars[E_rec_mom_kick]= pTrackData->getMomenta(ALG_KICK);
eNtupleVars[E_rec_mom_spline]= pTrackData->getMomenta(ALG_SPLINE);
eNtupleVars[E_rec_mom_rk]= pTrackData->getMomenta(ALG_RUNGEKUTTA);
eNtupleVars[E_inner_chi2] = pHitData->fInnerMdcChiSquare;
eNtupleVars[E_outer_chi2] = pHitData->fOuterMdcChiSquare;
eNtupleVars[E_rich_corr_kick] = pHitData->hasRingCorrelation[ALG_KICK];
eNtupleVars[E_rich_corr_spline] = pHitData->hasRingCorrelation[ALG_SPLINE];
eNtupleVars[E_rich_corr_rk] = pHitData->hasRingCorrelation[ALG_RUNGEKUTTA];
eNtupleVars[E_rich_flag] = pHitData->getFlagRICH();
eNtupleVars[E_rec_id]= pPart->getPid();
if(pKine->getID()==pPart->getPid() &&
pHitData->hasRingCorrelation[selectedMomAlg])
{
eNtupleVars[E_isIdentified] = kTRUE;
eNtupleVars[E_isRemovedByCPCandidateCutL]=isRemovedByCPCandidateCut(pPart,kTRUE);
eNtupleVars[E_isRemovedByCPCandidateCutH]=isRemovedByCPCandidateCut(pPart,kFALSE);
break;
}
}
}
}
}
if(isRemovedByActivePairCut(pKine))
{
eNtupleVars[E_isRemovedByActivePairCut]=kTRUE;
}
pItPairs->Reset();
Int_t cc = 0;
while((pPairSim = (HPairSim*)pItPairs->Next()) != NULL)
{
cc++;
Int_t correlationFlag = isPairedWithBackgroundTrack(pPairSim,pKine);
if(correlationFlag>0)
{
}
if(correlationFlag>=76)
{
eNtupleVars[E_isLegOfPair]=kTRUE;
if(pPairSim->getIsCutFlag()==0)
{
if(correlationFlag>=76 && eNtupleVars[E_isLegOfSurvivingPair]==kFALSE)
{
eNtupleVars[E_isLegOfSurvivingPair]=kTRUE;
}
if(correlationFlag<76)
{
}
}
}
}
effNtuple->Fill(eNtupleVars);
}
return 0;
}
Bool_t HPidEfficiencyCalculation::finalize(void) {
return writeAndCloseOutFile();
}
void HPidEfficiencyCalculation::resetNtuple() {
for(Int_t n=0;n<MAX_VARS;n++)
{
eNtupleVars[n]=-99;
}
eNtupleVars[E_isPartOfATrack] = kFALSE;
eNtupleVars[E_isFullyReconstructed] = kFALSE;
eNtupleVars[E_isIdentified]=kFALSE;
eNtupleVars[E_isAccepted]=kFALSE;
eNtupleVars[E_isCloseToOtherParticle]=kFALSE;
eNtupleVars[E_isEmbeddedTrack]=kFALSE;
eNtupleVars[E_isCloseToEmbeddedParticle]=kFALSE;
eNtupleVars[E_isRemovedByActivePairCut]=kFALSE;
eNtupleVars[E_isLegOfPair]=kFALSE;
eNtupleVars[E_isLegOfSurvivingPair]=kFALSE;
eNtupleVars[E_isRemovedByCPCandidateCutL]=kFALSE;
eNtupleVars[E_isRemovedByCPCandidateCutH]=kFALSE;
}
Bool_t HPidEfficiencyCalculation::isInAcceptance(HGeantKine* pKine)
{
Int_t nStatMDC1, nStatMDC2, nStatMDC3, nStatMDC4;
Int_t nStatTof, nStatShower;
Int_t lTrack,lId;
pKine->getParticle(lTrack,lId);
nStatMDC1 = nStatMDC2 = nStatMDC3 = nStatMDC4 = 0;
nStatTof = nStatShower = 0;
HGeantMdc *pMdc;
pItMdcGeant->Reset();
while((pMdc = (HGeantMdc*) pItMdcGeant->Next()) != NULL){
if (pMdc->getTrack() == lTrack)
{
switch (((Int_t)pMdc->getModule()))
{
case 0: nStatMDC1 = 1;
break;
case 1: nStatMDC2 = 1;
break;
case 2: nStatMDC3 = 1;
break;
case 3: nStatMDC4 = 1;
break;
default: cerr << "WRONG MDC module number!" << endl;
}
}
}
HGeantTof *pTof;
pItTOFGeant->Reset();
while((pTof = (HGeantTof*) pItTOFGeant->Next()) != NULL){
if (pTof->getTrack() == lTrack)
{
nStatTof = 1;
}
}
HGeantShower *pShower;
pItShowerGeant->Reset();
while((pShower = (HGeantShower*) pItShowerGeant->Next()) != NULL){
if (pShower->getTrack() == lTrack)
{
nStatShower = 1;
}
}
if (nStatMDC1 && nStatMDC2 && nStatMDC3 && (nStatTof || nStatShower))
{
return kTRUE;
}
return kFALSE;
}
Int_t HPidEfficiencyCalculation::findCorrespondingMDCSegment(Int_t geantTrackNumber,Bool_t bInnerSeg)
{
pItMdcSeg->Reset();
HMdcSegSim* pSeg=NULL;
Int_t idx=0;
while((pSeg=(HMdcSegSim*)pItMdcSeg->Next()))
{
for(Int_t i = 0; i < pSeg->getNTracks(); i++)
{
if(pSeg->getTrack(i) == geantTrackNumber)
{
if(pSeg->getIOSeg()==0 && bInnerSeg)
{
return idx;
}
if(pSeg->getIOSeg()==1 && !bInnerSeg)
{
return idx;
}
}
}
idx++;
}
return -99;
}
Int_t HPidEfficiencyCalculation::findCorrespondingMETAHit(Int_t geantTrackNumber)
{
pItTOFHit->Reset();
HTofHitSim* pTof=NULL;
Int_t idx=0;
while((pTof=(HTofHitSim*)pItTOFHit->Next()))
{
if(pTof->getNTrack1()== geantTrackNumber)
{
eNtupleVars[E_system]=1;
return idx;
}
if(pTof->getNTrack2()== geantTrackNumber)
{
eNtupleVars[E_system]=1;
return idx;
}
idx++;
}
idx=0;
pItShowerHitTof->Reset();
HShowerHitTofTrack* pShw=NULL;
while((pShw=(HShowerHitTofTrack*)pItShowerHitTof->Next()))
{
if(pShw->getTrack()== geantTrackNumber)
{
eNtupleVars[E_system]=0;
return idx;
}
idx++;
}
return -99;
}
Int_t HPidEfficiencyCalculation::findCorrespondingRichRing(Int_t geantTrackNumber)
{
pItRichHit->Reset();
HRichHitSim* pRich=NULL;
Int_t idx=0;
while((pRich=(HRichHitSim*)pItRichHit->Next()))
{
if(pRich->track1 == geantTrackNumber ||
pRich->track2 == geantTrackNumber ||
pRich->track3 == geantTrackNumber )
{
return idx;
}
idx++;
}
return -99;
}
Bool_t HPidEfficiencyCalculation::isCloseToOtherParticle(HGeantKine* pKine,Bool_t b_UseSourceParticlesOnly)
{
HGeantKine *pKine1=NULL;
Int_t trackNumber,particleID,trackNumber1,particleID1;
Int_t parentTrackNumber,medium,mechanism;
Int_t parentTrackNumber1,medium1,mechanism1;
parentTrackNumber1 = medium1 = mechanism1 = -1;
pKine->getParticle(trackNumber,particleID);
pKine->getCreator(parentTrackNumber,medium,mechanism);
pItGeant1->Reset();
while((pKine1 = (HGeantKine*) pItGeant1->Next()) != NULL){
pKine1->getParticle(trackNumber1,particleID1);
if(trackNumber!=trackNumber1) {
Double_t diff_Opang=calcOpeningAngleKine(pKine,pKine1);
if (diff_Opang<9.5)
{
if(b_UseSourceParticlesOnly)
{
if(isEmbeddedTrack(pKine1))
{
return kTRUE;
}
}
else
{
return kTRUE;
}
}
}
}
return kFALSE;
}
Bool_t HPidEfficiencyCalculation::isCloseToEmbeddedParticle(HGeantKine* pKine)
{
return isCloseToOtherParticle(pKine,kTRUE);
}
Bool_t HPidEfficiencyCalculation::isEmbeddedTrack(HGeantKine* pKine)
{
Int_t parentTrackNumber,medium,mechanism;
pKine->getCreator(parentTrackNumber,medium,mechanism);
Int_t pid = pKine->getID();
if(sourceID==URQMD_PURE)
{
if(medium==targetmedium && (pid==2 || pid==3))
{
return kTRUE;
}
else
{
return kFALSE;
}
}
if(sourceID==URQMD_EMBEDDED)
{
if(parentTrackNumber==0 && (pid==2 || pid==3))
{
return kTRUE;
}
else
{
return kFALSE;
}
}
return kFALSE;
}
Bool_t HPidEfficiencyCalculation::isHiddenLepton(const HPidGeantTrackSet* pGeantTrackSet)
{
Int_t itsPrimaryParticleID = pGeantTrackSet->getGeantPID(0);
Int_t itsPrimaryCorrelationFlag = pGeantTrackSet->getCorrelationFlag(0);
Int_t itsPrimaryProcessID;
itsPrimaryProcessID = pGeantTrackSet->getGeantProcessID(0);
Int_t itsSecondaryParticleID = -99;
Int_t itsSecondaryCorrelationFlag = -99;
Int_t itsSecondaryProcessID = -99;
Int_t itsTertiaryParticleID = -99;
Int_t itsTertiaryCorrelationFlag = -99;
Int_t itsTertiaryProcessID = -99;
if(pGeantTrackSet->getNCorrelatedTrackIds()>1)
{
itsSecondaryParticleID = pGeantTrackSet->getGeantPID(1);
itsSecondaryCorrelationFlag = pGeantTrackSet->getCorrelationFlag(1);
itsSecondaryProcessID = pGeantTrackSet->getGeantProcessID(1);
}
if(pGeantTrackSet->getNCorrelatedTrackIds()>2)
{
itsTertiaryParticleID = pGeantTrackSet->getGeantPID(2);
itsTertiaryCorrelationFlag = pGeantTrackSet->getCorrelationFlag(2);
itsTertiaryProcessID = pGeantTrackSet->getGeantProcessID(2);
}
if( itsPrimaryParticleID == 2 || itsPrimaryParticleID == 3)
{
if(itsPrimaryCorrelationFlag == 14 || itsPrimaryCorrelationFlag == 12)
{
if((itsSecondaryProcessID==10 || itsSecondaryProcessID == 6 || itsSecondaryProcessID == 7) && pGeantTrackSet->getNCorrelatedTrackIds()>1)
{
if(itsSecondaryCorrelationFlag == 64 && (itsSecondaryParticleID==2 || itsSecondaryParticleID==3))
{
return kTRUE;
}
}
if((itsTertiaryProcessID==10 || itsTertiaryProcessID == 6 || itsTertiaryProcessID == 7)&& pGeantTrackSet->getNCorrelatedTrackIds()>2)
{
if(itsTertiaryCorrelationFlag == 64 && (itsTertiaryParticleID==2 || itsTertiaryParticleID==3))
{
return kTRUE;
}
}
}
}
return kFALSE;
}
Bool_t HPidEfficiencyCalculation::isEmbeddedTrack(Int_t geantTrackNumber)
{
HGeantKine* pKine=HPidFL::getGeantKineObjFromTrkNbFast(geantTrackNumber);
return isEmbeddedTrack(pKine);
}
Bool_t HPidEfficiencyCalculation::isBackgroundTrack(Int_t geantTrackNumber)
{
return !isEmbeddedTrack(geantTrackNumber);
}
Float_t HPidEfficiencyCalculation::calcOpeningAngleKine(HGeantKine *kine1,HGeantKine *kine2)
{
if(kine1==kine2)
{
}
TVector3 vec1;
if (kine1)
{
Float_t xMom1,yMom1,zMom1;
kine1->getMomentum(xMom1,yMom1,zMom1);
vec1.SetX(xMom1);
vec1.SetY(yMom1);
vec1.SetZ(zMom1);
}
else
{
::Error("HPidEfficiencyCalcultion::calcOpeningAngleKine","first pointer is 0");
}
TVector3 vec2;
if (kine2)
{
Float_t xMom2,yMom2,zMom2;
kine2->getMomentum(xMom2,yMom2,zMom2);
vec2.SetX(xMom2);
vec2.SetY(yMom2);
vec2.SetZ(zMom2);
}
else
{
::Error("HPidEfficiencyCalcultion::calcOpeningAngleKine","second pointer is 0");
}
return TMath::RadToDeg()*vec1.Angle(vec2);
}
Bool_t HPidEfficiencyCalculation::isRemovedByCPCandidateCut(HPidParticle* pPart,Bool_t bSelectLeptons)
{
if(bSelectLeptons)
{
}
if (pPart->getTrackData()->getAngleWithClosestCandidate(-1,0)<9 && !bSelectLeptons)
{
return kTRUE;
}
if (pPart->getTrackData()->getAngleWithClosestCandidate(1,0)<9 && bSelectLeptons)
{
return kTRUE;
}
return kFALSE;
}
Bool_t HPidEfficiencyCalculation::isRemovedByActivePairCut(HGeantKine* pKine)
{
if(!pKine)
{
::Error("HPidEfficiencyCalculation::isRemovedByActivePairCut","kine-pointer==NULL");
}
Int_t trackNumber,particleID;
pKine->getParticle(trackNumber,particleID);
HPairSim * pPair=NULL;
pItPairs->Reset();
Int_t comDet=0;
Int_t paircnt;
paircnt = 0;
while(( pPair= (HPairSim *)pItPairs->Next()) != NULL)
{
comDet = isPairedWithBackgroundTrack(pPair,pKine);
if(comDet>=76){
if(pPair->getPid1()==3 ||
pPair->getPid1()==2 ||
pPair->getPid2()==3 ||
pPair->getPid2()==2)
{
if(pPair->getIsCutFlag()!=0){
return kTRUE;
}
}
}
}
return kFALSE;
}
Bool_t HPidEfficiencyCalculation::isDileptonPair(HPairSim *pair){
if((pair->getPid1()==2 || pair->getPid1()==3) && (pair->getPid2()==2 || pair->getPid2()==3))
{
return kTRUE;
}
else
{
return kFALSE;
}
}
Int_t HPidEfficiencyCalculation::isPairedWithBackgroundLepton(HPairSim* pPair, HGeantKine* pKine)
{
if(isDileptonPair(pPair))
{
return isPairedWithBackgroundTrack(pPair, pKine);
}
return kFALSE;
}
Int_t HPidEfficiencyCalculation::isPairedWithBackgroundTrack(HPairSim* pPair, HGeantKine* pKine)
{
Int_t geantTrackNumber = pKine->getTrack();
HPairGeantData pg(pPair);
if(pg.getTrackNumber1()==geantTrackNumber && isEmbeddedTrack(pKine))
{
if(isBackgroundTrack(pg.getTrackNumber2()))
{
return pg.getCommonDetectors1();
}
else return 0;
}
else if(pg.getTrackNumber2()==geantTrackNumber && isEmbeddedTrack(pKine))
{
if(isBackgroundTrack(pg.getTrackNumber1()))
{
return pg.getCommonDetectors2();
}
else return 0;
}
return 0;
}
Bool_t HPidEfficiencyCalculation::isSeenInSingleDetector(Int_t correlationFlag)
{
if(correlationFlag==1 ||
correlationFlag==2 ||
correlationFlag==4 ||
correlationFlag==8 ||
correlationFlag==64 )
{
return kTRUE;
}
return kFALSE;
}
Last change: Sat May 22 13:07:02 2010
Last generated: 2010-05-22 13:07
This page has been automatically generated. If you have any comments or suggestions about the page layout send a mail to ROOT support, or contact the developers with any questions or problems regarding ROOT.