Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Bug fixes to PU jet ID variables #15119

Merged
merged 7 commits into from Jul 15, 2016
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Jump to
Jump to file
Failed to load files.
Diff view
Diff view
25 changes: 25 additions & 0 deletions RecoJets/JetProducers/python/PileupJetIDCutParams_cfi.py
@@ -1,5 +1,30 @@
import FWCore.ParameterSet.Config as cms

###########################################################
## Working points for the 81X training (completed in 80X with variable fixes)
###########################################################
full_81x_chs_wp = cms.PSet(
#4 Eta Categories 0-2.5 2.5-2.75 2.75-3.0 3.0-5.0

#Tight Id
Pt010_Tight = cms.vdouble( 0.69, -0.35, -0.26, -0.21),
Pt1020_Tight = cms.vdouble( 0.69, -0.35, -0.26, -0.21),
Pt2030_Tight = cms.vdouble( 0.69, -0.35, -0.26, -0.21),
Pt3050_Tight = cms.vdouble( 0.86, -0.10, -0.05, -0.01),

#Medium Id
Pt010_Medium = cms.vdouble( 0.18, -0.55, -0.42, -0.36),
Pt1020_Medium = cms.vdouble( 0.18, -0.55, -0.42, -0.36),
Pt2030_Medium = cms.vdouble( 0.18, -0.55, -0.42, -0.36),
Pt3050_Medium = cms.vdouble( 0.61, -0.35, -0.23, -0.17),

#Loose Id
Pt010_Loose = cms.vdouble(-0.97, -0.68, -0.53, -0.47),
Pt1020_Loose = cms.vdouble(-0.97, -0.68, -0.53, -0.47),
Pt2030_Loose = cms.vdouble(-0.97, -0.68, -0.53, -0.47),
Pt3050_Loose = cms.vdouble(-0.89, -0.52, -0.38, -0.30)
)

###########################################################
## Working points for the 80X training
###########################################################
Expand Down
104 changes: 104 additions & 0 deletions RecoJets/JetProducers/python/PileupJetIDParams_cfi.py
@@ -1,6 +1,110 @@
import FWCore.ParameterSet.Config as cms
from RecoJets.JetProducers.PileupJetIDCutParams_cfi import *

####################################################################################################################
full_81x_chs = cms.PSet(
impactParTkThreshold = cms.double(1.),
cutBased = cms.bool(False),
etaBinnedWeights = cms.bool(True),
tmvaMethod = cms.string("JetIDMVAHighPt"),
version = cms.int32(-1),
nEtaBins = cms.int32(4),
trainings = cms.VPSet(
cms.PSet(
jEtaMin = cms.double(0.),
jEtaMax = cms.double(2.5),
tmvaWeights = cms.string("RecoJets/JetProducers/data/pileupJetId_80XvarFix_Eta0to2p5_BDT.weights.xml.gz"),
tmvaVariables = cms.vstring(
"nvtx",
"dR2Mean" ,
"nParticles" ,
"nCharged" ,
"majW" ,
"minW",
"frac01" ,
"frac02" ,
"frac03" ,
"frac04" ,
"ptD" ,
"beta" ,
"pull" ,
"jetR" ,
"jetRchg" ,
)
),
cms.PSet(
jEtaMin = cms.double(2.5),
jEtaMax = cms.double(2.75),
tmvaWeights = cms.string("RecoJets/JetProducers/data/pileupJetId_80XvarFix_Eta2p5to2p75_BDT.weights.xml.gz"),
tmvaVariables = cms.vstring(
"nvtx",
"dR2Mean" ,
"nParticles" ,
"nCharged" ,
"majW" ,
"minW",
"frac01" ,
"frac02" ,
"frac03" ,
"frac04" ,
"ptD" ,
"beta" ,
"pull" ,
"jetR" ,
"jetRchg" ,
)
),
cms.PSet(
jEtaMin = cms.double(2.75),
jEtaMax = cms.double(3.),
tmvaWeights = cms.string("RecoJets/JetProducers/data/pileupJetId_80XvarFix_Eta2p75to3_BDT.weights.xml.gz"),
tmvaVariables = cms.vstring(
"nvtx",
"dR2Mean" ,
"nParticles" ,
"nCharged" ,
"majW" ,
"minW",
"frac01" ,
"frac02" ,
"frac03" ,
"frac04" ,
"ptD" ,
"beta" ,
"pull" ,
"jetR" ,
"jetRchg" ,
)
),
cms.PSet(
jEtaMin = cms.double(3.),
jEtaMax = cms.double(5.),
tmvaWeights = cms.string("RecoJets/JetProducers/data/pileupJetId_80XvarFix_Eta3to5_BDT.weights.xml.gz"),
tmvaVariables = cms.vstring(
"nvtx",
"dR2Mean" ,
"nParticles" ,
"majW" ,
"minW",
"frac01" ,
"frac02" ,
"frac03" ,
"frac04" ,
"ptD" ,
"pull" ,
"jetR" ,
)
),
),
tmvaSpectators = cms.vstring(
"jetPt" ,
"jetEta" ,
),
JetIdParams = full_81x_chs_wp,
label = cms.string("full")
)


####################################################################################################################
full_80x_chs = cms.PSet(
impactParTkThreshold = cms.double(1.),
Expand Down
3 changes: 2 additions & 1 deletion RecoJets/JetProducers/python/PileupJetID_cfi.py
Expand Up @@ -10,8 +10,9 @@
_chsalgos_74x = cms.VPSet(full_74x_chs,cutbased)
_chsalgos_76x = cms.VPSet(full_76x_chs,cutbased)
_chsalgos_80x = cms.VPSet(full_80x_chs,cutbased)
_chsalgos_81x = cms.VPSet(full_81x_chs,cutbased)

_stdalgos = _chsalgos_80x
_stdalgos = _chsalgos_81x

# Calculate+store variables and run MVAs
pileupJetId = cms.EDProducer('PileupJetIdProducer',
Expand Down
40 changes: 30 additions & 10 deletions RecoJets/JetProducers/src/PileupJetIdAlgo.cc
Expand Up @@ -338,7 +338,7 @@ PileupJetIdentifier PileupJetIdAlgo::computeIdVariables(const reco::Jet * jet, f
float candPt = icand->pt();
float candPtFrac = candPt/jetPt;
float candDr = reco::deltaR(*icand,*jet);
float candDeta = std::abs( icand->eta() - jet->eta() );
float candDeta = icand->eta() - jet->eta();
float candDphi = reco::deltaPhi(*icand,*jet);
float candPtDr = candPt * candDr;
size_t icone = std::lower_bound(&cones[0],&cones[ncones],candDr) - &cones[0];
Expand Down Expand Up @@ -428,7 +428,7 @@ PileupJetIdentifier PileupJetIdAlgo::computeIdVariables(const reco::Jet * jet, f
sumTkPt += tkpt;
// 'classic' beta definition based on track-vertex association
bool inVtx0 = vtx->trackWeight ( lPF->trackRef()) > 0 ;

bool inAnyOther = false;
// alternative beta definition based on track-vertex distance of closest approach
double dZ0 = std::abs(lPF->trackRef()->dz(vtx->position()));
Expand Down Expand Up @@ -459,21 +459,41 @@ PileupJetIdentifier PileupJetIdAlgo::computeIdVariables(const reco::Jet * jet, f
internalId_.betaStar_ += tkpt;
}
}
}
}
else{
// setting classic and alternative to be the same for now
float tkpt = candPt;
sumTkPt += tkpt;
bool inVtx0 = false;
bool inVtxOther = false;
if (lPack->fromPV() == pat::PackedCandidate::PVUsedInFit) inVtx0 = true;
if (lPack->fromPV() == 0) inVtxOther = true;
double dZ0=9999.;
double dZ_tmp = 9999.;
for (unsigned vtx_i = 0 ; vtx_i < allvtx.size() ; vtx_i++ ) {
auto iv = allvtx[vtx_i];

if (iv.isFake())
continue;

// Match to vertex in case of copy as above
bool isVtx0 = (iv.position() - vtx->position()).r() < 0.02;

if (isVtx0) {
if (lPack->fromPV(vtx_i) == pat::PackedCandidate::PVUsedInFit) inVtx0 = true;
if (lPack->fromPV(vtx_i) == 0) inVtxOther = true;
dZ0 = lPack->dz(vtx_i);
}

if (fabs(lPack->dz(iv.position())) < fabs(dZ_tmp)) {
dZ_tmp = lPack->dz(iv.position());
}
}
if (inVtx0){
internalId_.betaClassic_ += tkpt;
internalId_.beta_ += tkpt;
}
if (inVtxOther){
} else if (inVtxOther){
internalId_.betaStarClassic_ += tkpt;
}
if (fabs(dZ0) < 0.2){
internalId_.beta_ += tkpt;
} else if (fabs(dZ_tmp) < 0.2){
internalId_.betaStar_ += tkpt;
}
}
Expand Down Expand Up @@ -617,7 +637,7 @@ PileupJetIdentifier PileupJetIdAlgo::computeIdVariables(const reco::Jet * jet, f
internalId_.sumNePt_ = sumPtNe;

internalId_.jetR_ = lLead->pt()/sumPt;
internalId_.jetRchg_ = lLeadEm->pt()/sumPt;
internalId_.jetRchg_ = lLeadCh->pt()/sumPt;
internalId_.dRMatch_ = dRmin;

if( sumTkPt != 0. ) {
Expand Down