CMS 3D CMS Logo

leptonTimeLifeInfo_common_cff.py
Go to the documentation of this file.
1 #
2 # Common definition of time-life variables for pat-leptons produced
3 # with {Electron,Muon,Tau}TimeLifeInfoTableProducer
4 #
5 import FWCore.ParameterSet.Config as cms
8 from PhysicsTools.PatAlgos.patRefitVertexProducer_cfi import patRefitVertexProducer
9 from PhysicsTools.NanoAOD.simpleVertexFlatTableProducer_cfi import simpleVertexFlatTableProducer
10 from PhysicsTools.PatAlgos.patElectronTimeLifeInfoProducer_cfi import patElectronTimeLifeInfoProducer
11 from PhysicsTools.PatAlgos.patMuonTimeLifeInfoProducer_cfi import patMuonTimeLifeInfoProducer
12 from PhysicsTools.PatAlgos.patTauTimeLifeInfoProducer_cfi import patTauTimeLifeInfoProducer
13 from PhysicsTools.NanoAOD.simpleCandidate2TrackTimeLifeInfoFlatTableProducer_cfi import simpleCandidate2TrackTimeLifeInfoFlatTableProducer
16 
17 # common settings of lepton life-time info producer
18 prod_common = cms.PSet(
19  pvSource = cms.InputTag("offlineSlimmedPrimaryVerticesWithBS"),
20  pvChoice = cms.int32(0) #0: PV[0], 1: smallest dz
21 )
22 
23 # impact parameter
24 ipVars = cms.PSet(
25  ipLength = Var("ipLength().value()", float, doc="lenght of impact parameter (3d)", precision=10),
26  ipLengthSig = Var("ipLength().significance()", float, doc="significance of impact parameter", precision=10),
27  IPx = Var("ipVector().x()", float, doc="x coordinate of impact parameter vector", precision=10),
28  IPy = Var("ipVector().y()", float, doc="y coordinate of impact parameter vector", precision=10),
29  IPz = Var("ipVector().z()", float, doc="z coordinate of impact parameter vector", precision=10)
30 )
31 
32 # track parameters and covariance at ref. point
33 trackVars = cms.PSet(
34  track_qoverp = Var("?hasTrack()?track().parameter(0):0", float, doc="track q/p", precision=10),
35  track_lambda = Var("?hasTrack()?track().parameter(1):0", float, doc="track lambda", precision=10),
36  track_phi = Var("?hasTrack()?track().parameter(2):0", float, doc="track phi", precision=10),
37  #track_deltaPhi = Var("?hasTrack()?deltaPhi(track().parameter(2), phi):0", float, doc="track phi minus lepton phi", precision=10),
38  track_dxy = Var("?hasTrack()?track().parameter(3):0", float, doc="track dxy", precision=10),
39  track_dsz = Var("?hasTrack()?track().parameter(4):0", float, doc="track dsz", precision=10),
40  bField_z = Var("?hasTrack()?bField_z:0", float, doc="z coordinate of magnetic field at track ref. point", precision=10),
41 )
42 # track covariance elements (adding to trackVars)
43 for i in range(0,5):
44  for j in range(i,5):
45  jistr = str(j)+str(i)
46  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))
47 
48 # secondary vertex
49 svVars = cms.PSet(
50  # SV
51  hasRefitSV = Var("hasSV()", bool, doc="has SV refit using miniAOD quantities"),
52  refitSVx = Var("?hasSV()?sv().x():0", float, doc="x coordinate of SV", precision=10),
53  refitSVy = Var("?hasSV()?sv().y():0", float, doc="y coordinate of SV", precision=10),
54  refitSVz = Var("?hasSV()?sv().z():0", float, doc="z coordinate of SV", precision=10),
55  refitSVchi2 = Var("?hasSV()?sv().chi2():0", float, doc="chi2 of SV fit", precision=8),
56  refitSVndof = Var("?hasSV()?sv().ndof():0", float, doc="ndof of SV fit", precision=8),
57  # flight-length
58  #refitFlightLength = Var("?hasSV()?flightLength().value():0", float, doc="flight-length,i.e. the PV to SV distance", precision=10),
59  #refitFlightLengthSig = Var("?hasSV()?flightLength().significance():0", float, doc="Significance of flight-length", precision=10)
60 )
61 # secondary vertex covariance elements (adding to svVars)
62 for i in range(0,3):
63  for j in range(i,3):
64  jistr = str(j)+str(i)
65  setattr(svVars, 'refitSVcov'+jistr, Var("?hasSV()?sv().covariance("+str(j)+","+str(i)+"):0", float, doc="Covariance of SV ("+str(j)+","+str(i)+")", precision=10))
66 
67 # primary vertex covariance elements
68 pvCovVars = cms.PSet()
69 for i in range(0,3):
70  for j in range(i,3):
71  jistr = str(j)+str(i)
72  setattr(pvCovVars, 'cov'+jistr, Var("covariance("+str(j)+","+str(i)+")", float, doc="vertex covariance ("+str(j)+","+str(i)+")", precision=10))
73 
74 # Module to refit PV with beam-spot constraint that is not present in Run-2 samples
75 refittedPV = patRefitVertexProducer.clone(
76  srcVertices = "offlineSlimmedPrimaryVertices",
77 )
78 run2_nanoAOD_ANY.toModify(
79  prod_common, pvSource = "refittedPV")
80 
81 # Definition of DQM plots
82 ipVarsPlots = cms.VPSet(
83  Plot1D('ipLength', 'ipLength', 25, -0.25, 0.25, 'signed lenght of impact parameter (3d)'),
84  Plot1D('ipLengthSig', 'ipLengthSig', 60, -5, 10, 'signed significance of impact parameter'),
85  Plot1D('IPx', 'IPx', 40, -0.02, 0.02, 'x coordinate of impact parameter vector'),
86  Plot1D('IPy', 'IPy', 40, -0.02, 0.02, 'y coordinate of impact parameter vector'),
87  Plot1D('IPz', 'IPz', 40, -0.02, 0.02, 'z coordinate of impact parameter vector')
88 )
89 trackVarsPlots = cms.VPSet(
90  Plot1D('track_qoverp', 'track_qoverp', 40, -0.2, 0.2, 'track q/p'),
91  Plot1D('track_lambda', 'track_lambda', 30, -1.5, 1.5, 'track lambda'),
92  Plot1D('track_phi', 'track_phi', 20, -3.14159, 3.14159, 'track phi'),
93  Plot1D('track_dxy', 'track_dxy', 20, -0.1, 0.1, 'track dxy'),
94  Plot1D('track_dsz', 'track_dsz', 20, -10, 10, 'track dsz'),
95  NoPlot('bField_z')
96 )
97 #no plots for track covariance elements, but store placeholders
98 for i in range(0,5):
99  for j in range(i,5):
100  trackVarsPlots.append(NoPlot('track_cov'+str(j)+str(i)))
101 svVarsPlots = cms.VPSet(
102  Plot1D('hasRefitSV', 'hasRefitSV', 2, 0, 2, 'has SV refit using miniAOD quantities'),
103  Plot1D('refitSVx', 'refitSVx', 20, -0.1, 0.1, 'x coordinate of refitted SV'),
104  Plot1D('refitSVy', 'refitSVy', 20, -0.1, 0.1, 'y coordinate of refitted SV'),
105  Plot1D('refitSVz', 'refitSVz', 20, -20, 20, 'z coordinate of refitted SV'),
106  Plot1D('refitSVchi2', 'refitSVchi2', 20, 0, 100, 'chi2 of SV fit'),
107  Plot1D('refitSVndof', 'refitSVndof', 10, 0, 10, 'ndof of SV fit')
108 )
109 #no plots for SV covariance elements, but store placeholders
110 for i in range(0,3):
111  for j in range(i,3):
112  svVarsPlots.append(NoPlot('refitSVcov'+str(j)+str(i)))
113 
114 #
115 # Customization sequences and functions
116 #
117 # electrons
118 electronVars = cms.PSet(
119  ipVars,
120  trackVars
121 )
122 for var in electronVars.parameters_():
123  setattr(getattr(electronVars, var), "src", cms.InputTag("electronTimeLifeInfos"))
125  process.electronTimeLifeInfos = patElectronTimeLifeInfoProducer.clone(
126  src = process.electronTable.src,
127  selection = 'pt > 15',
128  pvSource = prod_common.pvSource,
129  pvChoice = prod_common.pvChoice
130  )
131  process.electronTimeLifeInfoTable = simpleCandidate2TrackTimeLifeInfoFlatTableProducer.clone(
132  name = process.electronTable.name,
133  src = process.electronTable.src,
134  doc = cms.string("Additional time-life info for non-prompt electrons"),
135  extension = True,
136  externalTypedVariables = electronVars
137  )
138  process.electronTimeLifeInfoTask = cms.Task(
139  process.electronTimeLifeInfos,
140  process.electronTimeLifeInfoTable
141  )
142  # refit PV with beam-spot constraint that is not present in Run-2 samples
143  if not hasattr(process,'refittedPV'):
144  setattr(process,'refittedPV',refittedPV)
145  _electronTimeLifeInfoTaskRun2 = process.electronTimeLifeInfoTask.copy()
146  _electronTimeLifeInfoTaskRun2.add(process.refittedPV)
147  run2_nanoAOD_ANY.toReplaceWith(process.electronTimeLifeInfoTask,
148  _electronTimeLifeInfoTaskRun2)
149  process.electronTablesTask.add(process.electronTimeLifeInfoTask)
150  # add DQM plots if needed
151  if hasattr(process,'nanoDQM'):
152  process.nanoDQM.vplots.Electron.plots.extend(ipVarsPlots)
153  process.nanoDQM.vplots.Electron.plots.extend(trackVarsPlots)
154  return process
155 
156 # muons
157 muonVars = cms.PSet(
158  ipVars,
159  trackVars
160 )
161 for var in muonVars.parameters_():
162  setattr(getattr(muonVars, var), "src", cms.InputTag("muonTimeLifeInfos"))
164  process.muonTimeLifeInfos = patMuonTimeLifeInfoProducer.clone(
165  src = process.muonTable.src,
166  selection = 'pt > 15',
167  pvSource = prod_common.pvSource,
168  pvChoice = prod_common.pvChoice
169  )
170  process.muonTimeLifeInfoTable = simpleCandidate2TrackTimeLifeInfoFlatTableProducer.clone(
171  name = process.muonTable.name,
172  src = process.muonTable.src,
173  doc = cms.string("Additional time-life info for non-prompt muon"),
174  extension = True,
175  externalTypedVariables = muonVars
176  )
177  process.muonTimeLifeInfoTask = cms.Task(
178  process.muonTimeLifeInfos,
179  process.muonTimeLifeInfoTable
180  )
181  # refit PV with beam-spot constraint that is not present in Run-2 samples
182  if not hasattr(process,'refittedPV'):
183  setattr(process,'refittedPV',refittedPV)
184  _muonTimeLifeInfoTaskRun2 = process.muonTimeLifeInfoTask.copy()
185  _muonTimeLifeInfoTaskRun2.add(process.refittedPV)
186  run2_nanoAOD_ANY.toReplaceWith(process.muonTimeLifeInfoTask,
187  _muonTimeLifeInfoTaskRun2)
188  process.muonTablesTask.add(process.muonTimeLifeInfoTask)
189  # add DQM plots if needed
190  if hasattr(process,'nanoDQM'):
191  process.nanoDQM.vplots.Muon.plots.extend(ipVarsPlots)
192  process.nanoDQM.vplots.Muon.plots.extend(trackVarsPlots)
193  return process
194 
195 # taus
196 tauVars = cms.PSet(
197  svVars,
198  ipVars,
199  trackVars
200 )
201 for var in tauVars.parameters_():
202  setattr(getattr(tauVars, var), "src", cms.InputTag("tauTimeLifeInfos"))
204  process.tauTimeLifeInfos = patTauTimeLifeInfoProducer.clone(
205  src = process.tauTable.src,
206  pvSource = prod_common.pvSource,
207  pvChoice = prod_common.pvChoice
208  )
209  process.tauTimeLifeInfoTable = simpleCandidate2TrackTimeLifeInfoFlatTableProducer.clone(
210  name = process.tauTable.name,
211  src = process.tauTable.src,
212  doc = cms.string("Additional tau time-life info"),
213  extension = True,
214  externalTypedVariables = tauVars
215  )
216  process.tauTimeLifeInfoTask = cms.Task(
217  process.tauTimeLifeInfos,
218  process.tauTimeLifeInfoTable
219  )
220  # refit PV with beam-spot constraint that is not present in Run-2 samples
221  if not hasattr(process,'refittedPV'):
222  setattr(process,'refittedPV',refittedPV)
223  _tauTimeLifeInfoTaskRun2 = process.tauTimeLifeInfoTask.copy()
224  _tauTimeLifeInfoTaskRun2.add(process.refittedPV)
225  run2_nanoAOD_ANY.toReplaceWith(process.tauTimeLifeInfoTask,
226  _tauTimeLifeInfoTaskRun2)
227  process.tauTablesTask.add(process.tauTimeLifeInfoTask)
228  # add DQM plots if needed
229  if hasattr(process,'nanoDQM'):
230  process.nanoDQM.vplots.Tau.plots.extend(ipVarsPlots)
231  process.nanoDQM.vplots.Tau.plots.extend(trackVarsPlots)
232  process.nanoDQM.vplots.Tau.plots.extend(svVarsPlots)
233  return process
234 
235 # Vertices
236 def addExtendVertexInfo(process):
237  process.pvbsTable = simpleVertexFlatTableProducer.clone(
238  src = prod_common.pvSource,
239  name = "PVBS",
240  doc = "main primary vertex with beam-spot",
241  maxLen = 1,
242  variables = cms.PSet(
243  pvCovVars,
244  x = Var("position().x()", float, doc = "position x coordinate, in cm", precision = 10),
245  y = Var("position().y()", float, doc = "position y coordinate, in cm", precision = 10),
246  z = Var("position().z()", float, doc = "position z coordinate, in cm", precision = 16),
247  ndof = Var("ndof()", float, doc = "number of degrees of freedom", precision = 8),
248  chi2 = Var("normalizedChi2()", float, doc = "reduced chi2, i.e. chi2/ndof", precision = 8),
249  ),
250  )
251  process.pvbsTableTask = cms.Task(process.pvbsTable)
252  # refit PV with beam-spot constraint that is not present in Run-2 samples
253  if not hasattr(process,'refittedPV'):
254  setattr(process,'refittedPV',refittedPV)
255  _pvbsTableTaskRun2 = process.pvbsTableTask.copy()
256  _pvbsTableTaskRun2.add(process.refittedPV)
257  run2_nanoAOD_ANY.toReplaceWith(process.pvbsTableTask,
258  _pvbsTableTaskRun2)
259  process.vertexTablesTask.add(process.pvbsTableTask)
260  return process
261 
262 # Full
263 def addTimeLifeInfo(process):
265  addTimeLifeInfoToMuons(process)
266  addTimeLifeInfoToTaus(process)
267  addExtendVertexInfo(process)
268  return process
def Var(expr, valtype, doc=None, precision=-1)
Definition: common_cff.py:16
#define str(s)