Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 #
0002 # Common definition of time-life variables for pat-leptons produced
0003 # with {Electron,Muon,Tau}TimeLifeInfoTableProducer
0004 #
0005 import FWCore.ParameterSet.Config as cms
0006 from PhysicsTools.NanoAOD.common_cff import *
0007 from PhysicsTools.NanoAOD.nano_eras_cff import *
0008 from PhysicsTools.PatAlgos.patRefitVertexProducer_cfi import patRefitVertexProducer
0009 from PhysicsTools.NanoAOD.simpleVertexFlatTableProducer_cfi import simpleVertexFlatTableProducer
0010 from PhysicsTools.PatAlgos.patElectronTimeLifeInfoProducer_cfi import patElectronTimeLifeInfoProducer
0011 from PhysicsTools.PatAlgos.patMuonTimeLifeInfoProducer_cfi import patMuonTimeLifeInfoProducer
0012 from PhysicsTools.PatAlgos.patTauTimeLifeInfoProducer_cfi import patTauTimeLifeInfoProducer
0013 from PhysicsTools.NanoAOD.simplePATElectron2TrackTimeLifeInfoFlatTableProducer_cfi import simplePATElectron2TrackTimeLifeInfoFlatTableProducer
0014 from PhysicsTools.NanoAOD.simplePATMuon2TrackTimeLifeInfoFlatTableProducer_cfi import simplePATMuon2TrackTimeLifeInfoFlatTableProducer
0015 from PhysicsTools.NanoAOD.simplePATTau2TrackTimeLifeInfoFlatTableProducer_cfi import simplePATTau2TrackTimeLifeInfoFlatTableProducer
0016 from TrackingTools.TransientTrack.TransientTrackBuilder_cfi import *
0017 from PhysicsTools.NanoAOD.nanoDQM_tools_cff import *
0018 
0019 # common settings of lepton life-time info producer
0020 prod_common = cms.PSet(
0021     pvSource = cms.InputTag("offlineSlimmedPrimaryVerticesWithBS"),
0022     pvChoice = cms.int32(0) #0: PV[0], 1: smallest dz
0023 )
0024 
0025 # impact parameter
0026 ipVars = cms.PSet(
0027     #ipLength = Var("ipLength().value()", float, doc="lenght of impact parameter (3d)", precision=10),
0028     ipLengthSig = Var("ipLength().significance()", float, doc="significance of impact parameter", precision=10),
0029     IPx = Var("ipVector().x()", float, doc="x coordinate of impact parameter vector", precision=10),
0030     IPy = Var("ipVector().y()", float, doc="y coordinate of impact parameter vector", precision=10),
0031     IPz = Var("ipVector().z()", float, doc="z coordinate of impact parameter vector", precision=10)
0032 )
0033 
0034 # track parameters and covariance at ref. point
0035 trackVars = cms.PSet(
0036     track_qoverp = Var("?hasTrack()?track().parameter(0):0", float, doc="track q/p", precision=10),
0037     track_lambda = Var("?hasTrack()?track().parameter(1):0", float, doc="track lambda", precision=10),
0038     track_phi = Var("?hasTrack()?track().parameter(2):0", float, doc="track phi", precision=10),
0039     #track_deltaPhi = Var("?hasTrack()?deltaPhi(track().parameter(2), phi):0", float, doc="track phi minus lepton phi", precision=10),
0040     track_dxy = Var("?hasTrack()?track().parameter(3):0", float, doc="track dxy", precision=10),
0041     track_dsz = Var("?hasTrack()?track().parameter(4):0", float, doc="track dsz", precision=10),
0042     bField_z = Var("?hasTrack()?bField_z:0", float, doc="z coordinate of magnetic field at track ref. point", precision=10),
0043 )
0044 # track covariance elements (adding to trackVars)
0045 for i in range(0,5):
0046     for j in range(i,5):
0047         jistr = str(j)+str(i)
0048         setattr(trackVars, 'track_cov'+jistr, Var("?hasTrack()?track().covariance("+str(j)+","+str(i)+"):0", float, doc="track covariance element ("+str(j)+","+str(i)+")", precision=10))
0049 
0050 # secondary vertex
0051 svVars = cms.PSet(
0052     # SV
0053     hasRefitSV = Var("hasSV()", bool, doc="has SV refit using miniAOD quantities"),
0054     refitSVx = Var("?hasSV()?sv().x():0", float, doc="x coordinate of SV", precision=10),
0055     refitSVy = Var("?hasSV()?sv().y():0", float, doc="y coordinate of SV", precision=10),
0056     refitSVz = Var("?hasSV()?sv().z():0", float, doc="z coordinate of SV", precision=10),
0057     refitSVchi2 = Var("?hasSV()?sv().normalizedChi2():0", float, doc="reduced chi2, i.e. chi2/ndof, of SV fit", precision=8),
0058     #refitSVndof = Var("?hasSV()?sv().ndof():0", float, doc="ndof of SV fit", precision=8),
0059     # flight-length
0060     #refitFlightLength = Var("?hasSV()?flightLength().value():0", float, doc="flight-length,i.e. the PV to SV distance", precision=10),
0061     #refitFlightLengthSig = Var("?hasSV()?flightLength().significance():0", float, doc="Significance of flight-length", precision=10)
0062 )
0063 # secondary vertex covariance elements (adding to svVars)
0064 for i in range(0,3):
0065     for j in range(i,3):
0066         jistr = str(j)+str(i)
0067         setattr(svVars, 'refitSVcov'+jistr, Var("?hasSV()?sv().covariance("+str(j)+","+str(i)+"):0", float, doc="Covariance of SV ("+str(j)+","+str(i)+")", precision=10))
0068 
0069 # primary vertex covariance elements
0070 pvCovVars = cms.PSet()
0071 for i in range(0,3):
0072     for j in range(i,3):
0073         jistr = str(j)+str(i)
0074         setattr(pvCovVars, 'cov'+jistr, Var("covariance("+str(j)+","+str(i)+")", float, doc="vertex covariance ("+str(j)+","+str(i)+")", precision=10))
0075 
0076 # Module to refit PV with beam-spot constraint that is not present in Run-2 samples
0077 refittedPV = patRefitVertexProducer.clone(
0078     srcVertices = "offlineSlimmedPrimaryVertices",
0079 )
0080 run2_nanoAOD_ANY.toModify(
0081     prod_common, pvSource = "refittedPV")
0082 
0083 # Definition of DQM plots
0084 ipVarsPlots = cms.VPSet(
0085     #Plot1D('ipLength', 'ipLength', 25, -0.25, 0.25, 'signed lenght of impact parameter (3d)'),
0086     Plot1D('ipLengthSig', 'ipLengthSig', 60, -5, 10, 'signed significance of impact parameter'),
0087     Plot1D('IPx', 'IPx', 40, -0.02, 0.02, 'x coordinate of impact parameter vector'),
0088     Plot1D('IPy', 'IPy', 40, -0.02, 0.02, 'y coordinate of impact parameter vector'),
0089     Plot1D('IPz', 'IPz', 40, -0.02, 0.02, 'z coordinate of impact parameter vector')
0090 )
0091 trackVarsPlots = cms.VPSet(
0092     Plot1D('track_qoverp', 'track_qoverp', 40, -0.2, 0.2, 'track q/p'),
0093     Plot1D('track_lambda', 'track_lambda', 30, -1.5, 1.5, 'track lambda'),
0094     Plot1D('track_phi', 'track_phi', 20, -3.14159, 3.14159, 'track phi'),
0095     Plot1D('track_dxy', 'track_dxy', 20, -0.1, 0.1, 'track dxy'),
0096     Plot1D('track_dsz', 'track_dsz', 20, -10, 10, 'track dsz'),
0097     NoPlot('bField_z')
0098 )
0099 #no plots for track covariance elements, but store placeholders
0100 for i in range(0,5):
0101     for j in range(i,5):
0102         trackVarsPlots.append(NoPlot('track_cov'+str(j)+str(i)))
0103 svVarsPlots = cms.VPSet(
0104     Plot1D('hasRefitSV', 'hasRefitSV', 2, 0, 2, 'has SV refit using miniAOD quantities'),
0105     Plot1D('refitSVx', 'refitSVx', 20, -0.1, 0.1, 'x coordinate of refitted SV'),
0106     Plot1D('refitSVy', 'refitSVy', 20, -0.1, 0.1, 'y coordinate of refitted SV'),
0107     Plot1D('refitSVz', 'refitSVz', 20, -20, 20, 'z coordinate of refitted SV'),
0108     Plot1D('refitSVchi2', 'refitSVchi2', 20, 0, 40, 'reduced chi2 of SV fit'),
0109     #Plot1D('refitSVndof', 'refitSVndof', 10, 0, 10, 'ndof of SV fit')
0110 )
0111 #no plots for SV covariance elements, but store placeholders
0112 for i in range(0,3):
0113     for j in range(i,3):
0114         svVarsPlots.append(NoPlot('refitSVcov'+str(j)+str(i)))
0115 
0116 #
0117 # Customization sequences and functions
0118 #
0119 # electrons
0120 def addElectronTimeLifeInfoTask(process):
0121     process.electronTimeLifeInfos = patElectronTimeLifeInfoProducer.clone(
0122         src = process.electronTable.src,
0123         selection = 'pt > 15',
0124         pvSource = prod_common.pvSource,
0125         pvChoice = prod_common.pvChoice
0126     )
0127     process.electronTimeLifeInfoTable = simplePATElectron2TrackTimeLifeInfoFlatTableProducer.clone(
0128         name = process.electronTable.name,
0129         src = process.electronTable.src,
0130         doc = cms.string("Additional time-life info for non-prompt electrons"),
0131         extension = True,
0132         externalTypedVariables = cms.PSet()
0133     )
0134     process.electronTimeLifeInfoTask = cms.Task(
0135         process.electronTimeLifeInfos,
0136         process.electronTimeLifeInfoTable
0137     )
0138     # refit PV with beam-spot constraint that is not present in Run-2 samples
0139     if not hasattr(process,'refittedPV'):
0140         setattr(process,'refittedPV',refittedPV)
0141     _electronTimeLifeInfoTaskRun2 = process.electronTimeLifeInfoTask.copy()
0142     _electronTimeLifeInfoTaskRun2.add(process.refittedPV)
0143     run2_nanoAOD_ANY.toReplaceWith(process.electronTimeLifeInfoTask,
0144                                    _electronTimeLifeInfoTaskRun2)
0145     process.electronTablesTask.add(process.electronTimeLifeInfoTask)
0146     return process
0147 #base vars
0148 electronVars = cms.PSet(
0149     ipVars
0150 )
0151 for var in electronVars.parameters_():
0152     setattr(getattr(electronVars, var), "src", cms.InputTag("electronTimeLifeInfos"))
0153 def addTimeLifeInfoToElectrons(process):
0154     if not hasattr(process,'electronTimeLifeInfoTask'):
0155         process = addElectronTimeLifeInfoTask(process)
0156     electronExtVars = cms.PSet(
0157         process.electronTimeLifeInfoTable.externalTypedVariables,
0158         electronVars
0159     )
0160     process.electronTimeLifeInfoTable.externalTypedVariables = electronExtVars
0161     # add DQM plots if needed
0162     if hasattr(process,'nanoDQM'):
0163         process.nanoDQM.vplots.Electron.plots.extend(ipVarsPlots)
0164     return process
0165 #track vars
0166 electronTrackVars = cms.PSet(
0167     trackVars
0168 )
0169 for var in electronTrackVars.parameters_():
0170     setattr(getattr(electronTrackVars, var), "src", cms.InputTag("electronTimeLifeInfos"))
0171 def addElectronTrackVarsToTimeLifeInfo(process):
0172     if not hasattr(process,'electronTimeLifeInfoTask'):
0173         process = addElectronTimeLifeInfoTask(process)
0174     electronExtVars = cms.PSet(
0175         process.electronTimeLifeInfoTable.externalTypedVariables,
0176         electronTrackVars
0177     )
0178     process.electronTimeLifeInfoTable.externalTypedVariables = electronExtVars
0179     # add DQM plots if needed
0180     if hasattr(process,'nanoDQM'):
0181         process.nanoDQM.vplots.Electron.plots.extend(trackVarsPlots)
0182     return process
0183 
0184 # muons
0185 def addMuonTimeLifeInfoTask(process):
0186     process.muonTimeLifeInfos = patMuonTimeLifeInfoProducer.clone(
0187         src = process.muonTable.src,
0188         selection = 'pt > 15',
0189         pvSource = prod_common.pvSource,
0190         pvChoice = prod_common.pvChoice
0191     )
0192     process.muonTimeLifeInfoTable = simplePATMuon2TrackTimeLifeInfoFlatTableProducer.clone(
0193         name = process.muonTable.name,
0194         src = process.muonTable.src,
0195         doc = cms.string("Additional time-life info for non-prompt muon"),
0196         extension = True,
0197         externalTypedVariables = cms.PSet()
0198     )
0199     process.muonTimeLifeInfoTask = cms.Task(
0200         process.muonTimeLifeInfos,
0201         process.muonTimeLifeInfoTable
0202     )
0203     # refit PV with beam-spot constraint that is not present in Run-2 samples
0204     if not hasattr(process,'refittedPV'):
0205         setattr(process,'refittedPV',refittedPV)
0206     _muonTimeLifeInfoTaskRun2 = process.muonTimeLifeInfoTask.copy()
0207     _muonTimeLifeInfoTaskRun2.add(process.refittedPV)
0208     run2_nanoAOD_ANY.toReplaceWith(process.muonTimeLifeInfoTask,
0209                                    _muonTimeLifeInfoTaskRun2)
0210     process.muonTablesTask.add(process.muonTimeLifeInfoTask)
0211     return process
0212 #base vars
0213 muonVars = cms.PSet(
0214     ipVars
0215 )
0216 for var in muonVars.parameters_():
0217     setattr(getattr(muonVars, var), "src", cms.InputTag("muonTimeLifeInfos"))
0218 def addTimeLifeInfoToMuons(process):
0219     if not hasattr(process,'muonTimeLifeInfoTask'):
0220         process = addMuonTimeLifeInfoTask(process)
0221     muonExtVars = cms.PSet(
0222         process.muonTimeLifeInfoTable.externalTypedVariables,
0223         muonVars
0224     )
0225     process.muonTimeLifeInfoTable.externalTypedVariables = muonExtVars
0226     # add DQM plots if needed
0227     if hasattr(process,'nanoDQM'):
0228         process.nanoDQM.vplots.Muon.plots.extend(ipVarsPlots)
0229     return process
0230 #track vars
0231 muonTrackVars = cms.PSet(
0232     trackVars
0233 )
0234 for var in muonTrackVars.parameters_():
0235     setattr(getattr(muonTrackVars, var), "src", cms.InputTag("muonTimeLifeInfos"))
0236 def addMuonTrackVarsToTimeLifeInfo(process):
0237     if not hasattr(process,'muonTimeLifeInfoTask'):
0238         process = addMuonTimeLifeInfoTask(process)
0239     muonExtVars = cms.PSet(
0240         process.muonTimeLifeInfoTable.externalTypedVariables,
0241         muonTrackVars
0242     )
0243     process.muonTimeLifeInfoTable.externalTypedVariables = muonExtVars
0244     # add DQM plots if needed
0245     if hasattr(process,'nanoDQM'):
0246         process.nanoDQM.vplots.Muon.plots.extend(trackVarsPlots)
0247     return process
0248 
0249 # taus
0250 def addTauTimeLifeInfoTask(process):
0251     process.tauTimeLifeInfos = patTauTimeLifeInfoProducer.clone(
0252         src = process.tauTable.src,
0253         pvSource = prod_common.pvSource,
0254         pvChoice = prod_common.pvChoice
0255     )
0256     process.tauTimeLifeInfoTable = simplePATTau2TrackTimeLifeInfoFlatTableProducer.clone(
0257         name = process.tauTable.name,
0258         src = process.tauTable.src,
0259         doc = cms.string("Additional tau time-life info"),
0260         extension = True,
0261         externalTypedVariables = cms.PSet()
0262     )
0263     process.tauTimeLifeInfoTask = cms.Task(
0264         process.tauTimeLifeInfos,
0265         process.tauTimeLifeInfoTable
0266     )
0267     # refit PV with beam-spot constraint that is not present in Run-2 samples
0268     if not hasattr(process,'refittedPV'):
0269         setattr(process,'refittedPV',refittedPV)
0270     _tauTimeLifeInfoTaskRun2 = process.tauTimeLifeInfoTask.copy()
0271     _tauTimeLifeInfoTaskRun2.add(process.refittedPV)
0272     run2_nanoAOD_ANY.toReplaceWith(process.tauTimeLifeInfoTask,
0273                                    _tauTimeLifeInfoTaskRun2)
0274     process.tauTablesTask.add(process.tauTimeLifeInfoTask)
0275     return process
0276 # base vars
0277 tauVars = cms.PSet(
0278     svVars,
0279     ipVars
0280 )
0281 for var in tauVars.parameters_():
0282     setattr(getattr(tauVars, var), "src", cms.InputTag("tauTimeLifeInfos"))
0283 def addTimeLifeInfoToTaus(process):
0284     if not hasattr(process,'tauTimeLifeInfoTask'):
0285         process = addTauTimeLifeInfoTask(process)
0286     tauExtVars = cms.PSet(
0287         process.tauTimeLifeInfoTable.externalTypedVariables,
0288         tauVars
0289     )
0290     process.tauTimeLifeInfoTable.externalTypedVariables = tauExtVars
0291     # add DQM plots if needed
0292     if hasattr(process,'nanoDQM'):
0293         process.nanoDQM.vplots.Tau.plots.extend(ipVarsPlots)
0294         process.nanoDQM.vplots.Tau.plots.extend(svVarsPlots)
0295     return process
0296 #track vars
0297 tauTrackVars = cms.PSet(
0298     trackVars
0299 )
0300 for var in tauTrackVars.parameters_():
0301     setattr(getattr(tauTrackVars, var), "src", cms.InputTag("tauTimeLifeInfos"))
0302 def addTauTrackVarsToTimeLifeInfo(process):
0303     if not hasattr(process,'tauTimeLifeInfoTask'):
0304         process = addTauTimeLifeInfoTask(process)
0305     tauExtVars = cms.PSet(
0306         process.tauTimeLifeInfoTable.externalTypedVariables,
0307         tauTrackVars
0308     )
0309     process.tauTimeLifeInfoTable.externalTypedVariables = tauExtVars
0310     # add DQM plots if needed
0311     if hasattr(process,'nanoDQM'):
0312         process.nanoDQM.vplots.Tau.plots.extend(trackVarsPlots)
0313     return process
0314 
0315 # Vertices
0316 def addExtendVertexInfo(process):
0317     process.pvbsTable = simpleVertexFlatTableProducer.clone(
0318         src = prod_common.pvSource,
0319         name = "PVBS",
0320         doc = "main primary vertex with beam-spot",
0321         maxLen = 1,
0322         variables = cms.PSet(
0323             pvCovVars,
0324             x = Var("position().x()", float, doc = "position x coordinate, in cm", precision = 10),
0325             y = Var("position().y()", float, doc = "position y coordinate, in cm", precision = 10),
0326             z = Var("position().z()", float, doc = "position z coordinate, in cm", precision = 16),
0327             #ndof = Var("ndof()", float, doc = "number of degrees of freedom", precision = 8),#MB: not important
0328             chi2 = Var("normalizedChi2()", float, doc = "reduced chi2, i.e. chi2/ndof", precision = 8),
0329         ),
0330     )
0331     process.pvbsTableTask = cms.Task(process.pvbsTable)
0332     # refit PV with beam-spot constraint that is not present in Run-2 samples
0333     if not hasattr(process,'refittedPV'):
0334         setattr(process,'refittedPV',refittedPV)
0335     _pvbsTableTaskRun2 = process.pvbsTableTask.copy()
0336     _pvbsTableTaskRun2.add(process.refittedPV)
0337     run2_nanoAOD_ANY.toReplaceWith(process.pvbsTableTask,
0338                                    _pvbsTableTaskRun2)
0339     process.vertexTablesTask.add(process.pvbsTableTask)
0340     return process
0341 
0342 # Time-life info without track parameters
0343 def addTimeLifeInfoBase(process):
0344     addTimeLifeInfoToElectrons(process)
0345     addTimeLifeInfoToMuons(process)
0346     addTimeLifeInfoToTaus(process)
0347     addExtendVertexInfo(process)
0348     return process
0349 # Add track parameters to time-life info
0350 def addTrackVarsToTimeLifeInfo(process):
0351     addElectronTrackVarsToTimeLifeInfo(process)
0352     addMuonTrackVarsToTimeLifeInfo(process)
0353     addTauTrackVarsToTimeLifeInfo(process)
0354     return process
0355 # Full
0356 def addTimeLifeInfo(process):
0357     addTimeLifeInfoBase(process)
0358     addTrackVarsToTimeLifeInfo(process)
0359     return process