CMS 3D CMS Logo

DetachedTripletStep_cff.py
Go to the documentation of this file.
1 import FWCore.ParameterSet.Config as cms
2 from Configuration.Eras.Modifier_tracker_apv_vfp30_2016_cff import tracker_apv_vfp30_2016 as _tracker_apv_vfp30_2016
3 import RecoTracker.IterativeTracking.iterativeTkConfig as _cfg
4 
5 # for fastsim
6 from Configuration.Eras.Modifier_fastSim_cff import fastSim
7 from FastSimulation.Tracking.SeedingMigration import _hitSetProducerToFactoryPSet
8 
9 # for dnn classifier
10 from Configuration.ProcessModifiers.trackdnn_cff import trackdnn
11 from RecoTracker.IterativeTracking.dnnQualityCuts import qualityCutDictionary
12 
13 # for no-loopers
14 from Configuration.ProcessModifiers.trackingNoLoopers_cff import trackingNoLoopers
15 
16 
19 
20 # REMOVE HITS ASSIGNED TO GOOD TRACKS FROM PREVIOUS ITERATIONS
21 detachedTripletStepClusters = _cfg.clusterRemoverForIter('DetachedTripletStep')
22 for _eraName, _postfix, _era in _cfg.nonDefaultEras():
23  _era.toReplaceWith(detachedTripletStepClusters, _cfg.clusterRemoverForIter('DetachedTripletStep', _eraName, _postfix))
24 
25 # SEEDING LAYERS
27 detachedTripletStepSeedLayers = RecoTracker.TkSeedingLayers.PixelLayerTriplets_cfi.PixelLayerTriplets.clone(
28  BPix = dict(skipClusters = cms.InputTag('detachedTripletStepClusters')),
29  FPix = dict(skipClusters = cms.InputTag('detachedTripletStepClusters'))
30 )
31 _phase1LayerList = [
32  'BPix1+BPix2+BPix3',
33  'BPix2+BPix3+BPix4',
34 # 'BPix1+BPix3+BPix4', # has 'hole', not tested
35 # 'BPix1+BPix2+BPix4', # has 'hole', not tested
36  'BPix2+BPix3+FPix1_pos', 'BPix2+BPix3+FPix1_neg',
37 # 'BPix1+BPix2+FPix1_pos', 'BPix1+BPix2+FPix1_neg', # mostly fake tracks, lots of seeds
38 # 'BPix1+BPix3+FPix1_pos', 'BPix1+BPix3+FPix1_neg', # has 'hole', not tested
39  'BPix2+FPix1_pos+FPix2_pos', 'BPix2+FPix1_neg+FPix2_neg',
40 # 'BPix1+FPix1_pos+FPix2_pos', 'BPix1+FPix1_neg+FPix2_neg', # mostly fake tracks, lots of seeds
41 # 'BPix1+BPix2+FPix2_pos', 'BPix1+BPix2+FPix2_neg', # has 'hole', not tested
42  'FPix1_pos+FPix2_pos+FPix3_pos', 'FPix1_neg+FPix2_neg+FPix3_neg',
43 # 'BPix1+FPix2_pos+FPix3_pos', 'BPix1+FPix2_neg+FPix3_neg', # has 'hole', not tested
44 # 'BPix1+FPix1_pos+FPix3_pos', 'BPix1+FPix1_neg+FPix3_neg' # has 'hole', not tested
45  ]
46 from Configuration.Eras.Modifier_trackingPhase1_cff import trackingPhase1
47 trackingPhase1.toModify(detachedTripletStepSeedLayers, layerList=_phase1LayerList)
48 
49 # TrackingRegion
50 from RecoTracker.TkTrackingRegions.globalTrackingRegionFromBeamSpotFixedZ_cfi import globalTrackingRegionFromBeamSpotFixedZ as _globalTrackingRegionFromBeamSpotFixedZ
51 detachedTripletStepTrackingRegions = _globalTrackingRegionFromBeamSpotFixedZ.clone(RegionPSet = dict(
52  ptMin = 0.3,
53  originHalfLength = 15.0,
54  originRadius = 1.5
55 ))
56 trackingPhase1.toModify(detachedTripletStepTrackingRegions, RegionPSet = dict(ptMin = 0.25))
57 
58 from Configuration.Eras.Modifier_pp_on_XeXe_2017_cff import pp_on_XeXe_2017
59 from Configuration.ProcessModifiers.pp_on_AA_cff import pp_on_AA
60 from RecoTracker.TkTrackingRegions.globalTrackingRegionWithVertices_cff import globalTrackingRegionWithVertices as _globalTrackingRegionWithVertices
61 (pp_on_XeXe_2017 | pp_on_AA).toReplaceWith(detachedTripletStepTrackingRegions,
62  _globalTrackingRegionWithVertices.clone(RegionPSet=dict(
63  fixedError = 2.5,
64  ptMin = 0.9,
65  originRadius = 1.5)
66  )
67 )
68 from Configuration.Eras.Modifier_highBetaStar_2018_cff import highBetaStar_2018
69 highBetaStar_2018.toModify(detachedTripletStepTrackingRegions, RegionPSet = dict(ptMin = 0.05))
70 
71 # seeding
72 from RecoTracker.TkHitPairs.hitPairEDProducer_cfi import hitPairEDProducer as _hitPairEDProducer
73 detachedTripletStepHitDoublets = _hitPairEDProducer.clone(
74  seedingLayers = 'detachedTripletStepSeedLayers',
75  trackingRegions = 'detachedTripletStepTrackingRegions',
76  maxElement = 50000000,
77  produceIntermediateHitDoublets = True,
78 )
79 from RecoTracker.PixelSeeding.pixelTripletLargeTipEDProducer_cfi import pixelTripletLargeTipEDProducer as _pixelTripletLargeTipEDProducer
81 detachedTripletStepHitTriplets = _pixelTripletLargeTipEDProducer.clone(
82  doublets = 'detachedTripletStepHitDoublets',
83  produceSeedingHitSets = True,
84 )
85 from RecoTracker.TkSeedGenerator.seedCreatorFromRegionConsecutiveHitsTripletOnlyEDProducer_cff import seedCreatorFromRegionConsecutiveHitsTripletOnlyEDProducer as _seedCreatorFromRegionConsecutiveHitsTripletOnlyEDProducer
86 detachedTripletStepSeeds = _seedCreatorFromRegionConsecutiveHitsTripletOnlyEDProducer.clone(
87  seedingHitSets = 'detachedTripletStepHitTriplets',
88  SeedComparitorPSet = dict(# FIXME: is this defined in any cfi that could be imported instead of copy-paste?
89  ComponentName = 'PixelClusterShapeSeedComparitor',
90  FilterAtHelixStage = cms.bool(False),
91  FilterPixelHits = cms.bool(True),
92  FilterStripHits = cms.bool(False),
93  ClusterShapeHitFilterName = cms.string('ClusterShapeHitFilter'),
94  ClusterShapeCacheSrc = cms.InputTag('siPixelClusterShapeCache')
95  ),
96 )
97 
98 from RecoTracker.PixelSeeding.caHitTripletEDProducer_cfi import caHitTripletEDProducer as _caHitTripletEDProducer
99 trackingPhase1.toModify(detachedTripletStepHitDoublets, layerPairs = [0,1]) # layer pairs (0,1), (1,2)
100 trackingPhase1.toReplaceWith(detachedTripletStepHitTriplets, _caHitTripletEDProducer.clone(
101  doublets = 'detachedTripletStepHitDoublets',
102  extraHitRPhitolerance = detachedTripletStepHitTriplets.extraHitRPhitolerance,
103  maxChi2 = dict(
104  pt1 = 0.8, pt2 = 2,
105  value1 = 300, value2 = 10,
106  ),
107  useBendingCorrection = True,
108  CAThetaCut = 0.001,
109  CAPhiCut = 0,
110  CAHardPtCut = 0.2,
111 ))
112 highBetaStar_2018.toModify(detachedTripletStepHitTriplets,CAThetaCut = 0.002,CAPhiCut = 0.1,CAHardPtCut = 0)
113 
114 import FastSimulation.Tracking.TrajectorySeedProducer_cfi
115 _fastSim_detachedTripletStepSeeds = FastSimulation.Tracking.TrajectorySeedProducer_cfi.trajectorySeedProducer.clone(
116  trackingRegions = 'detachedTripletStepTrackingRegions',
117  hitMasks = cms.InputTag('detachedTripletStepMasks'),
118  seedFinderSelector = dict( pixelTripletGeneratorFactory = _hitSetProducerToFactoryPSet(detachedTripletStepHitTriplets),
119  layerList = detachedTripletStepSeedLayers.layerList.value())
120 )#new for phase1
121 trackingPhase1.toModify(_fastSim_detachedTripletStepSeeds, seedFinderSelector = dict(
122  pixelTripletGeneratorFactory = None,
123  CAHitTripletGeneratorFactory = _hitSetProducerToFactoryPSet(detachedTripletStepHitTriplets),
124  #new parameters required for phase1 seeding
125  BPix = dict(
126  TTRHBuilder = 'WithoutRefit',
127  HitProducer = 'TrackingRecHitProducer',
128  ),
129  FPix = dict(
130  TTRHBuilder = 'WithoutRefit',
131  HitProducer = 'TrackingRecHitProducer',
132  ),
133  layerPairs = detachedTripletStepHitDoublets.layerPairs.value()
134  )
135 )
136 fastSim.toReplaceWith(detachedTripletStepSeeds,_fastSim_detachedTripletStepSeeds)
137 
138 # QUALITY CUTS DURING TRACK BUILDING
140 _detachedTripletStepTrajectoryFilterBase = TrackingTools.TrajectoryFiltering.TrajectoryFilter_cff.CkfBaseTrajectoryFilter_block.clone(
141 # maxLostHitsFraction = cms.double(1./10.),
142 # constantValueForLostHitsFractionFilter = cms.double(0.701),
143  minimumNumberOfHits = 3,
144  minPt = 0.075,
145 )
146 detachedTripletStepTrajectoryFilterBase = _detachedTripletStepTrajectoryFilterBase.clone(
147  maxCCCLostHits = 0,
148  minGoodStripCharge = cms.PSet(refToPSet_ = cms.string('SiStripClusterChargeCutLoose'))
149 )
150 from Configuration.Eras.Modifier_tracker_apv_vfp30_2016_cff import tracker_apv_vfp30_2016
151 _tracker_apv_vfp30_2016.toModify(detachedTripletStepTrajectoryFilterBase, maxCCCLostHits = 2)
152 from Configuration.Eras.Modifier_trackingLowPU_cff import trackingLowPU
153 trackingLowPU.toReplaceWith(detachedTripletStepTrajectoryFilterBase, _detachedTripletStepTrajectoryFilterBase.clone(
154  maxLostHitsFraction = 1./10.,
155  constantValueForLostHitsFractionFilter = 0.701,
156 ))
157 
158 (pp_on_XeXe_2017 | pp_on_AA).toModify(detachedTripletStepTrajectoryFilterBase, minPt=0.9)
159 
160 import RecoTracker.PixelLowPtUtilities.StripSubClusterShapeTrajectoryFilter_cfi
161 detachedTripletStepTrajectoryFilterShape = RecoTracker.PixelLowPtUtilities.StripSubClusterShapeTrajectoryFilter_cfi.StripSubClusterShapeTrajectoryFilterTIX12.clone()
162 detachedTripletStepTrajectoryFilter = cms.PSet(
163  ComponentType = cms.string('CompositeTrajectoryFilter'),
164  filters = cms.VPSet(
165  cms.PSet( refToPSet_ = cms.string('detachedTripletStepTrajectoryFilterBase')),
166 # cms.PSet( refToPSet_ = cms.string('detachedTripletStepTrajectoryFilterShape'))
167  ),
168 )
169 
170 
171 import RecoTracker.MeasurementDet.Chi2ChargeMeasurementEstimator_cfi
172 detachedTripletStepChi2Est = RecoTracker.MeasurementDet.Chi2ChargeMeasurementEstimator_cfi.Chi2ChargeMeasurementEstimator.clone(
173  ComponentName = 'detachedTripletStepChi2Est',
174  nSigma = 3.0,
175  MaxChi2 = 9.0,
176  clusterChargeCut = cms.PSet(refToPSet_ = cms.string('SiStripClusterChargeCutTight')),
177 )
178 _tracker_apv_vfp30_2016.toModify(detachedTripletStepChi2Est,
179  clusterChargeCut = dict(refToPSet_ = 'SiStripClusterChargeCutTiny')
180 )
181 
182 # TRACK BUILDING
184 detachedTripletStepTrajectoryBuilder = RecoTracker.CkfPattern.GroupedCkfTrajectoryBuilder_cfi.GroupedCkfTrajectoryBuilderIterativeDefault.clone(
185  trajectoryFilter = dict(refToPSet_ = 'detachedTripletStepTrajectoryFilter'),
186  maxCand = 3,
187  alwaysUseInvalidHits = True,
188  estimator = 'detachedTripletStepChi2Est',
189  maxDPhiForLooperReconstruction = 2.0,
190  maxPtForLooperReconstruction = 0.7,
191 )
192 trackingNoLoopers.toModify(detachedTripletStepTrajectoryBuilder,
193  maxPtForLooperReconstruction = 0.0)
194 trackingLowPU.toModify(detachedTripletStepTrajectoryBuilder,
195  maxCand = 2,
196  alwaysUseInvalidHits = False,
197 )
198 
199 # MAKING OF TRACK CANDIDATES
201 # Give handle for CKF for HI
202 _detachedTripletStepTrackCandidatesCkf = RecoTracker.CkfPattern.CkfTrackCandidates_cfi.ckfTrackCandidatesIterativeDefault.clone(
203  src = 'detachedTripletStepSeeds',
204  clustersToSkip = 'detachedTripletStepClusters',
205 
206  numHitsForSeedCleaner = 50,
207  onlyPixelHitsForSeedCleaner = True,
208  TrajectoryBuilderPSet = dict(refToPSet_ = 'detachedTripletStepTrajectoryBuilder'),
209  doSeedingRegionRebuilding = True,
210  useHitsSplitting = True,
211 )
212 detachedTripletStepTrackCandidates = _detachedTripletStepTrackCandidatesCkf.clone()
213 
214 from TrackingTools.TrajectoryCleaning.TrajectoryCleanerBySharedHits_cfi import trajectoryCleanerBySharedHits
215 detachedTripletStepTrajectoryCleanerBySharedHits = trajectoryCleanerBySharedHits.clone(
216  ComponentName = 'detachedTripletStepTrajectoryCleanerBySharedHits',
217  fractionShared = 0.13,
218  allowSharedFirstHit = True
219 )
220 detachedTripletStepTrackCandidates.TrajectoryCleaner = 'detachedTripletStepTrajectoryCleanerBySharedHits'
221 trackingLowPU.toModify(detachedTripletStepTrajectoryCleanerBySharedHits, fractionShared = 0.19)
222 
223 from Configuration.ProcessModifiers.trackingMkFitDetachedTripletStep_cff import trackingMkFitDetachedTripletStep
224 import RecoTracker.MkFit.mkFitSeedConverter_cfi as mkFitSeedConverter_cfi
225 import RecoTracker.MkFit.mkFitIterationConfigESProducer_cfi as mkFitIterationConfigESProducer_cfi
226 import RecoTracker.MkFit.mkFitProducer_cfi as mkFitProducer_cfi
227 import RecoTracker.MkFit.mkFitOutputConverter_cfi as mkFitOutputConverter_cfi
228 detachedTripletStepTrackCandidatesMkFitSeeds = mkFitSeedConverter_cfi.mkFitSeedConverter.clone(
229  seeds = 'detachedTripletStepSeeds',
230 )
231 detachedTripletStepTrackCandidatesMkFitConfig = mkFitIterationConfigESProducer_cfi.mkFitIterationConfigESProducer.clone(
232  ComponentName = 'detachedTripletStepTrackCandidatesMkFitConfig',
233  config = 'RecoTracker/MkFit/data/mkfit-phase1-detachedTripletStep.json',
234 )
235 detachedTripletStepTrackCandidatesMkFit = mkFitProducer_cfi.mkFitProducer.clone(
236  seeds = 'detachedTripletStepTrackCandidatesMkFitSeeds',
237  config = ('', 'detachedTripletStepTrackCandidatesMkFitConfig'),
238  clustersToSkip = 'detachedTripletStepClusters',
239 )
240 trackingMkFitDetachedTripletStep.toReplaceWith(detachedTripletStepTrackCandidates, mkFitOutputConverter_cfi.mkFitOutputConverter.clone(
241  seeds = 'detachedTripletStepSeeds',
242  mkFitSeeds = 'detachedTripletStepTrackCandidatesMkFitSeeds',
243  tracks = 'detachedTripletStepTrackCandidatesMkFit',
244 ))
245 (pp_on_XeXe_2017 | pp_on_AA).toModify(detachedTripletStepTrackCandidatesMkFitConfig, minPt=0.9)
246 
247 import FastSimulation.Tracking.TrackCandidateProducer_cfi
248 _fastSim_detachedTripletStepTrackCandidates = FastSimulation.Tracking.TrackCandidateProducer_cfi.trackCandidateProducer.clone(
249  src = 'detachedTripletStepSeeds',
250  MinNumberOfCrossedLayers = 3,
251  hitMasks = cms.InputTag('detachedTripletStepMasks')
252 )
253 fastSim.toReplaceWith(detachedTripletStepTrackCandidates,_fastSim_detachedTripletStepTrackCandidates)
254 
255 
256 # TRACK FITTING
257 import RecoTracker.TrackProducer.TrackProducerIterativeDefault_cfi
258 detachedTripletStepTracks = RecoTracker.TrackProducer.TrackProducerIterativeDefault_cfi.TrackProducerIterativeDefault.clone(
259  AlgorithmName = 'detachedTripletStep',
260  src = 'detachedTripletStepTrackCandidates',
261  Fitter = 'FlexibleKFFittingSmoother'
262 )
263 fastSim.toModify(detachedTripletStepTracks,TTRHBuilder = 'WithoutRefit')
264 
265 # TRACK SELECTION AND QUALITY FLAG SETTING.
268 detachedTripletStepClassifier1 = TrackMVAClassifierDetached.clone(
269  src = 'detachedTripletStepTracks',
270  mva = dict(GBRForestLabel = 'MVASelectorIter3_13TeV'),
271  qualityCuts = [-0.5,0.0,0.5]
272 )
273 fastSim.toModify(detachedTripletStepClassifier1,vertices = 'firstStepPrimaryVerticesBeforeMixing')
274 
275 detachedTripletStepClassifier2 = TrackMVAClassifierPrompt.clone(
276  src = 'detachedTripletStepTracks',
277  mva = dict(GBRForestLabel = 'MVASelectorIter0_13TeV'),
278  qualityCuts = [-0.2,0.0,0.4]
279 )
280 fastSim.toModify(detachedTripletStepClassifier2,vertices = 'firstStepPrimaryVerticesBeforeMixing')
281 
283 detachedTripletStep = ClassifierMerger.clone(
284  inputClassifiers=['detachedTripletStepClassifier1','detachedTripletStepClassifier2']
285 )
286 trackingPhase1.toReplaceWith(detachedTripletStep, detachedTripletStepClassifier1.clone(
287  mva = dict(GBRForestLabel = 'MVASelectorDetachedTripletStep_Phase1'),
288  qualityCuts = [-0.2,0.3,0.8]
289 ))
290 pp_on_AA.toModify(detachedTripletStep,
291  mva = dict(GBRForestLabel = 'HIMVASelectorDetachedTripletStep_Phase1'),
292  qualityCuts = [-0.2, 0.4, 0.85],
293 )
294 
298 trackdnn.toReplaceWith(detachedTripletStep, trackTfClassifier.clone(
299  src = 'detachedTripletStepTracks',
300  qualityCuts = qualityCutDictionary.DetachedTripletStep.value()
301 ))
302 (trackdnn & fastSim).toModify(detachedTripletStep,vertices = 'firstStepPrimaryVerticesBeforeMixing')
303 
304 (pp_on_AA & trackdnn).toModify(detachedTripletStep, qualityCuts = [-0.32, 0.5, 0.98] )
305 
306 highBetaStar_2018.toModify(detachedTripletStep,qualityCuts = [-0.5,0.0,0.5])
307 
308 # For LowPU
309 import RecoTracker.FinalTrackSelectors.multiTrackSelector_cfi
310 detachedTripletStepSelector = RecoTracker.FinalTrackSelectors.multiTrackSelector_cfi.multiTrackSelector.clone(
311  src = 'detachedTripletStepTracks',
312  useAnyMVA = cms.bool(False),
313  GBRForestLabel = cms.string('MVASelectorIter3'),
314  trackSelectors = [
315  RecoTracker.FinalTrackSelectors.multiTrackSelector_cfi.looseMTS.clone(
316  name = 'detachedTripletStepVtxLoose',
317  chi2n_par = 1.6,
318  res_par = ( 0.003, 0.001 ),
319  minNumberLayers = 3,
320  d0_par1 = ( 1.2, 3.0 ),
321  dz_par1 = ( 1.2, 3.0 ),
322  d0_par2 = ( 1.3, 3.0 ),
323  dz_par2 = ( 1.3, 3.0 )
324  ),
325  RecoTracker.FinalTrackSelectors.multiTrackSelector_cfi.looseMTS.clone(
326  name = 'detachedTripletStepTrkLoose',
327  chi2n_par = 0.7,
328  res_par = ( 0.003, 0.001 ),
329  minNumberLayers = 3,
330  d0_par1 = ( 1.6, 4.0 ),
331  dz_par1 = ( 1.6, 4.0 ),
332  d0_par2 = ( 1.6, 4.0 ),
333  dz_par2 = ( 1.6, 4.0 )
334  ),
335  RecoTracker.FinalTrackSelectors.multiTrackSelector_cfi.tightMTS.clone(
336  name = 'detachedTripletStepVtxTight',
337  preFilterName = 'detachedTripletStepVtxLoose',
338  chi2n_par = 0.7,
339  res_par = ( 0.003, 0.001 ),
340  minNumberLayers = 3,
341  maxNumberLostLayers = 1,
342  minNumber3DLayers = 3,
343  d0_par1 = ( 0.95, 3.0 ),
344  dz_par1 = ( 0.9, 3.0 ),
345  d0_par2 = ( 1.0, 3.0 ),
346  dz_par2 = ( 1.0, 3.0 )
347  ),
348  RecoTracker.FinalTrackSelectors.multiTrackSelector_cfi.tightMTS.clone(
349  name = 'detachedTripletStepTrkTight',
350  preFilterName = 'detachedTripletStepTrkLoose',
351  chi2n_par = 0.5,
352  res_par = ( 0.003, 0.001 ),
353  minNumberLayers = 5,
354  maxNumberLostLayers = 1,
355  minNumber3DLayers = 3,
356  d0_par1 = ( 1.1, 4.0 ),
357  dz_par1 = ( 1.1, 4.0 ),
358  d0_par2 = ( 1.1, 4.0 ),
359  dz_par2 = ( 1.1, 4.0 )
360  ),
361  RecoTracker.FinalTrackSelectors.multiTrackSelector_cfi.highpurityMTS.clone(
362  name = 'detachedTripletStepVtx',
363  preFilterName = 'detachedTripletStepVtxTight',
364  chi2n_par = 0.7,
365  res_par = ( 0.003, 0.001 ),
366  minNumberLayers = 3,
367  maxNumberLostLayers = 1,
368  minNumber3DLayers = 3,
369  d0_par1 = ( 0.85, 3.0 ),
370  dz_par1 = ( 0.8, 3.0 ),
371  d0_par2 = ( 0.9, 3.0 ),
372  dz_par2 = ( 0.9, 3.0 )
373  ),
374  RecoTracker.FinalTrackSelectors.multiTrackSelector_cfi.highpurityMTS.clone(
375  name = 'detachedTripletStepTrk',
376  preFilterName = 'detachedTripletStepTrkTight',
377  chi2n_par = 0.4,
378  res_par = ( 0.003, 0.001 ),
379  minNumberLayers = 5,
380  maxNumberLostLayers = 1,
381  minNumber3DLayers = 4,
382  d0_par1 = ( 1.0, 4.0 ),
383  dz_par1 = ( 1.0, 4.0 ),
384  d0_par2 = ( 1.0, 4.0 ),
385  dz_par2 = ( 1.0, 4.0 )
386  )
387  ] #end of vpset
388 ) #end of clone
389 
390 from RecoTracker.FinalTrackSelectors.trackAlgoPriorityOrder_cfi import trackAlgoPriorityOrder
392 trackingLowPU.toReplaceWith(detachedTripletStep, RecoTracker.FinalTrackSelectors.trackListMerger_cfi.trackListMerger.clone(
393  TrackProducers = ['detachedTripletStepTracks',
394  'detachedTripletStepTracks'],
395  hasSelector = [1,1],
396  selectedTrackQuals = ['detachedTripletStepSelector:detachedTripletStepVtx',
397  'detachedTripletStepSelector:detachedTripletStepTrk'],
398  setsToMerge = [cms.PSet( tLists=cms.vint32(0,1), pQual=cms.bool(True) )],
399  writeOnlyTrkQuals = True
400 ))
401 
402 DetachedTripletStepTask = cms.Task(detachedTripletStepClusters,
403  detachedTripletStepSeedLayers,
404  detachedTripletStepTrackingRegions,
405  detachedTripletStepHitDoublets,
406  detachedTripletStepHitTriplets,
407  detachedTripletStepSeeds,
408  detachedTripletStepTrackCandidates,
409  detachedTripletStepTracks,
410  detachedTripletStepClassifier1,detachedTripletStepClassifier2,
411  detachedTripletStep)
412 DetachedTripletStep = cms.Sequence(DetachedTripletStepTask)
413 
414 _DetachedTripletStepTask_trackingMkFit = DetachedTripletStepTask.copy()
415 _DetachedTripletStepTask_trackingMkFit.add(detachedTripletStepTrackCandidatesMkFitSeeds, detachedTripletStepTrackCandidatesMkFit, detachedTripletStepTrackCandidatesMkFitConfig)
416 trackingMkFitDetachedTripletStep.toReplaceWith(DetachedTripletStepTask, _DetachedTripletStepTask_trackingMkFit)
417 
418 _DetachedTripletStepTask_LowPU = DetachedTripletStepTask.copyAndExclude([detachedTripletStepClassifier2])
419 _DetachedTripletStepTask_LowPU.replace(detachedTripletStepClassifier1, detachedTripletStepSelector)
420 trackingLowPU.toReplaceWith(DetachedTripletStepTask, _DetachedTripletStepTask_LowPU)
421 
422 # fast tracking mask producer
423 from FastSimulation.Tracking.FastTrackerRecHitMaskProducer_cfi import maskProducerFromClusterRemover
424 detachedTripletStepMasks = maskProducerFromClusterRemover(detachedTripletStepClusters)
425 fastSim.toReplaceWith(DetachedTripletStepTask,
426  cms.Task(detachedTripletStepMasks
427  ,detachedTripletStepTrackingRegions
428  ,detachedTripletStepSeeds
429  ,detachedTripletStepTrackCandidates
430  ,detachedTripletStepTracks
431  ,detachedTripletStepClassifier1,detachedTripletStepClassifier2
432  ,detachedTripletStep
433  ) )
TRIGGER SELECTION #####.
def _hitSetProducerToFactoryPSet(producer)