CMS 3D CMS Logo

TrackValidation_cff.py
Go to the documentation of this file.
1 import FWCore.ParameterSet.Config as cms
2 
7 from Validation.RecoTrack.trajectorySeedTracks_cfi import trajectorySeedTracks as _trajectorySeedTracks
11 import cutsRecoTracks_cfi
12 
15 from CommonTools.RecoAlgos.trackingParticleRefSelector_cfi import trackingParticleRefSelector as _trackingParticleRefSelector
16 from CommonTools.RecoAlgos.trackingParticleConversionRefSelector_cfi import trackingParticleConversionRefSelector as _trackingParticleConversionRefSelector
17 from SimTracker.TrackHistory.trackingParticleBHadronRefSelector_cfi import trackingParticleBHadronRefSelector as _trackingParticleBHadronRefSelector
19 from CommonTools.RecoAlgos.recoChargedRefCandidateToTrackRefProducer_cfi import recoChargedRefCandidateToTrackRefProducer as _recoChargedRefCandidateToTrackRefProducer
20 
21 import RecoTracker.IterativeTracking.iterativeTkConfig as _cfg
22 import RecoTracker.IterativeTracking.iterativeTkUtils as _utils
23 from Configuration.Eras.Modifier_fastSim_cff import fastSim
24 import six
25 
26 ### First define the stuff for the standard validation sequence
27 ## Track selectors
28 for _eraName, _postfix, _era in _cfg.allEras():
29  _seedProd = ["initialStepSeedsPreSplitting"]
30  _trackProd = ["initialStepTracksPreSplitting"]
31  if _eraName in ["trackingLowPU", "trackingPhase2PU140"]: # these don't have preSplitting
32  _seedProd = []
33  _trackProd = []
34 
35  locals()["_algos"+_postfix] = ["generalTracks"] + _cfg.iterationAlgos(_postfix) + ["duplicateMerge"]
36  locals()["_seedProducers"+_postfix] = _seedProd + _cfg.seedProducers(_postfix)
37  locals()["_trackProducers"+_postfix] = _trackProd + _cfg.trackProducers(_postfix)
38 
39  if _eraName != "trackingPhase2PU140":
40  locals()["_electronSeedProducers"+_postfix] = ["tripletElectronSeeds", "pixelPairElectronSeeds", "stripPairElectronSeeds"]
41  else:
42  locals()["_electronSeedProducers"+_postfix] = ["tripletElectronSeeds"]
43 
44 _removeForFastSimSeedProducers =["initialStepSeedsPreSplitting",
45  "jetCoreRegionalStepSeeds",
46  "muonSeededSeedsInOut",
47  "muonSeededSeedsOutIn"]
48 _seedProducers_fastSim = [ x for x in _seedProducers if x not in _removeForFastSimSeedProducers]
49 
50 _removeForFastTrackProducers = ["initialStepTracksPreSplitting",
51  "jetCoreRegionalStepTracks",
52  "muonSeededTracksInOut",
53  "muonSeededTracksOutIn"]
54 _trackProducers_fastSim = [ x for x in _trackProducers if x not in _removeForFastTrackProducers]
55 
56 def _algoToSelector(algo):
57  sel = ""
58  if algo != "generalTracks":
59  sel = algo[0].upper()+algo[1:]
60  return "cutsRecoTracks"+sel
61 
62 def _addSelectorsByAlgo(algos, modDict):
63  names = []
64  seq = cms.Sequence()
65  for algo in algos:
66  if algo == "generalTracks":
67  continue
68  modName = _algoToSelector(algo)
69  if modName not in modDict:
70  mod = cutsRecoTracks_cfi.cutsRecoTracks.clone(algorithm=[algo])
71  modDict[modName] = mod
72  else:
73  mod = modDict[modName]
74  names.append(modName)
75  seq += mod
76  return (names, seq)
77 def _addSelectorsByHp(algos, modDict):
78  seq = cms.Sequence()
79  names = []
80  for algo in algos:
81  modName = _algoToSelector(algo)
82  modNameHp = modName+"Hp"
83  if modNameHp not in modDict:
84  if algo == "generalTracks":
85  mod = cutsRecoTracks_cfi.cutsRecoTracks.clone(quality=["highPurity"])
86  else:
87  mod = modDict[modName].clone(quality=["highPurity"])
88  modDict[modNameHp] = mod
89  else:
90  mod = modDict[modNameHp]
91  names.append(modNameHp)
92  seq += mod
93  return (names, seq)
94 def _addSelectorsBySrc(modules, midfix, src, modDict):
95  seq = cms.Sequence()
96  names = []
97  for modName in modules:
98  modNameNew = modName.replace("cutsRecoTracks", "cutsRecoTracks"+midfix)
99  if modNameNew not in modDict:
100  mod = modDict[modName].clone(src=src)
101  modDict[modNameNew] = mod
102  else:
103  mod = modDict[modNameNew]
104  names.append(modNameNew)
105  seq += mod
106  return (names, seq)
107 def _addSelectorsByOriginalAlgoMask(modules, midfix, algoParam,modDict):
108  seq = cms.Sequence()
109  names = []
110  for modName in modules:
111  if modName[-2:] == "Hp":
112  modNameNew = modName[:-2] + midfix + "Hp"
113  else:
114  modNameNew = modName + midfix
115  if modNameNew not in modDict:
116  mod = modDict[modName].clone()
117  setattr(mod, algoParam, mod.algorithm.value())
118  mod.algorithm = []
119  modDict[modNameNew] = mod
120  else:
121  mod = modDict[modNameNew]
122  names.append(modNameNew)
123  seq += mod
124  return (names, seq)
125 def _addSeedToTrackProducers(seedProducers,modDict):
126  names = []
127  seq = cms.Sequence()
128  for seed in seedProducers:
129  modName = "seedTracks"+seed
130  if modName not in modDict:
131  mod = _trajectorySeedTracks.clone(src=seed)
132  modDict[modName] = mod
133  else:
134  mod = modDict[modName]
135  names.append(modName)
136  seq += mod
137  return (names, seq)
138 
139 _relevantEras = _cfg.allEras()
140 _relevantErasAndFastSim = _relevantEras + [("fastSim", "_fastSim", fastSim)]
141 def _translateArgs(args, postfix, modDict):
142  ret = []
143  for arg in args:
144  if isinstance(arg, list):
145  ret.append(_translateArgs(arg, postfix, modDict))
146  else:
147  ret.append(modDict[arg+postfix])
148  return ret
149 def _sequenceForEachEra(function, args, names, sequence, modDict, plainArgs=[], modifySequence=None, includeFastSim=False):
150  if sequence[0] != "_":
151  raise Exception("Sequence name is expected to begin with _")
152 
153  _eras = _relevantErasAndFastSim if includeFastSim else _relevantEras
154  for eraName, postfix, _era in _eras:
155  _args = _translateArgs(args, postfix, modDict)
156  _args.extend(plainArgs)
157  ret = function(*_args, modDict=modDict)
158  if len(ret) != 2:
159  raise Exception("_sequenceForEachEra is expected to return 2 values, but function returned %d" % len(ret))
160  modDict[names+postfix] = ret[0]
161  modDict[sequence+postfix] = ret[1]
162 
163  # The sequence of the first era will be the default one
164  defaultSequenceName = sequence+_eras[0][0]
165  defaultSequence = modDict[defaultSequenceName]
166  modDict[defaultSequenceName[1:]] = defaultSequence # remove leading underscore
167 
168  # Optionally modify sequences before applying the era
169  if modifySequence is not None:
170  for eraName, postfix, _era in _eras:
171  modifySequence(modDict[sequence+postfix])
172 
173  # Apply eras
174  for _eraName, _postfix, _era in _eras[1:]:
175  _era.toReplaceWith(defaultSequence, modDict[sequence+_postfix])
176 def _setForEra(module, eraName, era, **kwargs):
177  if eraName == "":
178  for key, value in six.iteritems(kwargs):
179  setattr(module, key, value)
180  else:
181  era.toModify(module, **kwargs)
182 
183 # Seeding layer sets
184 def _getSeedingLayers(seedProducers, config):
185  def _findSeedingLayers(name):
186  prod = getattr(config, name)
187  if hasattr(prod, "triplets"):
188  if hasattr(prod, "layerList"): # merger
189  return prod.layerList.refToPSet_.value()
190  return _findSeedingLayers(prod.triplets.getModuleLabel())
191  elif hasattr(prod, "doublets"):
192  return _findSeedingLayers(prod.doublets.getModuleLabel())
193  label = prod.trackingRegionsSeedingLayers.getModuleLabel()
194  if label != "":
195  return label
196  return prod.seedingLayers.getModuleLabel()
197 
198  seedingLayersMerged = []
199  for seedName in seedProducers:
200  seedProd = getattr(config, seedName)
201  seedingLayersName = None
202  seedingLayers = None
203  if hasattr(seedProd, "OrderedHitsFactoryPSet"): # old seeding framework
204  seedingLayersName = seedProd.OrderedHitsFactoryPSet.SeedingLayers.getModuleLabel()
205  elif hasattr(seedProd, "seedingHitSets"): # new seeding framework
206  seedingLayersName = _findSeedingLayers(seedProd.seedingHitSets.getModuleLabel())
207  elif hasattr(seedProd, "layerList"): # FastSim:
208  seedingLayers = seedProd.layerList.value()
209  else:
210  continue
211 
212  if seedingLayersName is not None:
213  seedingLayers = getattr(config, seedingLayersName).layerList.value()
214  for layerSet in seedingLayers:
215  if layerSet not in seedingLayersMerged:
216  seedingLayersMerged.append(layerSet)
217 
218  return seedingLayersMerged
219 import RecoTracker.IterativeTracking.iterativeTk_cff as _iterativeTk_cff
220 import RecoTracker.IterativeTracking.ElectronSeeds_cff as _ElectronSeeds_cff
221 for _eraName, _postfix, _era in _relevantErasAndFastSim:
222  _stdLayers = _getSeedingLayers(locals()["_seedProducers"+_postfix], _iterativeTk_cff)
223  _eleLayers = []
224  if "_electronSeedProducers"+_postfix in locals(): # doesn't exist for FastSim
225  for _layer in _getSeedingLayers(locals()["_electronSeedProducers"+_postfix], _ElectronSeeds_cff):
226  if _layer not in _stdLayers:
227  _eleLayers.append(_layer)
228 
229  locals()["_seedingLayerSets"+_postfix] = _stdLayers
230  locals()["_seedingLayerSetsForElectrons"+_postfix] = _eleLayers
231 
232 
233 # MVA selectors
234 def _getMVASelectors(postfix):
235  mvaSel = _utils.getMVASelectors(postfix)
236 
237  pset = cms.untracked.PSet()
238  for iteration, (trackProducer, classifiers) in six.iteritems(mvaSel):
239  setattr(pset, trackProducer, cms.untracked.vstring(classifiers))
240  return pset
241 for _eraName, _postfix, _era in _relevantEras:
242  locals()["_mvaSelectors"+_postfix] = _getMVASelectors(_postfix)
243 
244 # Validation iterative steps
245 _sequenceForEachEra(_addSelectorsByAlgo, args=["_algos"], names="_selectorsByAlgo", sequence="_tracksValidationSelectorsByAlgo", modDict=globals())
246 
247 # high purity
248 _sequenceForEachEra(_addSelectorsByHp, args=["_algos"], names="_selectorsByAlgoHp", sequence="_tracksValidationSelectorsByAlgoHp", modDict=globals())
249 
250 # by originalAlgo
251 for _eraName, _postfix, _era in _relevantEras:
252  locals()["_selectorsByAlgoAndHp"+_postfix] = locals()["_selectorsByAlgo"+_postfix] + locals()["_selectorsByAlgoHp"+_postfix]
253  # For ByAlgoMask
254  locals()["_selectorsByAlgoAndHpNoGenTk"+_postfix] = [n for n in locals()["_selectorsByAlgoAndHp"+_postfix] if n not in ["generalTracks", "cutsRecoTracksHp"]]
255  # For ByOriginalAlgo
256  locals()["_selectorsByAlgoAndHpNoGenTkDupMerge"+_postfix] = [n for n in locals()["_selectorsByAlgoAndHpNoGenTk"+_postfix] if n not in ["cutsRecoTracksDuplicateMerge", "cutsRecoTracksDuplicateMergeHp"]]
257 _sequenceForEachEra(_addSelectorsByOriginalAlgoMask, modDict = globals(),
258  args = ["_selectorsByAlgoAndHpNoGenTkDupMerge"], plainArgs = ["ByOriginalAlgo", "originalAlgorithm"],
259  names = "_selectorsByOriginalAlgo", sequence = "_tracksValidationSelectorsByOriginalAlgo")
260 
261 
262 for _eraName, _postfix, _era in _relevantEras:
263  selectors = locals()["_selectorsByAlgoHp"+_postfix]
264  locals()["_generalTracksHp"+_postfix] = selectors[0]
265  locals()["_selectorsByAlgoHp"+_postfix] = selectors[1:]
266 
267 # BTV-like selection
268 import PhysicsTools.RecoAlgos.btvTracks_cfi as btvTracks_cfi
269 cutsRecoTracksBtvLike = btvTracks_cfi.btvTrackRefs.clone()
270 
271 # Select tracks associated to AK4 jets
273 ak4JetTracksAssociatorExplicitAll = ak4JTA_cff.ak4JetTracksAssociatorExplicit.clone(
274  jets = "ak4PFJets"
275 )
277 import CommonTools.RecoAlgos.jetTracksAssociationToTrackRefs_cfi as jetTracksAssociationToTrackRefs_cfi
278 cutsRecoTracksAK4PFJets = jetTracksAssociationToTrackRefs_cfi.jetTracksAssociationToTrackRefs.clone(
279  association = "ak4JetTracksAssociatorExplicitAll",
280  jets = "ak4PFJets",
281  correctedPtMin = 10,
282 )
283 
284 
285 ## Select signal TrackingParticles, and do the corresponding associations
286 trackingParticlesSignal = _trackingParticleRefSelector.clone(
287  signalOnly = True,
288  chargedOnly = False,
289  tip = 1e5,
290  lip = 1e5,
291  minRapidity = -10,
292  maxRapidity = 10,
293  ptMin = 0,
294 )
295 
296 # select tracks with pT > 0.9 GeV (for upgrade fake rates)
297 generalTracksPt09 = cutsRecoTracks_cfi.cutsRecoTracks.clone(ptMin=0.9)
298 # and then the selectors
299 _sequenceForEachEra(_addSelectorsBySrc, modDict=globals(),
300  args=[["_generalTracksHp"]],
301  plainArgs=["Pt09", "generalTracksPt09"],
302  names="_selectorsPt09", sequence="_tracksValidationSelectorsPt09",
303  modifySequence=lambda seq:seq.insert(0, generalTracksPt09))
304 
305 # select tracks from the PV
306 from CommonTools.RecoAlgos.TrackWithVertexRefSelector_cfi import trackWithVertexRefSelector as _trackWithVertexRefSelector
307 generalTracksFromPV = _trackWithVertexRefSelector.clone(
308  src = "generalTracks",
309  ptMin = 0,
310  ptMax = 1e10,
311  ptErrorCut = 1e10,
312  quality = "loose",
313  vertexTag = "offlinePrimaryVertices",
314  nVertices = 1,
315  vtxFallback = False,
316  zetaVtx = 0.1, # 1 mm
317  rhoVtx = 1e10, # intentionally no dxy cut
318 )
319 # and then the selectors
320 _sequenceForEachEra(_addSelectorsBySrc, modDict=globals(),
321  args=[["_generalTracksHp"]],
322  plainArgs=["FromPV", "generalTracksFromPV"],
323  names="_selectorsFromPV", sequence="_tracksValidationSelectorsFromPV",
324  modifySequence=lambda seq: seq.insert(0, generalTracksFromPV))
325 
326 # select tracks with pT > 0.9 GeV from the PV
327 generalTracksFromPVPt09 = generalTracksPt09.clone(src="generalTracksFromPV")
328 # and then the selectors
329 _sequenceForEachEra(_addSelectorsBySrc, modDict=globals(),
330  args=[["_generalTracksHp"]],
331  plainArgs=["FromPVPt09", "generalTracksFromPVPt09"],
332  names="_selectorsFromPVPt09", sequence="_tracksValidationSelectorsFromPVPt09",
333  modifySequence=lambda seq: seq.insert(0, generalTracksFromPVPt09))
334 
335 ## Select conversion TrackingParticles, and define the corresponding associator
336 trackingParticlesConversion = _trackingParticleConversionRefSelector.clone()
337 
338 ## Select electron TPs
339 trackingParticlesElectron = _trackingParticleRefSelector.clone(
340  pdgId = [-11, 11],
341  signalOnly = False,
342  tip = 1e5,
343  lip = 1e5,
344  minRapidity = -10,
345  maxRapidity = 10,
346  ptMin = 0,
347 )
348 
349 # Select B-hadron TPs
350 trackingParticlesBHadron = _trackingParticleBHadronRefSelector.clone()
351 
352 ## MTV instances
353 trackValidator = Validation.RecoTrack.MultiTrackValidator_cfi.multiTrackValidator.clone(
354  useLogPt = cms.untracked.bool(True),
355  dodEdxPlots = True,
356  doPVAssociationPlots = True
357  #,minpT = cms.double(-1)
358  #,maxpT = cms.double(3)
359  #,nintpT = cms.int32(40)
360 )
361 fastSim.toModify(trackValidator,
362  dodEdxPlots = False)
363 
364 for _eraName, _postfix, _era in _relevantEras:
365  _setForEra(trackValidator, _eraName, _era,
366  label = ["generalTracks", locals()["_generalTracksHp"+_postfix]] +
367  locals()["_selectorsByAlgo"+_postfix] + locals()["_selectorsByAlgoHp"+_postfix] +
368  locals()["_selectorsByOriginalAlgo"+_postfix] +
369  ["generalTracksPt09"] + locals()["_selectorsPt09"+_postfix] +
370  [
371  "cutsRecoTracksBtvLike",
372  "cutsRecoTracksAK4PFJets"
373  ],
374  doResolutionPlotsForLabels = [
375  "generalTracks",
376  locals()["_generalTracksHp"+_postfix],
377  "generalTracksPt09",
378  "cutsRecoTracksBtvLike",
379  ]
380  )
381  _setForEra(trackValidator.histoProducerAlgoBlock, _eraName, _era, seedingLayerSets=locals()["_seedingLayerSets"+_postfix])
382 
383 # for low-pT
384 trackValidatorTPPtLess09 = trackValidator.clone(
385  dirName = "Tracking/TrackTPPtLess09/",
386  label = [x for x in trackValidator.label.value() if ("Pt09" not in x) and ("BtvLike" not in x) and ("AK4PFJets" not in x)],
387  ptMaxTP = 0.9, # set maximum pT globally
388  histoProducerAlgoBlock = dict(
389  TpSelectorForEfficiencyVsEta = dict(ptMin=0.05), # enough to set min pT here
390  TpSelectorForEfficiencyVsPhi = dict(ptMin=0.05),
391  TpSelectorForEfficiencyVsVTXR = dict(ptMin=0.05),
392  TpSelectorForEfficiencyVsVTXZ = dict(ptMin=0.05),
393  ),
394  doSimPlots = False, # same as in trackValidator, no need to repeat here
395  doRecoTrackPlots = False, # fake rates are same as in trackValidator, no need to repeat here
396  doResolutionPlotsForLabels = ["disabled"], # resolutions are same as in trackValidator, no need to repeat here
397 )
398 
399 # For efficiency of signal TPs vs. signal tracks, and fake rate of
400 # signal tracks vs. signal TPs
401 trackValidatorFromPV = trackValidator.clone(
402  dirName = "Tracking/TrackFromPV/",
403  label_tp_effic = "trackingParticlesSignal",
404  label_tp_fake = "trackingParticlesSignal",
405  label_tp_effic_refvector = True,
406  label_tp_fake_refvector = True,
407  trackCollectionForDrCalculation = "generalTracksFromPV",
408  doPlotsOnlyForTruePV = True,
409  doPVAssociationPlots = False,
410  doResolutionPlotsForLabels = ["disabled"],
411 )
412 for _eraName, _postfix, _era in _relevantEras:
413  _setForEra(trackValidatorFromPV, _eraName, _era,
414  label = ["generalTracksFromPV"] + locals()["_selectorsFromPV"+_postfix] + ["generalTracksFromPVPt09"] + locals()["_selectorsFromPVPt09"+_postfix],
415  doResolutionPlotsForLabels = [] # for standard "FromPV" do resolution plots for all input collections as they are already limited
416  )
417 
418 # For fake rate of signal tracks vs. all TPs, and pileup rate of
419 # signal tracks vs. non-signal TPs
420 trackValidatorFromPVAllTP = trackValidatorFromPV.clone(
421  dirName = "Tracking/TrackFromPVAllTP/",
422  label_tp_effic = trackValidator.label_tp_effic.value(),
423  label_tp_fake = trackValidator.label_tp_fake.value(),
424  label_tp_effic_refvector = False,
425  label_tp_fake_refvector = False,
426  doSimPlots = False,
427  doSimTrackPlots = False,
428  doResolutionPlotsForLabels = ["disabled"], # resolution plots are the same as in "trackValidatorFromPV"
429 )
430 
431 # For efficiency of all TPs vs. all tracks
432 trackValidatorAllTPEffic = trackValidator.clone(
433  dirName = "Tracking/TrackAllTPEffic/",
434  label = [x for x in trackValidator.label.value() if "Pt09" not in x],
435  doSimPlots = False,
436  doRecoTrackPlots = True, # Fake rate of all tracks vs. all TPs is already included in trackValidator, but we want the reco plots for other reasons
437  doPVAssociationPlots = False,
438  doResolutionPlotsForLabels = ["disabled"], # resolution plots are the same as in "trackValidator"
439 )
440 trackValidatorAllTPEffic.histoProducerAlgoBlock.generalTpSelector.signalOnly = False
441 trackValidatorAllTPEffic.histoProducerAlgoBlock.TpSelectorForEfficiencyVsEta.signalOnly = False
442 trackValidatorAllTPEffic.histoProducerAlgoBlock.TpSelectorForEfficiencyVsPhi.signalOnly = False
443 trackValidatorAllTPEffic.histoProducerAlgoBlock.TpSelectorForEfficiencyVsPt.signalOnly = False
444 trackValidatorAllTPEffic.histoProducerAlgoBlock.TpSelectorForEfficiencyVsVTXR.signalOnly = False
445 trackValidatorAllTPEffic.histoProducerAlgoBlock.TpSelectorForEfficiencyVsVTXZ.signalOnly = False
446 for _eraName, _postfix, _era in _relevantEras:
447  _setForEra(trackValidatorAllTPEffic, _eraName, _era, label = ["generalTracks", locals()["_generalTracksHp"+_postfix]])
448 
449 # Built tracks, in the standard sequence mainly for monitoring the track selection MVA
450 _trackValidatorSeedingBuilding = trackValidator.clone( # common for built tracks and seeds (in trackingOnly)
451  associators = ["quickTrackAssociatorByHits"],
452  UseAssociators = True,
453  dodEdxPlots = False,
454  doPVAssociationPlots = False,
455  doSimPlots = False,
456  doResolutionPlotsForLabels = ["disabled"],
457 )
458 trackValidatorBuilding = _trackValidatorSeedingBuilding.clone(
459  dirName = "Tracking/TrackBuilding/",
460  doMVAPlots = True,
461 )
462 for _eraName, _postfix, _era in _relevantErasAndFastSim:
463  _setForEra(trackValidatorBuilding, _eraName, _era, label = locals()["_trackProducers"+_postfix])
464 fastSim.toModify(trackValidatorBuilding, doMVAPlots=False)
465 for _eraName, _postfix, _era in _relevantEras:
466  _setForEra(trackValidatorBuilding, _eraName, _era, mvaLabels = locals()["_mvaSelectors"+_postfix])
467 
468 
469 # For conversions
470 trackValidatorConversion = trackValidator.clone(
471  dirName = "Tracking/TrackConversion/",
472  label = [
473  "convStepTracks",
474  "conversionStepTracks",
475  "ckfInOutTracksFromConversions",
476  "ckfOutInTracksFromConversions",
477  ],
478  label_tp_effic = "trackingParticlesConversion",
479  label_tp_effic_refvector = True,
480  associators = ["quickTrackAssociatorByHits"],
481  UseAssociators = True,
482  doSimPlots = True,
483  dodEdxPlots = False,
484  doPVAssociationPlots = False,
485  calculateDrSingleCollection = False,
486 )
487 from RecoTracker.ConversionSeedGenerators.ConversionStep_cff import convLayerPairs as _convLayerPairs
488 def _uniqueFirstLayers(layerList):
489  firstLayers = [layerSet.split("+")[0] for layerSet in layerList]
490  ret = []
491  for l in firstLayers:
492  if not l in ret:
493  ret.append(l)
494  # For conversions add also the mono-TEC to the list as 'TEC'
495  # is used for both matched and unmatched rphi/stereo hits
496  if l.startswith("TEC"):
497  ret.append("M"+l)
498  return ret
499 # PhotonConversionTrajectorySeedProducerFromSingleLeg keeps only the
500 # first hit of the pairs in the seed, bookkeeping those is the best we
501 # can do without major further development
502 trackValidatorConversion.histoProducerAlgoBlock.seedingLayerSets = _uniqueFirstLayers(_convLayerPairs.layerList.value())
503 # relax lip and tip
504 for n in ["Eta", "Phi", "Pt", "VTXR", "VTXZ"]:
505  pset = getattr(trackValidatorConversion.histoProducerAlgoBlock, "TpSelectorForEfficiencyVs"+n)
506  pset.lip = trackValidatorConversion.lipTP.value()
507  pset.tip = trackValidatorConversion.tipTP.value()
508 
509 # For electrons
510 trackValidatorGsfTracks = trackValidatorConversion.clone(
511  dirName = "Tracking/TrackGsf/",
512  label = ["electronGsfTracks"],
513  label_tp_effic = "trackingParticlesElectron",
514 )
515 # add the additional seeding layers from ElectronSeeds
516 for _eraName, _postfix, _era in _relevantEras:
517  _setForEra(trackValidatorGsfTracks.histoProducerAlgoBlock, _eraName, _era, seedingLayerSets=trackValidator.histoProducerAlgoBlock.seedingLayerSets.value()+locals()["_seedingLayerSetsForElectrons"+_postfix])
518 
519 
520 
521 # for B-hadrons
522 trackValidatorBHadron = trackValidator.clone(
523  dirName = "Tracking/TrackBHadron/",
524  label_tp_effic = "trackingParticlesBHadron",
525  label_tp_effic_refvector = True,
526  doSimPlots = True,
527  doRecoTrackPlots = False, # Fake rate is defined wrt. all TPs, and that is already included in trackValidator
528  dodEdxPlots = False,
529 )
530 for _eraName, _postfix, _era in _relevantEras:
531  _setForEra(trackValidatorBHadron, _eraName, _era,
532  label = ["generalTracks", locals()["_generalTracksHp"+_postfix], "cutsRecoTracksBtvLike"]
533  )
534 
535 
536 # the track selectors
537 tracksValidationSelectors = cms.Sequence(
538  tracksValidationSelectorsByAlgo +
539  tracksValidationSelectorsByAlgoHp +
540  tracksValidationSelectorsByOriginalAlgo +
541  cutsRecoTracksBtvLike +
542  ak4JetTracksAssociatorExplicitAll +
543  cutsRecoTracksAK4PFJets
544 )
545 tracksValidationTruth = cms.Sequence(
546  tpClusterProducer +
547  quickTrackAssociatorByHits +
548  trackingParticleRecoTrackAsssociation +
549  VertexAssociatorByPositionAndTracks +
550  trackingParticleNumberOfLayersProducer
551 )
552 fastSim.toModify(tracksValidationTruth, lambda x: x.remove(tpClusterProducer))
553 
554 tracksPreValidation = cms.Sequence(
555  tracksValidationSelectors +
556  tracksValidationSelectorsPt09 +
557  tracksValidationSelectorsFromPV +
558  tracksValidationSelectorsFromPVPt09 +
559  tracksValidationTruth +
560  cms.ignore(trackingParticlesSignal) +
561  cms.ignore(trackingParticlesElectron) +
562  trackingParticlesConversion
563 )
564 fastSim.toReplaceWith(tracksPreValidation, tracksPreValidation.copyAndExclude([
565  trackingParticlesElectron,
566  trackingParticlesConversion,
567 ]))
568 
569 tracksValidation = cms.Sequence(
570  tracksPreValidation +
571  trackValidator +
572  trackValidatorTPPtLess09 +
573  trackValidatorFromPV +
574  trackValidatorFromPVAllTP +
575  trackValidatorAllTPEffic +
576  trackValidatorBuilding +
577  trackValidatorConversion +
578  trackValidatorGsfTracks
579 )
580 fastSim.toReplaceWith(tracksValidation, tracksValidation.copyAndExclude([
581  trackValidatorConversion,
582  trackValidatorGsfTracks,
583 ]))
584 
585 ### Then define stuff for standalone mode (i.e. MTV with RECO+DIGI input)
586 
587 # Select by originalAlgo and algoMask
588 _sequenceForEachEra(_addSelectorsByOriginalAlgoMask, modDict = globals(),
589  args = ["_selectorsByAlgoAndHpNoGenTk"], plainArgs = ["ByAlgoMask", "algorithmMaskContains"],
590  names = "_selectorsByAlgoMask", sequence = "_tracksValidationSelectorsByAlgoMaskStandalone")
591 
592 # Select pT>0.9 by iteration
593 # Need to avoid generalTracks+HP because those are already included in the standard validator
594 _sequenceForEachEra(_addSelectorsBySrc, modDict = globals(),
595  args = ["_selectorsByAlgoAndHpNoGenTk"], plainArgs = ["Pt09", "generalTracksPt09"],
596  names = "_selectorsPt09Standalone", sequence = "_tracksValidationSelectorsPt09Standalone")
597 
598 # Select fromPV by iteration
599 # Need to avoid generalTracks+HP because those are already included in the standard validator
600 _sequenceForEachEra(_addSelectorsBySrc, modDict = globals(),
601  args = ["_selectorsByAlgoAndHpNoGenTk"], plainArgs = ["FromPV", "generalTracksFromPV"],
602  names = "_selectorsFromPVStandalone", sequence = "_tracksValidationSelectorsFromPVStandalone")
603 
604 # Select pt>0.9 and fromPV by iteration
605 # Need to avoid generalTracks+HP because those are already included in the standard validator
606 _sequenceForEachEra(_addSelectorsBySrc, modDict = globals(),
607  args = ["_selectorsByAlgoAndHpNoGenTk"], plainArgs = ["FromPVPt09", "generalTracksFromPVPt09"],
608  names = "_selectorsFromPVPt09Standalone", sequence = "_tracksValidationSelectorsFromPVPt09Standalone")
609 
610 # MTV instances
611 trackValidatorStandalone = trackValidator.clone()
612 trackValidatorTPPtLess09Standalone = trackValidatorTPPtLess09.clone()
613 for _eraName, _postfix, _era in _relevantEras:
614  _setForEra(trackValidatorStandalone, _eraName, _era, label = trackValidator.label + locals()["_selectorsByAlgoMask"+_postfix] + locals()["_selectorsPt09Standalone"+_postfix])
615  _setForEra(trackValidatorTPPtLess09Standalone, _eraName, _era, label = trackValidatorTPPtLess09.label + locals()["_selectorsByAlgoMask"+_postfix] + locals()["_selectorsPt09Standalone"+_postfix])
616 
617 trackValidatorFromPVStandalone = trackValidatorFromPV.clone()
618 for _eraName, _postfix, _era in _relevantEras:
619  _setForEra(trackValidatorFromPVStandalone, _eraName, _era, label = trackValidatorFromPV.label + locals()["_selectorsFromPVStandalone"+_postfix] + locals()["_selectorsFromPVPt09Standalone"+_postfix])
620 # do resolutions as in the standard version
621 
622 trackValidatorFromPVAllTPStandalone = trackValidatorFromPVAllTP.clone(
623  label = trackValidatorFromPVStandalone.label.value()
624 )
625 trackValidatorAllTPEfficStandalone = trackValidatorAllTPEffic.clone(
626  label = [ x for x in trackValidator.label.value() if x not in ["cutsRecoTracksBtvLike", "cutsRecoTracksAK4PFJets"] and "Pt09" not in x]
627 )
628 
629 trackValidatorConversionStandalone = trackValidatorConversion.clone( label = [x for x in trackValidatorConversion.label if x != "convStepTracks"])
630 
631 trackValidatorBHadronStandalone = trackValidatorBHadron.clone(label = [x for x in trackValidatorStandalone.label if "Pt09" not in x])
632 
633 # sequences
634 tracksPreValidationStandalone = tracksPreValidation.copy()
635 tracksPreValidationStandalone += trackingParticlesBHadron
636 fastSim.toReplaceWith(tracksPreValidationStandalone, tracksPreValidation)
637 
638 tracksValidationSelectorsStandalone = cms.Sequence(
639  tracksValidationSelectorsByAlgoMaskStandalone +
640  tracksValidationSelectorsPt09Standalone +
641  tracksValidationSelectorsFromPVStandalone +
642  tracksValidationSelectorsFromPVPt09Standalone
643 )
644 
645 # we copy this for both Standalone and TrackingOnly
646 # and later make modifications from it which change based on era
647 _trackValidatorsBase = cms.Sequence(
648  trackValidatorStandalone +
649  trackValidatorTPPtLess09Standalone +
650  trackValidatorFromPVStandalone +
651  trackValidatorFromPVAllTPStandalone +
652  trackValidatorAllTPEfficStandalone +
653  trackValidatorConversionStandalone +
654  trackValidatorGsfTracks +
655  trackValidatorBHadronStandalone
656 )
657 trackValidatorsStandalone = _trackValidatorsBase.copy()
658 fastSim.toModify(trackValidatorsStandalone, lambda x: x.remove(trackValidatorConversionStandalone) )
659 
660 tracksValidationStandalone = cms.Sequence(
661  ak4PFL1FastL2L3CorrectorChain +
662  tracksPreValidationStandalone +
663  tracksValidationSelectorsStandalone +
664  trackValidatorsStandalone
665 )
666 
667 ### TrackingOnly mode (i.e. MTV with DIGI input + tracking-only reconstruction)
668 
669 # selectors
670 tracksValidationSelectorsTrackingOnly = tracksValidationSelectors.copyAndExclude([ak4JetTracksAssociatorExplicitAll,cutsRecoTracksAK4PFJets]) # selectors using track information only (i.e. no PF)
671 _sequenceForEachEra(_addSeedToTrackProducers, args=["_seedProducers"], names="_seedSelectors", sequence="_tracksValidationSeedSelectorsTrackingOnly", includeFastSim=True, modDict=globals())
672 
673 # MTV instances
674 trackValidatorTrackingOnly = trackValidatorStandalone.clone(label = [ x for x in trackValidatorStandalone.label if x != "cutsRecoTracksAK4PFJets"] )
675 
676 trackValidatorSeedingTrackingOnly = _trackValidatorSeedingBuilding.clone(
677  dirName = "Tracking/TrackSeeding/",
678  label = _seedSelectors,
679  doSeedPlots = True,
680 )
681 for _eraName, _postfix, _era in _relevantErasAndFastSim:
682  _setForEra(trackValidatorSeedingTrackingOnly, _eraName, _era, label = locals()["_seedSelectors"+_postfix])
683 
684 
685 trackValidatorConversionTrackingOnly = trackValidatorConversion.clone(label = [x for x in trackValidatorConversion.label if x not in ["ckfInOutTracksFromConversions", "ckfOutInTracksFromConversions"]])
686 
687 trackValidatorBHadronTrackingOnly = trackValidatorBHadron.clone(label = [x for x in trackValidatorTrackingOnly.label if "Pt09" not in x])
688 
689 # sequences
690 tracksPreValidationTrackingOnly = tracksPreValidationStandalone.copy()
691 tracksPreValidationTrackingOnly.replace(tracksValidationSelectors, tracksValidationSelectorsTrackingOnly)
692 
693 trackValidatorsTrackingOnly = _trackValidatorsBase.copy()
694 trackValidatorsTrackingOnly.replace(trackValidatorStandalone, trackValidatorTrackingOnly)
695 trackValidatorsTrackingOnly += trackValidatorSeedingTrackingOnly
696 trackValidatorsTrackingOnly += trackValidatorBuilding
697 trackValidatorsTrackingOnly.replace(trackValidatorConversionStandalone, trackValidatorConversionTrackingOnly)
698 trackValidatorsTrackingOnly.remove(trackValidatorGsfTracks)
699 trackValidatorsTrackingOnly.replace(trackValidatorBHadronStandalone, trackValidatorBHadronTrackingOnly)
700 fastSim.toModify(trackValidatorsTrackingOnly, lambda x: x.remove(trackValidatorConversionTrackingOnly))
701 fastSim.toModify(trackValidatorsTrackingOnly, lambda x: x.remove(trackValidatorBHadronTrackingOnly))
702 
703 
704 tracksValidationTrackingOnly = cms.Sequence(
705  tracksPreValidationTrackingOnly +
706  tracksValidationSelectorsStandalone +
707  tracksValidationSeedSelectorsTrackingOnly +
708  trackValidatorsTrackingOnly
709 )
710 
711 
712 ### Pixel tracking only mode (placeholder for now)
713 tpClusterProducerPixelTrackingOnly = tpClusterProducer.clone(
714  pixelClusterSrc = "siPixelClustersPreSplitting"
715 )
716 quickTrackAssociatorByHitsPixelTrackingOnly = quickTrackAssociatorByHits.clone(
717  cluster2TPSrc = "tpClusterProducerPixelTrackingOnly"
718 )
719 trackingParticlePixelTrackAsssociation = trackingParticleRecoTrackAsssociation.clone(
720  label_tr = "pixelTracks",
721  associator = "quickTrackAssociatorByHitsPixelTrackingOnly",
722 )
723 PixelVertexAssociatorByPositionAndTracks = VertexAssociatorByPositionAndTracks.clone(
724  trackAssociation = "trackingParticlePixelTrackAsssociation"
725 )
726 
727 trackValidatorPixelTrackingOnly = trackValidator.clone(
728  dirName = "Tracking/PixelTrack/",
729  label = ["pixelTracks"],
730  doResolutionPlotsForLabels = [],
731  trackCollectionForDrCalculation = "pixelTracks",
732  associators = ["trackingParticlePixelTrackAsssociation"],
733  label_vertex = "pixelVertices",
734  vertexAssociator = "PixelVertexAssociatorByPositionAndTracks",
735  dodEdxPlots = False,
736 )
737 
738 tracksValidationTruthPixelTrackingOnly = tracksValidationTruth.copy()
739 tracksValidationTruthPixelTrackingOnly.replace(tpClusterProducer, tpClusterProducerPixelTrackingOnly)
740 tracksValidationTruthPixelTrackingOnly.replace(quickTrackAssociatorByHits, quickTrackAssociatorByHitsPixelTrackingOnly)
741 tracksValidationTruthPixelTrackingOnly.replace(trackingParticleRecoTrackAsssociation, trackingParticlePixelTrackAsssociation)
742 tracksValidationTruthPixelTrackingOnly.replace(VertexAssociatorByPositionAndTracks, PixelVertexAssociatorByPositionAndTracks)
743 tracksValidationPixelTrackingOnly = cms.Sequence(
744  tracksValidationTruthPixelTrackingOnly +
745  trackValidatorPixelTrackingOnly
746 )
747 
748 
749 ### Lite mode (only generalTracks and HP)
750 trackValidatorLite = trackValidator.clone(
751  label = ["generalTracks", "cutsRecoTracksHp"]
752 )
753 tracksValidationLite = cms.Sequence(
754  cutsRecoTracksHp +
755  tracksValidationTruth +
756  trackValidatorLite
757 )
758 
759 ## customization for timing
760 from Configuration.Eras.Modifier_phase2_timing_layer_cff import phase2_timing_layer
761 phase2_timing_layer.toModify( generalTracksFromPV,
762  vertexTag = cms.InputTag('offlinePrimaryVertices4D'),
763  timesTag = cms.InputTag('trackTimeValueMapProducer:generalTracksConfigurableFlatResolutionModel'),
764  timeResosTag = cms.InputTag('trackTimeValueMapProducer:generalTracksConfigurableFlatResolutionModelResolution'),
765  nSigmaDtVertex = cms.double(3) )
766 phase2_timing_layer.toModify( trackValidatorStandalone,
767  label_vertex = cms.untracked.InputTag('offlinePrimaryVertices4D') )
768 phase2_timing_layer.toModify( trackValidatorFromPVStandalone,
769  label_vertex = cms.untracked.InputTag('offlinePrimaryVertices4D') )
770 phase2_timing_layer.toModify( trackValidatorFromPVAllTPStandalone,
771  label_vertex = cms.untracked.InputTag('offlinePrimaryVertices4D') )
772 phase2_timing_layer.toModify( trackValidatorConversionStandalone,
773  label_vertex = cms.untracked.InputTag('offlinePrimaryVertices4D') )
774 phase2_timing_layer.toModify( trackValidatorGsfTracks,
775  label_vertex = cms.untracked.InputTag('offlinePrimaryVertices4D') )
Definition: vlib.h:256
def _sequenceForEachEra(function, args, names, sequence, modDict, plainArgs=[], modifySequence=None, includeFastSim=False)
def _setForEra(module, eraName, era, kwargs)
def _translateArgs(args, postfix, modDict)
def _uniqueFirstLayers(layerList)
_addSelectorsByOriginalAlgoMask
Then define stuff for standalone mode (i.e.
def _getSeedingLayers(seedProducers, config)
TEveGeoShape * clone(const TEveElement *element, TEveElement *parent)
Definition: eve_macros.cc:135