CMS 3D CMS Logo

jetsAK8_cff.py
Go to the documentation of this file.
1 import FWCore.ParameterSet.Config as cms
2 
5 from PhysicsTools.NanoAOD.simpleCandidateFlatTableProducer_cfi import simpleCandidateFlatTableProducer
6 
8 # Note: Safe to always add 'L2L3Residual' as MC contains dummy L2L3Residual corrections (always set to 1)
9 # (cf. https://twiki.cern.ch/twiki/bin/view/CMSPublic/WorkBookJetEnergyCorrections#CMSSW_7_6_4_and_above )
10 jetCorrFactorsAK8 = patJetCorrFactors.clone(src='slimmedJetsAK8',
11  levels = cms.vstring('L1FastJet',
12  'L2Relative',
13  'L3Absolute',
14  'L2L3Residual'),
15  payload = cms.string('AK8PFPuppi'),
16  primaryVertices = cms.InputTag("offlineSlimmedPrimaryVertices"),
17 )
18 
20 updatedJetsAK8 = updatedPatJets.clone(
21  addBTagInfo=False,
22  jetSource='slimmedJetsAK8',
23  jetCorrFactorsSource=cms.VInputTag(cms.InputTag("jetCorrFactorsAK8") ),
24 )
25 
26 #
27 # JetID
28 #
29 looseJetIdAK8 = cms.EDProducer("PatJetIDValueMapProducer",
30  filterParams=cms.PSet(
31  version = cms.string('WINTER16'),
32  quality = cms.string('LOOSE'),
33  ),
34  src = cms.InputTag("updatedJetsAK8")
35 )
36 tightJetIdAK8 = cms.EDProducer("PatJetIDValueMapProducer",
37  filterParams=cms.PSet(
38  version = cms.string('RUN3WINTER22PUPPI'),
39  quality = cms.string('TIGHT'),
40  ),
41  src = cms.InputTag("updatedJetsAK8")
42 )
43 tightJetIdLepVetoAK8 = cms.EDProducer("PatJetIDValueMapProducer",
44  filterParams=cms.PSet(
45  version = cms.string('RUN3WINTER22PUPPI'),
46  quality = cms.string('TIGHTLEPVETO'),
47  ),
48  src = cms.InputTag("updatedJetsAK8")
49 )
50 
51 run2_jme_2016.toModify(
52  tightJetIdAK8.filterParams, version = "RUN2UL16PUPPI"
53 ).toModify(
54  tightJetIdLepVetoAK8.filterParams, version = "RUN2UL16PUPPI"
55 )
56 
57 (run2_jme_2017 | run2_jme_2018 | run3_nanoAOD_122 | run3_nanoAOD_124).toModify(
58  tightJetIdAK8.filterParams, version = "RUN2ULPUPPI"
59 ).toModify(
60  tightJetIdLepVetoAK8.filterParams, version = "RUN2ULPUPPI"
61 )
62 
63 run3_jme_Winter22runsBCDEprompt.toModify(
64  tightJetIdAK8.filterParams, version = "RUN3WINTER22PUPPIrunsBCDEprompt"
65 ).toModify(
66  tightJetIdLepVetoAK8.filterParams, version = "RUN3WINTER22PUPPIrunsBCDEprompt"
67 )
68 
69 updatedJetsAK8WithUserData = cms.EDProducer("PATJetUserDataEmbedder",
70  src = cms.InputTag("updatedJetsAK8"),
71  userFloats = cms.PSet(),
72  userInts = cms.PSet(
73  tightId = cms.InputTag("tightJetIdAK8"),
74  tightIdLepVeto = cms.InputTag("tightJetIdLepVetoAK8"),
75  ),
76 )
77 
78 finalJetsAK8 = cms.EDFilter("PATJetRefSelector",
79  src = cms.InputTag("updatedJetsAK8WithUserData"),
80  cut = cms.string("pt > 170")
81 )
82 
83 
84 lepInAK8JetVars = cms.EDProducer("LepInJetProducer",
85  src = cms.InputTag("updatedJetsAK8WithUserData"),
86  srcEle = cms.InputTag("finalElectrons"),
87  srcMu = cms.InputTag("finalMuons")
88 )
89 
90 fatJetTable = simpleCandidateFlatTableProducer.clone(
91  src = cms.InputTag("finalJetsAK8"),
92  cut = cms.string(" pt > 170"), #probably already applied in miniaod
93  name = cms.string("FatJet"),
94  doc = cms.string("slimmedJetsAK8, i.e. ak8 fat jets for boosted analysis"),
95  variables = cms.PSet(P4Vars,
96  jetId = Var("userInt('tightId')*2+4*userInt('tightIdLepVeto')", "uint8",doc="Jet ID flags bit1 is loose (always false in 2017 since it does not exist), bit2 is tight, bit3 is tightLepVeto"),
97  area = Var("jetArea()", float, doc="jet catchment area, for JECs",precision=10),
98  rawFactor = Var("1.-jecFactor('Uncorrected')",float,doc="1 - Factor to get back to raw pT",precision=6),
99  tau1 = Var("userFloat('NjettinessAK8Puppi:tau1')",float, doc="Nsubjettiness (1 axis)",precision=10),
100  tau2 = Var("userFloat('NjettinessAK8Puppi:tau2')",float, doc="Nsubjettiness (2 axis)",precision=10),
101  tau3 = Var("userFloat('NjettinessAK8Puppi:tau3')",float, doc="Nsubjettiness (3 axis)",precision=10),
102  tau4 = Var("userFloat('NjettinessAK8Puppi:tau4')",float, doc="Nsubjettiness (4 axis)",precision=10),
103  n2b1 = Var("?hasUserFloat('nb1AK8PuppiSoftDrop:ecfN2')?userFloat('nb1AK8PuppiSoftDrop:ecfN2'):-99999.", float, doc="N2 with beta=1 (for jets with raw pT>250 GeV)", precision=10),
104  n3b1 = Var("?hasUserFloat('nb1AK8PuppiSoftDrop:ecfN3')?userFloat('nb1AK8PuppiSoftDrop:ecfN3'):-99999.", float, doc="N3 with beta=1 (for jets with raw pT>250 GeV)", precision=10),
105  msoftdrop = Var("groomedMass('SoftDropPuppi')",float, doc="Corrected soft drop mass with PUPPI",precision=10),
106  particleNetWithMass_QCD = Var("bDiscriminator('pfParticleNetJetTags:probQCDbb')+bDiscriminator('pfParticleNetJetTags:probQCDcc')+bDiscriminator('pfParticleNetJetTags:probQCDb')+bDiscriminator('pfParticleNetJetTags:probQCDc')+bDiscriminator('pfParticleNetJetTags:probQCDothers')",float,doc="ParticleNet tagger (w/ mass) QCD(bb,cc,b,c,others) sum",precision=10),
107  particleNetWithMass_TvsQCD = Var("bDiscriminator('pfParticleNetDiscriminatorsJetTags:TvsQCD')",float,doc="ParticleNet tagger (w/ mass) top vs QCD discriminator",precision=10),
108  particleNetWithMass_WvsQCD = Var("bDiscriminator('pfParticleNetDiscriminatorsJetTags:WvsQCD')",float,doc="ParticleNet tagger (w/ mass) W vs QCD discriminator",precision=10),
109  particleNetWithMass_ZvsQCD = Var("bDiscriminator('pfParticleNetDiscriminatorsJetTags:ZvsQCD')",float,doc="ParticleNet tagger (w/ mass) Z vs QCD discriminator",precision=10),
110  particleNetWithMass_H4qvsQCD = Var("bDiscriminator('pfParticleNetDiscriminatorsJetTags:H4qvsQCD')",float,doc="ParticleNet tagger (w/ mass) H(->VV->qqqq) vs QCD discriminator",precision=10),
111  particleNetWithMass_HbbvsQCD = Var("bDiscriminator('pfParticleNetDiscriminatorsJetTags:HbbvsQCD')",float,doc="ParticleNet tagger (w/mass) H(->bb) vs QCD discriminator",precision=10),
112  particleNetWithMass_HccvsQCD = Var("bDiscriminator('pfParticleNetDiscriminatorsJetTags:HccvsQCD')",float,doc="ParticleNet tagger (w/mass) H(->cc) vs QCD discriminator",precision=10),
113  particleNet_QCD = Var("bDiscriminator('pfParticleNetFromMiniAODAK8JetTags:probQCD2hf')+bDiscriminator('pfParticleNetFromMiniAODAK8JetTags:probQCD1hf')+bDiscriminator('pfParticleNetFromMiniAODAK8JetTags:probQCD0hf')",float,doc="ParticleNet tagger QCD(0+1+2HF) sum",precision=10),
114  particleNet_QCD2HF = Var("bDiscriminator('pfParticleNetFromMiniAODAK8JetTags:probQCD2hf')",float,doc="ParticleNet tagger QCD 2 HF (b/c) score",precision=10),
115  particleNet_QCD1HF = Var("bDiscriminator('pfParticleNetFromMiniAODAK8JetTags:probQCD1hf')",float,doc="ParticleNet tagger QCD 1 HF (b/c) score",precision=10),
116  particleNet_QCD0HF = Var("bDiscriminator('pfParticleNetFromMiniAODAK8JetTags:probQCD0hf')",float,doc="ParticleNet tagger QCD 0 HF (b/c) score",precision=10),
117  particleNet_massCorr = Var("bDiscriminator('pfParticleNetFromMiniAODAK8JetTags:masscorr')",float,doc="ParticleNet mass regression, relative correction to JEC-corrected jet mass (no softdrop)",precision=10),
118  particleNet_XbbVsQCD = Var("bDiscriminator('pfParticleNetFromMiniAODAK8DiscriminatorsJetTags:HbbvsQCD')",float,doc="ParticleNet X->bb vs. QCD score: Xbb/(Xbb+QCD)",precision=10),
119  particleNet_XccVsQCD = Var("bDiscriminator('pfParticleNetFromMiniAODAK8DiscriminatorsJetTags:HccvsQCD')",float,doc="ParticleNet X->cc vs. QCD score: Xcc/(Xcc+QCD)",precision=10),
120  particleNet_XqqVsQCD = Var("bDiscriminator('pfParticleNetFromMiniAODAK8DiscriminatorsJetTags:HqqvsQCD')",float,doc="ParticleNet X->qq (uds) vs. QCD score: Xqq/(Xqq+QCD)",precision=10),
121  particleNet_XggVsQCD = Var("bDiscriminator('pfParticleNetFromMiniAODAK8DiscriminatorsJetTags:HggvsQCD')",float,doc="ParticleNet X->gg vs. QCD score: Xgg/(Xgg+QCD)",precision=10),
122  particleNet_XttVsQCD = Var("bDiscriminator('pfParticleNetFromMiniAODAK8DiscriminatorsJetTags:HttvsQCD')",float,doc="ParticleNet X->tau_h tau_h vs. QCD score: Xtt/(Xtt+QCD)",precision=10),
123  particleNet_XtmVsQCD = Var("bDiscriminator('pfParticleNetFromMiniAODAK8DiscriminatorsJetTags:HtmvsQCD')",float,doc="ParticleNet X->mu tau_h vs. QCD score: Xtm/(Xtm+QCD)",precision=10),
124  particleNet_XteVsQCD = Var("bDiscriminator('pfParticleNetFromMiniAODAK8DiscriminatorsJetTags:HtevsQCD')",float,doc="ParticleNet X->e tau_h vs. QCD score: Xte/(Xte+QCD)",precision=10),
125  subJetIdx1 = Var("?nSubjetCollections()>0 && subjets('SoftDropPuppi').size()>0?subjets('SoftDropPuppi')[0].key():-1", "int16",
126  doc="index of first subjet"),
127  subJetIdx2 = Var("?nSubjetCollections()>0 && subjets('SoftDropPuppi').size()>1?subjets('SoftDropPuppi')[1].key():-1", "int16",
128  doc="index of second subjet"),
129  nConstituents = Var("numberOfDaughters()","uint8",doc="Number of particles in the jet"),
130  chMultiplicity = Var("?isPFJet()?chargedMultiplicity():-1","int16",doc="(Puppi-weighted) Number of charged particles in the jet"),
131  neMultiplicity = Var("?isPFJet()?neutralMultiplicity():-1","int16",doc="(Puppi-weighted) Number of neutral particles in the jet"),
132  chHEF = Var("?isPFJet()?chargedHadronEnergyFraction():-1", float, doc="charged Hadron Energy Fraction", precision=6),
133  neHEF = Var("?isPFJet()?neutralHadronEnergyFraction():-1", float, doc="neutral Hadron Energy Fraction", precision=6),
134  chEmEF = Var("?isPFJet()?chargedEmEnergyFraction():-1", float, doc="charged Electromagnetic Energy Fraction", precision=6),
135  neEmEF = Var("?isPFJet()?neutralEmEnergyFraction():-1", float, doc="neutral Electromagnetic Energy Fraction", precision=6),
136  muEF = Var("?isPFJet()?muonEnergyFraction():-1", float, doc="muon Energy Fraction", precision=6),
137  ),
138  externalVariables = cms.PSet(
139  lsf3 = ExtVar(cms.InputTag("lepInAK8JetVars:lsf3"),float, doc="Lepton Subjet Fraction (3 subjets)",precision=10),
140  muonIdx3SJ = ExtVar(cms.InputTag("lepInAK8JetVars:muIdx3SJ"),"int16", doc="index of muon matched to jet"),
141  electronIdx3SJ = ExtVar(cms.InputTag("lepInAK8JetVars:eleIdx3SJ"),"int16",doc="index of electron matched to jet"),
142  )
143 )
144 
145 run2_nanoAOD_ANY.toModify(
146  fatJetTable.variables,
147  btagCSVV2 = Var("bDiscriminator('pfCombinedInclusiveSecondaryVertexV2BJetTags')",float,doc=" pfCombinedInclusiveSecondaryVertexV2 b-tag discriminator (aka CSVV2)",precision=10),
148  # Remove for V9
149  chMultiplicity = None,
150  neMultiplicity = None,
151  chHEF = None,
152  neHEF = None,
153  chEmEF = None,
154  neEmEF = None,
155  muEF = None
156 )
157 (run3_nanoAOD_122 | run3_nanoAOD_124).toModify(
158  fatJetTable.variables,
159  # New ParticleNet trainings are not available in MiniAOD until Run3 13X
160  particleNet_QCD = None,
161  particleNet_QCD2HF = None,
162  particleNet_QCD1HF = None,
163  particleNet_QCD0HF = None,
164  particleNet_massCorr = None,
165  particleNet_XbbVsQCD = None,
166  particleNet_XccVsQCD = None,
167  particleNet_XqqVsQCD = None,
168  particleNet_XggVsQCD = None,
169  particleNet_XttVsQCD = None,
170  particleNet_XtmVsQCD = None,
171  particleNet_XteVsQCD = None,
172  # Remove for V11 and earlier versions
173  chMultiplicity = None,
174  neMultiplicity = None,
175  chHEF = None,
176  neHEF = None,
177  chEmEF = None,
178  neEmEF = None,
179  muEF = None
180 )
181 
182 (run2_nanoAOD_106Xv2 | run3_nanoAOD_122 | run3_nanoAOD_124).toModify(
183  fatJetTable.variables,
184  # Restore taggers that were decommisionned for Run-3
185  btagDeepB = Var("?(bDiscriminator('pfDeepCSVJetTags:probb')+bDiscriminator('pfDeepCSVJetTags:probbb'))>=0?bDiscriminator('pfDeepCSVJetTags:probb')+bDiscriminator('pfDeepCSVJetTags:probbb'):-1",float,doc="DeepCSV b+bb tag discriminator",precision=10),
186  btagHbb = Var("bDiscriminator('pfBoostedDoubleSecondaryVertexAK8BJetTags')",float,doc="Higgs to BB tagger discriminator",precision=10),
187  btagDDBvLV2 = Var("bDiscriminator('pfMassIndependentDeepDoubleBvLV2JetTags:probHbb')",float,doc="DeepDoubleX V2(mass-decorrelated) discriminator for H(Z)->bb vs QCD",precision=10),
188  btagDDCvLV2 = Var("bDiscriminator('pfMassIndependentDeepDoubleCvLV2JetTags:probHcc')",float,doc="DeepDoubleX V2 (mass-decorrelated) discriminator for H(Z)->cc vs QCD",precision=10),
189  btagDDCvBV2 = Var("bDiscriminator('pfMassIndependentDeepDoubleCvBV2JetTags:probHcc')",float,doc="DeepDoubleX V2 (mass-decorrelated) discriminator for H(Z)->cc vs H(Z)->bb",precision=10),
190  deepTag_TvsQCD = Var("bDiscriminator('pfDeepBoostedDiscriminatorsJetTags:TvsQCD')",float,doc="DeepBoostedJet tagger top vs QCD discriminator",precision=10),
191  deepTag_WvsQCD = Var("bDiscriminator('pfDeepBoostedDiscriminatorsJetTags:WvsQCD')",float,doc="DeepBoostedJet tagger W vs QCD discriminator",precision=10),
192  deepTag_ZvsQCD = Var("bDiscriminator('pfDeepBoostedDiscriminatorsJetTags:ZvsQCD')",float,doc="DeepBoostedJet tagger Z vs QCD discriminator",precision=10),
193  deepTag_H = Var("bDiscriminator('pfDeepBoostedJetTags:probHbb')+bDiscriminator('pfDeepBoostedJetTags:probHcc')+bDiscriminator('pfDeepBoostedJetTags:probHqqqq')",float,doc="DeepBoostedJet tagger H(bb,cc,4q) sum",precision=10),
194  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),
195  deepTag_QCDothers = Var("bDiscriminator('pfDeepBoostedJetTags:probQCDothers')",float,doc="DeepBoostedJet tagger QCDothers value",precision=10),
196  deepTagMD_TvsQCD = Var("bDiscriminator('pfMassDecorrelatedDeepBoostedDiscriminatorsJetTags:TvsQCD')",float,doc="Mass-decorrelated DeepBoostedJet tagger top vs QCD discriminator",precision=10),
197  deepTagMD_WvsQCD = Var("bDiscriminator('pfMassDecorrelatedDeepBoostedDiscriminatorsJetTags:WvsQCD')",float,doc="Mass-decorrelated DeepBoostedJet tagger W vs QCD discriminator",precision=10),
198  deepTagMD_ZvsQCD = Var("bDiscriminator('pfMassDecorrelatedDeepBoostedDiscriminatorsJetTags:ZvsQCD')",float,doc="Mass-decorrelated DeepBoostedJet tagger Z vs QCD discriminator",precision=10),
199  deepTagMD_ZHbbvsQCD = Var("bDiscriminator('pfMassDecorrelatedDeepBoostedDiscriminatorsJetTags:ZHbbvsQCD')",float,doc="Mass-decorrelated DeepBoostedJet tagger Z/H->bb vs QCD discriminator",precision=10),
200  deepTagMD_ZbbvsQCD = Var("bDiscriminator('pfMassDecorrelatedDeepBoostedDiscriminatorsJetTags:ZbbvsQCD')",float,doc="Mass-decorrelated DeepBoostedJet tagger Z->bb vs QCD discriminator",precision=10),
201  deepTagMD_HbbvsQCD = Var("bDiscriminator('pfMassDecorrelatedDeepBoostedDiscriminatorsJetTags:HbbvsQCD')",float,doc="Mass-decorrelated DeepBoostedJet tagger H->bb vs QCD discriminator",precision=10),
202  deepTagMD_ZHccvsQCD = Var("bDiscriminator('pfMassDecorrelatedDeepBoostedDiscriminatorsJetTags:ZHccvsQCD')",float,doc="Mass-decorrelated DeepBoostedJet tagger Z/H->cc vs QCD discriminator",precision=10),
203  deepTagMD_H4qvsQCD = Var("bDiscriminator('pfMassDecorrelatedDeepBoostedDiscriminatorsJetTags:H4qvsQCD')",float,doc="Mass-decorrelated DeepBoostedJet tagger H->4q vs QCD discriminator",precision=10),
204  deepTagMD_bbvsLight = Var("bDiscriminator('pfMassDecorrelatedDeepBoostedDiscriminatorsJetTags:bbvsLight')",float,doc="Mass-decorrelated DeepBoostedJet tagger Z/H/gluon->bb vs light flavour discriminator",precision=10),
205  deepTagMD_ccvsLight = Var("bDiscriminator('pfMassDecorrelatedDeepBoostedDiscriminatorsJetTags:ccvsLight')",float,doc="Mass-decorrelated DeepBoostedJet tagger Z/H/gluon->cc vs light flavour discriminator",precision=10),
206  particleNetLegacy_mass = Var("bDiscriminator('pfParticleNetMassRegressionJetTags:mass')",float,doc="ParticleNet Legacy Run-2 mass regression",precision=10),
207  particleNetLegacy_Xbb = Var("bDiscriminator('pfMassDecorrelatedParticleNetJetTags:probXbb')",float,doc="Mass-decorrelated ParticleNet Legacy Run-2 tagger raw X->bb score. For X->bb vs QCD tagging, use Xbb/(Xbb+QCD)",precision=10),
208  particleNetLegacy_Xcc = Var("bDiscriminator('pfMassDecorrelatedParticleNetJetTags:probXcc')",float,doc="Mass-decorrelated ParticleNet Legacy Run-2 tagger raw X->cc score. For X->cc vs QCD tagging, use Xcc/(Xcc+QCD)",precision=10),
209  particleNetLegacy_Xqq = Var("bDiscriminator('pfMassDecorrelatedParticleNetJetTags:probXqq')",float,doc="Mass-decorrelated ParticleNet Legacy Run-2 tagger raw X->qq (uds) score. For X->qq vs QCD tagging, use Xqq/(Xqq+QCD). For W vs QCD tagging, use (Xcc+Xqq)/(Xcc+Xqq+QCD)",precision=10),
210  particleNetLegacy_QCD = Var("bDiscriminator('pfMassDecorrelatedParticleNetJetTags:probQCDbb')+bDiscriminator('pfMassDecorrelatedParticleNetJetTags:probQCDcc')+bDiscriminator('pfMassDecorrelatedParticleNetJetTags:probQCDb')+bDiscriminator('pfMassDecorrelatedParticleNetJetTags:probQCDc')+bDiscriminator('pfMassDecorrelatedParticleNetJetTags:probQCDothers')",float,doc="Mass-decorrelated ParticleNet Legacy Run-2 tagger raw QCD score",precision=10),
211 )
212 
213 
217 from PhysicsTools.PatAlgos.tools.jetTools import updateJetCollection
218 def nanoAOD_addDeepInfoAK8(process, addDeepBTag, addDeepBoostedJet, addDeepDoubleX, addDeepDoubleXV2, addParticleNetMassLegacy, addParticleNet, jecPayload):
219  _btagDiscriminators=[]
220  if addDeepBTag:
221  print("Updating process to run DeepCSV btag to AK8 jets")
222  _btagDiscriminators += ['pfDeepCSVJetTags:probb','pfDeepCSVJetTags:probbb']
223  if addDeepBoostedJet:
224  print("Updating process to run DeepBoostedJet on datasets before 103X")
225  from RecoBTag.ONNXRuntime.pfDeepBoostedJet_cff import _pfDeepBoostedJetTagsAll as pfDeepBoostedJetTagsAll
226  _btagDiscriminators += pfDeepBoostedJetTagsAll
227  if addParticleNet:
228  print("Updating process to run ParticleNet joint classification and mass regression")
229  from RecoBTag.ONNXRuntime.pfParticleNetFromMiniAODAK8_cff import _pfParticleNetFromMiniAODAK8JetTagsAll as pfParticleNetFromMiniAODAK8JetTagsAll
230  _btagDiscriminators += pfParticleNetFromMiniAODAK8JetTagsAll
231  if addParticleNetMassLegacy:
232  from RecoBTag.ONNXRuntime.pfParticleNet_cff import _pfParticleNetMassRegressionOutputs
233  _btagDiscriminators += _pfParticleNetMassRegressionOutputs
234  if addDeepDoubleX:
235  print("Updating process to run DeepDoubleX on datasets before 104X")
236  _btagDiscriminators += ['pfDeepDoubleBvLJetTags:probHbb', \
237  'pfDeepDoubleCvLJetTags:probHcc', \
238  'pfDeepDoubleCvBJetTags:probHcc', \
239  'pfMassIndependentDeepDoubleBvLJetTags:probHbb', 'pfMassIndependentDeepDoubleCvLJetTags:probHcc', 'pfMassIndependentDeepDoubleCvBJetTags:probHcc']
240  if addDeepDoubleXV2:
241  print("Updating process to run DeepDoubleXv2 on datasets before 11X")
242  _btagDiscriminators += [
243  'pfMassIndependentDeepDoubleBvLV2JetTags:probHbb',
244  'pfMassIndependentDeepDoubleCvLV2JetTags:probHcc',
245  'pfMassIndependentDeepDoubleCvBV2JetTags:probHcc'
246  ]
247  if len(_btagDiscriminators)==0: return process
248  print("Will recalculate the following discriminators on AK8 jets: "+", ".join(_btagDiscriminators))
249  updateJetCollection(
250  process,
251  jetSource = cms.InputTag('slimmedJetsAK8'),
252  pvSource = cms.InputTag('offlineSlimmedPrimaryVertices'),
253  svSource = cms.InputTag('slimmedSecondaryVertices'),
254  rParam = 0.8,
255  jetCorrections = (jecPayload.value(), cms.vstring(['L1FastJet', 'L2Relative', 'L3Absolute', 'L2L3Residual']), 'None'),
256  btagDiscriminators = _btagDiscriminators,
257  postfix='AK8WithDeepInfo',
258  printWarning = False
259  )
260  process.jetCorrFactorsAK8.src="selectedUpdatedPatJetsAK8WithDeepInfo"
261  process.updatedJetsAK8.jetSource="selectedUpdatedPatJetsAK8WithDeepInfo"
262  return process
263 
264 nanoAOD_addDeepInfoAK8_switch = cms.PSet(
265  nanoAOD_addDeepBTag_switch = cms.untracked.bool(False),
266  nanoAOD_addDeepBoostedJet_switch = cms.untracked.bool(False),
267  nanoAOD_addDeepDoubleX_switch = cms.untracked.bool(False),
268  nanoAOD_addDeepDoubleXV2_switch = cms.untracked.bool(False),
269  nanoAOD_addParticleNet_switch = cms.untracked.bool(False),
270  nanoAOD_addParticleNetMassLegacy_switch = cms.untracked.bool(False),
271  jecPayload = cms.untracked.string('AK8PFPuppi')
272 )
273 
274 
275 # ParticleNet legacy jet tagger is already in 106Xv2 MINIAOD,
276 # add PartlceNet legacy mass regression and new combined tagger + mass regression
277 run2_nanoAOD_106Xv2.toModify(
278  nanoAOD_addDeepInfoAK8_switch,
279  nanoAOD_addParticleNet_switch = True,
280  nanoAOD_addParticleNetMassLegacy_switch = True
281 )
282 
283 
286 
287 subJetTable = simpleCandidateFlatTableProducer.clone(
288  src = cms.InputTag("slimmedJetsAK8PFPuppiSoftDropPacked","SubJets"),
289  name = cms.string("SubJet"),
290  doc = cms.string("slimmedJetsAK8PFPuppiSoftDropPacked::SubJets, i.e. soft-drop subjets for ak8 fat jets for boosted analysis"),
291  variables = cms.PSet(P4Vars,
292  btagDeepB = Var("bDiscriminator('pfDeepCSVJetTags:probb')+bDiscriminator('pfDeepCSVJetTags:probbb')",float,doc="DeepCSV b+bb tag discriminator",precision=10),
293  rawFactor = Var("1.-jecFactor('Uncorrected')",float,doc="1 - Factor to get back to raw pT",precision=6),
294  area = Var("jetArea()", float, doc="jet catchment area, for JECs",precision=10),
295  tau1 = Var("userFloat('NjettinessAK8Subjets:tau1')",float, doc="Nsubjettiness (1 axis)",precision=10),
296  tau2 = Var("userFloat('NjettinessAK8Subjets:tau2')",float, doc="Nsubjettiness (2 axis)",precision=10),
297  tau3 = Var("userFloat('NjettinessAK8Subjets:tau3')",float, doc="Nsubjettiness (3 axis)",precision=10),
298  tau4 = Var("userFloat('NjettinessAK8Subjets:tau4')",float, doc="Nsubjettiness (4 axis)",precision=10),
299  n2b1 = Var("userFloat('nb1AK8PuppiSoftDropSubjets:ecfN2')", float, doc="N2 with beta=1", precision=10),
300  n3b1 = Var("userFloat('nb1AK8PuppiSoftDropSubjets:ecfN3')", float, doc="N3 with beta=1", precision=10),
301  )
302 )
303 
304 run2_nanoAOD_ANY.toModify(
305  subJetTable.variables,
306  btagCSVV2 = Var("bDiscriminator('pfCombinedInclusiveSecondaryVertexV2BJetTags')",float,doc=" pfCombinedInclusiveSecondaryVertexV2 b-tag discriminator (aka CSVV2)",precision=10)
307 )
308 
309 (run2_nanoAOD_106Xv2 | run3_nanoAOD_122 | run3_nanoAOD_124).toModify(
310  subJetTable.variables,
311  area = None,
312 )
313 
314 #jets are not as precise as muons
315 fatJetTable.variables.pt.precision=10
316 subJetTable.variables.pt.precision=10
317 
318 jetAK8UserDataTask = cms.Task(tightJetIdAK8,tightJetIdLepVetoAK8)
319 jetAK8Task = cms.Task(jetCorrFactorsAK8,updatedJetsAK8,jetAK8UserDataTask,updatedJetsAK8WithUserData,finalJetsAK8)
320 
321 #after lepton collections have been run
322 jetAK8LepTask = cms.Task(lepInAK8JetVars)
323 
324 jetAK8TablesTask = cms.Task(fatJetTable,subJetTable)
def ExtVar(tag, valtype, doc=None, precision=-1)
Definition: common_cff.py:27
def Var(expr, valtype, doc=None, precision=-1)
Definition: common_cff.py:16
User floats producers, selectors ##########################.
void print(TMatrixD &m, const char *label=nullptr, bool mathematicaFormat=false)
Definition: Utilities.cc:47
def nanoAOD_addDeepInfoAK8(process, addDeepBTag, addDeepBoostedJet, addDeepDoubleX, addDeepDoubleXV2, addParticleNetMassLegacy, addParticleNet, jecPayload)
Definition: jetsAK8_cff.py:218
static std::string join(char **cmd)
Definition: RemoteFile.cc:19