1 import FWCore.ParameterSet.Config
as cms
2 from Configuration.Eras.Modifier_run2_miniAOD_80XLegacy_cff
import run2_miniAOD_80XLegacy
3 from Configuration.Eras.Modifier_run2_nanoAOD_94X2016_cff
import run2_nanoAOD_94X2016
4 from Configuration.Eras.Modifier_run2_nanoAOD_94XMiniAODv1_cff
import run2_nanoAOD_94XMiniAODv1
5 from Configuration.Eras.Modifier_run2_nanoAOD_94XMiniAODv2_cff
import run2_nanoAOD_94XMiniAODv2
14 chsForSATkJets = cms.EDFilter(
"CandPtrSelector", src = cms.InputTag(
"packedPFCandidates"), cut = cms.string(
'charge()!=0 && pvAssociationQuality()>=5 && vertexRef().key()==0'))
15 softActivityJets = ak4PFJets.clone(src =
'chsForSATkJets', doAreaFastjet =
False, jetPtMin=1)
16 softActivityJets10 = cms.EDFilter(
"CandPtrSelector", src = cms.InputTag(
"softActivityJets"), cut = cms.string(
'pt>10'))
17 softActivityJets5 = cms.EDFilter(
"CandPtrSelector", src = cms.InputTag(
"softActivityJets"), cut = cms.string(
'pt>5'))
18 softActivityJets2 = cms.EDFilter(
"CandPtrSelector", src = cms.InputTag(
"softActivityJets"), cut = cms.string(
'pt>2'))
23 jetCorrFactorsNano = patJetCorrFactors.clone(src=
'slimmedJets',
24 levels = cms.vstring(
'L1FastJet',
28 primaryVertices = cms.InputTag(
"offlineSlimmedPrimaryVertices"),
30 jetCorrFactorsAK8 = patJetCorrFactors.clone(src=
'slimmedJetsAK8',
31 levels = cms.vstring(
'L1FastJet',
35 payload = cms.string(
'AK8PFPuppi'),
36 primaryVertices = cms.InputTag(
"offlineSlimmedPrimaryVertices"),
38 run2_miniAOD_80XLegacy.toModify(jetCorrFactorsAK8, payload = cms.string(
'AK8PFchs'))
42 updatedJets = updatedPatJets.clone(
44 jetSource=
'slimmedJets',
45 jetCorrFactorsSource=cms.VInputTag(cms.InputTag(
"jetCorrFactorsNano") ),
48 updatedJetsAK8 = updatedPatJets.clone(
50 jetSource=
'slimmedJetsAK8',
51 jetCorrFactorsSource=cms.VInputTag(cms.InputTag(
"jetCorrFactorsAK8") ),
55 looseJetId = cms.EDProducer(
"PatJetIDValueMapProducer",
56 filterParams=cms.PSet(
57 version = cms.string(
'WINTER16'),
58 quality = cms.string(
'LOOSE'),
60 src = cms.InputTag(
"updatedJets")
62 tightJetId = cms.EDProducer(
"PatJetIDValueMapProducer",
63 filterParams=cms.PSet(
64 version = cms.string(
'SUMMER18'),
65 quality = cms.string(
'TIGHT'),
67 src = cms.InputTag(
"updatedJets")
69 tightJetIdLepVeto = cms.EDProducer(
"PatJetIDValueMapProducer",
70 filterParams=cms.PSet(
71 version = cms.string(
'SUMMER18'),
72 quality = cms.string(
'TIGHTLEPVETO'),
74 src = cms.InputTag(
"updatedJets")
76 for modifier
in run2_miniAOD_80XLegacy, run2_nanoAOD_94X2016:
77 modifier.toModify( tightJetId.filterParams, version =
"WINTER16" )
78 modifier.toModify( tightJetIdLepVeto.filterParams, version =
"WINTER16" )
79 for modifier
in run2_nanoAOD_94XMiniAODv1, run2_nanoAOD_94XMiniAODv2:
80 modifier.toModify( tightJetId.filterParams, version =
"WINTER17" )
81 modifier.toModify( tightJetIdLepVeto.filterParams, version =
"WINTER17" )
84 looseJetIdAK8 = cms.EDProducer(
"PatJetIDValueMapProducer",
85 filterParams=cms.PSet(
86 version = cms.string(
'WINTER16'),
87 quality = cms.string(
'LOOSE'),
89 src = cms.InputTag(
"updatedJetsAK8")
91 tightJetIdAK8 = cms.EDProducer(
"PatJetIDValueMapProducer",
92 filterParams=cms.PSet(
93 version = cms.string(
'SUMMER18PUPPI'),
94 quality = cms.string(
'TIGHT'),
96 src = cms.InputTag(
"updatedJetsAK8")
98 tightJetIdLepVetoAK8 = cms.EDProducer(
"PatJetIDValueMapProducer",
99 filterParams=cms.PSet(
100 version = cms.string(
'SUMMER18PUPPI'),
101 quality = cms.string(
'TIGHTLEPVETO'),
103 src = cms.InputTag(
"updatedJetsAK8")
105 for modifier
in run2_miniAOD_80XLegacy, run2_nanoAOD_94X2016:
106 modifier.toModify( tightJetIdAK8.filterParams, version =
"WINTER16" )
107 modifier.toModify( tightJetIdLepVetoAK8.filterParams, version =
"WINTER16" )
108 for modifier
in run2_nanoAOD_94XMiniAODv1, run2_nanoAOD_94XMiniAODv2:
109 modifier.toModify( tightJetIdAK8.filterParams, version =
"WINTER17PUPPI" )
110 modifier.toModify( tightJetIdLepVetoAK8.filterParams, version =
"WINTER17PUPPI" )
113 bJetVars = cms.EDProducer(
"JetRegressionVarProducer",
114 pvsrc = cms.InputTag(
"offlineSlimmedPrimaryVertices"),
115 src = cms.InputTag(
"updatedJets"),
116 svsrc = cms.InputTag(
"slimmedSecondaryVertices"),
117 gpsrc = cms.InputTag(
"prunedGenParticles"),
122 jercVars = cms.EDProducer(
"BetaStarPackedCandidateVarProducer",
123 srcJet = cms.InputTag(
"updatedJets"),
124 srcPF = cms.InputTag(
"packedPFCandidates"),
125 maxDR = cms.double(0.4)
129 updatedJetsWithUserData = cms.EDProducer(
"PATJetUserDataEmbedder",
130 src = cms.InputTag(
"updatedJets"),
131 userFloats = cms.PSet(
132 leadTrackPt = cms.InputTag(
"bJetVars:leadTrackPt"),
133 leptonPtRel = cms.InputTag(
"bJetVars:leptonPtRel"),
134 leptonPtRatio = cms.InputTag(
"bJetVars:leptonPtRatio"),
135 leptonPtRelInv = cms.InputTag(
"bJetVars:leptonPtRelInv"),
136 leptonPtRelv0 = cms.InputTag(
"bJetVars:leptonPtRelv0"),
137 leptonPtRatiov0 = cms.InputTag(
"bJetVars:leptonPtRatiov0"),
138 leptonPtRelInvv0 = cms.InputTag(
"bJetVars:leptonPtRelInvv0"),
139 leptonDeltaR = cms.InputTag(
"bJetVars:leptonDeltaR"),
140 leptonPt = cms.InputTag(
"bJetVars:leptonPt"),
141 vtxPt = cms.InputTag(
"bJetVars:vtxPt"),
142 vtxMass = cms.InputTag(
"bJetVars:vtxMass"),
143 vtx3dL = cms.InputTag(
"bJetVars:vtx3dL"),
144 vtx3deL = cms.InputTag(
"bJetVars:vtx3deL"),
145 ptD = cms.InputTag(
"bJetVars:ptD"),
146 genPtwNu = cms.InputTag(
"bJetVars:genPtwNu"),
147 qgl = cms.InputTag(
'qgtagger:qgLikelihood'),
148 jercCHPUF = cms.InputTag(
"jercVars:chargedHadronPUEnergyFraction"),
149 jercCHF = cms.InputTag(
"jercVars:chargedHadronCHSEnergyFraction"),
152 tightId = cms.InputTag(
"tightJetId"),
153 tightIdLepVeto = cms.InputTag(
"tightJetIdLepVeto"),
154 vtxNtrk = cms.InputTag(
"bJetVars:vtxNtrk"),
155 leptonPdgId = cms.InputTag(
"bJetVars:leptonPdgId"),
159 for modifier
in run2_miniAOD_80XLegacy, run2_nanoAOD_94X2016:
160 modifier.toModify( updatedJetsWithUserData.userInts,
161 looseId = cms.InputTag(
"looseJetId"),
164 updatedJetsAK8WithUserData = cms.EDProducer(
"PATJetUserDataEmbedder",
165 src = cms.InputTag(
"updatedJetsAK8"),
166 userFloats = cms.PSet(),
168 tightId = cms.InputTag(
"tightJetIdAK8"),
169 tightIdLepVeto = cms.InputTag(
"tightJetIdLepVetoAK8"),
172 for modifier
in run2_miniAOD_80XLegacy, run2_nanoAOD_94X2016:
173 modifier.toModify( updatedJetsAK8WithUserData.userInts,
174 looseId = cms.InputTag(
"looseJetIdAK8"),
178 finalJets = cms.EDFilter(
"PATJetRefSelector",
179 src = cms.InputTag(
"updatedJetsWithUserData"),
180 cut = cms.string(
"pt > 15")
183 finalJetsAK8 = cms.EDFilter(
"PATJetRefSelector",
184 src = cms.InputTag(
"updatedJetsAK8WithUserData"),
185 cut = cms.string(
"pt > 170")
196 jetTable = cms.EDProducer(
"SimpleCandidateFlatTableProducer",
197 src = cms.InputTag(
"linkedObjects",
"jets"),
198 cut = cms.string(
""),
199 name = cms.string(
"Jet"),
200 doc = cms.string(
"slimmedJets, i.e. ak4 PFJets CHS with JECs applied, after basic selection (" + finalJets.cut.value()+
")"),
201 singleton = cms.bool(
False),
202 extension = cms.bool(
False),
203 externalVariables = cms.PSet(
205 bRegCorr =
ExtVar(cms.InputTag(
"bjetNN:corr"),float, doc=
"pt correction for b-jet energy regression",precision=12),
206 bRegRes =
ExtVar(cms.InputTag(
"bjetNN:res"),float, doc=
"res on pt corrected with b-jet regression",precision=8),
208 variables = cms.PSet(P4Vars,
209 area =
Var(
"jetArea()", float, doc=
"jet catchment area, for JECs",precision=10),
210 nMuons =
Var(
"?hasOverlaps('muons')?overlaps('muons').size():0", int, doc=
"number of muons in the jet"),
211 muonIdx1 =
Var(
"?overlaps('muons').size()>0?overlaps('muons')[0].key():-1", int, doc=
"index of first matching muon"),
212 muonIdx2 =
Var(
"?overlaps('muons').size()>1?overlaps('muons')[1].key():-1", int, doc=
"index of second matching muon"),
213 electronIdx1 =
Var(
"?overlaps('electrons').size()>0?overlaps('electrons')[0].key():-1", int, doc=
"index of first matching electron"),
214 electronIdx2 =
Var(
"?overlaps('electrons').size()>1?overlaps('electrons')[1].key():-1", int, doc=
"index of second matching electron"),
215 nElectrons =
Var(
"?hasOverlaps('electrons')?overlaps('electrons').size():0", int, doc=
"number of electrons in the jet"),
216 btagCMVA =
Var(
"bDiscriminator('pfCombinedMVAV2BJetTags')",float,doc=
"CMVA V2 btag discriminator",precision=10),
217 btagDeepB =
Var(
"bDiscriminator('pfDeepCSVJetTags:probb')+bDiscriminator('pfDeepCSVJetTags:probbb')",float,doc=
"DeepCSV b+bb tag discriminator",precision=10),
218 btagDeepFlavB =
Var(
"bDiscriminator('pfDeepFlavourJetTags:probb')+bDiscriminator('pfDeepFlavourJetTags:probbb')+bDiscriminator('pfDeepFlavourJetTags:problepb')",float,doc=
"DeepFlavour b+bb+lepb tag discriminator",precision=10),
219 btagCSVV2 =
Var(
"bDiscriminator('pfCombinedInclusiveSecondaryVertexV2BJetTags')",float,doc=
" pfCombinedInclusiveSecondaryVertexV2 b-tag discriminator (aka CSVV2)",precision=10),
220 btagDeepC =
Var(
"bDiscriminator('pfDeepCSVJetTags:probc')",float,doc=
"DeepCSV charm btag discriminator",precision=10),
221 btagDeepFlavC =
Var(
"bDiscriminator('pfDeepFlavourJetTags:probc')",float,doc=
"DeepFlavour charm tag discriminator",precision=10),
223 puId =
Var(
"userInt('pileupJetId:fullId')",int,doc=
"Pilup ID flags"),
224 jetId =
Var(
"userInt('tightId')*2+4*userInt('tightIdLepVeto')",int,doc=
"Jet ID flags bit1 is loose (always false in 2017 since it does not exist), bit2 is tight, bit3 is tightLepVeto"),
225 qgl =
Var(
"userFloat('qgl')",float,doc=
"Quark vs Gluon likelihood discriminator",precision=10),
226 nConstituents =
Var(
"numberOfDaughters()",int,doc=
"Number of particles in the jet"),
227 rawFactor =
Var(
"1.-jecFactor('Uncorrected')",float,doc=
"1 - Factor to get back to raw pT",precision=6),
228 chHEF =
Var(
"chargedHadronEnergyFraction()", float, doc=
"charged Hadron Energy Fraction", precision= 6),
229 neHEF =
Var(
"neutralHadronEnergyFraction()", float, doc=
"neutral Hadron Energy Fraction", precision= 6),
230 chEmEF =
Var(
"chargedEmEnergyFraction()", float, doc=
"charged Electromagnetic Energy Fraction", precision= 6),
231 neEmEF =
Var(
"neutralEmEnergyFraction()", float, doc=
"neutral Electromagnetic Energy Fraction", precision= 6),
232 muEF =
Var(
"muonEnergyFraction()", float, doc=
"muon Energy Fraction", precision= 6),
233 jercCHPUF =
Var(
"userFloat('jercCHPUF')", float, doc=
"Pileup Charged Hadron Energy Fraction with the JERC group definition", precision= 6),
234 jercCHF =
Var(
"userFloat('jercCHF')", float, doc=
"Charged Hadron Energy Fraction with the JERC group definition", precision= 6),
239 jetTable.variables.pt.precision=10
242 for modifier
in run2_miniAOD_80XLegacy, run2_nanoAOD_94X2016:
243 modifier.toModify( jetTable.variables, jetId =
Var(
"userInt('tightIdLepVeto')*4+userInt('tightId')*2+userInt('looseId')",int,doc=
"Jet ID flags bit1 is loose, bit2 is tight, bit3 is tightLepVeto"))
245 bjetMVA= cms.EDProducer(
"BJetEnergyRegressionMVA",
246 backend = cms.string(
"TMVA"),
247 src = cms.InputTag(
"linkedObjects",
"jets"),
248 pvsrc = cms.InputTag(
"offlineSlimmedPrimaryVertices"),
249 svsrc = cms.InputTag(
"slimmedSecondaryVertices"),
250 rhosrc = cms.InputTag(
"fixedGridRhoFastjetAll"),
251 weightFile = cms.FileInPath(
"PhysicsTools/NanoAOD/data/bjet-regression.xml"),
252 name = cms.string(
"JetReg"),
253 isClassifier = cms.bool(
False),
254 variablesOrder = cms.vstring([
"Jet_pt",
"nPVs",
"Jet_eta",
"Jet_mt",
"Jet_leadTrackPt",
"Jet_leptonPtRel",
"Jet_leptonPt",
"Jet_leptonDeltaR",
"Jet_neHEF",
"Jet_neEmEF",
"Jet_vtxPt",
"Jet_vtxMass",
"Jet_vtx3dL",
"Jet_vtxNtrk",
"Jet_vtx3deL"]),
255 variables = cms.PSet(
256 Jet_pt = cms.string(
"pt"),
257 Jet_eta = cms.string(
"eta"),
258 Jet_mt = cms.string(
"mt"),
259 Jet_leadTrackPt = cms.string(
"userFloat('leadTrackPt')"),
260 Jet_vtxNtrk = cms.string(
"userInt('vtxNtrk')"),
261 Jet_vtxMass = cms.string(
"userFloat('vtxMass')"),
262 Jet_vtx3dL = cms.string(
"userFloat('vtx3dL')"),
263 Jet_vtx3deL = cms.string(
"userFloat('vtx3deL')"),
264 Jet_vtxPt = cms.string(
"userFloat('vtxPt')"),
265 Jet_leptonPtRel = cms.string(
"userFloat('leptonPtRelv0')"),
266 Jet_leptonPt = cms.string(
"?overlaps('muons').size()>0?overlaps('muons')[0].pt():(?overlaps('electrons').size()>0?overlaps('electrons')[0].pt():0)"),
267 Jet_neHEF = cms.string(
"neutralHadronEnergy()/energy()"),
268 Jet_neEmEF = cms.string(
"neutralEmEnergy()/energy()"),
269 Jet_leptonDeltaR = cms.string(
'''?overlaps('muons').size()>0?deltaR(eta,phi,overlaps('muons')[0].eta,overlaps('muons')[0].phi): 270 (?overlaps('electrons').size()>0?deltaR(eta,phi,overlaps('electrons')[0].eta,overlaps('electrons')[0].phi): 276 bjetNN= cms.EDProducer(
"BJetEnergyRegressionMVA",
277 backend = cms.string(
"TF"),
278 src = cms.InputTag(
"linkedObjects",
"jets"),
279 pvsrc = cms.InputTag(
"offlineSlimmedPrimaryVertices"),
280 svsrc = cms.InputTag(
"slimmedSecondaryVertices"),
281 rhosrc = cms.InputTag(
"fixedGridRhoFastjetAll"),
283 weightFile = cms.FileInPath(
"PhysicsTools/NanoAOD/data/breg_training_2017.pb"),
284 name = cms.string(
"JetRegNN"),
285 isClassifier = cms.bool(
False),
286 variablesOrder = cms.vstring([
"Jet_pt",
"Jet_eta",
"rho",
"Jet_mt",
"Jet_leadTrackPt",
"Jet_leptonPtRel",
"Jet_leptonDeltaR",
"Jet_neHEF",
"Jet_neEmEF",
"Jet_vtxPt",
"Jet_vtxMass",
"Jet_vtx3dL",
"Jet_vtxNtrk",
"Jet_vtx3deL",
"Jet_numDaughters_pt03",
"Jet_energyRing_dR0_em_Jet_rawEnergy",
"Jet_energyRing_dR1_em_Jet_rawEnergy",
"Jet_energyRing_dR2_em_Jet_rawEnergy",
"Jet_energyRing_dR3_em_Jet_rawEnergy",
"Jet_energyRing_dR4_em_Jet_rawEnergy",
"Jet_energyRing_dR0_neut_Jet_rawEnergy",
"Jet_energyRing_dR1_neut_Jet_rawEnergy",
"Jet_energyRing_dR2_neut_Jet_rawEnergy",
"Jet_energyRing_dR3_neut_Jet_rawEnergy",
"Jet_energyRing_dR4_neut_Jet_rawEnergy",
"Jet_energyRing_dR0_ch_Jet_rawEnergy",
"Jet_energyRing_dR1_ch_Jet_rawEnergy",
"Jet_energyRing_dR2_ch_Jet_rawEnergy",
"Jet_energyRing_dR3_ch_Jet_rawEnergy",
"Jet_energyRing_dR4_ch_Jet_rawEnergy",
"Jet_energyRing_dR0_mu_Jet_rawEnergy",
"Jet_energyRing_dR1_mu_Jet_rawEnergy",
"Jet_energyRing_dR2_mu_Jet_rawEnergy",
"Jet_energyRing_dR3_mu_Jet_rawEnergy",
"Jet_energyRing_dR4_mu_Jet_rawEnergy",
"Jet_chHEF",
"Jet_chEmEF",
"Jet_leptonPtRelInv",
"isEle",
"isMu",
"isOther",
"Jet_mass",
"Jet_ptd"]),
287 variables = cms.PSet(
288 Jet_pt = cms.string(
"pt*jecFactor('Uncorrected')"),
289 Jet_mt = cms.string(
"mt*jecFactor('Uncorrected')"),
290 Jet_eta = cms.string(
"eta"),
291 Jet_mass = cms.string(
"mass*jecFactor('Uncorrected')"),
292 Jet_ptd = cms.string(
"userFloat('ptD')"),
293 Jet_leadTrackPt = cms.string(
"userFloat('leadTrackPt')"),
294 Jet_vtxNtrk = cms.string(
"userInt('vtxNtrk')"),
295 Jet_vtxMass = cms.string(
"userFloat('vtxMass')"),
296 Jet_vtx3dL = cms.string(
"userFloat('vtx3dL')"),
297 Jet_vtx3deL = cms.string(
"userFloat('vtx3deL')"),
298 Jet_vtxPt = cms.string(
"userFloat('vtxPt')"),
300 Jet_leptonPtRel = cms.string(
"userFloat('leptonPtRelv0')"),
301 Jet_leptonPtRelInv = cms.string(
"userFloat('leptonPtRelInvv0')*jecFactor('Uncorrected')"),
302 Jet_leptonDeltaR = cms.string(
"userFloat('leptonDeltaR')"),
304 Jet_neHEF = cms.string(
"neutralHadronEnergyFraction()"),
305 Jet_neEmEF = cms.string(
"neutralEmEnergyFraction()"),
306 Jet_chHEF = cms.string(
"chargedHadronEnergyFraction()"),
307 Jet_chEmEF = cms.string(
"chargedEmEnergyFraction()"),
308 isMu = cms.string(
"?abs(userInt('leptonPdgId'))==13?1:0"),
309 isEle = cms.string(
"?abs(userInt('leptonPdgId'))==11?1:0"),
310 isOther = cms.string(
"?userInt('leptonPdgId')==0?1:0"),
312 inputTensorName = cms.string(
"ffwd_inp"),
313 outputTensorName = cms.string(
"ffwd_out/BiasAdd"),
314 outputNames = cms.vstring([
"corr",
"res"]),
315 outputFormulas = cms.vstring([
"at(0)*0.28492164611816406+1.0596693754196167",
"0.5*(at(2)-at(1))*0.28492164611816406"]),
316 nThreads = cms.uint32(1),
317 singleThreadPool = cms.string(
"no_threads"),
323 saJetTable = cms.EDProducer(
"SimpleCandidateFlatTableProducer",
324 src = cms.InputTag(
"softActivityJets"),
325 cut = cms.string(
""),
326 maxLen = cms.uint32(6),
327 name = cms.string(
"SoftActivityJet"),
328 doc = cms.string(
"jets clustered from charged candidates compatible with primary vertex (" + chsForSATkJets.cut.value()+
")"),
329 singleton = cms.bool(
False),
330 extension = cms.bool(
False),
331 variables = cms.PSet(P3Vars,
335 saJetTable.variables.pt.precision=10
336 saJetTable.variables.eta.precision=8
337 saJetTable.variables.phi.precision=8
339 saTable = cms.EDProducer(
"GlobalVariablesTableProducer",
340 variables = cms.PSet(
341 SoftActivityJetHT =
ExtVar( cms.InputTag(
"softActivityJets"),
"candidatescalarsum", doc =
"scalar sum of soft activity jet pt, pt>1" ),
342 SoftActivityJetHT10 =
ExtVar( cms.InputTag(
"softActivityJets10"),
"candidatescalarsum", doc =
"scalar sum of soft activity jet pt , pt >10" ),
343 SoftActivityJetHT5 =
ExtVar( cms.InputTag(
"softActivityJets5"),
"candidatescalarsum", doc =
"scalar sum of soft activity jet pt, pt>5" ),
344 SoftActivityJetHT2 =
ExtVar( cms.InputTag(
"softActivityJets2"),
"candidatescalarsum", doc =
"scalar sum of soft activity jet pt, pt >2" ),
345 SoftActivityJetNjets10 =
ExtVar( cms.InputTag(
"softActivityJets10"),
"candidatesize", doc =
"number of soft activity jet pt, pt >2" ),
346 SoftActivityJetNjets5 =
ExtVar( cms.InputTag(
"softActivityJets5"),
"candidatesize", doc =
"number of soft activity jet pt, pt >5" ),
347 SoftActivityJetNjets2 =
ExtVar( cms.InputTag(
"softActivityJets2"),
"candidatesize", doc =
"number of soft activity jet pt, pt >10" ),
355 fatJetTable = cms.EDProducer(
"SimpleCandidateFlatTableProducer",
356 src = cms.InputTag(
"finalJetsAK8"),
357 cut = cms.string(
" pt > 170"),
358 name = cms.string(
"FatJet"),
359 doc = cms.string(
"slimmedJetsAK8, i.e. ak8 fat jets for boosted analysis"),
360 singleton = cms.bool(
False),
361 extension = cms.bool(
False),
362 variables = cms.PSet(P4Vars,
363 jetId =
Var(
"userInt('tightId')*2+4*userInt('tightIdLepVeto')",int,doc=
"Jet ID flags bit1 is loose (always false in 2017 since it does not exist), bit2 is tight, bit3 is tightLepVeto"),
364 area =
Var(
"jetArea()", float, doc=
"jet catchment area, for JECs",precision=10),
365 rawFactor =
Var(
"1.-jecFactor('Uncorrected')",float,doc=
"1 - Factor to get back to raw pT",precision=6),
366 tau1 =
Var(
"userFloat('NjettinessAK8Puppi:tau1')",float, doc=
"Nsubjettiness (1 axis)",precision=10),
367 tau2 =
Var(
"userFloat('NjettinessAK8Puppi:tau2')",float, doc=
"Nsubjettiness (2 axis)",precision=10),
368 tau3 =
Var(
"userFloat('NjettinessAK8Puppi:tau3')",float, doc=
"Nsubjettiness (3 axis)",precision=10),
369 tau4 =
Var(
"userFloat('NjettinessAK8Puppi:tau4')",float, doc=
"Nsubjettiness (4 axis)",precision=10),
370 n2b1 =
Var(
"userFloat('ak8PFJetsPuppiSoftDropValueMap:nb1AK8PuppiSoftDropN2')", float, doc=
"N2 with beta=1", precision=10),
371 n3b1 =
Var(
"userFloat('ak8PFJetsPuppiSoftDropValueMap:nb1AK8PuppiSoftDropN3')", float, doc=
"N3 with beta=1", precision=10),
372 msoftdrop =
Var(
"groomedMass('SoftDropPuppi')",float, doc=
"Corrected soft drop mass with PUPPI",precision=10),
373 btagCMVA =
Var(
"bDiscriminator('pfCombinedMVAV2BJetTags')",float,doc=
"CMVA V2 btag discriminator",precision=10),
374 btagDeepB =
Var(
"bDiscriminator('pfDeepCSVJetTags:probb')+bDiscriminator('pfDeepCSVJetTags:probbb')",float,doc=
"DeepCSV b+bb tag discriminator",precision=10),
375 btagCSVV2 =
Var(
"bDiscriminator('pfCombinedInclusiveSecondaryVertexV2BJetTags')",float,doc=
" pfCombinedInclusiveSecondaryVertexV2 b-tag discriminator (aka CSVV2)",precision=10),
376 btagHbb =
Var(
"bDiscriminator('pfBoostedDoubleSecondaryVertexAK8BJetTags')",float,doc=
"Higgs to BB tagger discriminator",precision=10),
377 btagDDBvL_noMD =
Var(
"bDiscriminator('pfDeepDoubleBvLJetTags:probHbb')",float,doc=
"DeepDoubleX discriminator (no mass-decorrelation) for H(Z)->bb vs QCD",precision=10),
378 btagDDCvL_noMD =
Var(
"bDiscriminator('pfDeepDoubleCvLJetTags:probHcc')",float,doc=
"DeepDoubleX discriminator (no mass-decorrelation) for H(Z)->cc vs QCD",precision=10),
379 btagDDCvB_noMD =
Var(
"bDiscriminator('pfDeepDoubleCvBJetTags:probHcc')",float,doc=
"DeepDoubleX discriminator (no mass-decorrelation) for H(Z)->cc vs H(Z)->bb",precision=10),
380 btagDDBvL =
Var(
"bDiscriminator('pfMassIndependentDeepDoubleBvLJetTags:probHbb')",float,doc=
"DeepDoubleX (mass-decorrelated) discriminator for H(Z)->bb vs QCD",precision=10),
381 btagDDCvL =
Var(
"bDiscriminator('pfMassIndependentDeepDoubleCvLJetTags:probHcc')",float,doc=
"DeepDoubleX (mass-decorrelated) discriminator for H(Z)->cc vs QCD",precision=10),
382 btagDDCvB =
Var(
"bDiscriminator('pfMassIndependentDeepDoubleCvBJetTags:probHcc')",float,doc=
"DeepDoubleX (mass-decorrelated) discriminator for H(Z)->cc vs H(Z)->bb",precision=10),
383 deepTag_TvsQCD =
Var(
"bDiscriminator('pfDeepBoostedDiscriminatorsJetTags:TvsQCD')",float,doc=
"DeepBoostedJet tagger top vs QCD discriminator",precision=10),
384 deepTag_WvsQCD =
Var(
"bDiscriminator('pfDeepBoostedDiscriminatorsJetTags:WvsQCD')",float,doc=
"DeepBoostedJet tagger W vs QCD discriminator",precision=10),
385 deepTag_ZvsQCD =
Var(
"bDiscriminator('pfDeepBoostedDiscriminatorsJetTags:ZvsQCD')",float,doc=
"DeepBoostedJet tagger Z vs QCD discriminator",precision=10),
386 deepTag_H =
Var(
"bDiscriminator('pfDeepBoostedJetTags:probHbb')+bDiscriminator('pfDeepBoostedJetTags:probHcc')+bDiscriminator('pfDeepBoostedJetTags:probHqqqq')",float,doc=
"DeepBoostedJet tagger H(bb,cc,4q) sum",precision=10),
387 deepTag_QCD =
Var(
"bDiscriminator('pfDeepBoostedJetTags:probQCDbb')+bDiscriminator('pfDeepBoostedJetTags:probQCDcc')+bDiscriminator('pfDeepBoostedJetTags:probQCDb')+bDiscriminator('pfDeepBoostedJetTags:probQCDc')+bDiscriminator('pfDeepBoostedJetTags:probQCDothers')",float,doc=
"DeepBoostedJet tagger QCD(bb,cc,b,c,others) sum",precision=10),
388 deepTag_QCDothers =
Var(
"bDiscriminator('pfDeepBoostedJetTags:probQCDothers')",float,doc=
"DeepBoostedJet tagger QCDothers value",precision=10),
389 deepTagMD_TvsQCD =
Var(
"bDiscriminator('pfMassDecorrelatedDeepBoostedDiscriminatorsJetTags:TvsQCD')",float,doc=
"Mass-decorrelated DeepBoostedJet tagger top vs QCD discriminator",precision=10),
390 deepTagMD_WvsQCD =
Var(
"bDiscriminator('pfMassDecorrelatedDeepBoostedDiscriminatorsJetTags:WvsQCD')",float,doc=
"Mass-decorrelated DeepBoostedJet tagger W vs QCD discriminator",precision=10),
391 deepTagMD_ZvsQCD =
Var(
"bDiscriminator('pfMassDecorrelatedDeepBoostedDiscriminatorsJetTags:ZvsQCD')",float,doc=
"Mass-decorrelated DeepBoostedJet tagger Z vs QCD discriminator",precision=10),
392 deepTagMD_ZHbbvsQCD =
Var(
"bDiscriminator('pfMassDecorrelatedDeepBoostedDiscriminatorsJetTags:ZHbbvsQCD')",float,doc=
"Mass-decorrelated DeepBoostedJet tagger Z/H->bb vs QCD discriminator",precision=10),
393 deepTagMD_ZbbvsQCD =
Var(
"bDiscriminator('pfMassDecorrelatedDeepBoostedDiscriminatorsJetTags:ZbbvsQCD')",float,doc=
"Mass-decorrelated DeepBoostedJet tagger Z->bb vs QCD discriminator",precision=10),
394 deepTagMD_HbbvsQCD =
Var(
"bDiscriminator('pfMassDecorrelatedDeepBoostedDiscriminatorsJetTags:HbbvsQCD')",float,doc=
"Mass-decorrelated DeepBoostedJet tagger H->bb vs QCD discriminator",precision=10),
395 deepTagMD_ZHccvsQCD =
Var(
"bDiscriminator('pfMassDecorrelatedDeepBoostedDiscriminatorsJetTags:ZHccvsQCD')",float,doc=
"Mass-decorrelated DeepBoostedJet tagger Z/H->cc vs QCD discriminator",precision=10),
396 deepTagMD_H4qvsQCD =
Var(
"bDiscriminator('pfMassDecorrelatedDeepBoostedDiscriminatorsJetTags:H4qvsQCD')",float,doc=
"Mass-decorrelated DeepBoostedJet tagger H->4q vs QCD discriminator",precision=10),
397 deepTagMD_bbvsLight =
Var(
"bDiscriminator('pfMassDecorrelatedDeepBoostedDiscriminatorsJetTags:bbvsLight')",float,doc=
"Mass-decorrelated DeepBoostedJet tagger Z/H/gluon->bb vs light flavour discriminator",precision=10),
398 deepTagMD_ccvsLight =
Var(
"bDiscriminator('pfMassDecorrelatedDeepBoostedDiscriminatorsJetTags:ccvsLight')",float,doc=
"Mass-decorrelated DeepBoostedJet tagger Z/H/gluon->cc vs light flavour discriminator",precision=10),
399 subJetIdx1 =
Var(
"?nSubjetCollections()>0 && subjets('SoftDropPuppi').size()>0?subjets('SoftDropPuppi')[0].key():-1", int,
400 doc=
"index of first subjet"),
401 subJetIdx2 =
Var(
"?nSubjetCollections()>0 && subjets('SoftDropPuppi').size()>1?subjets('SoftDropPuppi')[1].key():-1", int,
402 doc=
"index of second subjet"),
411 run2_miniAOD_80XLegacy.toModify( bjetNN,outputFormulas = cms.vstring([
"at(0)*0.31628304719924927+1.0454729795455933",
"0.5*(at(2)-at(1))*0.31628304719924927"]))
412 run2_miniAOD_80XLegacy.toModify( fatJetTable.variables, msoftdrop_chs =
Var(
"userFloat('ak8PFJetsCHSSoftDropMass')",float, doc=
"Legacy uncorrected soft drop mass with CHS",precision=10))
413 run2_miniAOD_80XLegacy.toModify( fatJetTable.variables.tau1, expr = cms.string(
"userFloat(\'ak8PFJetsPuppiValueMap:NjettinessAK8PuppiTau1\')"),)
414 run2_miniAOD_80XLegacy.toModify( fatJetTable.variables.tau2, expr = cms.string(
"userFloat(\'ak8PFJetsPuppiValueMap:NjettinessAK8PuppiTau2\')"),)
415 run2_miniAOD_80XLegacy.toModify( fatJetTable.variables.tau3, expr = cms.string(
"userFloat(\'ak8PFJetsPuppiValueMap:NjettinessAK8PuppiTau3\')"),)
416 run2_miniAOD_80XLegacy.toModify( fatJetTable.variables, tau4 =
None)
417 run2_miniAOD_80XLegacy.toModify( fatJetTable.variables, n2b1 =
None)
418 run2_miniAOD_80XLegacy.toModify( fatJetTable.variables, n3b1 =
None)
419 for modifier
in run2_miniAOD_80XLegacy, run2_nanoAOD_94X2016:
420 modifier.toModify( fatJetTable.variables, jetId =
Var(
"userInt('tightId')*2+userInt('looseId')",int,doc=
"Jet ID flags bit1 is loose, bit2 is tight"))
425 subJetTable = cms.EDProducer(
"SimpleCandidateFlatTableProducer",
426 src = cms.InputTag(
"slimmedJetsAK8PFPuppiSoftDropPacked",
"SubJets"),
427 cut = cms.string(
""),
428 name = cms.string(
"SubJet"),
429 doc = cms.string(
"slimmedJetsAK8, i.e. ak8 fat jets for boosted analysis"),
430 singleton = cms.bool(
False),
431 extension = cms.bool(
False),
432 variables = cms.PSet(P4Vars,
433 btagCMVA =
Var(
"bDiscriminator('pfCombinedMVAV2BJetTags')",float,doc=
"CMVA V2 btag discriminator",precision=10),
434 btagDeepB =
Var(
"bDiscriminator('pfDeepCSVJetTags:probb')+bDiscriminator('pfDeepCSVJetTags:probbb')",float,doc=
"DeepCSV b+bb tag discriminator",precision=10),
435 btagCSVV2 =
Var(
"bDiscriminator('pfCombinedInclusiveSecondaryVertexV2BJetTags')",float,doc=
" pfCombinedInclusiveSecondaryVertexV2 b-tag discriminator (aka CSVV2)",precision=10),
436 rawFactor =
Var(
"1.-jecFactor('Uncorrected')",float,doc=
"1 - Factor to get back to raw pT",precision=6),
437 tau1 =
Var(
"userFloat('NjettinessAK8Subjets:tau1')",float, doc=
"Nsubjettiness (1 axis)",precision=10),
438 tau2 =
Var(
"userFloat('NjettinessAK8Subjets:tau2')",float, doc=
"Nsubjettiness (2 axis)",precision=10),
439 tau3 =
Var(
"userFloat('NjettinessAK8Subjets:tau3')",float, doc=
"Nsubjettiness (3 axis)",precision=10),
440 tau4 =
Var(
"userFloat('NjettinessAK8Subjets:tau4')",float, doc=
"Nsubjettiness (4 axis)",precision=10),
441 n2b1 =
Var(
"userFloat('nb1AK8PuppiSoftDropSubjets:ecfN2')", float, doc=
"N2 with beta=1", precision=10),
442 n3b1 =
Var(
"userFloat('nb1AK8PuppiSoftDropSubjets:ecfN3')", float, doc=
"N3 with beta=1", precision=10),
447 fatJetTable.variables.pt.precision=10
448 subJetTable.variables.pt.precision=10
450 run2_miniAOD_80XLegacy.toModify( subJetTable.variables, tau1 =
None)
451 run2_miniAOD_80XLegacy.toModify( subJetTable.variables, tau2 =
None)
452 run2_miniAOD_80XLegacy.toModify( subJetTable.variables, tau3 =
None)
453 run2_miniAOD_80XLegacy.toModify( subJetTable.variables, tau4 =
None)
454 run2_miniAOD_80XLegacy.toModify( subJetTable.variables, n2b1 =
None)
455 run2_miniAOD_80XLegacy.toModify( subJetTable.variables, n3b1 =
None)
456 run2_miniAOD_80XLegacy.toModify( subJetTable.variables, btagCMVA =
None, btagDeepB =
None)
459 corrT1METJetTable = cms.EDProducer(
"SimpleCandidateFlatTableProducer",
460 src = cms.InputTag(
"corrT1METJets"),
461 cut = cms.string(
""),
462 name = cms.string(
"CorrT1METJet"),
463 doc = cms.string(
"Additional low-pt jets for Type-1 MET re-correction"),
464 singleton = cms.bool(
False),
465 extension = cms.bool(
False),
466 variables = cms.PSet(
467 rawPt =
Var(
"pt()*jecFactor('Uncorrected')",float,precision=10),
468 eta =
Var(
"eta", float,precision=12),
469 phi =
Var(
"phi", float, precision=12),
470 area =
Var(
"jetArea()", float, doc=
"jet catchment area, for JECs",precision=10),
477 jetMCTable = cms.EDProducer(
"SimpleCandidateFlatTableProducer",
478 src = cms.InputTag(
"linkedObjects",
"jets"),
479 cut = cms.string(
""),
480 name = cms.string(
"Jet"),
481 singleton = cms.bool(
False),
482 extension = cms.bool(
True),
483 variables = cms.PSet(
484 partonFlavour =
Var(
"partonFlavour()", int, doc=
"flavour from parton matching"),
485 hadronFlavour =
Var(
"hadronFlavour()", int, doc=
"flavour from hadron ghost clustering"),
486 genJetIdx =
Var(
"?genJetFwdRef().backRef().isNonnull()?genJetFwdRef().backRef().key():-1", int, doc=
"index of matched gen jet"),
489 genJetTable = cms.EDProducer(
"SimpleCandidateFlatTableProducer",
490 src = cms.InputTag(
"slimmedGenJets"),
491 cut = cms.string(
"pt > 10"),
492 name = cms.string(
"GenJet"),
493 doc = cms.string(
"slimmedGenJets, i.e. ak4 Jets made with visible genparticles"),
494 singleton = cms.bool(
False),
495 extension = cms.bool(
False),
496 variables = cms.PSet(P4Vars,
500 patJetPartons = cms.EDProducer(
'HadronAndPartonSelector',
501 src = cms.InputTag(
"generator"),
502 particles = cms.InputTag(
"prunedGenParticles"),
503 partonMode = cms.string(
"Auto"),
504 fullChainPhysPartons = cms.bool(
True)
506 genJetFlavourAssociation = cms.EDProducer(
"JetFlavourClustering",
507 jets = genJetTable.src,
508 bHadrons = cms.InputTag(
"patJetPartons",
"bHadrons"),
509 cHadrons = cms.InputTag(
"patJetPartons",
"cHadrons"),
510 partons = cms.InputTag(
"patJetPartons",
"physicsPartons"),
511 leptons = cms.InputTag(
"patJetPartons",
"leptons"),
512 jetAlgorithm = cms.string(
"AntiKt"),
513 rParam = cms.double(0.4),
514 ghostRescaling = cms.double(1e-18),
515 hadronFlavourHasPriority = cms.bool(
False)
517 genJetFlavourTable = cms.EDProducer(
"GenJetFlavourTableProducer",
518 name = genJetTable.name,
519 src = genJetTable.src,
520 cut = genJetTable.cut,
521 deltaR = cms.double(0.1),
522 jetFlavourInfos = cms.InputTag(
"slimmedGenJetsFlavourInfos"),
525 genJetAK8Table = cms.EDProducer(
"SimpleCandidateFlatTableProducer",
526 src = cms.InputTag(
"slimmedGenJetsAK8"),
527 cut = cms.string(
"pt > 100."),
528 name = cms.string(
"GenJetAK8"),
529 doc = cms.string(
"slimmedGenJetsAK8, i.e. ak8 Jets made with visible genparticles"),
530 singleton = cms.bool(
False),
531 extension = cms.bool(
False),
532 variables = cms.PSet(P4Vars,
536 genJetAK8FlavourAssociation = cms.EDProducer(
"JetFlavourClustering",
537 jets = genJetAK8Table.src,
538 bHadrons = cms.InputTag(
"patJetPartons",
"bHadrons"),
539 cHadrons = cms.InputTag(
"patJetPartons",
"cHadrons"),
540 partons = cms.InputTag(
"patJetPartons",
"physicsPartons"),
541 leptons = cms.InputTag(
"patJetPartons",
"leptons"),
542 jetAlgorithm = cms.string(
"AntiKt"),
543 rParam = cms.double(0.8),
544 ghostRescaling = cms.double(1e-18),
545 hadronFlavourHasPriority = cms.bool(
False)
547 genJetAK8FlavourTable = cms.EDProducer(
"GenJetFlavourTableProducer",
548 name = genJetAK8Table.name,
549 src = genJetAK8Table.src,
550 cut = genJetAK8Table.cut,
551 deltaR = cms.double(0.1),
552 jetFlavourInfos = cms.InputTag(
"genJetAK8FlavourAssociation"),
554 genSubJetAK8Table = cms.EDProducer(
"SimpleCandidateFlatTableProducer",
555 src = cms.InputTag(
"slimmedGenJetsAK8SoftDropSubJets"),
556 cut = cms.string(
""),
557 name = cms.string(
"SubGenJetAK8"),
558 doc = cms.string(
"slimmedGenJetsAK8SoftDropSubJets, i.e. subjets of ak8 Jets made with visible genparticles"),
559 singleton = cms.bool(
False),
560 extension = cms.bool(
False),
561 variables = cms.PSet(P4Vars,
566 run2_miniAOD_80XLegacy.toModify( genJetFlavourTable, jetFlavourInfos = cms.InputTag(
"genJetFlavourAssociation"),)
569 qgtagger=QGTagger.clone(srcJets=
"updatedJets",srcVertexCollection=
"offlineSlimmedPrimaryVertices")
572 jetSequence = cms.Sequence(jetCorrFactorsNano+updatedJets+tightJetId+tightJetIdLepVeto+bJetVars+jercVars+qgtagger+updatedJetsWithUserData+jetCorrFactorsAK8+updatedJetsAK8+tightJetIdAK8+tightJetIdLepVetoAK8+updatedJetsAK8WithUserData+chsForSATkJets+softActivityJets+softActivityJets2+softActivityJets5+softActivityJets10+finalJets+finalJetsAK8)
574 _jetSequence_2016 = jetSequence.copy()
575 _jetSequence_2016.insert(_jetSequence_2016.index(tightJetId), looseJetId)
576 _jetSequence_2016.insert(_jetSequence_2016.index(tightJetIdAK8), looseJetIdAK8)
577 for modifier
in run2_miniAOD_80XLegacy, run2_nanoAOD_94X2016:
578 modifier.toReplaceWith(jetSequence, _jetSequence_2016)
581 jetTables = cms.Sequence(bjetMVA+bjetNN+jetTable+fatJetTable+subJetTable+saJetTable+saTable)
584 jetMC = cms.Sequence(jetMCTable+genJetTable+patJetPartons+genJetFlavourTable+genJetAK8Table+genJetAK8FlavourAssociation+genJetAK8FlavourTable+genSubJetAK8Table)
585 _jetMC_pre94X = jetMC.copy()
586 _jetMC_pre94X.insert(_jetMC_pre94X.index(genJetFlavourTable),genJetFlavourAssociation)
587 _jetMC_pre94X.remove(genSubJetAK8Table)
588 run2_miniAOD_80XLegacy.toReplaceWith(jetMC, _jetMC_pre94X)
def ExtVar(tag, valtype, compression=None, doc=None, mcOnly=False, precision=-1)
def Var(expr, valtype, compression=None, doc=None, mcOnly=False, precision=-1)