CMS 3D CMS Logo

DisplacedRegionalStep_cff.py
Go to the documentation of this file.
1 import FWCore.ParameterSet.Config as cms
2 import RecoTracker.IterativeTracking.iterativeTkConfig as _cfg
3 
4 #for dnn classifier
5 from Configuration.ProcessModifiers.trackdnn_cff import trackdnn
6 from RecoTracker.IterativeTracking.dnnQualityCuts import qualityCutDictionary
7 
8 # for no-loopers
9 from Configuration.ProcessModifiers.trackingNoLoopers_cff import trackingNoLoopers
10 
11 
14 
15 displacedRegionalStepClusters = _cfg.clusterRemoverForIter("DisplacedRegionalStep")
16 for _eraName, _postfix, _era in _cfg.nonDefaultEras():
17  _era.toReplaceWith(displacedRegionalStepClusters, _cfg.clusterRemoverForIter("DisplacedRegionalStep", _eraName, _postfix))
18 
19 # TRIPLET SEEDING LAYERS
22 
23 displacedRegionalStepSeedLayersTripl = _mod.seedingLayersEDProducer.clone(
24  layerList = [
25  #TOB
26  'TOB1+TOB2+MTOB3',#'TOB1+TOB2+MTOB4',
27  #TOB+MTEC
28  'TOB1+TOB2+MTEC1_pos','TOB1+TOB2+MTEC1_neg',
29  ],
30  TOB = dict(
31  TTRHBuilder = cms.string('WithTrackAngle'),
32  clusterChargeCut = cms.PSet(refToPSet_ = cms.string('SiStripClusterChargeCutTight')),
33  matchedRecHits = cms.InputTag('siStripMatchedRecHits','matchedRecHit'),
34  skipClusters = cms.InputTag('displacedRegionalStepClusters')
35  ),
36  MTOB = dict(
37  TTRHBuilder = cms.string('WithTrackAngle'),
38  clusterChargeCut = cms.PSet(refToPSet_ = cms.string('SiStripClusterChargeCutTight')),
39  skipClusters = cms.InputTag('displacedRegionalStepClusters'),
40  rphiRecHits = cms.InputTag('siStripMatchedRecHits','rphiRecHit')
41  ),
42  MTEC = dict(
43  rphiRecHits = cms.InputTag('siStripMatchedRecHits','rphiRecHit'),
44  skipClusters = cms.InputTag('displacedRegionalStepClusters'),
45  useRingSlector = cms.bool(True),
46  TTRHBuilder = cms.string('WithTrackAngle'),
47  clusterChargeCut = cms.PSet(refToPSet_ = cms.string('SiStripClusterChargeCutTight')),
48  minRing = cms.int32(6),
49  maxRing = cms.int32(7)
50  )
51 )
52 
53 # Triplet TrackingRegion
54 from RecoTracker.FinalTrackSelectors.displacedRegionalStepInputTracks_cfi import displacedRegionalStepInputTracks
55 from RecoVertex.V0Producer.generalV0Candidates_cfi import generalV0Candidates as _generalV0Candidates
56 from RecoTracker.DisplacedRegionalTracking.displacedRegionProducer_cfi import displacedRegionProducer as _displacedRegionProducer
57 displacedRegionalStepSeedingV0Candidates = _generalV0Candidates.clone(
58  trackRecoAlgorithm = "displacedRegionalStepInputTracks",
59  doLambdas = False,
60  doFit = False,
61  useRefTracks = False,
62  vtxDecayXYCut = 1.,
63  ssVtxDecayXYCut = 5.,
64  allowSS = True,
65  innerTkDCACut = 0.2,
66  allowWideAngleVtx = True,
67  mPiPiCut = 13000.,
68  cosThetaXYCut = 0.,
69  kShortMassCut = 13000.,
70 )
71 displacedRegionalStepSeedingVertices = _displacedRegionProducer.clone(
72  minRadius = 2.0,
73  nearThreshold = 20.0,
74  farThreshold = -1.0,
75  discriminatorCut = 0.5,
76  trackClusters = ["displacedRegionalStepSeedingV0Candidates", "Kshort"],
77 )
78 
79 from RecoTracker.TkTrackingRegions.globalTrackingRegionWithVertices_cfi import globalTrackingRegionWithVertices as _globalTrackingRegionWithVertices
80 displacedRegionalStepFarTrackingRegions = _globalTrackingRegionWithVertices.clone(RegionPSet = dict(
81  originRadius = 1.0,
82  fixedError = 1.0,
83  VertexCollection = ["displacedRegionalStepSeedingVertices", "farRegionsOfInterest"],
84  useFakeVertices = True,
85  ptMin = 0.6,
86  allowEmpty = True
87 ))
88 
89 # Triplet seeding
90 from RecoTracker.PixelLowPtUtilities.ClusterShapeHitFilterESProducer_cfi import ClusterShapeHitFilterESProducer as _ClusterShapeHitFilterESProducer
91 displacedRegionalStepClusterShapeHitFilter = _ClusterShapeHitFilterESProducer.clone(
92  ComponentName = 'displacedRegionalStepClusterShapeHitFilter',
93  doStripShapeCut = False,
94  clusterChargeCut = dict(refToPSet_ = 'SiStripClusterChargeCutTight')
95 )
96 
97 from RecoTracker.TkHitPairs.hitPairEDProducer_cfi import hitPairEDProducer as _hitPairEDProducer
98 displacedRegionalStepHitDoubletsTripl = _hitPairEDProducer.clone(
99  seedingLayers = "displacedRegionalStepSeedLayersTripl",
100  trackingRegions = "displacedRegionalStepFarTrackingRegions",
101  maxElement = 50000000,
102  produceIntermediateHitDoublets = True,
103 )
104 from RecoTracker.TkSeedGenerator.multiHitFromChi2EDProducer_cfi import multiHitFromChi2EDProducer as _multiHitFromChi2EDProducer
105 displacedRegionalStepHitTripletsTripl = _multiHitFromChi2EDProducer.clone(
106  doublets = "displacedRegionalStepHitDoubletsTripl",
107  extraPhiKDBox = 0.01,
108 )
109 from RecoTracker.TkSeedGenerator.seedCreatorFromRegionConsecutiveHitsEDProducer_cff import seedCreatorFromRegionConsecutiveHitsEDProducer as _seedCreatorFromRegionConsecutiveHitsTripletOnlyEDProducer
110 _displacedRegionalStepSeedComparitorPSet = dict(
111  ComponentName = 'CombinedSeedComparitor',
112  mode = cms.string("and"),
113  comparitors = cms.VPSet(
114  cms.PSet(# FIXME: is this defined in any cfi that could be imported instead of copy-paste?
115  ComponentName = cms.string('PixelClusterShapeSeedComparitor'),
116  FilterAtHelixStage = cms.bool(True),
117  FilterPixelHits = cms.bool(False),
118  FilterStripHits = cms.bool(True),
119  ClusterShapeHitFilterName = cms.string('displacedRegionalStepClusterShapeHitFilter'),
120  ClusterShapeCacheSrc = cms.InputTag("siPixelClusterShapeCache") # not really needed here since FilterPixelHits=False
121  )
122  )
123 )
124 
125 displacedRegionalStepSeedsTripl = _seedCreatorFromRegionConsecutiveHitsTripletOnlyEDProducer.clone(#empirically better than 'SeedFromConsecutiveHitsTripletOnlyCreator'
126  seedingHitSets = "displacedRegionalStepHitTripletsTripl",
127  SeedComparitorPSet = _displacedRegionalStepSeedComparitorPSet,
128 )
129 
130 from RecoTracker.PixelLowPtUtilities.StripSubClusterShapeSeedFilter_cfi import StripSubClusterShapeSeedFilter as _StripSubClusterShapeSeedFilter
131 from Configuration.ProcessModifiers.approxSiStripClusters_cff import approxSiStripClusters
132 (~approxSiStripClusters).toModify(displacedRegionalStepSeedsTripl.SeedComparitorPSet.comparitors, func = lambda list: list.append(_StripSubClusterShapeSeedFilter.clone()) )
133 
134 # PAIR SEEDING LAYERS
135 displacedRegionalStepSeedLayersPair = _mod.seedingLayersEDProducer.clone(
136  layerList = ['TOB1+TEC1_pos','TOB1+TEC1_neg',
137  'TEC1_pos+TEC2_pos','TEC1_neg+TEC2_neg',
138  'TEC2_pos+TEC3_pos','TEC2_neg+TEC3_neg',
139  'TEC3_pos+TEC4_pos','TEC3_neg+TEC4_neg',
140  'TEC4_pos+TEC5_pos','TEC4_neg+TEC5_neg',
141  'TEC5_pos+TEC6_pos','TEC5_neg+TEC6_neg',
142  'TEC6_pos+TEC7_pos','TEC6_neg+TEC7_neg'],
143  TOB = dict(
144  TTRHBuilder = cms.string('WithTrackAngle'),
145  clusterChargeCut = cms.PSet(refToPSet_ = cms.string('SiStripClusterChargeCutTight')),
146  matchedRecHits = cms.InputTag('siStripMatchedRecHits','matchedRecHit'),
147  skipClusters = cms.InputTag('displacedRegionalStepClusters')
148  ),
149  TEC = dict(
150  matchedRecHits = cms.InputTag('siStripMatchedRecHits','matchedRecHit'),
151  skipClusters = cms.InputTag('displacedRegionalStepClusters'),
152  useRingSlector = cms.bool(True),
153  TTRHBuilder = cms.string('WithTrackAngle'),
154  clusterChargeCut = cms.PSet(refToPSet_ = cms.string('SiStripClusterChargeCutTight')),
155  minRing = cms.int32(5),
156  maxRing = cms.int32(5)
157  )
158 )
159 
160 
161 # Pair seeds
162 displacedRegionalStepHitDoubletsPair = _hitPairEDProducer.clone(
163  seedingLayers = "displacedRegionalStepSeedLayersPair",
164  trackingRegions = "displacedRegionalStepFarTrackingRegions",
165  produceSeedingHitSets = True,
166  maxElementTotal = 12000000,
167 )
168 from RecoTracker.TkSeedGenerator.seedCreatorFromRegionConsecutiveHitsEDProducer_cff import seedCreatorFromRegionConsecutiveHitsEDProducer as _seedCreatorFromRegionConsecutiveHitsEDProducer
169 displacedRegionalStepSeedsPair = _seedCreatorFromRegionConsecutiveHitsEDProducer.clone(
170  seedingHitSets = "displacedRegionalStepHitDoubletsPair",
171  SeedComparitorPSet = _displacedRegionalStepSeedComparitorPSet,
172 )
173 
174 # FROM PIXELLESS
175 displacedRegionalStepPLSeedLayersTripl = _mod.seedingLayersEDProducer.clone(
176  layerList = [
177  #TIB
178  'TIB1+TIB2+MTIB3',#'TIB1+TIB2+MTIB4',
179  #TIB+TID
180  'TIB1+TIB2+MTID1_pos','TIB1+TIB2+MTID1_neg',
181  #TID
182  'TID1_pos+TID2_pos+TID3_pos','TID1_neg+TID2_neg+TID3_neg',#ring 1-2 (matched)
183  'TID1_pos+TID2_pos+MTID3_pos','TID1_neg+TID2_neg+MTID3_neg',#ring 3 (mono)
184  'TID1_pos+TID2_pos+MTEC1_pos','TID1_neg+TID2_neg+MTEC1_neg',
185  #TID+TEC RING 1-3
186  'TID2_pos+TID3_pos+TEC1_pos','TID2_neg+TID3_neg+TEC1_neg',#ring 1-2 (matched)
187  'TID2_pos+TID3_pos+MTEC1_pos','TID2_neg+TID3_neg+MTEC1_neg',#ring 3 (mono)
188  #TEC RING 1-3
189  #'TEC1_pos+TEC2_pos+TEC3_pos', 'TEC1_neg+TEC2_neg+TEC3_neg',
190  #'TEC1_pos+TEC2_pos+MTEC3_pos','TEC1_neg+TEC2_neg+MTEC3_neg',
191  #'TEC1_pos+TEC2_pos+TEC4_pos', 'TEC1_neg+TEC2_neg+TEC4_neg',
192  #'TEC1_pos+TEC2_pos+MTEC4_pos','TEC1_neg+TEC2_neg+MTEC4_neg',
193  #'TEC2_pos+TEC3_pos+TEC4_pos', 'TEC2_neg+TEC3_neg+TEC4_neg',
194  #'TEC2_pos+TEC3_pos+MTEC4_pos','TEC2_neg+TEC3_neg+MTEC4_neg',
195  #'TEC2_pos+TEC3_pos+TEC5_pos', 'TEC2_neg+TEC3_neg+TEC5_neg',
196  #'TEC2_pos+TEC3_pos+TEC6_pos', 'TEC2_neg+TEC3_neg+TEC6_neg',
197  #'TEC3_pos+TEC4_pos+TEC5_pos', 'TEC3_neg+TEC4_neg+TEC5_neg',
198  #'TEC3_pos+TEC4_pos+MTEC5_pos','TEC3_neg+TEC4_neg+MTEC5_neg',
199  #'TEC3_pos+TEC5_pos+TEC6_pos', 'TEC3_neg+TEC5_neg+TEC6_neg',
200  #'TEC4_pos+TEC5_pos+TEC6_pos', 'TEC4_neg+TEC5_neg+TEC6_neg'
201  ],
202  TIB = dict(
203  TTRHBuilder = cms.string('WithTrackAngle'),
204  clusterChargeCut = cms.PSet(refToPSet_ = cms.string('SiStripClusterChargeCutTight')),
205  matchedRecHits = cms.InputTag('siStripMatchedRecHits','matchedRecHit'),
206  skipClusters = cms.InputTag('displacedRegionalStepClusters')
207  ),
208  MTIB = dict(
209  TTRHBuilder = cms.string('WithTrackAngle'),
210  clusterChargeCut = cms.PSet(refToPSet_ = cms.string('SiStripClusterChargeCutTight')),
211  skipClusters = cms.InputTag('displacedRegionalStepClusters'),
212  rphiRecHits = cms.InputTag('siStripMatchedRecHits','rphiRecHit')
213  ),
214  TID = dict(
215  matchedRecHits = cms.InputTag('siStripMatchedRecHits','matchedRecHit'),
216  skipClusters = cms.InputTag('displacedRegionalStepClusters'),
217  useRingSlector = cms.bool(True),
218  TTRHBuilder = cms.string('WithTrackAngle'),
219  clusterChargeCut = cms.PSet(refToPSet_ = cms.string('SiStripClusterChargeCutTight')),
220  minRing = cms.int32(1),
221  maxRing = cms.int32(2)
222  ),
223  MTID = dict(
224  rphiRecHits = cms.InputTag('siStripMatchedRecHits','rphiRecHit'),
225  skipClusters = cms.InputTag('displacedRegionalStepClusters'),
226  useRingSlector = cms.bool(True),
227  TTRHBuilder = cms.string('WithTrackAngle'),
228  clusterChargeCut = cms.PSet(refToPSet_ = cms.string('SiStripClusterChargeCutTight')),
229  minRing = cms.int32(3),
230  maxRing = cms.int32(3)
231  ),
232  TEC = dict(
233  matchedRecHits = cms.InputTag('siStripMatchedRecHits','matchedRecHit'),
234  skipClusters = cms.InputTag('displacedRegionalStepClusters'),
235  useRingSlector = cms.bool(True),
236  TTRHBuilder = cms.string('WithTrackAngle'),
237  clusterChargeCut = cms.PSet(refToPSet_ = cms.string('SiStripClusterChargeCutTight')),
238  minRing = cms.int32(1),
239  maxRing = cms.int32(2)
240  ),
241  MTEC = dict(
242  rphiRecHits = cms.InputTag('siStripMatchedRecHits','rphiRecHit'),
243  skipClusters = cms.InputTag('displacedRegionalStepClusters'),
244  useRingSlector = cms.bool(True),
245  TTRHBuilder = cms.string('WithTrackAngle'),
246  clusterChargeCut = cms.PSet(refToPSet_ = cms.string('SiStripClusterChargeCutTight')),
247  minRing = cms.int32(3),
248  maxRing = cms.int32(3)
249  )
250 )
251 
252 # Pair TrackingRegion
253 displacedRegionalStepNearTrackingRegions = _globalTrackingRegionWithVertices.clone(RegionPSet = dict(
254  originRadius = 1.0,
255  fixedError = 1.0,
256  VertexCollection = ["displacedRegionalStepSeedingVertices", "nearRegionsOfInterest"],
257  useFakeVertices = True,
258  ptMin = 0.6,
259  allowEmpty = True
260 ))
261 
262 displacedRegionalStepPLHitDoubletsTripl = _hitPairEDProducer.clone(
263  seedingLayers = 'displacedRegionalStepPLSeedLayersTripl',
264  trackingRegions = 'displacedRegionalStepNearTrackingRegions',
265  maxElement = 50000000,
266  produceIntermediateHitDoublets = True,
267 )
268 
269 displacedRegionalStepPLHitTripletsTripl = _multiHitFromChi2EDProducer.clone(
270  doublets = 'displacedRegionalStepPLHitDoubletsTripl',
271 )
272 
273 displacedRegionalStepPLSeedsTripl = _seedCreatorFromRegionConsecutiveHitsTripletOnlyEDProducer.clone(
274  seedingHitSets = 'displacedRegionalStepPLHitTripletsTripl',
275  SeedComparitorPSet = _displacedRegionalStepSeedComparitorPSet,
276 )
277 
278 # Combined seeds
280 displacedRegionalStepSeeds = RecoTracker.TkSeedGenerator.GlobalCombinedSeeds_cfi.globalCombinedSeeds.clone(
281  seedCollections = ['displacedRegionalStepSeedsTripl', 'displacedRegionalStepSeedsPair', 'displacedRegionalStepPLSeedsTripl']
282 )
283 
284 # QUALITY CUTS DURING TRACK BUILDING (for inwardss and outwards track building steps)
286 _displacedRegionalStepTrajectoryFilterBase = TrackingTools.TrajectoryFiltering.TrajectoryFilter_cff.CkfBaseTrajectoryFilter_block.clone(
287  maxLostHits = 0,
288  minimumNumberOfHits = 5,
289  minPt = 0.1,
290  minHitsMinPt = 3
291  )
292 displacedRegionalStepTrajectoryFilter = _displacedRegionalStepTrajectoryFilterBase.clone(
293  seedPairPenalty = 1,
294 )
295 
296 displacedRegionalStepInOutTrajectoryFilter = displacedRegionalStepTrajectoryFilter.clone(
297  minimumNumberOfHits = 4,
298 )
299 
300 import RecoTracker.MeasurementDet.Chi2ChargeMeasurementEstimator_cfi
301 displacedRegionalStepChi2Est = RecoTracker.MeasurementDet.Chi2ChargeMeasurementEstimator_cfi.Chi2ChargeMeasurementEstimator.clone(
302  ComponentName = 'displacedRegionalStepChi2Est',
303  nSigma = 3.0,
304  MaxChi2 = 16.0,
305  clusterChargeCut = cms.PSet(refToPSet_ = cms.string('SiStripClusterChargeCutTight'))
306 )
307 
308 # TRACK BUILDING
310 displacedRegionalStepTrajectoryBuilder = RecoTracker.CkfPattern.GroupedCkfTrajectoryBuilder_cfi.GroupedCkfTrajectoryBuilderIterativeDefault.clone(
311  trajectoryFilter = dict(refToPSet_ = 'displacedRegionalStepTrajectoryFilter'),
312  inOutTrajectoryFilter = dict(refToPSet_ = 'displacedRegionalStepInOutTrajectoryFilter'),
313  useSameTrajFilter = False,
314  minNrOfHitsForRebuild = 4,
315  alwaysUseInvalidHits = False,
316  maxCand = 2,
317  estimator = 'displacedRegionalStepChi2Est',
318  #startSeedHitsInRebuild = True,
319  maxDPhiForLooperReconstruction = 2.0,
320  maxPtForLooperReconstruction = 0.7,
321 )
322 trackingNoLoopers.toModify(displacedRegionalStepTrajectoryBuilder,
323  maxPtForLooperReconstruction = 0.0)
324 
325 # MAKING OF TRACK CANDIDATES
327 # Give handle for CKF for HI
328 _displacedRegionalStepTrackCandidatesCkf = RecoTracker.CkfPattern.CkfTrackCandidates_cfi.ckfTrackCandidatesIterativeDefault.clone(
329  src = 'displacedRegionalStepSeeds',
330  clustersToSkip = 'displacedRegionalStepClusters',
331 
332  numHitsForSeedCleaner = 50,
333  onlyPixelHitsForSeedCleaner = False,
334  TrajectoryBuilderPSet = dict(refToPSet_ = 'displacedRegionalStepTrajectoryBuilder'),
335  doSeedingRegionRebuilding = True,
336  useHitsSplitting = True,
337  cleanTrajectoryAfterInOut = True,
338  TrajectoryCleaner = 'displacedRegionalStepTrajectoryCleanerBySharedHits',
339 )
340 displacedRegionalStepTrackCandidates = _displacedRegionalStepTrackCandidatesCkf.clone()
341 
342 from Configuration.ProcessModifiers.trackingMkFitDisplacedRegionalStep_cff import trackingMkFitDisplacedRegionalStep
343 import RecoTracker.MkFit.mkFitSeedConverter_cfi as mkFitSeedConverter_cfi
344 import RecoTracker.MkFit.mkFitIterationConfigESProducer_cfi as mkFitIterationConfigESProducer_cfi
345 import RecoTracker.MkFit.mkFitProducer_cfi as mkFitProducer_cfi
346 import RecoTracker.MkFit.mkFitOutputConverter_cfi as mkFitOutputConverter_cfi
347 displacedRegionalStepTrackCandidatesMkFitSeeds = mkFitSeedConverter_cfi.mkFitSeedConverter.clone(
348  seeds = 'displacedRegionalStepSeeds',
349 )
350 displacedRegionalStepTrackCandidatesMkFitConfig = mkFitIterationConfigESProducer_cfi.mkFitIterationConfigESProducer.clone(
351  ComponentName = 'displacedRegionalStepTrackCandidatesMkFitConfig',
352  config = 'RecoTracker/MkFit/data/mkfit-phase1-tobTecStep.json',
353 )
354 displacedRegionalStepTrackCandidatesMkFit = mkFitProducer_cfi.mkFitProducer.clone(
355  seeds = 'displacedRegionalStepTrackCandidatesMkFitSeeds',
356  config = ('', 'displacedRegionalStepTrackCandidatesMkFitConfig'),
357  clustersToSkip = 'displacedRegionalStepClusters',
358 )
359 trackingMkFitDisplacedRegionalStep.toReplaceWith(displacedRegionalStepTrackCandidates, mkFitOutputConverter_cfi.mkFitOutputConverter.clone(
360  seeds = 'displacedRegionalStepSeeds',
361  mkFitSeeds = 'displacedRegionalStepTrackCandidatesMkFitSeeds',
362  tracks = 'displacedRegionalStepTrackCandidatesMkFit',
363 ))
364 
365 from TrackingTools.TrajectoryCleaning.TrajectoryCleanerBySharedHits_cfi import trajectoryCleanerBySharedHits
366 displacedRegionalStepTrajectoryCleanerBySharedHits = trajectoryCleanerBySharedHits.clone(
367  ComponentName = 'displacedRegionalStepTrajectoryCleanerBySharedHits',
368  fractionShared = 0.09,
369  allowSharedFirstHit = True
370  )
371 
372 # TRACK FITTING AND SMOOTHING OPTIONS
374 displacedRegionalStepFitterSmoother = TrackingTools.TrackFitters.RungeKuttaFitters_cff.KFFittingSmootherWithOutliersRejectionAndRK.clone(
375  ComponentName = 'displacedRegionalStepFitterSmoother',
376  EstimateCut = 30,
377  MinNumberOfHits = 7,
378  Fitter = 'displacedRegionalStepRKFitter',
379  Smoother = 'displacedRegionalStepRKSmoother'
380  )
381 
382 displacedRegionalStepFitterSmootherForLoopers = displacedRegionalStepFitterSmoother.clone(
383  ComponentName = 'displacedRegionalStepFitterSmootherForLoopers',
384  Fitter = 'displacedRegionalStepRKFitterForLoopers',
385  Smoother = 'displacedRegionalStepRKSmootherForLoopers'
386 )
387 
388 # Also necessary to specify minimum number of hits after final track fit
389 displacedRegionalStepRKTrajectoryFitter = TrackingTools.TrackFitters.RungeKuttaFitters_cff.RKTrajectoryFitter.clone(
390  ComponentName = 'displacedRegionalStepRKFitter',
391  minHits = 7
392 )
393 
394 displacedRegionalStepRKTrajectoryFitterForLoopers = displacedRegionalStepRKTrajectoryFitter.clone(
395  ComponentName = 'displacedRegionalStepRKFitterForLoopers',
396  Propagator = 'PropagatorWithMaterialForLoopers',
397 )
398 
399 displacedRegionalStepRKTrajectorySmoother = TrackingTools.TrackFitters.RungeKuttaFitters_cff.RKTrajectorySmoother.clone(
400  ComponentName = 'displacedRegionalStepRKSmoother',
401  errorRescaling = 10.0,
402  minHits = 7
403 )
404 
405 displacedRegionalStepRKTrajectorySmootherForLoopers = displacedRegionalStepRKTrajectorySmoother.clone(
406  ComponentName = 'displacedRegionalStepRKSmootherForLoopers',
407  Propagator = 'PropagatorWithMaterialForLoopers',
408 )
409 
411 displacedRegionalFlexibleKFFittingSmoother = TrackingTools.TrackFitters.FlexibleKFFittingSmoother_cfi.FlexibleKFFittingSmoother.clone(
412  ComponentName = 'displacedRegionalFlexibleKFFittingSmoother',
413  standardFitter = 'displacedRegionalStepFitterSmoother',
414  looperFitter = 'displacedRegionalStepFitterSmootherForLoopers',
415 )
416 
417 # TRACK FITTING
418 import RecoTracker.TrackProducer.TrackProducerIterativeDefault_cfi
419 displacedRegionalStepTracks = RecoTracker.TrackProducer.TrackProducerIterativeDefault_cfi.TrackProducerIterativeDefault.clone(
420  src = 'displacedRegionalStepTrackCandidates',
421  AlgorithmName = 'displacedRegionalStep',
422  #Fitter = 'displacedRegionalStepFitterSmoother',
423  Fitter = 'displacedRegionalFlexibleKFFittingSmoother',
424 )
425 
426 # TRACK SELECTION AND QUALITY FLAG SETTING.
429 displacedRegionalStepClassifier1 = TrackMVAClassifierDetached.clone(
430  src = 'displacedRegionalStepTracks',
431  mva = dict(GBRForestLabel = 'MVASelectorIter6_13TeV'),
432  qualityCuts = [-0.6,-0.45,-0.3]
433 )
434 
435 displacedRegionalStepClassifier2 = TrackMVAClassifierPrompt.clone(
436  src = 'displacedRegionalStepTracks',
437  mva = dict(GBRForestLabel = 'MVASelectorIter0_13TeV'),
438  qualityCuts = [0.0,0.0,0.0]
439 )
440 
442 displacedRegionalStep = ClassifierMerger.clone(
443  inputClassifiers=['displacedRegionalStepClassifier1','displacedRegionalStepClassifier2']
444 )
445 
446 from Configuration.Eras.Modifier_trackingPhase1_cff import trackingPhase1
447 trackingPhase1.toReplaceWith(displacedRegionalStep, displacedRegionalStepClassifier1.clone(
448  mva = dict(GBRForestLabel = 'MVASelectorTobTecStep_Phase1'),
449  qualityCuts = [-0.6,-0.45,-0.3]
450 ))
451 
455 trackdnn.toReplaceWith(displacedRegionalStep, trackTfClassifier.clone(
456  src = 'displacedRegionalStepTracks',
457  qualityCuts = qualityCutDictionary.DisplacedRegionalStep.value()
458 ))
459 
460 DisplacedRegionalStepTask = cms.Task(displacedRegionalStepClusters,
461  displacedRegionalStepSeedLayersTripl,
462  displacedRegionalStepInputTracks,
463  displacedRegionalStepSeedingV0Candidates,
464  displacedRegionalStepSeedingVertices,
465  displacedRegionalStepFarTrackingRegions,
466  displacedRegionalStepHitDoubletsTripl,
467  displacedRegionalStepHitTripletsTripl,
468  displacedRegionalStepSeedsTripl,
469  displacedRegionalStepSeedLayersPair,
470  displacedRegionalStepHitDoubletsPair,
471  displacedRegionalStepSeedsPair,
472  displacedRegionalStepPLSeedLayersTripl,
473  displacedRegionalStepNearTrackingRegions,
474  displacedRegionalStepPLHitDoubletsTripl,
475  displacedRegionalStepPLHitTripletsTripl,
476  displacedRegionalStepPLSeedsTripl,
477  displacedRegionalStepSeeds,
478  displacedRegionalStepTrackCandidates,
479  displacedRegionalStepTracks,
480  displacedRegionalStepClassifier1,displacedRegionalStepClassifier2,
481  displacedRegionalStep)
482 DisplacedRegionalStep = cms.Sequence(DisplacedRegionalStepTask)
483 
484 _DisplacedRegionalStepTask_trackingMkFit = DisplacedRegionalStepTask.copy()
485 _DisplacedRegionalStepTask_trackingMkFit.add(displacedRegionalStepTrackCandidatesMkFitSeeds, displacedRegionalStepTrackCandidatesMkFit, displacedRegionalStepTrackCandidatesMkFitConfig)
486 trackingMkFitDisplacedRegionalStep.toReplaceWith(DisplacedRegionalStepTask, _DisplacedRegionalStepTask_trackingMkFit)