File indexing completed on 2023-05-10 03:54:41
0001 import FWCore.ParameterSet.Config as cms
0002 import RecoTracker.IterativeTracking.iterativeTkConfig as _cfg
0003
0004
0005 from Configuration.ProcessModifiers.trackdnn_cff import trackdnn
0006 from RecoTracker.IterativeTracking.dnnQualityCuts import qualityCutDictionary
0007
0008
0009 from Configuration.ProcessModifiers.trackingNoLoopers_cff import trackingNoLoopers
0010
0011
0012
0013
0014
0015 displacedRegionalStepClusters = _cfg.clusterRemoverForIter("DisplacedRegionalStep")
0016 for _eraName, _postfix, _era in _cfg.nonDefaultEras():
0017 _era.toReplaceWith(displacedRegionalStepClusters, _cfg.clusterRemoverForIter("DisplacedRegionalStep", _eraName, _postfix))
0018
0019
0020 from RecoLocalTracker.SiStripClusterizer.SiStripClusterChargeCut_cfi import *
0021 import RecoTracker.TkSeedingLayers.seedingLayersEDProducer_cfi as _mod
0022
0023 displacedRegionalStepSeedLayersTripl = _mod.seedingLayersEDProducer.clone(
0024 layerList = [
0025
0026 'TOB1+TOB2+MTOB3',
0027
0028 'TOB1+TOB2+MTEC1_pos','TOB1+TOB2+MTEC1_neg',
0029 ],
0030 TOB = dict(
0031 TTRHBuilder = cms.string('WithTrackAngle'),
0032 clusterChargeCut = cms.PSet(refToPSet_ = cms.string('SiStripClusterChargeCutTight')),
0033 matchedRecHits = cms.InputTag('siStripMatchedRecHits','matchedRecHit'),
0034 skipClusters = cms.InputTag('displacedRegionalStepClusters')
0035 ),
0036 MTOB = dict(
0037 TTRHBuilder = cms.string('WithTrackAngle'),
0038 clusterChargeCut = cms.PSet(refToPSet_ = cms.string('SiStripClusterChargeCutTight')),
0039 skipClusters = cms.InputTag('displacedRegionalStepClusters'),
0040 rphiRecHits = cms.InputTag('siStripMatchedRecHits','rphiRecHit')
0041 ),
0042 MTEC = dict(
0043 rphiRecHits = cms.InputTag('siStripMatchedRecHits','rphiRecHit'),
0044 skipClusters = cms.InputTag('displacedRegionalStepClusters'),
0045 useRingSlector = cms.bool(True),
0046 TTRHBuilder = cms.string('WithTrackAngle'),
0047 clusterChargeCut = cms.PSet(refToPSet_ = cms.string('SiStripClusterChargeCutTight')),
0048 minRing = cms.int32(6),
0049 maxRing = cms.int32(7)
0050 )
0051 )
0052
0053
0054 from RecoTracker.FinalTrackSelectors.displacedRegionalStepInputTracks_cfi import displacedRegionalStepInputTracks
0055 from RecoVertex.V0Producer.generalV0Candidates_cfi import generalV0Candidates as _generalV0Candidates
0056 from RecoTracker.DisplacedRegionalTracking.displacedRegionProducer_cfi import displacedRegionProducer as _displacedRegionProducer
0057 displacedRegionalStepSeedingV0Candidates = _generalV0Candidates.clone(
0058 trackRecoAlgorithm = "displacedRegionalStepInputTracks",
0059 doLambdas = False,
0060 doFit = False,
0061 useRefTracks = False,
0062 vtxDecayXYCut = 1.,
0063 ssVtxDecayXYCut = 5.,
0064 allowSS = True,
0065 innerTkDCACut = 0.2,
0066 allowWideAngleVtx = True,
0067 mPiPiCut = 13000.,
0068 cosThetaXYCut = 0.,
0069 kShortMassCut = 13000.,
0070 )
0071 displacedRegionalStepSeedingVertices = _displacedRegionProducer.clone(
0072 minRadius = 2.0,
0073 nearThreshold = 20.0,
0074 farThreshold = -1.0,
0075 discriminatorCut = 0.5,
0076 trackClusters = ["displacedRegionalStepSeedingV0Candidates", "Kshort"],
0077 )
0078
0079 from RecoTracker.TkTrackingRegions.globalTrackingRegionWithVertices_cfi import globalTrackingRegionWithVertices as _globalTrackingRegionWithVertices
0080 displacedRegionalStepFarTrackingRegions = _globalTrackingRegionWithVertices.clone(RegionPSet = dict(
0081 originRadius = 1.0,
0082 fixedError = 1.0,
0083 VertexCollection = ["displacedRegionalStepSeedingVertices", "farRegionsOfInterest"],
0084 useFakeVertices = True,
0085 ptMin = 0.6,
0086 allowEmpty = True
0087 ))
0088
0089
0090 from RecoTracker.PixelLowPtUtilities.ClusterShapeHitFilterESProducer_cfi import ClusterShapeHitFilterESProducer as _ClusterShapeHitFilterESProducer
0091 displacedRegionalStepClusterShapeHitFilter = _ClusterShapeHitFilterESProducer.clone(
0092 ComponentName = 'displacedRegionalStepClusterShapeHitFilter',
0093 doStripShapeCut = False,
0094 clusterChargeCut = dict(refToPSet_ = 'SiStripClusterChargeCutTight')
0095 )
0096
0097 from RecoTracker.TkHitPairs.hitPairEDProducer_cfi import hitPairEDProducer as _hitPairEDProducer
0098 displacedRegionalStepHitDoubletsTripl = _hitPairEDProducer.clone(
0099 seedingLayers = "displacedRegionalStepSeedLayersTripl",
0100 trackingRegions = "displacedRegionalStepFarTrackingRegions",
0101 maxElement = 50000000,
0102 produceIntermediateHitDoublets = True,
0103 )
0104 from RecoTracker.TkSeedGenerator.multiHitFromChi2EDProducer_cfi import multiHitFromChi2EDProducer as _multiHitFromChi2EDProducer
0105 displacedRegionalStepHitTripletsTripl = _multiHitFromChi2EDProducer.clone(
0106 doublets = "displacedRegionalStepHitDoubletsTripl",
0107 extraPhiKDBox = 0.01,
0108 )
0109 from RecoTracker.TkSeedGenerator.seedCreatorFromRegionConsecutiveHitsEDProducer_cff import seedCreatorFromRegionConsecutiveHitsEDProducer as _seedCreatorFromRegionConsecutiveHitsTripletOnlyEDProducer
0110 _displacedRegionalStepSeedComparitorPSet = dict(
0111 ComponentName = 'CombinedSeedComparitor',
0112 mode = cms.string("and"),
0113 comparitors = cms.VPSet(
0114 cms.PSet(
0115 ComponentName = cms.string('PixelClusterShapeSeedComparitor'),
0116 FilterAtHelixStage = cms.bool(True),
0117 FilterPixelHits = cms.bool(False),
0118 FilterStripHits = cms.bool(True),
0119 ClusterShapeHitFilterName = cms.string('displacedRegionalStepClusterShapeHitFilter'),
0120 ClusterShapeCacheSrc = cms.InputTag("siPixelClusterShapeCache")
0121 )
0122 )
0123 )
0124
0125 displacedRegionalStepSeedsTripl = _seedCreatorFromRegionConsecutiveHitsTripletOnlyEDProducer.clone(
0126 seedingHitSets = "displacedRegionalStepHitTripletsTripl",
0127 SeedComparitorPSet = _displacedRegionalStepSeedComparitorPSet,
0128 )
0129
0130 from RecoTracker.PixelLowPtUtilities.StripSubClusterShapeSeedFilter_cfi import StripSubClusterShapeSeedFilter as _StripSubClusterShapeSeedFilter
0131 from Configuration.ProcessModifiers.approxSiStripClusters_cff import approxSiStripClusters
0132 (~approxSiStripClusters).toModify(displacedRegionalStepSeedsTripl.SeedComparitorPSet.comparitors, func = lambda list: list.append(_StripSubClusterShapeSeedFilter.clone()) )
0133
0134
0135 displacedRegionalStepSeedLayersPair = _mod.seedingLayersEDProducer.clone(
0136 layerList = ['TOB1+TEC1_pos','TOB1+TEC1_neg',
0137 'TEC1_pos+TEC2_pos','TEC1_neg+TEC2_neg',
0138 'TEC2_pos+TEC3_pos','TEC2_neg+TEC3_neg',
0139 'TEC3_pos+TEC4_pos','TEC3_neg+TEC4_neg',
0140 'TEC4_pos+TEC5_pos','TEC4_neg+TEC5_neg',
0141 'TEC5_pos+TEC6_pos','TEC5_neg+TEC6_neg',
0142 'TEC6_pos+TEC7_pos','TEC6_neg+TEC7_neg'],
0143 TOB = dict(
0144 TTRHBuilder = cms.string('WithTrackAngle'),
0145 clusterChargeCut = cms.PSet(refToPSet_ = cms.string('SiStripClusterChargeCutTight')),
0146 matchedRecHits = cms.InputTag('siStripMatchedRecHits','matchedRecHit'),
0147 skipClusters = cms.InputTag('displacedRegionalStepClusters')
0148 ),
0149 TEC = dict(
0150 matchedRecHits = cms.InputTag('siStripMatchedRecHits','matchedRecHit'),
0151 skipClusters = cms.InputTag('displacedRegionalStepClusters'),
0152 useRingSlector = cms.bool(True),
0153 TTRHBuilder = cms.string('WithTrackAngle'),
0154 clusterChargeCut = cms.PSet(refToPSet_ = cms.string('SiStripClusterChargeCutTight')),
0155 minRing = cms.int32(5),
0156 maxRing = cms.int32(5)
0157 )
0158 )
0159
0160
0161
0162 displacedRegionalStepHitDoubletsPair = _hitPairEDProducer.clone(
0163 seedingLayers = "displacedRegionalStepSeedLayersPair",
0164 trackingRegions = "displacedRegionalStepFarTrackingRegions",
0165 produceSeedingHitSets = True,
0166 maxElementTotal = 12000000,
0167 )
0168 from RecoTracker.TkSeedGenerator.seedCreatorFromRegionConsecutiveHitsEDProducer_cff import seedCreatorFromRegionConsecutiveHitsEDProducer as _seedCreatorFromRegionConsecutiveHitsEDProducer
0169 displacedRegionalStepSeedsPair = _seedCreatorFromRegionConsecutiveHitsEDProducer.clone(
0170 seedingHitSets = "displacedRegionalStepHitDoubletsPair",
0171 SeedComparitorPSet = _displacedRegionalStepSeedComparitorPSet,
0172 )
0173
0174
0175 displacedRegionalStepPLSeedLayersTripl = _mod.seedingLayersEDProducer.clone(
0176 layerList = [
0177
0178 'TIB1+TIB2+MTIB3',
0179
0180 'TIB1+TIB2+MTID1_pos','TIB1+TIB2+MTID1_neg',
0181
0182 'TID1_pos+TID2_pos+TID3_pos','TID1_neg+TID2_neg+TID3_neg',
0183 'TID1_pos+TID2_pos+MTID3_pos','TID1_neg+TID2_neg+MTID3_neg',
0184 'TID1_pos+TID2_pos+MTEC1_pos','TID1_neg+TID2_neg+MTEC1_neg',
0185
0186 'TID2_pos+TID3_pos+TEC1_pos','TID2_neg+TID3_neg+TEC1_neg',
0187 'TID2_pos+TID3_pos+MTEC1_pos','TID2_neg+TID3_neg+MTEC1_neg',
0188
0189
0190
0191
0192
0193
0194
0195
0196
0197
0198
0199
0200
0201 ],
0202 TIB = dict(
0203 TTRHBuilder = cms.string('WithTrackAngle'),
0204 clusterChargeCut = cms.PSet(refToPSet_ = cms.string('SiStripClusterChargeCutTight')),
0205 matchedRecHits = cms.InputTag('siStripMatchedRecHits','matchedRecHit'),
0206 skipClusters = cms.InputTag('displacedRegionalStepClusters')
0207 ),
0208 MTIB = dict(
0209 TTRHBuilder = cms.string('WithTrackAngle'),
0210 clusterChargeCut = cms.PSet(refToPSet_ = cms.string('SiStripClusterChargeCutTight')),
0211 skipClusters = cms.InputTag('displacedRegionalStepClusters'),
0212 rphiRecHits = cms.InputTag('siStripMatchedRecHits','rphiRecHit')
0213 ),
0214 TID = dict(
0215 matchedRecHits = cms.InputTag('siStripMatchedRecHits','matchedRecHit'),
0216 skipClusters = cms.InputTag('displacedRegionalStepClusters'),
0217 useRingSlector = cms.bool(True),
0218 TTRHBuilder = cms.string('WithTrackAngle'),
0219 clusterChargeCut = cms.PSet(refToPSet_ = cms.string('SiStripClusterChargeCutTight')),
0220 minRing = cms.int32(1),
0221 maxRing = cms.int32(2)
0222 ),
0223 MTID = dict(
0224 rphiRecHits = cms.InputTag('siStripMatchedRecHits','rphiRecHit'),
0225 skipClusters = cms.InputTag('displacedRegionalStepClusters'),
0226 useRingSlector = cms.bool(True),
0227 TTRHBuilder = cms.string('WithTrackAngle'),
0228 clusterChargeCut = cms.PSet(refToPSet_ = cms.string('SiStripClusterChargeCutTight')),
0229 minRing = cms.int32(3),
0230 maxRing = cms.int32(3)
0231 ),
0232 TEC = dict(
0233 matchedRecHits = cms.InputTag('siStripMatchedRecHits','matchedRecHit'),
0234 skipClusters = cms.InputTag('displacedRegionalStepClusters'),
0235 useRingSlector = cms.bool(True),
0236 TTRHBuilder = cms.string('WithTrackAngle'),
0237 clusterChargeCut = cms.PSet(refToPSet_ = cms.string('SiStripClusterChargeCutTight')),
0238 minRing = cms.int32(1),
0239 maxRing = cms.int32(2)
0240 ),
0241 MTEC = dict(
0242 rphiRecHits = cms.InputTag('siStripMatchedRecHits','rphiRecHit'),
0243 skipClusters = cms.InputTag('displacedRegionalStepClusters'),
0244 useRingSlector = cms.bool(True),
0245 TTRHBuilder = cms.string('WithTrackAngle'),
0246 clusterChargeCut = cms.PSet(refToPSet_ = cms.string('SiStripClusterChargeCutTight')),
0247 minRing = cms.int32(3),
0248 maxRing = cms.int32(3)
0249 )
0250 )
0251
0252
0253 displacedRegionalStepNearTrackingRegions = _globalTrackingRegionWithVertices.clone(RegionPSet = dict(
0254 originRadius = 1.0,
0255 fixedError = 1.0,
0256 VertexCollection = ["displacedRegionalStepSeedingVertices", "nearRegionsOfInterest"],
0257 useFakeVertices = True,
0258 ptMin = 0.6,
0259 allowEmpty = True
0260 ))
0261
0262 displacedRegionalStepPLHitDoubletsTripl = _hitPairEDProducer.clone(
0263 seedingLayers = 'displacedRegionalStepPLSeedLayersTripl',
0264 trackingRegions = 'displacedRegionalStepNearTrackingRegions',
0265 maxElement = 50000000,
0266 produceIntermediateHitDoublets = True,
0267 )
0268
0269 displacedRegionalStepPLHitTripletsTripl = _multiHitFromChi2EDProducer.clone(
0270 doublets = 'displacedRegionalStepPLHitDoubletsTripl',
0271 )
0272
0273 displacedRegionalStepPLSeedsTripl = _seedCreatorFromRegionConsecutiveHitsTripletOnlyEDProducer.clone(
0274 seedingHitSets = 'displacedRegionalStepPLHitTripletsTripl',
0275 SeedComparitorPSet = _displacedRegionalStepSeedComparitorPSet,
0276 )
0277
0278
0279 import RecoTracker.TkSeedGenerator.GlobalCombinedSeeds_cfi
0280 displacedRegionalStepSeeds = RecoTracker.TkSeedGenerator.GlobalCombinedSeeds_cfi.globalCombinedSeeds.clone(
0281 seedCollections = ['displacedRegionalStepSeedsTripl', 'displacedRegionalStepSeedsPair', 'displacedRegionalStepPLSeedsTripl']
0282 )
0283
0284
0285 import TrackingTools.TrajectoryFiltering.TrajectoryFilter_cff
0286 _displacedRegionalStepTrajectoryFilterBase = TrackingTools.TrajectoryFiltering.TrajectoryFilter_cff.CkfBaseTrajectoryFilter_block.clone(
0287 maxLostHits = 0,
0288 minimumNumberOfHits = 5,
0289 minPt = 0.1,
0290 minHitsMinPt = 3
0291 )
0292 displacedRegionalStepTrajectoryFilter = _displacedRegionalStepTrajectoryFilterBase.clone(
0293 seedPairPenalty = 1,
0294 )
0295
0296 displacedRegionalStepInOutTrajectoryFilter = displacedRegionalStepTrajectoryFilter.clone(
0297 minimumNumberOfHits = 4,
0298 )
0299
0300 import RecoTracker.MeasurementDet.Chi2ChargeMeasurementEstimator_cfi
0301 displacedRegionalStepChi2Est = RecoTracker.MeasurementDet.Chi2ChargeMeasurementEstimator_cfi.Chi2ChargeMeasurementEstimator.clone(
0302 ComponentName = 'displacedRegionalStepChi2Est',
0303 nSigma = 3.0,
0304 MaxChi2 = 16.0,
0305 clusterChargeCut = cms.PSet(refToPSet_ = cms.string('SiStripClusterChargeCutTight'))
0306 )
0307
0308
0309 import RecoTracker.CkfPattern.GroupedCkfTrajectoryBuilder_cfi
0310 displacedRegionalStepTrajectoryBuilder = RecoTracker.CkfPattern.GroupedCkfTrajectoryBuilder_cfi.GroupedCkfTrajectoryBuilderIterativeDefault.clone(
0311 trajectoryFilter = dict(refToPSet_ = 'displacedRegionalStepTrajectoryFilter'),
0312 inOutTrajectoryFilter = dict(refToPSet_ = 'displacedRegionalStepInOutTrajectoryFilter'),
0313 useSameTrajFilter = False,
0314 minNrOfHitsForRebuild = 4,
0315 alwaysUseInvalidHits = False,
0316 maxCand = 2,
0317 estimator = 'displacedRegionalStepChi2Est',
0318
0319 maxDPhiForLooperReconstruction = 2.0,
0320 maxPtForLooperReconstruction = 0.7,
0321 )
0322 trackingNoLoopers.toModify(displacedRegionalStepTrajectoryBuilder,
0323 maxPtForLooperReconstruction = 0.0)
0324
0325
0326 import RecoTracker.CkfPattern.CkfTrackCandidates_cfi
0327
0328 _displacedRegionalStepTrackCandidatesCkf = RecoTracker.CkfPattern.CkfTrackCandidates_cfi.ckfTrackCandidatesIterativeDefault.clone(
0329 src = 'displacedRegionalStepSeeds',
0330 clustersToSkip = 'displacedRegionalStepClusters',
0331
0332 numHitsForSeedCleaner = 50,
0333 onlyPixelHitsForSeedCleaner = False,
0334 TrajectoryBuilderPSet = dict(refToPSet_ = 'displacedRegionalStepTrajectoryBuilder'),
0335 doSeedingRegionRebuilding = True,
0336 useHitsSplitting = True,
0337 cleanTrajectoryAfterInOut = True,
0338 TrajectoryCleaner = 'displacedRegionalStepTrajectoryCleanerBySharedHits',
0339 )
0340 displacedRegionalStepTrackCandidates = _displacedRegionalStepTrackCandidatesCkf.clone()
0341
0342 from Configuration.ProcessModifiers.trackingMkFitDisplacedRegionalStep_cff import trackingMkFitDisplacedRegionalStep
0343 import RecoTracker.MkFit.mkFitSeedConverter_cfi as mkFitSeedConverter_cfi
0344 import RecoTracker.MkFit.mkFitIterationConfigESProducer_cfi as mkFitIterationConfigESProducer_cfi
0345 import RecoTracker.MkFit.mkFitProducer_cfi as mkFitProducer_cfi
0346 import RecoTracker.MkFit.mkFitOutputConverter_cfi as mkFitOutputConverter_cfi
0347 displacedRegionalStepTrackCandidatesMkFitSeeds = mkFitSeedConverter_cfi.mkFitSeedConverter.clone(
0348 seeds = 'displacedRegionalStepSeeds',
0349 )
0350 displacedRegionalStepTrackCandidatesMkFitConfig = mkFitIterationConfigESProducer_cfi.mkFitIterationConfigESProducer.clone(
0351 ComponentName = 'displacedRegionalStepTrackCandidatesMkFitConfig',
0352 config = 'RecoTracker/MkFit/data/mkfit-phase1-tobTecStep.json',
0353 )
0354 displacedRegionalStepTrackCandidatesMkFit = mkFitProducer_cfi.mkFitProducer.clone(
0355 seeds = 'displacedRegionalStepTrackCandidatesMkFitSeeds',
0356 config = ('', 'displacedRegionalStepTrackCandidatesMkFitConfig'),
0357 clustersToSkip = 'displacedRegionalStepClusters',
0358 )
0359 trackingMkFitDisplacedRegionalStep.toReplaceWith(displacedRegionalStepTrackCandidates, mkFitOutputConverter_cfi.mkFitOutputConverter.clone(
0360 seeds = 'displacedRegionalStepSeeds',
0361 mkFitSeeds = 'displacedRegionalStepTrackCandidatesMkFitSeeds',
0362 tracks = 'displacedRegionalStepTrackCandidatesMkFit',
0363 ))
0364
0365 from TrackingTools.TrajectoryCleaning.TrajectoryCleanerBySharedHits_cfi import trajectoryCleanerBySharedHits
0366 displacedRegionalStepTrajectoryCleanerBySharedHits = trajectoryCleanerBySharedHits.clone(
0367 ComponentName = 'displacedRegionalStepTrajectoryCleanerBySharedHits',
0368 fractionShared = 0.09,
0369 allowSharedFirstHit = True
0370 )
0371
0372
0373 import TrackingTools.TrackFitters.RungeKuttaFitters_cff
0374 displacedRegionalStepFitterSmoother = TrackingTools.TrackFitters.RungeKuttaFitters_cff.KFFittingSmootherWithOutliersRejectionAndRK.clone(
0375 ComponentName = 'displacedRegionalStepFitterSmoother',
0376 EstimateCut = 30,
0377 MinNumberOfHits = 7,
0378 Fitter = 'displacedRegionalStepRKFitter',
0379 Smoother = 'displacedRegionalStepRKSmoother'
0380 )
0381
0382 displacedRegionalStepFitterSmootherForLoopers = displacedRegionalStepFitterSmoother.clone(
0383 ComponentName = 'displacedRegionalStepFitterSmootherForLoopers',
0384 Fitter = 'displacedRegionalStepRKFitterForLoopers',
0385 Smoother = 'displacedRegionalStepRKSmootherForLoopers'
0386 )
0387
0388
0389 displacedRegionalStepRKTrajectoryFitter = TrackingTools.TrackFitters.RungeKuttaFitters_cff.RKTrajectoryFitter.clone(
0390 ComponentName = 'displacedRegionalStepRKFitter',
0391 minHits = 7
0392 )
0393
0394 displacedRegionalStepRKTrajectoryFitterForLoopers = displacedRegionalStepRKTrajectoryFitter.clone(
0395 ComponentName = 'displacedRegionalStepRKFitterForLoopers',
0396 Propagator = 'PropagatorWithMaterialForLoopers',
0397 )
0398
0399 displacedRegionalStepRKTrajectorySmoother = TrackingTools.TrackFitters.RungeKuttaFitters_cff.RKTrajectorySmoother.clone(
0400 ComponentName = 'displacedRegionalStepRKSmoother',
0401 errorRescaling = 10.0,
0402 minHits = 7
0403 )
0404
0405 displacedRegionalStepRKTrajectorySmootherForLoopers = displacedRegionalStepRKTrajectorySmoother.clone(
0406 ComponentName = 'displacedRegionalStepRKSmootherForLoopers',
0407 Propagator = 'PropagatorWithMaterialForLoopers',
0408 )
0409
0410 import TrackingTools.TrackFitters.FlexibleKFFittingSmoother_cfi
0411 displacedRegionalFlexibleKFFittingSmoother = TrackingTools.TrackFitters.FlexibleKFFittingSmoother_cfi.FlexibleKFFittingSmoother.clone(
0412 ComponentName = 'displacedRegionalFlexibleKFFittingSmoother',
0413 standardFitter = 'displacedRegionalStepFitterSmoother',
0414 looperFitter = 'displacedRegionalStepFitterSmootherForLoopers',
0415 )
0416
0417
0418 import RecoTracker.TrackProducer.TrackProducerIterativeDefault_cfi
0419 displacedRegionalStepTracks = RecoTracker.TrackProducer.TrackProducerIterativeDefault_cfi.TrackProducerIterativeDefault.clone(
0420 src = 'displacedRegionalStepTrackCandidates',
0421 AlgorithmName = 'displacedRegionalStep',
0422
0423 Fitter = 'displacedRegionalFlexibleKFFittingSmoother',
0424 )
0425
0426
0427 from RecoTracker.FinalTrackSelectors.TrackMVAClassifierPrompt_cfi import *
0428 from RecoTracker.FinalTrackSelectors.TrackMVAClassifierDetached_cfi import *
0429 displacedRegionalStepClassifier1 = TrackMVAClassifierDetached.clone(
0430 src = 'displacedRegionalStepTracks',
0431 mva = dict(GBRForestLabel = 'MVASelectorIter6_13TeV'),
0432 qualityCuts = [-0.6,-0.45,-0.3]
0433 )
0434
0435 displacedRegionalStepClassifier2 = TrackMVAClassifierPrompt.clone(
0436 src = 'displacedRegionalStepTracks',
0437 mva = dict(GBRForestLabel = 'MVASelectorIter0_13TeV'),
0438 qualityCuts = [0.0,0.0,0.0]
0439 )
0440
0441 from RecoTracker.FinalTrackSelectors.ClassifierMerger_cfi import *
0442 displacedRegionalStep = ClassifierMerger.clone(
0443 inputClassifiers=['displacedRegionalStepClassifier1','displacedRegionalStepClassifier2']
0444 )
0445
0446 from Configuration.Eras.Modifier_trackingPhase1_cff import trackingPhase1
0447 trackingPhase1.toReplaceWith(displacedRegionalStep, displacedRegionalStepClassifier1.clone(
0448 mva = dict(GBRForestLabel = 'MVASelectorTobTecStep_Phase1'),
0449 qualityCuts = [-0.6,-0.45,-0.3]
0450 ))
0451
0452 from RecoTracker.FinalTrackSelectors.trackTfClassifier_cfi import *
0453 from RecoTracker.FinalTrackSelectors.trackSelectionTf_cfi import *
0454 from RecoTracker.FinalTrackSelectors.trackSelectionTf_CKF_cfi import *
0455 trackdnn.toReplaceWith(displacedRegionalStep, trackTfClassifier.clone(
0456 src = 'displacedRegionalStepTracks',
0457 qualityCuts = qualityCutDictionary.DisplacedRegionalStep.value()
0458 ))
0459
0460 DisplacedRegionalStepTask = cms.Task(displacedRegionalStepClusters,
0461 displacedRegionalStepSeedLayersTripl,
0462 displacedRegionalStepInputTracks,
0463 displacedRegionalStepSeedingV0Candidates,
0464 displacedRegionalStepSeedingVertices,
0465 displacedRegionalStepFarTrackingRegions,
0466 displacedRegionalStepHitDoubletsTripl,
0467 displacedRegionalStepHitTripletsTripl,
0468 displacedRegionalStepSeedsTripl,
0469 displacedRegionalStepSeedLayersPair,
0470 displacedRegionalStepHitDoubletsPair,
0471 displacedRegionalStepSeedsPair,
0472 displacedRegionalStepPLSeedLayersTripl,
0473 displacedRegionalStepNearTrackingRegions,
0474 displacedRegionalStepPLHitDoubletsTripl,
0475 displacedRegionalStepPLHitTripletsTripl,
0476 displacedRegionalStepPLSeedsTripl,
0477 displacedRegionalStepSeeds,
0478 displacedRegionalStepTrackCandidates,
0479 displacedRegionalStepTracks,
0480 displacedRegionalStepClassifier1,displacedRegionalStepClassifier2,
0481 displacedRegionalStep)
0482 DisplacedRegionalStep = cms.Sequence(DisplacedRegionalStepTask)
0483
0484 _DisplacedRegionalStepTask_trackingMkFit = DisplacedRegionalStepTask.copy()
0485 _DisplacedRegionalStepTask_trackingMkFit.add(displacedRegionalStepTrackCandidatesMkFitSeeds, displacedRegionalStepTrackCandidatesMkFit, displacedRegionalStepTrackCandidatesMkFitConfig)
0486 trackingMkFitDisplacedRegionalStep.toReplaceWith(DisplacedRegionalStepTask, _DisplacedRegionalStepTask_trackingMkFit)