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