Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-25 02:14:03

0001 import FWCore.ParameterSet.Config as cms
0002 from PhysicsTools.NanoAOD.common_cff import Var
0003 from PhysicsTools.NanoAOD.jetsAK4_Puppi_cff import jetPuppiTable, jetPuppiCorrFactorsNano, updatedJetsPuppi, updatedJetsPuppiWithUserData
0004 from PhysicsTools.NanoAOD.jetsAK8_cff import fatJetTable, subJetTable
0005 from PhysicsTools.PatAlgos.tools.jetTools import updateJetCollection
0006 from PhysicsTools.PatAlgos.tools.helpers import addToProcessAndTask, getPatAlgosToolsTask
0007 from PhysicsTools.NanoAOD.common_cff import Var, CandVars
0008 from PhysicsTools.NanoAOD.simpleCandidateFlatTableProducer_cfi import simpleCandidateFlatTableProducer
0009 from PhysicsTools.NanoAOD.btvMC_cff import allPFPFCandsMCSequence,ak4ak8PFCandsMCSequence,ak8onlyPFCandsMCSequence,ak4onlyPFCandsMCSequence
0010 ## Move PFNano (https://github.com/cms-jet/PFNano/) to NanoAOD
0011 
0012 ## From: https://github.com/cms-jet/PFNano/blob/13_0_7_from124MiniAOD/python/addBTV.py
0013 def update_jets_AK4(process):
0014     # Based on ``nanoAOD_addDeepInfo``
0015     # in https://github.com/cms-sw/cmssw/blob/master/PhysicsTools/NanoAOD/python/nano_cff.py
0016     # DeepJet flav_names as found in
0017     # https://github.com/cms-sw/cmssw/blob/master/RecoBTag/ONNXRuntime/plugins/DeepFlavourONNXJetTagsProducer.cc#L86
0018     # and https://twiki.cern.ch/twiki/bin/view/CMS/DeepJet
0019     from RecoBTag.ONNXRuntime.pfParticleTransformerAK4_cff import _pfParticleTransformerAK4JetTagsAll as pfParticleTransformerAK4JetTagsAll
0020     from RecoBTag.ONNXRuntime.pfParticleNetFromMiniAODAK4_cff import _pfParticleNetFromMiniAODAK4PuppiCentralJetTagsAll as pfParticleNetFromMiniAODAK4PuppiCentralJetTagsAll
0021     from RecoBTag.ONNXRuntime.pfParticleTransformerAK4_cff import _pfNegativeParticleTransformerAK4JetTagsProbs as pfNegativeParticleTransformerAK4JetTagsProbs
0022     from RecoBTag.ONNXRuntime.pfParticleNetFromMiniAODAK4_cff import _pfNegativeParticleNetFromMiniAODAK4PuppiCentralJetTagsProbs as pfNegativeParticleNetFromMiniAODAK4PuppiCentralJetTagsProbs
0023 
0024     _btagDiscriminators = [
0025         'pfJetProbabilityBJetTags',
0026         'pfJetBProbabilityBJetTags',
0027         'pfNegativeOnlyJetProbabilityBJetTags',
0028         'pfNegativeOnlyJetBProbabilityBJetTags',
0029         'pfDeepFlavourJetTags:probb',
0030         'pfDeepFlavourJetTags:probbb',
0031         'pfDeepFlavourJetTags:problepb',
0032         'pfDeepFlavourJetTags:probc',
0033         'pfDeepFlavourJetTags:probuds',
0034         'pfDeepFlavourJetTags:probg',
0035         'pfNegativeDeepFlavourJetTags:probb',
0036         'pfNegativeDeepFlavourJetTags:probbb',
0037         'pfNegativeDeepFlavourJetTags:problepb',
0038         'pfNegativeDeepFlavourJetTags:probc',
0039         'pfNegativeDeepFlavourJetTags:probuds',
0040         'pfNegativeDeepFlavourJetTags:probg',
0041     ]        + pfParticleNetFromMiniAODAK4PuppiCentralJetTagsAll + pfNegativeParticleNetFromMiniAODAK4PuppiCentralJetTagsProbs
0042     # \ #+ pfParticleTransformerAK4JetTagsAll + pfNegativeParticleTransformerAK4JetTagsProbs \
0043     updateJetCollection(
0044         process,
0045         jetSource=cms.InputTag('slimmedJetsPuppi'),
0046         jetCorrections=('AK4PFPuppi',
0047                         cms.vstring(
0048                             ['L1FastJet', 'L2Relative', 'L3Absolute',
0049                              'L2L3Residual']), 'None'),
0050         btagDiscriminators=_btagDiscriminators,
0051         postfix='PuppiWithDeepInfo',
0052     )
0053     process.load("Configuration.StandardSequences.MagneticField_cff")
0054     process.jetPuppiCorrFactorsNano.src = "selectedUpdatedPatJetsPuppiWithDeepInfo"
0055     process.updatedJetsPuppi.jetSource = "selectedUpdatedPatJetsPuppiWithDeepInfo"
0056     
0057     
0058     process.updatedPatJetsTransientCorrectedPuppiWithDeepInfo.tagInfoSources.append(cms.InputTag("pfDeepFlavourTagInfosPuppiWithDeepInfo"))
0059     process.updatedPatJetsTransientCorrectedPuppiWithDeepInfo.addTagInfos = cms.bool(True)
0060 
0061     
0062     
0063     return process
0064 
0065 def update_jets_AK8(process):
0066     # Based on ``nanoAOD_addDeepInfoAK8``
0067     # in https://github.com/cms-sw/cmssw/blob/master/PhysicsTools/NanoAOD/python/nano_cff.py
0068     # Care needs to be taken to make sure no discriminators from stock Nano are excluded -> would results in unfilled vars
0069     _btagDiscriminators = [
0070         'pfMassIndependentDeepDoubleBvLV2JetTags:probHbb',
0071         'pfMassIndependentDeepDoubleCvLV2JetTags:probHcc',
0072         'pfMassIndependentDeepDoubleCvBV2JetTags:probHcc',
0073         ]
0074     from RecoBTag.ONNXRuntime.pfParticleNet_cff import _pfParticleNetJetTagsAll as pfParticleNetJetTagsAll
0075     _btagDiscriminators += pfParticleNetJetTagsAll
0076     updateJetCollection(
0077         process,
0078         jetSource=cms.InputTag('slimmedJetsAK8'),
0079         pvSource=cms.InputTag('offlineSlimmedPrimaryVertices'),
0080         svSource=cms.InputTag('slimmedSecondaryVertices'),
0081         rParam=0.8,
0082         jetCorrections=('AK8PFPuppi',
0083                         cms.vstring([
0084                             'L1FastJet', 'L2Relative', 'L3Absolute',
0085                             'L2L3Residual'
0086                         ]), 'None'),
0087         btagDiscriminators=_btagDiscriminators,
0088         postfix='AK8WithDeepInfo',
0089         # this should work but doesn't seem to enable the tag info with addTagInfos
0090         # btagInfos=['pfDeepDoubleXTagInfos'],
0091         printWarning=False)
0092     process.jetCorrFactorsAK8.src = "selectedUpdatedPatJetsAK8WithDeepInfo"
0093     process.updatedJetsAK8.jetSource = "selectedUpdatedPatJetsAK8WithDeepInfo"
0094     # add DeepDoubleX taginfos
0095     process.updatedPatJetsTransientCorrectedAK8WithDeepInfo.tagInfoSources.append(cms.InputTag("pfDeepDoubleXTagInfosAK8WithDeepInfo"))
0096     process.updatedPatJetsTransientCorrectedAK8WithDeepInfo.addTagInfos = cms.bool(True)
0097     return process
0098 
0099 def update_jets_AK8_subjet(process):
0100     # Based on ``nanoAOD_addDeepInfoAK8``
0101     # in https://github.com/cms-sw/cmssw/blob/master/PhysicsTools/NanoAOD/python/nano_cff.py
0102     # and https://github.com/alefisico/RecoBTag-PerformanceMeasurements/blob/10_2_X_boostedCommissioning/test/runBTagAnalyzer_cfg.py
0103     _btagDiscriminators = [
0104         'pfJetProbabilityBJetTags',
0105         'pfDeepCSVJetTags:probb',
0106         'pfDeepCSVJetTags:probc',
0107         'pfDeepCSVJetTags:probbb',
0108         'pfDeepCSVJetTags:probudsg',
0109         ]
0110     updateJetCollection(
0111         process,
0112         labelName='SoftDropSubjetsPF',
0113         jetSource=cms.InputTag("slimmedJetsAK8PFPuppiSoftDropPacked", "SubJets"),
0114         jetCorrections=('AK4PFPuppi',
0115                         ['L2Relative', 'L3Absolute'], 'None'),
0116         btagDiscriminators=list(_btagDiscriminators),
0117         explicitJTA=True,  # needed for subjet b tagging
0118         svClustering=False,  # needed for subjet b tagging (IMPORTANT: Needs to be set to False to disable ghost-association which does not work with slimmed jets)
0119         fatJets=cms.InputTag('slimmedJetsAK8'),  # needed for subjet b tagging
0120         rParam=0.8,  # needed for subjet b tagging
0121         sortByPt=False, # Don't change order (would mess with subJetIdx for FatJets)
0122         postfix='AK8SubjetsWithDeepInfo')
0123 
0124     process.subJetTable.src = 'selectedUpdatedPatJetsSoftDropSubjetsPFAK8SubjetsWithDeepInfo' 
0125     
0126 
0127     return process
0128 
0129 def get_DDX_vars():
0130     # retreive 27 jet-level features used in double-b and deep double-x taggers
0131     # defined in arXiv:1712.07158
0132     DDXVars = cms.PSet(
0133         DDX_jetNTracks = Var("tagInfo(\'pfDeepDoubleX\').features().tag_info_features.jetNTracks", int, doc="number of tracks associated with the jet"),
0134         DDX_jetNSecondaryVertices = Var("tagInfo(\'pfDeepDoubleX\').features().tag_info_features.jetNSecondaryVertices", int, doc="number of SVs associated with the jet"),
0135         DDX_tau1_trackEtaRel_0 = Var("tagInfo(\'pfDeepDoubleX\').features().tag_info_features.tau1_trackEtaRel_0", float, doc="1st smallest track pseudorapidity, relative to the jet axis, associated to the 1st N-subjettiness axis", precision=10),
0136         DDX_tau1_trackEtaRel_1 = Var("tagInfo(\'pfDeepDoubleX\').features().tag_info_features.tau1_trackEtaRel_1", float, doc="2nd smallest track pseudorapidity, relative to the jet axis, associated to the 1st N-subjettiness axis", precision=10),
0137         DDX_tau1_trackEtaRel_2 = Var("tagInfo(\'pfDeepDoubleX\').features().tag_info_features.tau1_trackEtaRel_2", float, doc="3rd smallest track pseudorapidity, relative to the jet axis, associated to the 1st N-subjettiness axis", precision=10),
0138         DDX_tau2_trackEtaRel_0 = Var("tagInfo(\'pfDeepDoubleX\').features().tag_info_features.tau2_trackEtaRel_0", float, doc="1st smallest track pseudorapidity, relative to the jet axis, associated to the 2nd N-subjettiness axis", precision=10),
0139         DDX_tau2_trackEtaRel_1 = Var("tagInfo(\'pfDeepDoubleX\').features().tag_info_features.tau2_trackEtaRel_1", float, doc="2nd smallest track pseudorapidity, relative to the jet axis, associated to the 2nd N-subjettiness axis", precision=10),
0140         DDX_tau2_trackEtaRel_3 = Var("tagInfo(\'pfDeepDoubleX\').features().tag_info_features.tau2_trackEtaRel_2", float, doc="3rd smallest track pseudorapidity, relative to the jet axis, associated to the 2nd N-subjettiness axis", precision=10),
0141         DDX_tau1_flightDistance2dSig = Var("tagInfo(\'pfDeepDoubleX\').features().tag_info_features.tau1_flightDistance2dSig", float, doc="transverse distance significance between primary and secondary vertex associated to the 1st N-subjettiness axis", precision=10),
0142         DDX_tau2_flightDistance2dSig = Var("tagInfo(\'pfDeepDoubleX\').features().tag_info_features.tau2_flightDistance2dSig", float, doc="transverse distance significance between primary and secondary vertex associated to the 2nd N-subjettiness axis", precision=10),
0143         DDX_tau1_vertexDeltaR = Var("tagInfo(\'pfDeepDoubleX\').features().tag_info_features.tau1_vertexDeltaR", float, doc="deltaR between the 1st N-subjettiness axis and secondary vertex direction", precision=10),
0144         DDX_tau1_vertexEnergyRatio = Var("tagInfo(\'pfDeepDoubleX\').features().tag_info_features.tau1_vertexEnergyRatio", float, doc="ratio of energy at secondary vertex over total energy associated to the 1st N-subjettiness axis", precision=10),
0145         DDX_tau2_vertexEnergyRatio = Var("tagInfo(\'pfDeepDoubleX\').features().tag_info_features.tau2_vertexEnergyRatio", float, doc="ratio of energy at secondary vertex over total energy associated to the 2nd N-subjettiness axis", precision=10),
0146         DDX_tau1_vertexMass = Var("tagInfo(\'pfDeepDoubleX\').features().tag_info_features.tau1_vertexMass", float, doc="mass of track sum at secondary vertex associated to the 1st N-subjettiness axis", precision=10),
0147         DDX_tau2_vertexMass = Var("tagInfo(\'pfDeepDoubleX\').features().tag_info_features.tau2_vertexMass", float, doc="mass of track sum at secondary vertex associated to the 2nd N-subjettiness axis", precision=10),
0148         DDX_trackSip2dSigAboveBottom_0 = Var("tagInfo(\'pfDeepDoubleX\').features().tag_info_features.trackSip2dSigAboveBottom_0", float, doc="track 2D signed impact parameter significance of 1st track lifting mass above bottom", precision=10),
0149         DDX_trackSip2dSigAboveBottom_1 = Var("tagInfo(\'pfDeepDoubleX\').features().tag_info_features.trackSip2dSigAboveBottom_1", float, doc="track 2D signed impact parameter significance of 2nd track lifting mass above bottom", precision=10),
0150         DDX_trackSip2dSigAboveCharm = Var("tagInfo(\'pfDeepDoubleX\').features().tag_info_features.trackSip2dSigAboveCharm", float, doc="track 2D signed impact parameter significance of 1st track lifting mass above charm", precision=10),
0151         DDX_trackSip3dSig_0 = Var("tagInfo(\'pfDeepDoubleX\').features().tag_info_features.trackSip3dSig_0", float, doc="1st largest track 3D signed impact parameter significance", precision=10),
0152         DDX_tau1_trackSip3dSig_0 = Var("tagInfo(\'pfDeepDoubleX\').features().tag_info_features.tau1_trackSip3dSig_0", float, doc="1st largest track 3D signed impact parameter significance associated to the 1st N-subjettiness axis", precision=10),
0153         DDX_tau1_trackSip3dSig_1 = Var("tagInfo(\'pfDeepDoubleX\').features().tag_info_features.tau1_trackSip3dSig_1", float, doc="2nd largest track 3D signed impact parameter significance associated to the 1st N-subjettiness axis", precision=10),
0154         DDX_trackSip3dSig_1 = Var("tagInfo(\'pfDeepDoubleX\').features().tag_info_features.trackSip3dSig_1", float, doc="2nd largest track 3D signed impact parameter significance", precision=10),
0155         DDX_tau2_trackSip3dSig_0 = Var("tagInfo(\'pfDeepDoubleX\').features().tag_info_features.tau2_trackSip3dSig_0", float, doc="1st largest track 3D signed impact parameter significance associated to the 2nd N-subjettiness axis", precision=10),
0156         DDX_tau2_trackSip3dSig_1 = Var("tagInfo(\'pfDeepDoubleX\').features().tag_info_features.tau2_trackSip3dSig_1", float, doc="2nd largest track 3D signed impact parameter significance associated to the 2nd N-subjettiness axis", precision=10),
0157         DDX_trackSip3dSig_2 = Var("tagInfo(\'pfDeepDoubleX\').features().tag_info_features.trackSip3dSig_2", float, doc="3rd largest track 3D signed impact parameter significance", precision=10),
0158         DDX_trackSip3dSig_3 = Var("tagInfo(\'pfDeepDoubleX\').features().tag_info_features.trackSip3dSig_3", float, doc="4th largest track 3D signed impact parameter significance", precision=10),
0159         DDX_z_ratio = Var("tagInfo(\'pfDeepDoubleX\').features().tag_info_features.z_ratio", float, doc="z = deltaR(SV0,SV1)*pT(SV1)/m(SV0,SV1), defined in Eq. 7 of arXiv:1712.07158", precision=10)
0160     )
0161     return DDXVars
0162 
0163 def get_DeepCSV_vars():
0164     DeepCSVVars = cms.PSet(
0165         # Tagger inputs also include jet pt and eta
0166         # Track based (keep only jet-based features for DeepCSV from Run 3 commissioning)
0167         # DeepCSV_trackPtRel_0 = Var("?tagInfo(\'pfDeepCSV\').taggingVariables.checkTag(\'trackPtRel\')?tagInfo(\'pfDeepCSV\').taggingVariables.getList(\'trackPtRel\')[0]:-999", float, doc="track transverse momentum, relative to the jet axis", precision=10),
0168         # DeepCSV_trackPtRel_1 = Var("?tagInfo(\'pfDeepCSV\').taggingVariables.checkTag(\'trackPtRel\')?tagInfo(\'pfDeepCSV\').taggingVariables.getList(\'trackPtRel\')[1]:-999", float, doc="track transverse momentum, relative to the jet axis", precision=10),
0169         # DeepCSV_trackPtRel_2 = Var("?tagInfo(\'pfDeepCSV\').taggingVariables.checkTag(\'trackPtRel\')?tagInfo(\'pfDeepCSV\').taggingVariables.getList(\'trackPtRel\')[2]:-999", float, doc="track transverse momentum, relative to the jet axis", precision=10),
0170         # DeepCSV_trackPtRel_3 = Var("?tagInfo(\'pfDeepCSV\').taggingVariables.checkTag(\'trackPtRel\')?tagInfo(\'pfDeepCSV\').taggingVariables.getList(\'trackPtRel\')[3]:-999", float, doc="track transverse momentum, relative to the jet axis", precision=10),
0171         # DeepCSV_trackPtRel_4 = Var("?tagInfo(\'pfDeepCSV\').taggingVariables.checkTag(\'trackPtRel\')?tagInfo(\'pfDeepCSV\').taggingVariables.getList(\'trackPtRel\')[4]:-999", float, doc="track transverse momentum, relative to the jet axis", precision=10),
0172         # DeepCSV_trackPtRel_5 = Var("?tagInfo(\'pfDeepCSV\').taggingVariables.checkTag(\'trackPtRel\')?tagInfo(\'pfDeepCSV\').taggingVariables.getList(\'trackPtRel\')[5]:-999", float, doc="track transverse momentum, relative to the jet axis", precision=10),
0173         # DeepCSV_trackJetDistVal_0 = Var("?tagInfo(\'pfDeepCSV\').taggingVariables.checkTag(\'trackJetDistVal\')?tagInfo(\'pfDeepCSV\').taggingVariables.getList(\'trackJetDistVal\')[0]:-999", float, doc="minimum track approach distance to jet axis", precision=10),
0174         # DeepCSV_trackJetDistVal_1 = Var("?tagInfo(\'pfDeepCSV\').taggingVariables.checkTag(\'trackJetDistVal\')?tagInfo(\'pfDeepCSV\').taggingVariables.getList(\'trackJetDistVal\')[1]:-999", float, doc="minimum track approach distance to jet axis", precision=10),
0175         # DeepCSV_trackJetDistVal_2 = Var("?tagInfo(\'pfDeepCSV\').taggingVariables.checkTag(\'trackJetDistVal\')?tagInfo(\'pfDeepCSV\').taggingVariables.getList(\'trackJetDistVal\')[2]:-999", float, doc="minimum track approach distance to jet axis", precision=10),
0176         # DeepCSV_trackJetDistVal_3 = Var("?tagInfo(\'pfDeepCSV\').taggingVariables.checkTag(\'trackJetDistVal\')?tagInfo(\'pfDeepCSV\').taggingVariables.getList(\'trackJetDistVal\')[3]:-999", float, doc="minimum track approach distance to jet axis", precision=10),
0177         # DeepCSV_trackJetDistVal_4 = Var("?tagInfo(\'pfDeepCSV\').taggingVariables.checkTag(\'trackJetDistVal\')?tagInfo(\'pfDeepCSV\').taggingVariables.getList(\'trackJetDistVal\')[4]:-999", float, doc="minimum track approach distance to jet axis", precision=10),
0178         # DeepCSV_trackJetDistVal_5 = Var("?tagInfo(\'pfDeepCSV\').taggingVariables.checkTag(\'trackJetDistVal\')?tagInfo(\'pfDeepCSV\').taggingVariables.getList(\'trackJetDistVal\')[5]:-999", float, doc="minimum track approach distance to jet axis", precision=10),
0179         # DeepCSV_trackDeltaR_0 = Var("?tagInfo(\'pfDeepCSV\').taggingVariables.checkTag(\'trackDeltaR\')?tagInfo(\'pfDeepCSV\').taggingVariables.getList(\'trackDeltaR\')[0]:-999", float, doc="track pseudoangular distance from the jet axis", precision=10),
0180         # DeepCSV_trackDeltaR_1 = Var("?tagInfo(\'pfDeepCSV\').taggingVariables.checkTag(\'trackDeltaR\')?tagInfo(\'pfDeepCSV\').taggingVariables.getList(\'trackDeltaR\')[1]:-999", float, doc="track pseudoangular distance from the jet axis", precision=10),
0181         # DeepCSV_trackDeltaR_2 = Var("?tagInfo(\'pfDeepCSV\').taggingVariables.checkTag(\'trackDeltaR\')?tagInfo(\'pfDeepCSV\').taggingVariables.getList(\'trackDeltaR\')[2]:-999", float, doc="track pseudoangular distance from the jet axis", precision=10),
0182         # DeepCSV_trackDeltaR_3 = Var("?tagInfo(\'pfDeepCSV\').taggingVariables.checkTag(\'trackDeltaR\')?tagInfo(\'pfDeepCSV\').taggingVariables.getList(\'trackDeltaR\')[3]:-999", float, doc="track pseudoangular distance from the jet axis", precision=10),
0183         # DeepCSV_trackDeltaR_4 = Var("?tagInfo(\'pfDeepCSV\').taggingVariables.checkTag(\'trackDeltaR\')?tagInfo(\'pfDeepCSV\').taggingVariables.getList(\'trackDeltaR\')[4]:-999", float, doc="track pseudoangular distance from the jet axis", precision=10),
0184         # DeepCSV_trackDeltaR_5 = Var("?tagInfo(\'pfDeepCSV\').taggingVariables.checkTag(\'trackDeltaR\')?tagInfo(\'pfDeepCSV\').taggingVariables.getList(\'trackDeltaR\')[5]:-999", float, doc="track pseudoangular distance from the jet axis", precision=10),
0185         # DeepCSV_trackPtRatio_0 = Var("?tagInfo(\'pfDeepCSV\').taggingVariables.checkTag(\'trackPtRatio\')?tagInfo(\'pfDeepCSV\').taggingVariables.getList(\'trackPtRatio\')[0]:-999", float, doc="track transverse momentum, relative to the jet axis, normalized to its energy", precision=10),
0186         # DeepCSV_trackPtRatio_1 = Var("?tagInfo(\'pfDeepCSV\').taggingVariables.checkTag(\'trackPtRatio\')?tagInfo(\'pfDeepCSV\').taggingVariables.getList(\'trackPtRatio\')[1]:-999", float, doc="track transverse momentum, relative to the jet axis, normalized to its energy", precision=10),
0187         # DeepCSV_trackPtRatio_2 = Var("?tagInfo(\'pfDeepCSV\').taggingVariables.checkTag(\'trackPtRatio\')?tagInfo(\'pfDeepCSV\').taggingVariables.getList(\'trackPtRatio\')[2]:-999", float, doc="track transverse momentum, relative to the jet axis, normalized to its energy", precision=10),
0188         # DeepCSV_trackPtRatio_3 = Var("?tagInfo(\'pfDeepCSV\').taggingVariables.checkTag(\'trackPtRatio\')?tagInfo(\'pfDeepCSV\').taggingVariables.getList(\'trackPtRatio\')[3]:-999", float, doc="track transverse momentum, relative to the jet axis, normalized to its energy", precision=10),
0189         # DeepCSV_trackPtRatio_4 = Var("?tagInfo(\'pfDeepCSV\').taggingVariables.checkTag(\'trackPtRatio\')?tagInfo(\'pfDeepCSV\').taggingVariables.getList(\'trackPtRatio\')[4]:-999", float, doc="track transverse momentum, relative to the jet axis, normalized to its energy", precision=10),
0190         # DeepCSV_trackPtRatio_5 = Var("?tagInfo(\'pfDeepCSV\').taggingVariables.checkTag(\'trackPtRatio\')?tagInfo(\'pfDeepCSV\').taggingVariables.getList(\'trackPtRatio\')[5]:-999", float, doc="track transverse momentum, relative to the jet axis, normalized to its energy", precision=10),
0191         # DeepCSV_trackSip3dSig_0 = Var("?tagInfo(\'pfDeepCSV\').taggingVariables.checkTag(\'trackSip3dSig\')?tagInfo(\'pfDeepCSV\').taggingVariables.getList(\'trackSip3dSig\')[0]:-999", float, doc="track 3D signed impact parameter significance", precision=10),
0192         # DeepCSV_trackSip3dSig_1 = Var("?tagInfo(\'pfDeepCSV\').taggingVariables.checkTag(\'trackSip3dSig\')?tagInfo(\'pfDeepCSV\').taggingVariables.getList(\'trackSip3dSig\')[1]:-999", float, doc="track 3D signed impact parameter significance", precision=10),
0193         # DeepCSV_trackSip3dSig_2 = Var("?tagInfo(\'pfDeepCSV\').taggingVariables.checkTag(\'trackSip3dSig\')?tagInfo(\'pfDeepCSV\').taggingVariables.getList(\'trackSip3dSig\')[2]:-999", float, doc="track 3D signed impact parameter significance", precision=10),
0194         # DeepCSV_trackSip3dSig_3 = Var("?tagInfo(\'pfDeepCSV\').taggingVariables.checkTag(\'trackSip3dSig\')?tagInfo(\'pfDeepCSV\').taggingVariables.getList(\'trackSip3dSig\')[3]:-999", float, doc="track 3D signed impact parameter significance", precision=10),
0195         # DeepCSV_trackSip3dSig_4 = Var("?tagInfo(\'pfDeepCSV\').taggingVariables.checkTag(\'trackSip3dSig\')?tagInfo(\'pfDeepCSV\').taggingVariables.getList(\'trackSip3dSig\')[4]:-999", float, doc="track 3D signed impact parameter significance", precision=10),
0196         # DeepCSV_trackSip3dSig_5 = Var("?tagInfo(\'pfDeepCSV\').taggingVariables.checkTag(\'trackSip3dSig\')?tagInfo(\'pfDeepCSV\').taggingVariables.getList(\'trackSip3dSig\')[5]:-999", float, doc="track 3D signed impact parameter significance", precision=10),
0197         # DeepCSV_trackSip2dSig_0 = Var("?tagInfo(\'pfDeepCSV\').taggingVariables.checkTag(\'trackSip2dSig\')?tagInfo(\'pfDeepCSV\').taggingVariables.getList(\'trackSip2dSig\')[0]:-999", float, doc="track 2D signed impact parameter significance", precision=10),
0198         # DeepCSV_trackSip2dSig_1 = Var("?tagInfo(\'pfDeepCSV\').taggingVariables.checkTag(\'trackSip2dSig\')?tagInfo(\'pfDeepCSV\').taggingVariables.getList(\'trackSip2dSig\')[1]:-999", float, doc="track 2D signed impact parameter significance", precision=10),
0199         # DeepCSV_trackSip2dSig_2 = Var("?tagInfo(\'pfDeepCSV\').taggingVariables.checkTag(\'trackSip2dSig\')?tagInfo(\'pfDeepCSV\').taggingVariables.getList(\'trackSip2dSig\')[2]:-999", float, doc="track 2D signed impact parameter significance", precision=10),
0200         # DeepCSV_trackSip2dSig_3 = Var("?tagInfo(\'pfDeepCSV\').taggingVariables.checkTag(\'trackSip2dSig\')?tagInfo(\'pfDeepCSV\').taggingVariables.getList(\'trackSip2dSig\')[3]:-999", float, doc="track 2D signed impact parameter significance", precision=10),
0201         # DeepCSV_trackSip2dSig_4 = Var("?tagInfo(\'pfDeepCSV\').taggingVariables.checkTag(\'trackSip2dSig\')?tagInfo(\'pfDeepCSV\').taggingVariables.getList(\'trackSip2dSig\')[4]:-999", float, doc="track 2D signed impact parameter significance", precision=10),
0202         # DeepCSV_trackSip2dSig_5 = Var("?tagInfo(\'pfDeepCSV\').taggingVariables.checkTag(\'trackSip2dSig\')?tagInfo(\'pfDeepCSV\').taggingVariables.getList(\'trackSip2dSig\')[5]:-999", float, doc="track 2D signed impact parameter significance", precision=10),
0203         # DeepCSV_trackDecayLenVal_0 = Var("?tagInfo(\'pfDeepCSV\').taggingVariables.checkTag(\'trackDecayLenVal\')?tagInfo(\'pfDeepCSV\').taggingVariables.getList(\'trackDecayLenVal\')[0]:-999", float, doc="track decay length", precision=10),
0204         # DeepCSV_trackDecayLenVal_1 = Var("?tagInfo(\'pfDeepCSV\').taggingVariables.checkTag(\'trackDecayLenVal\')?tagInfo(\'pfDeepCSV\').taggingVariables.getList(\'trackDecayLenVal\')[1]:-999", float, doc="track decay length", precision=10),
0205         # DeepCSV_trackDecayLenVal_2 = Var("?tagInfo(\'pfDeepCSV\').taggingVariables.checkTag(\'trackDecayLenVal\')?tagInfo(\'pfDeepCSV\').taggingVariables.getList(\'trackDecayLenVal\')[2]:-999", float, doc="track decay length", precision=10),
0206         # DeepCSV_trackDecayLenVal_3 = Var("?tagInfo(\'pfDeepCSV\').taggingVariables.checkTag(\'trackDecayLenVal\')?tagInfo(\'pfDeepCSV\').taggingVariables.getList(\'trackDecayLenVal\')[3]:-999", float, doc="track decay length", precision=10),
0207         # DeepCSV_trackDecayLenVal_4 = Var("?tagInfo(\'pfDeepCSV\').taggingVariables.checkTag(\'trackDecayLenVal\')?tagInfo(\'pfDeepCSV\').taggingVariables.getList(\'trackDecayLenVal\')[4]:-999", float, doc="track decay length", precision=10),
0208         # DeepCSV_trackDecayLenVal_5 = Var("?tagInfo(\'pfDeepCSV\').taggingVariables.checkTag(\'trackDecayLenVal\')?tagInfo(\'pfDeepCSV\').taggingVariables.getList(\'trackDecayLenVal\')[5]:-999", float, doc="track decay length", precision=10),
0209         # DeepCSV_trackEtaRel_0 = Var("?tagInfo(\'pfDeepCSV\').taggingVariables.checkTag(\'trackEtaRel\')?tagInfo(\'pfDeepCSV\').taggingVariables.getList(\'trackEtaRel\')[0]:-999", float, doc="track pseudorapidity, relative to the jet axis", precision=10),
0210         # DeepCSV_trackEtaRel_1 = Var("?tagInfo(\'pfDeepCSV\').taggingVariables.checkTag(\'trackEtaRel\')?tagInfo(\'pfDeepCSV\').taggingVariables.getList(\'trackEtaRel\')[1]:-999", float, doc="track pseudorapidity, relative to the jet axis", precision=10),
0211         # DeepCSV_trackEtaRel_2 = Var("?tagInfo(\'pfDeepCSV\').taggingVariables.checkTag(\'trackEtaRel\')?tagInfo(\'pfDeepCSV\').taggingVariables.getList(\'trackEtaRel\')[2]:-999", float, doc="track pseudorapidity, relative to the jet axis", precision=10),
0212         # DeepCSV_trackEtaRel_3 = Var("?tagInfo(\'pfDeepCSV\').taggingVariables.checkTag(\'trackEtaRel\')?tagInfo(\'pfDeepCSV\').taggingVariables.getList(\'trackEtaRel\')[3]:-999", float, doc="track pseudorapidity, relative to the jet axis", precision=10),
0213         # Jet based
0214         DeepCSV_trackJetPt = Var("tagInfo(\'pfDeepCSV\').taggingVariables.get(\'trackJetPt\', -999)", float, doc="track-based jet transverse momentum", precision=10),
0215         DeepCSV_vertexCategory = Var("tagInfo(\'pfDeepCSV\').taggingVariables.get(\'vertexCategory\', -999)", float, doc="category of secondary vertex (Reco, Pseudo, No)", precision=10),
0216         DeepCSV_jetNSecondaryVertices = Var("tagInfo(\'pfDeepCSV\').taggingVariables.get(\'jetNSecondaryVertices\', -999)", int, doc="number of reconstructed possible secondary vertices in jet"),
0217         DeepCSV_jetNSelectedTracks = Var("tagInfo(\'pfDeepCSV\').taggingVariables.get(\'jetNSelectedTracks\', -999)", int, doc="selected tracks in the jet"), 
0218         DeepCSV_jetNTracksEtaRel = Var("tagInfo(\'pfDeepCSV\').taggingVariables.get(\'jetNTracksEtaRel\', -999)", int, doc="number of tracks for which etaRel is computed"), 
0219         DeepCSV_trackSumJetEtRatio = Var("tagInfo(\'pfDeepCSV\').taggingVariables.get(\'trackSumJetEtRatio\', -999)", float, doc="ratio of track sum transverse energy over jet energy", precision=10),
0220         DeepCSV_trackSumJetDeltaR = Var("tagInfo(\'pfDeepCSV\').taggingVariables.get(\'trackSumJetDeltaR\', -999)", float, doc="pseudoangular distance between jet axis and track fourvector sum", precision=10),
0221         DeepCSV_trackSip2dValAboveCharm = Var("tagInfo(\'pfDeepCSV\').taggingVariables.get(\'trackSip2dValAboveCharm\', -999)", float, doc="track 2D signed impact parameter of first track lifting mass above charm", precision=10),
0222         DeepCSV_trackSip2dSigAboveCharm = Var("tagInfo(\'pfDeepCSV\').taggingVariables.get(\'trackSip2dSigAboveCharm\', -999)", float, doc="track 2D signed impact parameter significance of first track lifting mass above charm", precision=10),
0223         DeepCSV_trackSip3dValAboveCharm = Var("tagInfo(\'pfDeepCSV\').taggingVariables.get(\'trackSip3dValAboveCharm\', -999)", float, doc="track 3D signed impact parameter of first track lifting mass above charm", precision=10),
0224         DeepCSV_trackSip3dSigAboveCharm = Var("tagInfo(\'pfDeepCSV\').taggingVariables.get(\'trackSip3dSigAboveCharm\', -999)", float, doc="track 3D signed impact parameter significance of first track lifting mass above charm", precision=10),
0225         DeepCSV_vertexMass = Var("tagInfo(\'pfDeepCSV\').taggingVariables.get(\'vertexMass\', -999)", float, doc="mass of track sum at secondary vertex", precision=10),
0226         DeepCSV_vertexNTracks = Var("tagInfo(\'pfDeepCSV\').taggingVariables.get(\'vertexNTracks\', -999)", int, doc="number of tracks at secondary vertex"),
0227         DeepCSV_vertexEnergyRatio = Var("tagInfo(\'pfDeepCSV\').taggingVariables.get(\'vertexEnergyRatio\', -999)", float, doc="ratio of energy at secondary vertex over total energy", precision=10),
0228         DeepCSV_vertexJetDeltaR = Var("tagInfo(\'pfDeepCSV\').taggingVariables.get(\'vertexJetDeltaR\', -999)", float, doc="pseudoangular distance between jet axis and secondary vertex direction", precision=10),
0229         DeepCSV_flightDistance2dVal = Var("tagInfo(\'pfDeepCSV\').taggingVariables.get(\'flightDistance2dVal\', -999)", float, doc="transverse distance between primary and secondary vertex", precision=10),
0230         DeepCSV_flightDistance2dSig = Var("tagInfo(\'pfDeepCSV\').taggingVariables.get(\'flightDistance2dSig\', -999)", float, doc="transverse distance significance between primary and secondary vertex", precision=10),
0231         DeepCSV_flightDistance3dVal = Var("tagInfo(\'pfDeepCSV\').taggingVariables.get(\'flightDistance3dVal\', -999)", float, doc="distance between primary and secondary vertex", precision=10),
0232         DeepCSV_flightDistance3dSig = Var("tagInfo(\'pfDeepCSV\').taggingVariables.get(\'flightDistance3dSig\', -999)", float, doc="distance significance between primary and secondary vertex", precision=10),
0233     )
0234     return DeepCSVVars
0235 
0236 ## Store all output nodes, negative tagger for SF 
0237 def get_DeepJet_outputs():
0238     DeepJetOutputVars = cms.PSet(
0239         btagDeepFlavB_b=Var("bDiscriminator('pfDeepFlavourJetTags:probb')",
0240                             float,
0241                             doc="DeepJet b tag probability",
0242                             precision=10),
0243         btagDeepFlavB_bb=Var("bDiscriminator('pfDeepFlavourJetTags:probbb')",
0244                             float,
0245                             doc="DeepJet bb tag probability",
0246                             precision=10),
0247         btagDeepFlavB_lepb=Var("bDiscriminator('pfDeepFlavourJetTags:problepb')",
0248                             float,
0249                             doc="DeepJet lepb tag probability",
0250                             precision=10),
0251         btagDeepFlavC=Var("bDiscriminator('pfDeepFlavourJetTags:probc')",
0252                             float,
0253                             doc="DeepJet c tag probability",
0254                             precision=10),
0255         btagDeepFlavUDS=Var("bDiscriminator('pfDeepFlavourJetTags:probuds')",
0256                             float,
0257                             doc="DeepJet uds tag probability",
0258                             precision=10),
0259         btagDeepFlavG=Var("bDiscriminator('pfDeepFlavourJetTags:probg')",
0260                             float,
0261                             doc="DeepJet gluon tag probability",
0262                             precision=10),
0263         # discriminators are already part of jets_cff.py from NanoAOD and therefore not added here     
0264 
0265         # negative taggers
0266         btagNegDeepFlavB = Var("bDiscriminator('pfNegativeDeepFlavourJetTags:probb')+bDiscriminator('pfNegativeDeepFlavourJetTags:probbb')+bDiscriminator('pfNegativeDeepFlavourJetTags:problepb')",
0267                             float,
0268                             doc="Negative DeepJet b+bb+lepb tag discriminator",
0269                             precision=10),
0270         btagNegDeepFlavCvL = Var("?(bDiscriminator('pfNegativeDeepFlavourJetTags:probc')+bDiscriminator('pfNegativeDeepFlavourJetTags:probuds')+bDiscriminator('pfNegativeDeepFlavourJetTags:probg'))>0?bDiscriminator('pfNegativeDeepFlavourJetTags:probc')/(bDiscriminator('pfNegativeDeepFlavourJetTags:probc')+bDiscriminator('pfNegativeDeepFlavourJetTags:probuds')+bDiscriminator('pfNegativeDeepFlavourJetTags:probg')):-1",
0271                             float,
0272                             doc="Negative DeepJet c vs uds+g discriminator",
0273                             precision=10),
0274         btagNegDeepFlavCvB = Var("?(bDiscriminator('pfNegativeDeepFlavourJetTags:probc')+bDiscriminator('pfNegativeDeepFlavourJetTags:probb')+bDiscriminator('pfNegativeDeepFlavourJetTags:probbb')+bDiscriminator('pfNegativeDeepFlavourJetTags:problepb'))>0?bDiscriminator('pfNegativeDeepFlavourJetTags:probc')/(bDiscriminator('pfNegativeDeepFlavourJetTags:probc')+bDiscriminator('pfNegativeDeepFlavourJetTags:probb')+bDiscriminator('pfNegativeDeepFlavourJetTags:probbb')+bDiscriminator('pfNegativeDeepFlavourJetTags:problepb')):-1",
0275                             float,
0276                             doc="Negative DeepJet c vs b+bb+lepb discriminator",
0277                             precision=10),
0278         btagNegDeepFlavQG = Var("?(bDiscriminator('pfNegativeDeepFlavourJetTags:probg')+bDiscriminator('pfNegativeDeepFlavourJetTags:probuds'))>0?bDiscriminator('pfNegativeDeepFlavourJetTags:probg')/(bDiscriminator('pfNegativeDeepFlavourJetTags:probg')+bDiscriminator('pfNegativeDeepFlavourJetTags:probuds')):-1",
0279                             float,
0280                             doc="Negative DeepJet g vs uds discriminator",
0281                             precision=10),
0282         btagNegDeepFlavB_b = Var("bDiscriminator('pfNegativeDeepFlavourJetTags:probb')",
0283                             float,
0284                             doc="Negative DeepJet b tag probability",
0285                             precision=10),
0286         btagNegDeepFlavB_bb = Var("bDiscriminator('pfNegativeDeepFlavourJetTags:probbb')",
0287                             float,
0288                             doc="Negative DeepJet bb tag probability",
0289                             precision=10),
0290         btagNegDeepFlavB_lepb = Var("bDiscriminator('pfNegativeDeepFlavourJetTags:problepb')",
0291                             float,
0292                             doc="Negative DeepJet lepb tag probability",
0293                             precision=10),
0294         btagNegDeepFlavC = Var("bDiscriminator('pfNegativeDeepFlavourJetTags:probc')",
0295                             float,
0296                             doc="Negative DeepJet c tag probability",
0297                             precision=10),
0298         btagNegDeepFlavUDS = Var("bDiscriminator('pfNegativeDeepFlavourJetTags:probuds')",
0299                             float,
0300                             doc="Negative DeepJet uds tag probability",
0301                             precision=10),
0302         btagNegDeepFlavG = Var("bDiscriminator('pfNegativeDeepFlavourJetTags:probg')",
0303                             float,
0304                             doc="Negative DeepJet gluon tag probability",
0305                             precision=10),
0306     )
0307     return DeepJetOutputVars
0308 
0309 def get_ParticleNetAK4_outputs():
0310     ## default scores in jetAK4_Puppi_cff.py collections
0311     ParticleNetAK4OutputVars = cms.PSet(
0312         # raw scores
0313         btagPNetProbB = Var("?bDiscriminator('pfParticleNetFromMiniAODAK4PuppiCentralJetTags:probb')>0?bDiscriminator('pfParticleNetFromMiniAODAK4PuppiCentralJetTags:probb'):-1",
0314                             float,
0315                             doc="ParticleNet b tag probability",
0316                             precision=10),
0317         btagPNetProbC = Var("?bDiscriminator('pfParticleNetFromMiniAODAK4PuppiCentralJetTags:probc')>0?bDiscriminator('pfParticleNetFromMiniAODAK4PuppiCentralJetTags:probc'):-1",
0318                             float,
0319                             doc="ParticleNet c tag probability",
0320                             precision=10),
0321         btagPNetProbUDS = Var("?bDiscriminator('pfParticleNetFromMiniAODAK4PuppiCentralJetTags:probuds')>0?bDiscriminator('pfParticleNetFromMiniAODAK4PuppiCentralJetTags:probuds'):-1",
0322                             float,
0323                             doc="ParticleNet uds tag probability",
0324                             precision=10),
0325         btagPNetProbG = Var("?bDiscriminator('pfParticleNetFromMiniAODAK4PuppiCentralJetTags:probg')>0?bDiscriminator('pfParticleNetFromMiniAODAK4PuppiCentralJetTags:probg'):-1",
0326                             float,
0327                             doc="ParticleNet gluon tag probability",
0328                             precision=10),
0329         
0330         # negative taggers
0331         btagNegPNetB = Var("?(bDiscriminator('pfNegativeParticleNetFromMiniAODAK4PuppiCentralJetTags:probb')+bDiscriminator('pfNegativeParticleNetFromMiniAODAK4PuppiCentralJetTags:probc')+bDiscriminator('pfNegativeParticleNetFromMiniAODAK4PuppiCentralJetTags:probuds')+bDiscriminator('pfNegativeParticleNetFromMiniAODAK4PuppiCentralJetTags:probg'))>0?(bDiscriminator('pfNegativeParticleNetFromMiniAODAK4PuppiCentralJetTags:probb'))/(bDiscriminator('pfNegativeParticleNetFromMiniAODAK4PuppiCentralJetTags:probb')+bDiscriminator('pfNegativeParticleNetFromMiniAODAK4PuppiCentralJetTags:probc')+bDiscriminator('pfNegativeParticleNetFromMiniAODAK4PuppiCentralJetTags:probuds')+bDiscriminator('pfNegativeParticleNetFromMiniAODAK4PuppiCentralJetTags:probg')):-1",
0332                             float,
0333                             doc="Negative ParticleNet b vs. udscg",
0334                             precision=10),
0335         btagNegPNetCvL = Var("?(bDiscriminator('pfNegativeParticleNetFromMiniAODAK4PuppiCentralJetTags:probc')+bDiscriminator('pfNegativeParticleNetFromMiniAODAK4PuppiCentralJetTags:probuds')+bDiscriminator('pfNegativeParticleNetFromMiniAODAK4PuppiCentralJetTags:probg'))>0?(bDiscriminator('pfNegativeParticleNetFromMiniAODAK4PuppiCentralJetTags:probc'))/(bDiscriminator('pfNegativeParticleNetFromMiniAODAK4PuppiCentralJetTags:probc')+bDiscriminator('pfNegativeParticleNetFromMiniAODAK4PuppiCentralJetTags:probuds')+bDiscriminator('pfNegativeParticleNetFromMiniAODAK4PuppiCentralJetTags:probg')):-1",
0336                             float,
0337                             doc="Negative ParticleNet c vs. udsg",
0338                             precision=10),
0339         btagNegPNetCvB = Var("?(bDiscriminator('pfNegativeParticleNetFromMiniAODAK4PuppiCentralJetTags:probc')+bDiscriminator('pfNegativeParticleNetFromMiniAODAK4PuppiCentralJetTags:probb'))>0?(bDiscriminator('pfNegativeParticleNetFromMiniAODAK4PuppiCentralJetTags:probc'))/(bDiscriminator('pfNegativeParticleNetFromMiniAODAK4PuppiCentralJetTags:probc')+bDiscriminator('pfNegativeParticleNetFromMiniAODAK4PuppiCentralJetTags:probb')):-1",
0340                             float,
0341                             doc="Negative ParticleNet c vs. b",
0342                             precision=10),
0343         btagNegPNetProbB = Var("?bDiscriminator('pfNegativeParticleNetFromMiniAODAK4PuppiCentralJetTags:probb')>0?bDiscriminator('pfNegativeParticleNetFromMiniAODAK4PuppiCentralJetTags:probb'):-1",
0344                             float,
0345                             doc="Negative ParticleNet b tag probability",
0346                             precision=10),
0347         btagNegPNetProbC = Var("?bDiscriminator('pfNegativeParticleNetFromMiniAODAK4PuppiCentralJetTags:probc')>0?bDiscriminator('pfNegativeParticleNetFromMiniAODAK4PuppiCentralJetTags:probc'):-1",
0348                             float,
0349                             doc="Negative ParticleNet c tag probability",
0350                             precision=10),
0351         btagNegPNetProbUDS = Var("?bDiscriminator('pfNegativeParticleNetFromMiniAODAK4PuppiCentralJetTags:probuds')>0?bDiscriminator('pfNegativeParticleNetFromMiniAODAK4PuppiCentralJetTags:probuds'):-1",
0352                             float,
0353                             doc="Negative ParticleNet uds tag probability",
0354                             precision=10),
0355         btagNegPNetProbG = Var("?bDiscriminator('pfNegativeParticleNetFromMiniAODAK4PuppiCentralJetTags:probg')>0?bDiscriminator('pfNegativeParticleNetFromMiniAODAK4PuppiCentralJetTags:probg'):-1",
0356                             float,
0357                             doc="Negative ParticleNet gluon tag probability",
0358                             precision=10),
0359     )
0360 
0361     return ParticleNetAK4OutputVars
0362 
0363 def get_ParticleTransformerAK4_outputs():
0364     ParticleTransformerAK4OutputVars = cms.PSet(       
0365         btagRobustParTAK4B_b=Var("bDiscriminator('pfParticleTransformerAK4JetTags:probb')",
0366                             float,
0367                             doc="RobustParTAK4 b tag probability",
0368                             precision=10),
0369         btagRobustParTAK4B_bb=Var("bDiscriminator('pfParticleTransformerAK4JetTags:probbb')",
0370                             float,
0371                             doc="RobustParTAK4 bb tag probability",
0372                             precision=10),
0373         btagRobustParTAK4B_lepb=Var("bDiscriminator('pfParticleTransformerAK4JetTags:problepb')",
0374                             float,
0375                             doc="RobustParTAK4 lepb tag probability",
0376                             precision=10),
0377         btagRobustParTAK4C=Var("bDiscriminator('pfParticleTransformerAK4JetTags:probc')",
0378                             float,
0379                             doc="RobustParTAK4 c tag probability",
0380                             precision=10),
0381         btagRobustParTAK4UDS=Var("bDiscriminator('pfParticleTransformerAK4JetTags:probuds')",
0382                             float,
0383                             doc="RobustParTAK4 uds tag probability",
0384                             precision=10),
0385         btagRobustParTAK4G=Var("bDiscriminator('pfParticleTransformerAK4JetTags:probg')",
0386                             float,
0387                             doc="RobustParTAK4 gluon tag probability",
0388                             precision=10),
0389 
0390         # negative taggers
0391         btagNegRobustParTAK4B = Var("bDiscriminator('pfNegativeParticleTransformerAK4JetTags:probb')+bDiscriminator('pfNegativeParticleTransformerAK4JetTags:probbb')+bDiscriminator('pfNegativeParticleTransformerAK4JetTags:problepb')",
0392                             float,
0393                             doc="Negative RobustParTAK4 b+bb+lepb tag discriminator",
0394                             precision=10),
0395         btagNegRobustParTAK4CvL = Var("?(bDiscriminator('pfNegativeParticleTransformerAK4JetTags:probc')+bDiscriminator('pfNegativeParticleTransformerAK4JetTags:probuds')+bDiscriminator('pfNegativeParticleTransformerAK4JetTags:probg'))>0?bDiscriminator('pfNegativeParticleTransformerAK4JetTags:probc')/(bDiscriminator('pfNegativeParticleTransformerAK4JetTags:probc')+bDiscriminator('pfNegativeParticleTransformerAK4JetTags:probuds')+bDiscriminator('pfNegativeParticleTransformerAK4JetTags:probg')):-1",
0396                             float,
0397                             doc="Negative RobustParTAK4 c vs uds+g discriminator",
0398                             precision=10),
0399         btagNegRobustParTAK4CvB = Var("?(bDiscriminator('pfNegativeParticleTransformerAK4JetTags:probc')+bDiscriminator('pfNegativeParticleTransformerAK4JetTags:probb')+bDiscriminator('pfNegativeParticleTransformerAK4JetTags:probbb')+bDiscriminator('pfNegativeParticleTransformerAK4JetTags:problepb'))>0?bDiscriminator('pfNegativeParticleTransformerAK4JetTags:probc')/(bDiscriminator('pfNegativeParticleTransformerAK4JetTags:probc')+bDiscriminator('pfNegativeParticleTransformerAK4JetTags:probb')+bDiscriminator('pfNegativeParticleTransformerAK4JetTags:probbb')+bDiscriminator('pfNegativeParticleTransformerAK4JetTags:problepb')):-1",
0400                             float,
0401                             doc="Negative RobustParTAK4 c vs b+bb+lepb discriminator",
0402                             precision=10),
0403         btagNegRobustParTAK4QG = Var("?(bDiscriminator('pfNegativeParticleTransformerAK4JetTags:probg')+bDiscriminator('pfNegativeParticleTransformerAK4JetTags:probuds'))>0?bDiscriminator('pfNegativeParticleTransformerAK4JetTags:probg')/(bDiscriminator('pfNegativeParticleTransformerAK4JetTags:probg')+bDiscriminator('pfNegativeParticleTransformerAK4JetTags:probuds')):-1",
0404                             float,
0405                             doc="Negative RobustParTAK4 g vs uds discriminator",
0406                             precision=10),
0407         btagNegRobustParTAK4B_b = Var("bDiscriminator('pfNegativeParticleTransformerAK4JetTags:probb')",
0408                             float,
0409                             doc="Negative RobustParTAK4 b tag probability",
0410                             precision=10),
0411         btagNegRobustParTAK4B_bb = Var("bDiscriminator('pfNegativeParticleTransformerAK4JetTags:probbb')",
0412                             float,
0413                             doc="Negative RobustParTAK4 bb tag probability",
0414                             precision=10),
0415         btagNegRobustParTAK4B_lepb = Var("bDiscriminator('pfNegativeParticleTransformerAK4JetTags:problepb')",
0416                             float,
0417                             doc="Negative RobustParTAK4 lepb tag probability",
0418                             precision=10),
0419         btagNegRobustParTAK4C = Var("bDiscriminator('pfNegativeParticleTransformerAK4JetTags:probc')",
0420                             float,
0421                             doc="Negative RobustParTAK4 c tag probability",
0422                             precision=10),
0423         btagNegRobustParTAK4UDS = Var("bDiscriminator('pfNegativeParticleTransformerAK4JetTags:probuds')",
0424                             float,
0425                             doc="Negative RobustParTAK4 uds tag probability",
0426                             precision=10),
0427         btagNegRobustParTAK4G = Var("bDiscriminator('pfNegativeParticleTransformerAK4JetTags:probg')",
0428                             float,
0429                             doc="Negative RobustParTAK4 gluon tag probability",
0430                             precision=10),
0431     )
0432 
0433     return ParticleTransformerAK4OutputVars
0434 
0435 def add_BTV(process,  addAK4=False, addAK8=False, scheme="btvSF"):
0436     process.customizeJetTask = cms.Task()
0437     process.schedule.associate(process.customizeJetTask)
0438     
0439     CommonVars = cms.PSet(
0440         Proba=Var("bDiscriminator('pfJetProbabilityBJetTags')",
0441                   float,
0442                   doc="Jet Probability (Usage:BTV)",
0443                   precision=10),
0444         ProbaN=Var("bDiscriminator('pfNegativeOnlyJetProbabilityBJetTags')",
0445                   float,
0446                   doc="Negative-only Jet Probability (Usage:BTV)",
0447                   precision=10),
0448         Bprob=Var("bDiscriminator('pfJetBProbabilityBJetTags')",
0449                   float,
0450                   doc="Jet B Probability (Usage:BTV)",
0451                   precision=10),
0452         BprobN=Var("bDiscriminator('pfNegativeOnlyJetBProbabilityBJetTags')",
0453                   float,
0454                   doc="Negative-only Jet B Probability (Usage:BTV)",
0455                   precision=10),
0456     )
0457     
0458     # decouple these from CommonVars, not relevant for data
0459     HadronCountingVars = cms.PSet(
0460         nBHadrons=Var("jetFlavourInfo().getbHadrons().size()",
0461                       int,
0462                       doc="number of b-hadrons"),
0463         nCHadrons=Var("jetFlavourInfo().getcHadrons().size()",
0464                       int,
0465                       doc="number of c-hadrons")
0466     )
0467 
0468     # AK4
0469     if addAK4:
0470         if scheme == "btvSF":
0471             _n_cpf = 3 
0472             _n_npf = 3
0473             _n_sv = 4
0474         elif scheme == "DeepJet":
0475             _n_cpf = 25 
0476             _n_npf = 25
0477             _n_sv = 4
0478         elif scheme == "RobustParTAK4":
0479             _n_cpf = 25 
0480             _n_npf = 25
0481             _n_sv = 12
0482         process = update_jets_AK4(process)
0483             
0484         process.customJetExtTable = cms.EDProducer(
0485             "SimplePATJetFlatTableProducer",
0486             src=jetPuppiTable.src,
0487             cut=jetPuppiTable.cut,
0488             name=jetPuppiTable.name,
0489             doc=jetPuppiTable.doc,
0490             singleton=cms.bool(False),  # the number of entries is variable
0491             extension=cms.bool(True),  # this is the extension table for Jets
0492             variables=cms.PSet(
0493                 CommonVars,
0494                 get_DeepCSV_vars(),
0495                 get_DeepJet_outputs(),  # outputs are added in any case, inputs only if requested
0496                 get_ParticleNetAK4_outputs(),
0497                 #get_ParticleTransformerAK4_outputs(),# removed in 2024
0498             ))
0499 
0500         # disable the ParT branches in default jetPuppi table
0501         from PhysicsTools.NanoAOD.nano_eras_cff import run3_nanoAOD_122, run3_nanoAOD_124
0502         (run3_nanoAOD_122 | run3_nanoAOD_124).toModify(
0503             process.jetPuppiTable.variables,
0504             btagRobustParTAK4B = None,
0505             btagRobustParTAK4CvL = None,
0506             btagRobustParTAK4CvB = None,
0507             btagRobustParTAK4QG = None,
0508         )
0509     
0510     
0511          # from Run3 onwards, always set storeAK4Truth to True for MC
0512         process.customAK4ConstituentsForDeepJetTable = cms.EDProducer("PatJetDeepJetTableProducer",
0513                                                                         jets = cms.InputTag("linkedObjects","jets"),
0514                                                                         n_cpf=cms.uint32(_n_cpf),
0515                                                                         n_npf=cms.uint32(_n_npf),
0516                                                                         n_sv=cms.uint32(_n_sv)
0517                                                                       )
0518         process.customizeJetTask.add(process.customJetExtTable)
0519         process.customizeJetTask.add(process.customAK4ConstituentsForDeepJetTable)
0520     # AK8
0521     if addAK8:
0522         process = update_jets_AK8(process)
0523         process = update_jets_AK8_subjet(process)
0524         process.customFatJetExtTable = cms.EDProducer(
0525             "SimplePATJetFlatTableProducer",
0526             src=fatJetTable.src,
0527             cut=fatJetTable.cut,
0528             name=fatJetTable.name,
0529             doc=fatJetTable.doc,
0530             singleton=cms.bool(False),  # the number of entries is variable
0531             extension=cms.bool(True),  # this is the extension table for FatJets
0532             variables=cms.PSet(
0533                 CommonVars,
0534                 #HadronCountingVars if runOnMC else cms.PSet(), # only necessary before 106x
0535                 get_DDX_vars() ,
0536             ))
0537 
0538 
0539         # Subjets
0540         process.customSubJetExtTable = cms.EDProducer(
0541             "SimplePATJetFlatTableProducer",
0542             src=subJetTable.src,
0543             cut=subJetTable.cut,
0544             name=subJetTable.name,
0545             doc=subJetTable.doc,
0546             singleton=cms.bool(False),  # the number of entries is variable
0547             extension=cms.bool(True),  # this is the extension table for FatJets
0548             variables=cms.PSet(
0549                 CommonVars,
0550                 #HadronCountingVars if runOnMC else cms.PSet(), # only necessary before 106x
0551         ))
0552 
0553         process.customizeJetTask.add(process.customFatJetExtTable)
0554         process.customizeJetTask.add(process.customSubJetExtTable)
0555 
0556 
0557 
0558 #  From https://github.com/cms-jet/PFNano/blob/13_0_7_from124MiniAOD/python/addPFCands_cff.py
0559 def addPFCands(process, allPF = False, addAK4=False, addAK8=False):
0560     process.customizedPFCandsTask = cms.Task()
0561     process.schedule.associate(process.customizedPFCandsTask)
0562 
0563     process.finalJetsAK8Constituents = cms.EDProducer("PatJetConstituentPtrSelector",
0564                                             src = cms.InputTag("finalJetsAK8"),
0565                                             cut = cms.string("")
0566                                             )
0567     process.finalJetsAK4Constituents = cms.EDProducer("PatJetConstituentPtrSelector",
0568                                             src = cms.InputTag("finalJetsPuppi"),
0569                                             cut = cms.string("")
0570                                             )
0571     if allPF:
0572         candInput = cms.InputTag("packedPFCandidates")
0573     elif not addAK8:
0574         candList = cms.VInputTag(cms.InputTag("finalJetsAK4Constituents", "constituents"))
0575         process.customizedPFCandsTask.add(process.finalJetsAK4Constituents)
0576         process.finalJetsConstituentsTable = cms.EDProducer("PackedCandidatePtrMerger", src = candList, skipNulls = cms.bool(True), warnOnSkip = cms.bool(True))
0577         candInput = cms.InputTag("finalJetsConstituentsTable")
0578     elif not addAK4:
0579         candList = cms.VInputTag(cms.InputTag("finalJetsAK8Constituents", "constituents"))
0580         process.customizedPFCandsTask.add(process.finalJetsAK8Constituents)
0581         process.finalJetsConstituentsTable = cms.EDProducer("PackedCandidatePtrMerger", src = candList, skipNulls = cms.bool(True), warnOnSkip = cms.bool(True))
0582         candInput = cms.InputTag("finalJetsConstituentsTable")
0583     else:
0584         candList = cms.VInputTag(cms.InputTag("finalJetsAK4Constituents", "constituents"), cms.InputTag("finalJetsAK8Constituents", "constituents"))
0585         process.customizedPFCandsTask.add(process.finalJetsAK4Constituents)
0586         process.customizedPFCandsTask.add(process.finalJetsAK8Constituents)
0587         process.finalJetsConstituentsTable = cms.EDProducer("PackedCandidatePtrMerger", src = candList, skipNulls = cms.bool(True), warnOnSkip = cms.bool(True))
0588         candInput = cms.InputTag("finalJetsConstituentsTable")
0589     
0590     process.customConstituentsExtTable = cms.EDProducer("SimplePATCandidateFlatTableProducer",
0591                                                         src = candInput,
0592                                                         cut = cms.string(""), #we should not filter after pruning
0593                                                         name = cms.string("PFCands"),
0594                                                         doc = cms.string("interesting particles from AK4 and AK8 jets"),
0595                                                         singleton = cms.bool(False), # the number of entries is variable
0596                                                         extension = cms.bool(False), # this is the extension table for the AK8 constituents
0597                                                         variables = cms.PSet(CandVars,
0598                                                             puppiWeight = Var("puppiWeight()", float, doc="Puppi weight",precision=10),
0599                                                             puppiWeightNoLep = Var("puppiWeightNoLep()", float, doc="Puppi weight removing leptons",precision=10),
0600                                                             vtxChi2 = Var("?hasTrackDetails()?vertexChi2():-1", float, doc="vertex chi2",precision=10),
0601                                                             trkChi2 = Var("?hasTrackDetails()?pseudoTrack().normalizedChi2():-1", float, doc="normalized trk chi2", precision=10),
0602                                                             dz = Var("?hasTrackDetails()?dz():-1", float, doc="pf dz", precision=10),
0603                                                             dzErr = Var("?hasTrackDetails()?dzError():-1", float, doc="pf dz err", precision=10),
0604                                                             d0 = Var("?hasTrackDetails()?dxy():-1", float, doc="pf d0", precision=10),
0605                                                             d0Err = Var("?hasTrackDetails()?dxyError():-1", float, doc="pf d0 err", precision=10),
0606                                                             pvAssocQuality = Var("pvAssociationQuality()", int, doc="primary vertex association quality. 0: NotReconstructedPrimary, 1: OtherDeltaZ, 4: CompatibilityBTag, 5: CompatibilityDz, 6: UsedInFitLoose, 7: UsedInFitTight"),
0607                                                             lostInnerHits = Var("lostInnerHits()", int, doc="lost inner hits. -1: validHitInFirstPixelBarrelLayer, 0: noLostInnerHits, 1: oneLostInnerHit, 2: moreLostInnerHits"),
0608                                                             lostOuterHits = Var("?hasTrackDetails()?pseudoTrack().hitPattern().numberOfLostHits('MISSING_OUTER_HITS'):0", int, doc="lost outer hits"),
0609                                                             numberOfHits = Var("numberOfHits()", int, doc="number of hits"),
0610                                                             numberOfPixelHits = Var("numberOfPixelHits()", int, doc="number of pixel hits"),
0611                                                             trkQuality = Var("?hasTrackDetails()?pseudoTrack().qualityMask():0", int, doc="track quality mask"),
0612                                                             trkHighPurity = Var("?hasTrackDetails()?pseudoTrack().quality('highPurity'):0", bool, doc="track is high purity"),
0613                                                             trkAlgo = Var("?hasTrackDetails()?pseudoTrack().algo():-1", int, doc="track algorithm"),
0614                                                             trkP = Var("?hasTrackDetails()?pseudoTrack().p():-1", float, doc="track momemtum", precision=-1),
0615                                                             trkPt = Var("?hasTrackDetails()?pseudoTrack().pt():-1", float, doc="track pt", precision=-1),
0616                                                             trkEta = Var("?hasTrackDetails()?pseudoTrack().eta():-1", float, doc="track pt", precision=12),
0617                                                             trkPhi = Var("?hasTrackDetails()?pseudoTrack().phi():-1", float, doc="track phi", precision=12),
0618                                                          )
0619                                     )
0620     process.customAK8ConstituentsTable = cms.EDProducer("PatJetConstituentTableProducer",
0621                                                         candidates = candInput,
0622                                                         jets = cms.InputTag("finalJetsAK8"),
0623                                                         jet_radius = cms.double(0.8),
0624                                                         name = cms.string("FatJetPFCands"),
0625                                                         idx_name = cms.string("pFCandsIdx"),
0626                                                         nameSV = cms.string("FatJetSVs"),
0627                                                         idx_nameSV = cms.string("sVIdx"),
0628                                                         )
0629     process.customAK4ConstituentsTable = cms.EDProducer("PatJetConstituentTableProducer",
0630                                                         candidates = candInput,
0631                                                         jets = cms.InputTag("finalJetsPuppi"), # was finalJets before
0632                                                         jet_radius = cms.double(0.4),
0633                                                         name = cms.string("JetPFCands"),
0634                                                         idx_name = cms.string("pFCandsIdx"),
0635                                                         nameSV = cms.string("JetSVs"),
0636                                                         idx_nameSV = cms.string("sVIdx"),
0637                                                         )
0638     process.customizedPFCandsTask.add(process.customConstituentsExtTable)
0639 
0640     if not allPF:
0641         process.customizedPFCandsTask.add(process.finalJetsConstituentsTable)
0642     # linkedObjects are WIP for Run3
0643     if addAK8:
0644         process.customizedPFCandsTask.add(process.customAK8ConstituentsTable)
0645     if addAK4: 
0646         process.customizedPFCandsTask.add(process.customAK4ConstituentsTable)
0647     
0648         
0649     return process
0650 
0651 ## Switches for BTV nano
0652 # Default(store SFs PFCands+TaggerInputs) for both AK4 & AK8 jets
0653 # btvNano_addAK4_switch, btvNano_addAK8_switch True, btvNano_addPF_switch  False, TaggerInput = "btvSF"
0654 
0655 btvNano_switch = cms.PSet(
0656     btvNano_addAK4_switch = cms.untracked.bool(True),
0657     btvNano_addAK8_switch = cms.untracked.bool(False),
0658     btvNano_addallPF_switch = cms.untracked.bool(False),
0659     TaggerInput = cms.string("btvSF")
0660   )
0661 
0662 def BTVCustomNanoAOD(process):
0663     addPFCands(process,btvNano_switch.btvNano_addallPF_switch,btvNano_switch.btvNano_addAK4_switch,btvNano_switch.btvNano_addAK8_switch)
0664     add_BTV(process, btvNano_switch.btvNano_addAK4_switch,btvNano_switch.btvNano_addAK8_switch,btvNano_switch.TaggerInput)
0665     process.load("PhysicsTools.NanoAOD.btvMC_cff")
0666     ### for MC
0667     if btvNano_switch.btvNano_addallPF_switch:
0668         process.nanoSequenceMC+=allPFPFCandsMCSequence
0669     else:
0670         if btvNano_switch.btvNano_addAK4_switch and btvNano_switch.btvNano_addAK8_switch :
0671             process.nanoSequenceMC+=ak4ak8PFCandsMCSequence
0672         elif btvNano_switch.btvNano_addAK4_switch and not btvNano_switch.btvNano_addAK8_switch :
0673             process.nanoSequenceMC+=ak4onlyPFCandsMCSequence
0674         elif not btvNano_switch.btvNano_addAK4_switch and btvNano_switch.btvNano_addAK8_switch:
0675             process.nanoSequenceMC+=ak8onlyPFCandsMCSequence
0676     
0677     return process