File indexing completed on 2024-11-05 05:20:42
0001 import FWCore.ParameterSet.Config as cms
0002
0003 from PhysicsTools.NanoAOD.nano_eras_cff import *
0004 from PhysicsTools.NanoAOD.common_cff import *
0005 from PhysicsTools.NanoAOD.simplePATJetFlatTableProducer_cfi import simplePATJetFlatTableProducer
0006
0007 from PhysicsTools.PatAlgos.recoLayer0.jetCorrFactors_cfi import *
0008
0009
0010 jetCorrFactorsAK8 = patJetCorrFactors.clone(src='slimmedJetsAK8',
0011 levels = cms.vstring('L1FastJet',
0012 'L2Relative',
0013 'L3Absolute',
0014 'L2L3Residual'),
0015 payload = cms.string('AK8PFPuppi'),
0016 primaryVertices = cms.InputTag("offlineSlimmedPrimaryVertices"),
0017 )
0018
0019 from PhysicsTools.PatAlgos.producersLayer1.jetUpdater_cfi import *
0020 updatedJetsAK8 = updatedPatJets.clone(
0021 addBTagInfo=False,
0022 jetSource='slimmedJetsAK8',
0023 jetCorrFactorsSource=cms.VInputTag(cms.InputTag("jetCorrFactorsAK8") ),
0024 )
0025
0026 updatedJetsAK8WithUserData = cms.EDProducer("PATJetUserDataEmbedder",
0027 src = cms.InputTag("updatedJetsAK8"),
0028 userFloats = cms.PSet(),
0029 userInts = cms.PSet(),
0030 )
0031
0032 finalJetsAK8 = cms.EDFilter("PATJetRefSelector",
0033 src = cms.InputTag("updatedJetsAK8WithUserData"),
0034 cut = cms.string("pt > 170")
0035 )
0036
0037
0038 lepInAK8JetVars = cms.EDProducer("LepInJetProducer",
0039 src = cms.InputTag("updatedJetsAK8WithUserData"),
0040 srcEle = cms.InputTag("finalElectrons"),
0041 srcMu = cms.InputTag("finalMuons")
0042 )
0043
0044 fatJetTable = simplePATJetFlatTableProducer.clone(
0045 src = cms.InputTag("finalJetsAK8"),
0046 cut = cms.string(" pt > 170"),
0047 name = cms.string("FatJet"),
0048 doc = cms.string("slimmedJetsAK8, i.e. ak8 fat jets for boosted analysis"),
0049 variables = cms.PSet(P4Vars,
0050 area = Var("jetArea()", float, doc="jet catchment area, for JECs",precision=10),
0051 rawFactor = Var("1.-jecFactor('Uncorrected')",float,doc="1 - Factor to get back to raw pT",precision=6),
0052 tau1 = Var("userFloat('NjettinessAK8Puppi:tau1')",float, doc="Nsubjettiness (1 axis)",precision=10),
0053 tau2 = Var("userFloat('NjettinessAK8Puppi:tau2')",float, doc="Nsubjettiness (2 axis)",precision=10),
0054 tau3 = Var("userFloat('NjettinessAK8Puppi:tau3')",float, doc="Nsubjettiness (3 axis)",precision=10),
0055 tau4 = Var("userFloat('NjettinessAK8Puppi:tau4')",float, doc="Nsubjettiness (4 axis)",precision=10),
0056 n2b1 = Var("?hasUserFloat('nb1AK8PuppiSoftDrop:ecfN2')?userFloat('nb1AK8PuppiSoftDrop:ecfN2'):-99999.", float, doc="N2 with beta=1 (for jets with raw pT>250 GeV)", precision=10),
0057 n3b1 = Var("?hasUserFloat('nb1AK8PuppiSoftDrop:ecfN3')?userFloat('nb1AK8PuppiSoftDrop:ecfN3'):-99999.", float, doc="N3 with beta=1 (for jets with raw pT>250 GeV)", precision=10),
0058 msoftdrop = Var("groomedMass('SoftDropPuppi')",float, doc="Corrected soft drop mass with PUPPI",precision=10),
0059 globalParT3_Xbb = Var("bDiscriminator('pfGlobalParticleTransformerAK8JetTags:probXbb')",float,doc="Mass-decorrelated GlobalParT-3 X->bb score. Note: For sig vs bkg (e.g. bkg=QCD) tagging, use sig/(sig+bkg) to construct the discriminator",precision=10),
0060 globalParT3_Xcc = Var("bDiscriminator('pfGlobalParticleTransformerAK8JetTags:probXcc')",float,doc="Mass-decorrelated GlobalParT-3 X->cc score",precision=10),
0061 globalParT3_Xcs = Var("bDiscriminator('pfGlobalParticleTransformerAK8JetTags:probXcs')",float,doc="Mass-decorrelated GlobalParT-3 X->cs score",precision=10),
0062 globalParT3_Xqq = Var("bDiscriminator('pfGlobalParticleTransformerAK8JetTags:probXqq')",float,doc="Mass-decorrelated GlobalParT-3 X->qq (ss/dd/uu) score",precision=10),
0063 globalParT3_Xtauhtaue = Var("bDiscriminator('pfGlobalParticleTransformerAK8JetTags:probXtauhtaue')",float,doc="Mass-decorrelated GlobalParT-3 X->tauhtaue score",precision=10),
0064 globalParT3_Xtauhtaum = Var("bDiscriminator('pfGlobalParticleTransformerAK8JetTags:probXtauhtaum')",float,doc="Mass-decorrelated GlobalParT-3 X->tauhtaum score",precision=10),
0065 globalParT3_Xtauhtauh = Var("bDiscriminator('pfGlobalParticleTransformerAK8JetTags:probXtauhtauh')",float,doc="Mass-decorrelated GlobalParT-3 X->tauhtauh score",precision=10),
0066 globalParT3_XWW4q = Var("bDiscriminator('pfGlobalParticleTransformerAK8JetTags:probXWW4q')",float,doc="Mass-decorrelated GlobalParT-3 X->WW4q score",precision=10),
0067 globalParT3_XWW3q = Var("bDiscriminator('pfGlobalParticleTransformerAK8JetTags:probXWW3q')",float,doc="Mass-decorrelated GlobalParT-3 X->WW3q score",precision=10),
0068 globalParT3_XWWqqev = Var("bDiscriminator('pfGlobalParticleTransformerAK8JetTags:probXWWqqev')",float,doc="Mass-decorrelated GlobalParT-3 X->WWqqev score",precision=10),
0069 globalParT3_XWWqqmv = Var("bDiscriminator('pfGlobalParticleTransformerAK8JetTags:probXWWqqmv')",float,doc="Mass-decorrelated GlobalParT-3 X->WWqqmv score",precision=10),
0070 globalParT3_TopbWqq = Var("bDiscriminator('pfGlobalParticleTransformerAK8JetTags:probTopbWqq')",float,doc="Mass-decorrelated GlobalParT-3 Top->bWqq score",precision=10),
0071 globalParT3_TopbWq = Var("bDiscriminator('pfGlobalParticleTransformerAK8JetTags:probTopbWq')",float,doc="Mass-decorrelated GlobalParT-3 Top->bWq score",precision=10),
0072 globalParT3_TopbWev = Var("bDiscriminator('pfGlobalParticleTransformerAK8JetTags:probTopbWev')",float,doc="Mass-decorrelated GlobalParT-3 Top->bWev score",precision=10),
0073 globalParT3_TopbWmv = Var("bDiscriminator('pfGlobalParticleTransformerAK8JetTags:probTopbWmv')",float,doc="Mass-decorrelated GlobalParT-3 Top->bWmv score",precision=10),
0074 globalParT3_TopbWtauhv = Var("bDiscriminator('pfGlobalParticleTransformerAK8JetTags:probTopbWtauhv')",float,doc="Mass-decorrelated GlobalParT-3 Top->bWtauhv score",precision=10),
0075 globalParT3_QCD = Var("bDiscriminator('pfGlobalParticleTransformerAK8JetTags:probQCD')",float,doc="Mass-decorrelated GlobalParT-3 QCD score.",precision=10),
0076 globalParT3_massCorrX2p = Var("bDiscriminator('pfGlobalParticleTransformerAK8JetTags:massCorrX2p')",float,doc="GlobalParT-3 mass regression corrector with respect to the original jet mass, optimised for resonance 2-prong (bb/cc/cs/ss/qq) jets. Use (massCorrX2p * mass * (1 - rawFactor)) to get the regressed mass",precision=10),
0077 globalParT3_massCorrGeneric = Var("bDiscriminator('pfGlobalParticleTransformerAK8JetTags:massCorrGeneric')",float,doc="GlobalParT-3 mass regression corrector with respect to the original jet mass, optimised for generic jet cases. Use (massCorrGeneric * mass * (1 - rawFactor)) to get the regressed mass",precision=10),
0078 globalParT3_withMassTopvsQCD = Var("bDiscriminator('pfGlobalParticleTransformerAK8JetTags:probWithMassTopvsQCD')",float,doc="GlobalParT-3 tagger (w/mass) Top vs QCD discriminator",precision=10),
0079 globalParT3_withMassWvsQCD = Var("bDiscriminator('pfGlobalParticleTransformerAK8JetTags:probWithMassWvsQCD')",float,doc="GlobalParT-3 tagger (w/mass) W vs QCD discriminator",precision=10),
0080 globalParT3_withMassZvsQCD = Var("bDiscriminator('pfGlobalParticleTransformerAK8JetTags:probWithMassZvsQCD')",float,doc="GlobalParT-3 tagger (w/mass) Z vs QCD discriminator",precision=10),
0081 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),
0082 particleNetWithMass_TvsQCD = Var("bDiscriminator('pfParticleNetDiscriminatorsJetTags:TvsQCD')",float,doc="ParticleNet tagger (w/ mass) top vs QCD discriminator",precision=10),
0083 particleNetWithMass_WvsQCD = Var("bDiscriminator('pfParticleNetDiscriminatorsJetTags:WvsQCD')",float,doc="ParticleNet tagger (w/ mass) W vs QCD discriminator",precision=10),
0084 particleNetWithMass_ZvsQCD = Var("bDiscriminator('pfParticleNetDiscriminatorsJetTags:ZvsQCD')",float,doc="ParticleNet tagger (w/ mass) Z vs QCD discriminator",precision=10),
0085 particleNetWithMass_H4qvsQCD = Var("bDiscriminator('pfParticleNetDiscriminatorsJetTags:H4qvsQCD')",float,doc="ParticleNet tagger (w/ mass) H(->VV->qqqq) vs QCD discriminator",precision=10),
0086 particleNetWithMass_HbbvsQCD = Var("bDiscriminator('pfParticleNetDiscriminatorsJetTags:HbbvsQCD')",float,doc="ParticleNet tagger (w/mass) H(->bb) vs QCD discriminator",precision=10),
0087 particleNetWithMass_HccvsQCD = Var("bDiscriminator('pfParticleNetDiscriminatorsJetTags:HccvsQCD')",float,doc="ParticleNet tagger (w/mass) H(->cc) vs QCD discriminator",precision=10),
0088 particleNet_QCD = Var("bDiscriminator('pfParticleNetFromMiniAODAK8JetTags:probQCD2hf')+bDiscriminator('pfParticleNetFromMiniAODAK8JetTags:probQCD1hf')+bDiscriminator('pfParticleNetFromMiniAODAK8JetTags:probQCD0hf')",float,doc="ParticleNet tagger QCD(0+1+2HF) sum",precision=10),
0089 particleNet_QCD2HF = Var("bDiscriminator('pfParticleNetFromMiniAODAK8JetTags:probQCD2hf')",float,doc="ParticleNet tagger QCD 2 HF (b/c) score",precision=10),
0090 particleNet_QCD1HF = Var("bDiscriminator('pfParticleNetFromMiniAODAK8JetTags:probQCD1hf')",float,doc="ParticleNet tagger QCD 1 HF (b/c) score",precision=10),
0091 particleNet_QCD0HF = Var("bDiscriminator('pfParticleNetFromMiniAODAK8JetTags:probQCD0hf')",float,doc="ParticleNet tagger QCD 0 HF (b/c) score",precision=10),
0092 particleNet_massCorr = Var("bDiscriminator('pfParticleNetFromMiniAODAK8JetTags:masscorr')",float,doc="ParticleNet mass regression, relative correction to JEC-corrected jet mass (no softdrop)",precision=10),
0093 particleNet_XbbVsQCD = Var("bDiscriminator('pfParticleNetFromMiniAODAK8DiscriminatorsJetTags:HbbvsQCD')",float,doc="ParticleNet X->bb vs. QCD score: Xbb/(Xbb+QCD)",precision=10),
0094 particleNet_XccVsQCD = Var("bDiscriminator('pfParticleNetFromMiniAODAK8DiscriminatorsJetTags:HccvsQCD')",float,doc="ParticleNet X->cc vs. QCD score: Xcc/(Xcc+QCD)",precision=10),
0095 particleNet_XqqVsQCD = Var("bDiscriminator('pfParticleNetFromMiniAODAK8DiscriminatorsJetTags:HqqvsQCD')",float,doc="ParticleNet X->qq (uds) vs. QCD score: Xqq/(Xqq+QCD)",precision=10),
0096 particleNet_XggVsQCD = Var("bDiscriminator('pfParticleNetFromMiniAODAK8DiscriminatorsJetTags:HggvsQCD')",float,doc="ParticleNet X->gg vs. QCD score: Xgg/(Xgg+QCD)",precision=10),
0097 particleNet_XttVsQCD = Var("bDiscriminator('pfParticleNetFromMiniAODAK8DiscriminatorsJetTags:HttvsQCD')",float,doc="ParticleNet X->tau_h tau_h vs. QCD score: Xtt/(Xtt+QCD)",precision=10),
0098 particleNet_XtmVsQCD = Var("bDiscriminator('pfParticleNetFromMiniAODAK8DiscriminatorsJetTags:HtmvsQCD')",float,doc="ParticleNet X->mu tau_h vs. QCD score: Xtm/(Xtm+QCD)",precision=10),
0099 particleNet_XteVsQCD = Var("bDiscriminator('pfParticleNetFromMiniAODAK8DiscriminatorsJetTags:HtevsQCD')",float,doc="ParticleNet X->e tau_h vs. QCD score: Xte/(Xte+QCD)",precision=10),
0100 particleNet_WVsQCD = Var("bDiscriminator('pfParticleNetFromMiniAODAK8DiscriminatorsJetTags:WvsQCD')",float,doc="ParticleNet W->qq vs. QCD score: Xqq+Xcc/(Xqq+Xcc+QCD)",precision=10),
0101 particleNetLegacy_mass = Var("bDiscriminator('pfParticleNetMassRegressionJetTags:mass')",float,doc="ParticleNet Legacy Run-2 mass regression",precision=10),
0102 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),
0103 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),
0104 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),
0105 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),
0106 subJetIdx1 = Var("?nSubjetCollections()>0 && subjets('SoftDropPuppi').size()>0?subjets('SoftDropPuppi')[0].key():-1", "int16",
0107 doc="index of first subjet"),
0108 subJetIdx2 = Var("?nSubjetCollections()>0 && subjets('SoftDropPuppi').size()>1?subjets('SoftDropPuppi')[1].key():-1", "int16",
0109 doc="index of second subjet"),
0110 nConstituents = Var("numberOfDaughters()","uint8",doc="Number of particles in the jet"),
0111 chMultiplicity = Var("?isPFJet()?chargedMultiplicity():-1","int16",doc="(Puppi-weighted) Number of charged particles in the jet"),
0112 neMultiplicity = Var("?isPFJet()?neutralMultiplicity():-1","int16",doc="(Puppi-weighted) Number of neutral particles in the jet"),
0113 chHEF = Var("?isPFJet()?chargedHadronEnergyFraction():-1", float, doc="charged Hadron Energy Fraction", precision=10),
0114 neHEF = Var("?isPFJet()?neutralHadronEnergyFraction():-1", float, doc="neutral Hadron Energy Fraction", precision=10),
0115 chEmEF = Var("?isPFJet()?chargedEmEnergyFraction():-1", float, doc="charged Electromagnetic Energy Fraction", precision=10),
0116 neEmEF = Var("?isPFJet()?neutralEmEnergyFraction():-1", float, doc="neutral Electromagnetic Energy Fraction", precision=10),
0117 hfHEF = Var("?isPFJet()?HFHadronEnergyFraction():-1",float,doc="hadronic Energy Fraction in HF",precision=10),
0118 hfEmEF = Var("?isPFJet()?HFEMEnergyFraction():-1",float,doc="electromagnetic Energy Fraction in HF",precision=10),
0119 muEF = Var("?isPFJet()?muonEnergyFraction():-1", float, doc="muon Energy Fraction", precision=10),
0120 ),
0121 externalVariables = cms.PSet(
0122 lsf3 = ExtVar(cms.InputTag("lepInAK8JetVars:lsf3"),float, doc="Lepton Subjet Fraction (3 subjets)",precision=10),
0123 muonIdx3SJ = ExtVar(cms.InputTag("lepInAK8JetVars:muIdx3SJ"),"int16", doc="index of muon matched to jet"),
0124 electronIdx3SJ = ExtVar(cms.InputTag("lepInAK8JetVars:eleIdx3SJ"),"int16",doc="index of electron matched to jet"),
0125 )
0126 )
0127
0128 run2_nanoAOD_ANY.toModify(
0129 fatJetTable.variables,
0130 btagCSVV2 = Var("bDiscriminator('pfCombinedInclusiveSecondaryVertexV2BJetTags')",float,doc=" pfCombinedInclusiveSecondaryVertexV2 b-tag discriminator (aka CSVV2)",precision=10),
0131
0132 chMultiplicity = None,
0133 neMultiplicity = None,
0134 chHEF = None,
0135 neHEF = None,
0136 chEmEF = None,
0137 neEmEF = None,
0138 muEF = None
0139 )
0140
0141 (run2_nanoAOD_106Xv2).toModify(
0142 fatJetTable.variables,
0143
0144 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),
0145 btagHbb = Var("bDiscriminator('pfBoostedDoubleSecondaryVertexAK8BJetTags')",float,doc="Higgs to BB tagger discriminator",precision=10),
0146 btagDDBvLV2 = Var("bDiscriminator('pfMassIndependentDeepDoubleBvLV2JetTags:probHbb')",float,doc="DeepDoubleX V2(mass-decorrelated) discriminator for H(Z)->bb vs QCD",precision=10),
0147 btagDDCvLV2 = Var("bDiscriminator('pfMassIndependentDeepDoubleCvLV2JetTags:probHcc')",float,doc="DeepDoubleX V2 (mass-decorrelated) discriminator for H(Z)->cc vs QCD",precision=10),
0148 btagDDCvBV2 = Var("bDiscriminator('pfMassIndependentDeepDoubleCvBV2JetTags:probHcc')",float,doc="DeepDoubleX V2 (mass-decorrelated) discriminator for H(Z)->cc vs H(Z)->bb",precision=10),
0149 deepTag_TvsQCD = Var("bDiscriminator('pfDeepBoostedDiscriminatorsJetTags:TvsQCD')",float,doc="DeepBoostedJet tagger top vs QCD discriminator",precision=10),
0150 deepTag_WvsQCD = Var("bDiscriminator('pfDeepBoostedDiscriminatorsJetTags:WvsQCD')",float,doc="DeepBoostedJet tagger W vs QCD discriminator",precision=10),
0151 deepTag_ZvsQCD = Var("bDiscriminator('pfDeepBoostedDiscriminatorsJetTags:ZvsQCD')",float,doc="DeepBoostedJet tagger Z vs QCD discriminator",precision=10),
0152 deepTag_H = Var("bDiscriminator('pfDeepBoostedJetTags:probHbb')+bDiscriminator('pfDeepBoostedJetTags:probHcc')+bDiscriminator('pfDeepBoostedJetTags:probHqqqq')",float,doc="DeepBoostedJet tagger H(bb,cc,4q) sum",precision=10),
0153 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),
0154 deepTag_QCDothers = Var("bDiscriminator('pfDeepBoostedJetTags:probQCDothers')",float,doc="DeepBoostedJet tagger QCDothers value",precision=10),
0155 deepTagMD_TvsQCD = Var("bDiscriminator('pfMassDecorrelatedDeepBoostedDiscriminatorsJetTags:TvsQCD')",float,doc="Mass-decorrelated DeepBoostedJet tagger top vs QCD discriminator",precision=10),
0156 deepTagMD_WvsQCD = Var("bDiscriminator('pfMassDecorrelatedDeepBoostedDiscriminatorsJetTags:WvsQCD')",float,doc="Mass-decorrelated DeepBoostedJet tagger W vs QCD discriminator",precision=10),
0157 deepTagMD_ZvsQCD = Var("bDiscriminator('pfMassDecorrelatedDeepBoostedDiscriminatorsJetTags:ZvsQCD')",float,doc="Mass-decorrelated DeepBoostedJet tagger Z vs QCD discriminator",precision=10),
0158 deepTagMD_ZHbbvsQCD = Var("bDiscriminator('pfMassDecorrelatedDeepBoostedDiscriminatorsJetTags:ZHbbvsQCD')",float,doc="Mass-decorrelated DeepBoostedJet tagger Z/H->bb vs QCD discriminator",precision=10),
0159 deepTagMD_ZbbvsQCD = Var("bDiscriminator('pfMassDecorrelatedDeepBoostedDiscriminatorsJetTags:ZbbvsQCD')",float,doc="Mass-decorrelated DeepBoostedJet tagger Z->bb vs QCD discriminator",precision=10),
0160 deepTagMD_HbbvsQCD = Var("bDiscriminator('pfMassDecorrelatedDeepBoostedDiscriminatorsJetTags:HbbvsQCD')",float,doc="Mass-decorrelated DeepBoostedJet tagger H->bb vs QCD discriminator",precision=10),
0161 deepTagMD_ZHccvsQCD = Var("bDiscriminator('pfMassDecorrelatedDeepBoostedDiscriminatorsJetTags:ZHccvsQCD')",float,doc="Mass-decorrelated DeepBoostedJet tagger Z/H->cc vs QCD discriminator",precision=10),
0162 deepTagMD_H4qvsQCD = Var("bDiscriminator('pfMassDecorrelatedDeepBoostedDiscriminatorsJetTags:H4qvsQCD')",float,doc="Mass-decorrelated DeepBoostedJet tagger H->4q vs QCD discriminator",precision=10),
0163 deepTagMD_bbvsLight = Var("bDiscriminator('pfMassDecorrelatedDeepBoostedDiscriminatorsJetTags:bbvsLight')",float,doc="Mass-decorrelated DeepBoostedJet tagger Z/H/gluon->bb vs light flavour discriminator",precision=10),
0164 deepTagMD_ccvsLight = Var("bDiscriminator('pfMassDecorrelatedDeepBoostedDiscriminatorsJetTags:ccvsLight')",float,doc="Mass-decorrelated DeepBoostedJet tagger Z/H/gluon->cc vs light flavour discriminator",precision=10),
0165 )
0166
0167
0168
0169
0170
0171 from PhysicsTools.PatAlgos.tools.jetTools import updateJetCollection
0172 def nanoAOD_addDeepInfoAK8(process, addDeepBTag, addDeepBoostedJet, addDeepDoubleX, addDeepDoubleXV2, addParticleNetMassLegacy, addParticleNet, addGlobalParT, jecPayload):
0173 _btagDiscriminators=[]
0174 if addDeepBTag:
0175 print("Updating process to run DeepCSV btag to AK8 jets")
0176 _btagDiscriminators += ['pfDeepCSVJetTags:probb','pfDeepCSVJetTags:probbb']
0177 if addDeepBoostedJet:
0178 print("Updating process to run DeepBoostedJet on datasets before 103X")
0179 from RecoBTag.ONNXRuntime.pfDeepBoostedJet_cff import _pfDeepBoostedJetTagsAll as pfDeepBoostedJetTagsAll
0180 _btagDiscriminators += pfDeepBoostedJetTagsAll
0181 if addGlobalParT:
0182 print("Updating process to run GlobalParT")
0183 from RecoBTag.ONNXRuntime.pfGlobalParticleTransformerAK8_cff import _pfGlobalParticleTransformerAK8JetTagsAll as pfGlobalParticleTransformerAK8JetTagsAll
0184 _btagDiscriminators += pfGlobalParticleTransformerAK8JetTagsAll
0185 if addParticleNet:
0186 print("Updating process to run ParticleNet joint classification and mass regression")
0187 from RecoBTag.ONNXRuntime.pfParticleNetFromMiniAODAK8_cff import _pfParticleNetFromMiniAODAK8JetTagsAll as pfParticleNetFromMiniAODAK8JetTagsAll
0188 _btagDiscriminators += pfParticleNetFromMiniAODAK8JetTagsAll
0189 if addParticleNetMassLegacy:
0190 from RecoBTag.ONNXRuntime.pfParticleNet_cff import _pfParticleNetMassRegressionOutputs
0191 _btagDiscriminators += _pfParticleNetMassRegressionOutputs
0192 if addDeepDoubleX:
0193 print("Updating process to run DeepDoubleX on datasets before 104X")
0194 _btagDiscriminators += ['pfDeepDoubleBvLJetTags:probHbb', \
0195 'pfDeepDoubleCvLJetTags:probHcc', \
0196 'pfDeepDoubleCvBJetTags:probHcc', \
0197 'pfMassIndependentDeepDoubleBvLJetTags:probHbb', 'pfMassIndependentDeepDoubleCvLJetTags:probHcc', 'pfMassIndependentDeepDoubleCvBJetTags:probHcc']
0198 if addDeepDoubleXV2:
0199 print("Updating process to run DeepDoubleXv2 on datasets before 11X")
0200 _btagDiscriminators += [
0201 'pfMassIndependentDeepDoubleBvLV2JetTags:probHbb',
0202 'pfMassIndependentDeepDoubleCvLV2JetTags:probHcc',
0203 'pfMassIndependentDeepDoubleCvBV2JetTags:probHcc'
0204 ]
0205 if len(_btagDiscriminators)==0: return process
0206 print("Will recalculate the following discriminators on AK8 jets: "+", ".join(_btagDiscriminators))
0207 updateJetCollection(
0208 process,
0209 jetSource = cms.InputTag('slimmedJetsAK8'),
0210 pvSource = cms.InputTag('offlineSlimmedPrimaryVertices'),
0211 svSource = cms.InputTag('slimmedSecondaryVertices'),
0212 rParam = 0.8,
0213 jetCorrections = (jecPayload.value(), cms.vstring(['L1FastJet', 'L2Relative', 'L3Absolute', 'L2L3Residual']), 'None'),
0214 btagDiscriminators = _btagDiscriminators,
0215 postfix='AK8WithDeepInfo',
0216 printWarning = False
0217 )
0218 process.jetCorrFactorsAK8.src="selectedUpdatedPatJetsAK8WithDeepInfo"
0219 process.updatedJetsAK8.jetSource="selectedUpdatedPatJetsAK8WithDeepInfo"
0220 return process
0221
0222 nanoAOD_addDeepInfoAK8_switch = cms.PSet(
0223 nanoAOD_addDeepBTag_switch = cms.untracked.bool(False),
0224 nanoAOD_addDeepBoostedJet_switch = cms.untracked.bool(False),
0225 nanoAOD_addDeepDoubleX_switch = cms.untracked.bool(False),
0226 nanoAOD_addDeepDoubleXV2_switch = cms.untracked.bool(False),
0227 nanoAOD_addParticleNetMassLegacy_switch = cms.untracked.bool(False),
0228 nanoAOD_addParticleNet_switch = cms.untracked.bool(False),
0229 nanoAOD_addGlobalParT_switch = cms.untracked.bool(False),
0230 jecPayload = cms.untracked.string('AK8PFPuppi')
0231 )
0232
0233
0234
0235
0236 run2_nanoAOD_106Xv2.toModify(
0237 nanoAOD_addDeepInfoAK8_switch,
0238 nanoAOD_addParticleNetMassLegacy_switch = True,
0239 nanoAOD_addParticleNet_switch = True,
0240 nanoAOD_addGlobalParT_switch = True,
0241 )
0242
0243
0244
0245
0246
0247 subJetTable = simplePATJetFlatTableProducer.clone(
0248 src = cms.InputTag("slimmedJetsAK8PFPuppiSoftDropPacked","SubJets"),
0249 name = cms.string("SubJet"),
0250 doc = cms.string("slimmedJetsAK8PFPuppiSoftDropPacked::SubJets, i.e. soft-drop subjets for ak8 fat jets for boosted analysis"),
0251 variables = cms.PSet(P4Vars,
0252 btagDeepFlavB = Var("bDiscriminator('pfDeepFlavourJetTags:probb')+bDiscriminator('pfDeepFlavourJetTags:probbb')+bDiscriminator('pfDeepFlavourJetTags:problepb')",float,doc="DeepJet b+bb+lepb tag discriminator",precision=10),
0253 btagUParTAK4B = Var("?bDiscriminator('pfUnifiedParticleTransformerAK4DiscriminatorsJetTags:BvsAll')>0?bDiscriminator('pfUnifiedParticleTransformerAK4DiscriminatorsJetTags:BvsAll'):-1",float,precision=10,doc="UnifiedParT b vs. udscg"),
0254 rawFactor = Var("1.-jecFactor('Uncorrected')",float,doc="1 - Factor to get back to raw pT",precision=6),
0255 area = Var("jetArea()", float, doc="jet catchment area, for JECs",precision=10),
0256 tau1 = Var("userFloat('NjettinessAK8Subjets:tau1')",float, doc="Nsubjettiness (1 axis)",precision=10),
0257 tau2 = Var("userFloat('NjettinessAK8Subjets:tau2')",float, doc="Nsubjettiness (2 axis)",precision=10),
0258 tau3 = Var("userFloat('NjettinessAK8Subjets:tau3')",float, doc="Nsubjettiness (3 axis)",precision=10),
0259 tau4 = Var("userFloat('NjettinessAK8Subjets:tau4')",float, doc="Nsubjettiness (4 axis)",precision=10),
0260 n2b1 = Var("userFloat('nb1AK8PuppiSoftDropSubjets:ecfN2')", float, doc="N2 with beta=1", precision=10),
0261 n3b1 = Var("userFloat('nb1AK8PuppiSoftDropSubjets:ecfN3')", float, doc="N3 with beta=1", precision=10),
0262 )
0263 )
0264
0265 run2_nanoAOD_ANY.toModify(
0266 subJetTable.variables,
0267 btagCSVV2 = Var("bDiscriminator('pfCombinedInclusiveSecondaryVertexV2BJetTags')",float,doc=" pfCombinedInclusiveSecondaryVertexV2 b-tag discriminator (aka CSVV2)",precision=10)
0268 )
0269
0270 (run2_nanoAOD_106Xv2).toModify(
0271 subJetTable.variables,
0272 area = None,
0273 )
0274
0275 run3_nanoAOD_pre142X.toModify(
0276 subJetTable.variables,
0277 btagDeepFlavB = None,
0278 btagUParTAK4B = None,
0279 btagDeepB = Var("bDiscriminator('pfDeepCSVJetTags:probb')+bDiscriminator('pfDeepCSVJetTags:probbb')",float,doc="DeepCSV b+bb tag discriminator",precision=10),
0280 )
0281
0282
0283
0284 fatJetTable.variables.pt.precision=10
0285 subJetTable.variables.pt.precision=10
0286
0287
0288
0289
0290 finalJetsAK8PFConstituents = cms.EDProducer("PatJetConstituentPtrSelector",
0291 src = fatJetTable.src,
0292 cut = cms.string("abs(eta) <= 2.5")
0293 )
0294
0295 finalJetsPFConstituents = cms.EDProducer("PackedCandidatePtrMerger",
0296 src = cms.VInputTag(cms.InputTag("finalJetsAK8PFConstituents", "constituents")),
0297 skipNulls = cms.bool(True),
0298 warnOnSkip = cms.bool(True)
0299 )
0300
0301 pfCandidatesTable = cms.EDProducer("SimplePATCandidateFlatTableProducer",
0302 src = cms.InputTag("finalJetsPFConstituents"),
0303 cut = cms.string(""),
0304 name = cms.string("PFCand"),
0305 doc = cms.string("PF candidate constituents of AK8 puppi jets (FatJet) with |eta| <= 2.5."),
0306 singleton = cms.bool(False),
0307 extension = cms.bool(False),
0308 variables = cms.PSet(
0309 pt = Var("pt * puppiWeight()", float, doc="Puppi-weighted pt", precision=10),
0310 mass = Var("mass * puppiWeight()", float, doc="Puppi-weighted mass", precision=10),
0311 eta = Var("eta", float, precision=12),
0312 phi = Var("phi", float, precision=12),
0313 pdgId = Var("pdgId", int, doc="PF candidate type (+/-211 = ChgHad, 130 = NeuHad, 22 = Photon, +/-11 = Electron, +/-13 = Muon, 1 = HFHad, 2 = HFEM)")
0314 )
0315 )
0316
0317 finalJetsAK8ConstituentsTable = cms.EDProducer("SimplePatJetConstituentTableProducer",
0318 name = cms.string(fatJetTable.name.value()+"PFCand"),
0319 candIdxName = cms.string("PFCandIdx"),
0320 candIdxDoc = cms.string("Index in the PFCand table"),
0321 jets = fatJetTable.src,
0322 candidates = pfCandidatesTable.src,
0323 jetCut = fatJetTable.cut
0324 )
0325
0326 jetAK8UserDataTask = cms.Task()
0327 jetAK8Task = cms.Task(jetCorrFactorsAK8,updatedJetsAK8,jetAK8UserDataTask,updatedJetsAK8WithUserData,finalJetsAK8,finalJetsAK8PFConstituents,finalJetsPFConstituents)
0328
0329
0330 jetAK8LepTask = cms.Task(lepInAK8JetVars)
0331
0332 jetAK8TablesTask = cms.Task(fatJetTable,subJetTable,pfCandidatesTable,finalJetsAK8ConstituentsTable)
0333