Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2025-06-20 01:53:40

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