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 _removeForFastSimSeedProducers =["initialStepSeedsPreSplitting",
39  "jetCoreRegionalStepSeeds",
40  "muonSeededSeedsInOut",
41  "muonSeededSeedsOutIn"]
42 _seedProducers_fastSim = [ x for x in _seedProducers if x not in _removeForFastSimSeedProducers]
43 
44 _removeForFastTrackProducers = ["initialStepTracksPreSplitting",
45  "jetCoreRegionalStepTracks",
46  "muonSeededTracksInOut",
47  "muonSeededTracksOutIn"]
48 _trackProducers_fastSim = [ x for x in _trackProducers if x not in _removeForFastTrackProducers]
49 
50 def _algoToSelector(algo):
51  sel = ""
52  if algo != "generalTracks":
53  sel = algo[0].upper()+algo[1:]
54  return "cutsRecoTracks"+sel
55 
56 def _addSelectorsByAlgo(algos, modDict):
57  names = []
58  seq = cms.Sequence()
59  for algo in algos:
60  if algo == "generalTracks":
61  continue
62  modName = _algoToSelector(algo)
63  if modName not in modDict:
64  mod = cutsRecoTracks_cfi.cutsRecoTracks.clone(algorithm=[algo])
65  modDict[modName] = mod
66  else:
67  mod = modDict[modName]
68  names.append(modName)
69  seq += mod
70  return (names, seq)
71 def _addSelectorsByHp(algos, modDict):
72  seq = cms.Sequence()
73  names = []
74  for algo in algos:
75  modName = _algoToSelector(algo)
76  modNameHp = modName+"Hp"
77  if modNameHp not in modDict:
78  if algo == "generalTracks":
79  mod = cutsRecoTracks_cfi.cutsRecoTracks.clone(quality=["highPurity"])
80  else:
81  mod = modDict[modName].clone(quality=["highPurity"])
82  modDict[modNameHp] = mod
83  else:
84  mod = modDict[modNameHp]
85  names.append(modNameHp)
86  seq += mod
87  return (names, seq)
88 def _addSelectorsBySrc(modules, midfix, src, modDict):
89  seq = cms.Sequence()
90  names = []
91  for modName in modules:
92  modNameNew = modName.replace("cutsRecoTracks", "cutsRecoTracks"+midfix)
93  if modNameNew not in modDict:
94  mod = modDict[modName].clone(src=src)
95  modDict[modNameNew] = mod
96  else:
97  mod = modDict[modNameNew]
98  names.append(modNameNew)
99  seq += mod
100  return (names, seq)
101 def _addSelectorsByOriginalAlgoMask(modules, midfix, algoParam,modDict):
102  seq = cms.Sequence()
103  names = []
104  for modName in modules:
105  if modName[-2:] == "Hp":
106  modNameNew = modName[:-2] + midfix + "Hp"
107  else:
108  modNameNew = modName + midfix
109  if modNameNew not in modDict:
110  mod = modDict[modName].clone()
111  setattr(mod, algoParam, mod.algorithm.value())
112  mod.algorithm = []
113  modDict[modNameNew] = mod
114  else:
115  mod = modDict[modNameNew]
116  names.append(modNameNew)
117  seq += mod
118  return (names, seq)
119 def _addSeedToTrackProducers(seedProducers,modDict):
120  names = []
121  seq = cms.Sequence()
122  for seed in seedProducers:
123  modName = "seedTracks"+seed
124  if modName not in modDict:
125  mod = _trajectorySeedTracks.clone(src=seed)
126  modDict[modName] = mod
127  else:
128  mod = modDict[modName]
129  names.append(modName)
130  seq += mod
131  return (names, seq)
132 
133 _relevantEras = _cfg.allEras()
134 _relevantErasAndFastSim = _relevantEras + [("fastSim", "_fastSim", fastSim)]
135 def _translateArgs(args, postfix, modDict):
136  ret = []
137  for arg in args:
138  if isinstance(arg, list):
139  ret.append(_translateArgs(arg, postfix, modDict))
140  else:
141  ret.append(modDict[arg+postfix])
142  return ret
143 def _sequenceForEachEra(function, args, names, sequence, modDict, plainArgs=[], modifySequence=None, includeFastSim=False):
144  if sequence[0] != "_":
145  raise Exception("Sequence name is expected to begin with _")
146 
147  _eras = _relevantErasAndFastSim if includeFastSim else _relevantEras
148  for eraName, postfix, _era in _eras:
149  _args = _translateArgs(args, postfix, modDict)
150  _args.extend(plainArgs)
151  ret = function(*_args, modDict=modDict)
152  if len(ret) != 2:
153  raise Exception("_sequenceForEachEra is expected to return 2 values, but function returned %d" % len(ret))
154  modDict[names+postfix] = ret[0]
155  modDict[sequence+postfix] = ret[1]
156 
157  # The sequence of the first era will be the default one
158  defaultSequenceName = sequence+_eras[0][0]
159  defaultSequence = modDict[defaultSequenceName]
160  modDict[defaultSequenceName[1:]] = defaultSequence # remove leading underscore
161 
162  # Optionally modify sequences before applying the era
163  if modifySequence is not None:
164  for eraName, postfix, _era in _eras:
165  modifySequence(modDict[sequence+postfix])
166 
167  # Apply eras
168  for _eraName, _postfix, _era in _eras[1:]:
169  _era.toReplaceWith(defaultSequence, modDict[sequence+_postfix])
170 def _setForEra(module, eraName, era, **kwargs):
171  if eraName == "":
172  for key, value in kwargs.iteritems():
173  setattr(module, key, value)
174  else:
175  era.toModify(module, **kwargs)
176 
177 # Seeding layer sets
178 def _getSeedingLayers(seedProducers):
179  import RecoTracker.IterativeTracking.iterativeTk_cff as _iterativeTk_cff
180 
181  def _findSeedingLayers(name):
182  prod = getattr(_iterativeTk_cff, name)
183  if hasattr(prod, "triplets"):
184  if hasattr(prod, "layerList"): # merger
185  return prod.layerList.refToPSet_.value()
186  return _findSeedingLayers(prod.triplets.getModuleLabel())
187  elif hasattr(prod, "doublets"):
188  return _findSeedingLayers(prod.doublets.getModuleLabel())
189  return prod.seedingLayers.getModuleLabel()
190 
191  seedingLayersMerged = []
192  for seedName in seedProducers:
193  seedProd = getattr(_iterativeTk_cff, seedName)
194  if hasattr(seedProd, "OrderedHitsFactoryPSet"):
195  if hasattr(seedProd, "SeedMergerPSet"):
196  seedingLayersName = seedProd.SeedMergerPSet.layerList.refToPSet_.value()
197  else:
198  seedingLayersName = seedProd.OrderedHitsFactoryPSet.SeedingLayers.getModuleLabel()
199  elif hasattr(seedProd, "seedingHitSets"):
200  seedingLayersName = _findSeedingLayers(seedProd.seedingHitSets.getModuleLabel())
201  else:
202  continue
203 
204  seedingLayers = getattr(_iterativeTk_cff, seedingLayersName).layerList.value()
205  for layerSet in seedingLayers:
206  if layerSet not in seedingLayersMerged:
207  seedingLayersMerged.append(layerSet)
208  return seedingLayersMerged
209 for _eraName, _postfix, _era in _relevantEras:
210  locals()["_seedingLayerSets"+_postfix] = _getSeedingLayers(locals()["_seedProducers"+_postfix])
211 
212 # MVA selectors
213 def _getMVASelectors(postfix):
214  mvaSel = _utils.getMVASelectors(postfix)
215 
216  pset = cms.untracked.PSet()
217  for iteration, (trackProducer, classifiers) in mvaSel.iteritems():
218  setattr(pset, trackProducer, cms.untracked.vstring(classifiers))
219  return pset
220 for _eraName, _postfix, _era in _relevantEras:
221  locals()["_mvaSelectors"+_postfix] = _getMVASelectors(_postfix)
222 
223 # Validation iterative steps
224 _sequenceForEachEra(_addSelectorsByAlgo, args=["_algos"], names="_selectorsByAlgo", sequence="_tracksValidationSelectorsByAlgo", modDict=globals())
225 
226 # high purity
227 _sequenceForEachEra(_addSelectorsByHp, args=["_algos"], names="_selectorsByAlgoHp", sequence="_tracksValidationSelectorsByAlgoHp", modDict=globals())
228 
229 for _eraName, _postfix, _era in _relevantEras:
230  selectors = locals()["_selectorsByAlgoHp"+_postfix]
231  locals()["_generalTracksHp"+_postfix] = selectors[0]
232  locals()["_selectorsByAlgoHp"+_postfix] = selectors[1:]
233 
234 # BTV-like selection
235 import PhysicsTools.RecoAlgos.btvTracks_cfi as btvTracks_cfi
236 cutsRecoTracksBtvLike = btvTracks_cfi.btvTrackRefs.clone()
237 
238 # Select tracks associated to AK4 jets
240 ak4JetTracksAssociatorExplicitAll = ak4JTA_cff.ak4JetTracksAssociatorExplicit.clone(
241  jets = "ak4PFJets"
242 )
244 import CommonTools.RecoAlgos.jetTracksAssociationToTrackRefs_cfi as jetTracksAssociationToTrackRefs_cfi
245 cutsRecoTracksAK4PFJets = jetTracksAssociationToTrackRefs_cfi.jetTracksAssociationToTrackRefs.clone(
246  association = "ak4JetTracksAssociatorExplicitAll",
247  jets = "ak4PFJets",
248  correctedPtMin = 10,
249 )
250 
251 
252 ## Select signal TrackingParticles, and do the corresponding associations
253 trackingParticlesSignal = _trackingParticleRefSelector.clone(
254  signalOnly = True,
255  chargedOnly = False,
256  tip = 1e5,
257  lip = 1e5,
258  minRapidity = -10,
259  maxRapidity = 10,
260  ptMin = 0,
261 )
262 
263 # select tracks with pT > 0.9 GeV (for upgrade fake rates)
264 generalTracksPt09 = cutsRecoTracks_cfi.cutsRecoTracks.clone(ptMin=0.9)
265 # and then the selectors
266 _sequenceForEachEra(_addSelectorsBySrc, modDict=globals(),
267  args=[["_generalTracksHp"]],
268  plainArgs=["Pt09", "generalTracksPt09"],
269  names="_selectorsPt09", sequence="_tracksValidationSelectorsPt09",
270  modifySequence=lambda seq:seq.insert(0, generalTracksPt09))
271 
272 # select tracks from the PV
273 from CommonTools.RecoAlgos.TrackWithVertexRefSelector_cfi import trackWithVertexRefSelector as _trackWithVertexRefSelector
274 generalTracksFromPV = _trackWithVertexRefSelector.clone(
275  src = "generalTracks",
276  ptMin = 0,
277  ptMax = 1e10,
278  ptErrorCut = 1e10,
279  quality = "loose",
280  vertexTag = "offlinePrimaryVertices",
281  nVertices = 1,
282  vtxFallback = False,
283  zetaVtx = 0.1, # 1 mm
284  rhoVtx = 1e10, # intentionally no dxy cut
285 )
286 # and then the selectors
287 _sequenceForEachEra(_addSelectorsBySrc, modDict=globals(),
288  args=[["_generalTracksHp"]],
289  plainArgs=["FromPV", "generalTracksFromPV"],
290  names="_selectorsFromPV", sequence="_tracksValidationSelectorsFromPV",
291  modifySequence=lambda seq: seq.insert(0, generalTracksFromPV))
292 
293 # select tracks with pT > 0.9 GeV from the PV
294 generalTracksFromPVPt09 = generalTracksPt09.clone(src="generalTracksFromPV")
295 # and then the selectors
296 _sequenceForEachEra(_addSelectorsBySrc, modDict=globals(),
297  args=[["_generalTracksHp"]],
298  plainArgs=["FromPVPt09", "generalTracksFromPVPt09"],
299  names="_selectorsFromPVPt09", sequence="_tracksValidationSelectorsFromPVPt09",
300  modifySequence=lambda seq: seq.insert(0, generalTracksFromPVPt09))
301 
302 ## Select conversion TrackingParticles, and define the corresponding associator
303 trackingParticlesConversion = _trackingParticleConversionRefSelector.clone()
304 
305 ## Select electron TPs
306 trackingParticlesElectron = _trackingParticleRefSelector.clone(
307  pdgId = [-11, 11],
308  signalOnly = False,
309  tip = 1e5,
310  lip = 1e5,
311  minRapidity = -10,
312  maxRapidity = 10,
313  ptMin = 0,
314 )
315 
316 # Select B-hadron TPs
317 trackingParticlesBHadron = _trackingParticleBHadronRefSelector.clone()
318 
319 ## MTV instances
320 trackValidator = Validation.RecoTrack.MultiTrackValidator_cfi.multiTrackValidator.clone(
321  useLogPt = cms.untracked.bool(True),
322  dodEdxPlots = True,
323  doPVAssociationPlots = True
324  #,minpT = cms.double(-1)
325  #,maxpT = cms.double(3)
326  #,nintpT = cms.int32(40)
327 )
328 fastSim.toModify(trackValidator,
329  dodEdxPlots = False)
330 
331 for _eraName, _postfix, _era in _relevantEras:
332  _setForEra(trackValidator, _eraName, _era,
333  label = ["generalTracks", locals()["_generalTracksHp"+_postfix]] +
334  locals()["_selectorsByAlgo"+_postfix] + locals()["_selectorsByAlgoHp"+_postfix] +
335  ["generalTracksPt09"] + locals()["_selectorsPt09"+_postfix] +
336  [
337  "cutsRecoTracksBtvLike",
338  "cutsRecoTracksAK4PFJets"
339  ]
340  )
341  _setForEra(trackValidator.histoProducerAlgoBlock, _eraName, _era, seedingLayerSets=locals()["_seedingLayerSets"+_postfix])
342 
343 # For efficiency of signal TPs vs. signal tracks, and fake rate of
344 # signal tracks vs. signal TPs
345 trackValidatorFromPV = trackValidator.clone(
346  dirName = "Tracking/TrackFromPV/",
347  label_tp_effic = "trackingParticlesSignal",
348  label_tp_fake = "trackingParticlesSignal",
349  label_tp_effic_refvector = True,
350  label_tp_fake_refvector = True,
351  trackCollectionForDrCalculation = "generalTracksFromPV",
352  doPlotsOnlyForTruePV = True,
353  doPVAssociationPlots = False,
354 )
355 for _eraName, _postfix, _era in _relevantEras:
356  _setForEra(trackValidatorFromPV, _eraName, _era, label = ["generalTracksFromPV"] + locals()["_selectorsFromPV"+_postfix] + ["generalTracksFromPVPt09"] + locals()["_selectorsFromPVPt09"+_postfix])
357 
358 # For fake rate of signal tracks vs. all TPs, and pileup rate of
359 # signal tracks vs. non-signal TPs
360 trackValidatorFromPVAllTP = trackValidatorFromPV.clone(
361  dirName = "Tracking/TrackFromPVAllTP/",
362  label_tp_effic = trackValidator.label_tp_effic.value(),
363  label_tp_fake = trackValidator.label_tp_fake.value(),
364  label_tp_effic_refvector = False,
365  label_tp_fake_refvector = False,
366  doSimPlots = False,
367  doSimTrackPlots = False,
368 )
369 
370 # For efficiency of all TPs vs. all tracks
371 trackValidatorAllTPEffic = trackValidator.clone(
372  dirName = "Tracking/TrackAllTPEffic/",
373  label = [x for x in trackValidator.label.value() if "Pt09" not in x],
374  doSimPlots = False,
375  doRecoTrackPlots = True, # Fake rate of all tracks vs. all TPs is already included in trackValidator, but we want the reco plots for other reasons
376  doPVAssociationPlots = False,
377 )
378 trackValidatorAllTPEffic.histoProducerAlgoBlock.generalTpSelector.signalOnly = False
379 trackValidatorAllTPEffic.histoProducerAlgoBlock.TpSelectorForEfficiencyVsEta.signalOnly = False
380 trackValidatorAllTPEffic.histoProducerAlgoBlock.TpSelectorForEfficiencyVsPhi.signalOnly = False
381 trackValidatorAllTPEffic.histoProducerAlgoBlock.TpSelectorForEfficiencyVsPt.signalOnly = False
382 trackValidatorAllTPEffic.histoProducerAlgoBlock.TpSelectorForEfficiencyVsVTXR.signalOnly = False
383 trackValidatorAllTPEffic.histoProducerAlgoBlock.TpSelectorForEfficiencyVsVTXZ.signalOnly = False
384 for _eraName, _postfix, _era in _relevantEras:
385  _setForEra(trackValidatorAllTPEffic, _eraName, _era, label = ["generalTracks", locals()["_generalTracksHp"+_postfix]])
386 
387 # For conversions
388 trackValidatorConversion = trackValidator.clone(
389  dirName = "Tracking/TrackConversion/",
390  label = [
391  "convStepTracks",
392  "conversionStepTracks",
393  "ckfInOutTracksFromConversions",
394  "ckfOutInTracksFromConversions",
395  ],
396  label_tp_effic = "trackingParticlesConversion",
397  label_tp_effic_refvector = True,
398  associators = ["quickTrackAssociatorByHits"],
399  UseAssociators = True,
400  doSimPlots = True,
401  dodEdxPlots = False,
402  doPVAssociationPlots = False,
403  calculateDrSingleCollection = False,
404 )
405 # relax lip and tip
406 for n in ["Eta", "Phi", "Pt", "VTXR", "VTXZ"]:
407  pset = getattr(trackValidatorConversion.histoProducerAlgoBlock, "TpSelectorForEfficiencyVs"+n)
408  pset.lip = trackValidatorConversion.lipTP.value()
409  pset.tip = trackValidatorConversion.tipTP.value()
410 
411 # For electrons
412 trackValidatorGsfTracks = trackValidatorConversion.clone(
413  dirName = "Tracking/TrackGsf/",
414  label = ["electronGsfTracks"],
415  label_tp_effic = "trackingParticlesElectron",
416 )
417 
418 # for B-hadrons
419 trackValidatorBHadron = trackValidator.clone(
420  dirName = "Tracking/TrackBHadron/",
421  label_tp_effic = "trackingParticlesBHadron",
422  label_tp_effic_refvector = True,
423  doSimPlots = True,
424  doRecoTrackPlots = False, # Fake rate is defined wrt. all TPs, and that is already included in trackValidator
425  dodEdxPlots = False,
426 )
427 for _eraName, _postfix, _era in _relevantEras:
428  _setForEra(trackValidatorBHadron, _eraName, _era,
429  label = ["generalTracks", locals()["_generalTracksHp"+_postfix], "cutsRecoTracksBtvLike"]
430  )
431 
432 
433 # the track selectors
434 tracksValidationSelectors = cms.Sequence(
435  tracksValidationSelectorsByAlgo +
436  tracksValidationSelectorsByAlgoHp +
437  cutsRecoTracksBtvLike +
438  ak4JetTracksAssociatorExplicitAll +
439  cutsRecoTracksAK4PFJets
440 )
441 tracksValidationTruth = cms.Sequence(
442  tpClusterProducer +
443  quickTrackAssociatorByHits +
444  trackingParticleRecoTrackAsssociation +
445  VertexAssociatorByPositionAndTracks +
446  trackingParticleNumberOfLayersProducer
447 )
448 fastSim.toModify(tracksValidationTruth, lambda x: x.remove(tpClusterProducer))
449 
450 tracksPreValidation = cms.Sequence(
451  tracksValidationSelectors +
452  tracksValidationSelectorsPt09 +
453  tracksValidationSelectorsFromPV +
454  tracksValidationSelectorsFromPVPt09 +
455  tracksValidationTruth +
456  cms.ignore(trackingParticlesSignal) +
457  cms.ignore(trackingParticlesElectron) +
458  trackingParticlesConversion
459 )
460 fastSim.toReplaceWith(tracksPreValidation, tracksPreValidation.copyAndExclude([
461  trackingParticlesElectron,
462  trackingParticlesConversion,
463 ]))
464 
465 tracksValidation = cms.Sequence(
466  tracksPreValidation +
467  trackValidator +
468  trackValidatorFromPV +
469  trackValidatorFromPVAllTP +
470  trackValidatorAllTPEffic +
471  trackValidatorConversion +
472  trackValidatorGsfTracks
473 )
474 fastSim.toReplaceWith(tracksValidation, tracksValidation.copyAndExclude([
475  trackValidatorConversion,
476  trackValidatorGsfTracks,
477 ]))
478 
479 ### Then define stuff for standalone mode (i.e. MTV with RECO+DIGI input)
480 
481 # Select by originalAlgo and algoMask
482 for _eraName, _postfix, _era in _relevantEras:
483  locals()["_selectorsByAlgoAndHp"+_postfix] = locals()["_selectorsByAlgo"+_postfix] + locals()["_selectorsByAlgoHp"+_postfix]
484 _sequenceForEachEra(_addSelectorsByOriginalAlgoMask, modDict = globals(),
485  args = ["_selectorsByAlgoAndHp"], plainArgs = ["ByOriginalAlgo", "originalAlgorithm"],
486  names = "_selectorsByOriginalAlgo", sequence = "_tracksValidationSelectorsByOriginalAlgoStandalone")
487 _sequenceForEachEra(_addSelectorsByOriginalAlgoMask, modDict = globals(),
488  args = ["_selectorsByAlgoAndHp"], plainArgs = ["ByAlgoMask", "algorithmMaskContains"],
489  names = "_selectorsByAlgoMask", sequence = "_tracksValidationSelectorsByAlgoMaskStandalone")
490 
491 # Select pT>0.9 by iteration
492 _sequenceForEachEra(_addSelectorsBySrc, modDict = globals(),
493  args = ["_selectorsByAlgoAndHp"], plainArgs = ["Pt09", "generalTracksPt09"],
494  names = "_selectorsPt09Standalone", sequence = "_tracksValidationSelectorsPt09Standalone")
495 
496 # Select fromPV by iteration
497 _sequenceForEachEra(_addSelectorsBySrc, modDict = globals(),
498  args = ["_selectorsByAlgoAndHp"], plainArgs = ["FromPV", "generalTracksFromPV"],
499  names = "_selectorsFromPVStandalone", sequence = "_tracksValidationSelectorsFromPVStandalone")
500 
501 # Select pt>0.9 and fromPV by iteration
502 _sequenceForEachEra(_addSelectorsBySrc, modDict = globals(),
503  args = ["_selectorsByAlgoAndHp"], plainArgs = ["FromPVPt09", "generalTracksFromPVPt09"],
504  names = "_selectorsFromPVPt09Standalone", sequence = "_tracksValidationSelectorsFromPVPt09Standalone")
505 
506 # MTV instances
507 trackValidatorStandalone = trackValidator.clone()
508 for _eraName, _postfix, _era in _relevantEras:
509  _setForEra(trackValidatorStandalone, _eraName, _era, label = trackValidator.label + locals()["_selectorsByOriginalAlgo"+_postfix] + locals()["_selectorsByAlgoMask"+_postfix] + locals()["_selectorsPt09Standalone"+_postfix])
510 
511 trackValidatorFromPVStandalone = trackValidatorFromPV.clone()
512 for _eraName, _postfix, _era in _relevantEras:
513  _setForEra(trackValidatorFromPVStandalone, _eraName, _era, label = trackValidatorFromPV.label + locals()["_selectorsFromPVStandalone"+_postfix] + locals()["_selectorsFromPVPt09Standalone"+_postfix])
514 
515 trackValidatorFromPVAllTPStandalone = trackValidatorFromPVAllTP.clone(
516  label = trackValidatorFromPVStandalone.label.value()
517 )
518 trackValidatorAllTPEfficStandalone = trackValidatorAllTPEffic.clone(
519  label = [ x for x in trackValidator.label.value() if x not in ["cutsRecoTracksBtvLike", "cutsRecoTracksAK4PFJets"] and "Pt09" not in x]
520 )
521 
522 trackValidatorConversionStandalone = trackValidatorConversion.clone( label = [x for x in trackValidatorConversion.label if x != "convStepTracks"])
523 
524 trackValidatorBHadronStandalone = trackValidatorBHadron.clone(label = [x for x in trackValidatorStandalone.label if "Pt09" not in x])
525 
526 # sequences
527 tracksPreValidationStandalone = tracksPreValidation.copy()
528 tracksPreValidationStandalone += trackingParticlesBHadron
529 fastSim.toReplaceWith(tracksPreValidationStandalone, tracksPreValidation)
530 
531 tracksValidationSelectorsStandalone = cms.Sequence(
532  tracksValidationSelectorsByOriginalAlgoStandalone +
533  tracksValidationSelectorsByAlgoMaskStandalone +
534  tracksValidationSelectorsPt09Standalone +
535  tracksValidationSelectorsFromPVStandalone +
536  tracksValidationSelectorsFromPVPt09Standalone
537 )
538 
539 # we copy this for both Standalone and TrackingOnly
540 # and later make modifications from it which change based on era
541 _trackValidatorsBase = cms.Sequence(
542  trackValidatorStandalone +
543  trackValidatorFromPVStandalone +
544  trackValidatorFromPVAllTPStandalone +
545  trackValidatorAllTPEfficStandalone +
546  trackValidatorConversionStandalone +
547  trackValidatorGsfTracks +
548  trackValidatorBHadronStandalone
549 )
550 trackValidatorsStandalone = _trackValidatorsBase.copy()
551 fastSim.toModify(trackValidatorsStandalone, lambda x: x.remove(trackValidatorConversionStandalone) )
552 
553 tracksValidationStandalone = cms.Sequence(
554  ak4PFL1FastL2L3CorrectorChain +
555  tracksPreValidationStandalone +
556  tracksValidationSelectorsStandalone +
557  trackValidatorsStandalone
558 )
559 
560 ### TrackingOnly mode (i.e. MTV with DIGI input + tracking-only reconstruction)
561 
562 # selectors
563 tracksValidationSelectorsTrackingOnly = tracksValidationSelectors.copyAndExclude([ak4JetTracksAssociatorExplicitAll,cutsRecoTracksAK4PFJets]) # selectors using track information only (i.e. no PF)
564 _sequenceForEachEra(_addSeedToTrackProducers, args=["_seedProducers"], names="_seedSelectors", sequence="_tracksValidationSeedSelectorsTrackingOnly", includeFastSim=True, modDict=globals())
565 
566 # MTV instances
567 trackValidatorTrackingOnly = trackValidatorStandalone.clone(label = [ x for x in trackValidatorStandalone.label if x != "cutsRecoTracksAK4PFJets"] )
568 
569 _trackValidatorSeedingBuildingTrackingOnly = trackValidatorTrackingOnly.clone( # common for seeds and built tracks
570  associators = ["quickTrackAssociatorByHits"],
571  UseAssociators = True,
572  dodEdxPlots = False,
573  doPVAssociationPlots = False,
574  doSimPlots = False,
575 )
576 trackValidatorBuildingTrackingOnly = _trackValidatorSeedingBuildingTrackingOnly.clone(
577  dirName = "Tracking/TrackBuilding/",
578  doMVAPlots = True,
579 )
580 for _eraName, _postfix, _era in _relevantErasAndFastSim:
581  _setForEra(trackValidatorBuildingTrackingOnly, _eraName, _era, label = locals()["_trackProducers"+_postfix])
582 fastSim.toModify(trackValidatorBuildingTrackingOnly, doMVAPlots=False)
583 for _eraName, _postfix, _era in _relevantEras:
584  _setForEra(trackValidatorBuildingTrackingOnly, _eraName, _era, mvaLabels = locals()["_mvaSelectors"+_postfix])
585 
586 trackValidatorSeedingTrackingOnly = _trackValidatorSeedingBuildingTrackingOnly.clone(
587  dirName = "Tracking/TrackSeeding/",
588  label = _seedSelectors,
589  doSeedPlots = True,
590 )
591 for _eraName, _postfix, _era in _relevantErasAndFastSim:
592  _setForEra(trackValidatorSeedingTrackingOnly, _eraName, _era, label = locals()["_seedSelectors"+_postfix])
593 
594 
595 trackValidatorConversionTrackingOnly = trackValidatorConversion.clone(label = [x for x in trackValidatorConversion.label if x not in ["ckfInOutTracksFromConversions", "ckfOutInTracksFromConversions"]])
596 
597 trackValidatorBHadronTrackingOnly = trackValidatorBHadron.clone(label = [x for x in trackValidatorTrackingOnly.label if "Pt09" not in x])
598 
599 # sequences
600 tracksPreValidationTrackingOnly = tracksPreValidationStandalone.copy()
601 tracksPreValidationTrackingOnly.replace(tracksValidationSelectors, tracksValidationSelectorsTrackingOnly)
602 
603 trackValidatorsTrackingOnly = _trackValidatorsBase.copy()
604 trackValidatorsTrackingOnly.replace(trackValidatorStandalone, trackValidatorTrackingOnly)
605 trackValidatorsTrackingOnly += (
606  trackValidatorSeedingTrackingOnly +
607  trackValidatorBuildingTrackingOnly
608 )
609 trackValidatorsTrackingOnly.replace(trackValidatorConversionStandalone, trackValidatorConversionTrackingOnly)
610 trackValidatorsTrackingOnly.remove(trackValidatorGsfTracks)
611 trackValidatorsTrackingOnly.replace(trackValidatorBHadronStandalone, trackValidatorBHadronTrackingOnly)
612 fastSim.toModify(trackValidatorsTrackingOnly, lambda x: x.remove(trackValidatorConversionTrackingOnly))
613 fastSim.toModify(trackValidatorsTrackingOnly, lambda x: x.remove(trackValidatorBHadronTrackingOnly))
614 
615 
616 tracksValidationTrackingOnly = cms.Sequence(
617  tracksPreValidationTrackingOnly +
618  tracksValidationSelectorsStandalone +
619  tracksValidationSeedSelectorsTrackingOnly +
620  trackValidatorsTrackingOnly
621 )
622 
623 
624 
625 ### Lite mode (only generalTracks and HP)
626 trackValidatorLite = trackValidator.clone(
627  label = ["generalTracks", "cutsRecoTracksHp"]
628 )
629 tracksValidationLite = cms.Sequence(
630  cutsRecoTracksHp +
631  tracksValidationTruth +
632  trackValidatorLite
633 )
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)
_addSelectorsByOriginalAlgoMask
Then define stuff for standalone mode (i.e.
def _getSeedingLayers(seedProducers)
TEveGeoShape * clone(const TEveElement *element, TEveElement *parent)
Definition: eve_macros.cc:135