File indexing completed on 2024-10-23 22:48:12
0001 import FWCore.ParameterSet.Config as cms
0002
0003
0004 from Configuration.ProcessModifiers.trackdnn_cff import trackdnn
0005 from RecoTracker.IterativeTracking.dnnQualityCuts import qualityCutDictionary
0006
0007
0008 from Configuration.ProcessModifiers.trackingNoLoopers_cff import trackingNoLoopers
0009
0010
0011
0012
0013 jetsForCoreTracking = cms.EDFilter('CandPtrSelector', src = cms.InputTag('ak4CaloJetsForTrk'), cut = cms.string('pt > 100 && abs(eta) < 2.5'), filter = cms.bool(False))
0014
0015 jetsForCoreTrackingBarrel = jetsForCoreTracking.clone( cut = 'pt > 100 && abs(eta) < 1.4' )
0016 jetsForCoreTrackingEndcap = jetsForCoreTracking.clone( cut = 'pt > 100 && abs(eta) < 2.5' )
0017
0018
0019 firstStepGoodPrimaryVertices = cms.EDFilter('PrimaryVertexObjectFilter',
0020 filterParams = cms.PSet(
0021 minNdof = cms.double(25.0),
0022 maxZ = cms.double(15.0),
0023 maxRho = cms.double(2.0)
0024 ),
0025 src=cms.InputTag('firstStepPrimaryVertices')
0026 )
0027
0028 import RecoTracker.TkSeedingLayers.seedingLayersEDProducer_cfi as _mod
0029
0030
0031 jetCoreRegionalStepSeedLayers = _mod.seedingLayersEDProducer.clone(
0032 layerList = ['BPix1+BPix2', 'BPix1+BPix3', 'BPix2+BPix3',
0033 'BPix1+FPix1_pos', 'BPix1+FPix1_neg',
0034 'BPix2+FPix1_pos', 'BPix2+FPix1_neg',
0035 'FPix1_pos+FPix2_pos', 'FPix1_neg+FPix2_neg',
0036
0037 'BPix3+TIB1','BPix3+TIB2'],
0038 TIB = dict(
0039 matchedRecHits = cms.InputTag('siStripMatchedRecHits','matchedRecHit'),
0040 TTRHBuilder = cms.string('WithTrackAngle'),
0041 clusterChargeCut = cms.PSet(refToPSet_ = cms.string('SiStripClusterChargeCutNone'))
0042 ),
0043 BPix = dict(
0044 useErrorsFromParam = cms.bool(True),
0045 hitErrorRPhi = cms.double(0.0027),
0046 hitErrorRZ = cms.double(0.006),
0047 TTRHBuilder = cms.string('WithTrackAngle'),
0048 HitProducer = cms.string('siPixelRecHits'),
0049
0050 ),
0051 FPix = dict(
0052 useErrorsFromParam = cms.bool(True),
0053 hitErrorRPhi = cms.double(0.0051),
0054 hitErrorRZ = cms.double(0.0036),
0055 TTRHBuilder = cms.string('WithTrackAngle'),
0056 HitProducer = cms.string('siPixelRecHits'),
0057
0058 )
0059 )
0060 from Configuration.Eras.Modifier_trackingPhase1_cff import trackingPhase1
0061 _layerListForPhase1 = [
0062 'BPix1+BPix2', 'BPix1+BPix3', 'BPix1+BPix4',
0063 'BPix2+BPix3', 'BPix2+BPix4',
0064 'BPix3+BPix4',
0065 'BPix1+FPix1_pos', 'BPix1+FPix1_neg',
0066 'BPix2+FPix1_pos', 'BPix2+FPix1_neg',
0067 'FPix1_pos+FPix2_pos', 'FPix1_neg+FPix2_neg',
0068 'FPix1_pos+FPix3_pos', 'FPix1_neg+FPix3_neg',
0069 'FPix2_pos+FPix3_pos', 'FPix2_neg+FPix3_neg',
0070
0071 'BPix4+TIB1','BPix4+TIB2'
0072 ]
0073 trackingPhase1.toModify(jetCoreRegionalStepSeedLayers, layerList = _layerListForPhase1)
0074 from Configuration.ProcessModifiers.jetCoreInPhase2_cff import jetCoreInPhase2
0075 _layerListForPhase2 = [
0076 'BPix1+BPix2', 'BPix1+BPix3', 'BPix1+BPix4',
0077 'BPix2+BPix3', 'BPix2+BPix4',
0078 'BPix3+BPix4',
0079 'BPix1+FPix1_pos', 'BPix1+FPix1_neg',
0080 'BPix2+FPix1_pos', 'BPix2+FPix1_neg',
0081 'FPix1_pos+FPix2_pos', 'FPix1_neg+FPix2_neg',
0082 'FPix1_pos+FPix3_pos', 'FPix1_neg+FPix3_neg',
0083 'FPix2_pos+FPix3_pos', 'FPix2_neg+FPix3_neg',
0084
0085
0086 ]
0087 jetCoreInPhase2.toModify(jetCoreRegionalStepSeedLayers, layerList = _layerListForPhase2)
0088
0089
0090 from RecoTauTag.HLTProducers.tauRegionalPixelSeedTrackingRegions_cfi import tauRegionalPixelSeedTrackingRegions as _tauRegionalPixelSeedTrackingRegions
0091 jetCoreRegionalStepTrackingRegions = _tauRegionalPixelSeedTrackingRegions.clone(
0092 RegionPSet=dict(
0093 ptMin = 10,
0094 deltaPhiRegion = 0.20,
0095 deltaEtaRegion = 0.20,
0096 JetSrc = 'jetsForCoreTracking',
0097 vertexSrc = 'firstStepGoodPrimaryVertices',
0098 howToUseMeasurementTracker = 'Never')
0099 )
0100 jetCoreInPhase2.toModify(jetCoreRegionalStepTrackingRegions,
0101 RegionPSet=dict(
0102 deltaPhiRegion = 0.10,
0103 deltaEtaRegion = 0.10,
0104 ),
0105 )
0106 jetCoreRegionalStepEndcapTrackingRegions = jetCoreRegionalStepTrackingRegions.clone(
0107 RegionPSet=dict(
0108 JetSrc = 'jetsForCoreTrackingEndcap')
0109 )
0110
0111
0112 from RecoTracker.TkHitPairs.hitPairEDProducer_cfi import hitPairEDProducer as _hitPairEDProducer
0113 jetCoreRegionalStepHitDoublets = _hitPairEDProducer.clone(
0114 seedingLayers = 'jetCoreRegionalStepSeedLayers',
0115 trackingRegions = 'jetCoreRegionalStepTrackingRegions',
0116 produceSeedingHitSets = True,
0117 maxElementTotal = 12000000,
0118 )
0119 jetCoreRegionalStepEndcapHitDoublets = jetCoreRegionalStepHitDoublets.clone(
0120 trackingRegions = 'jetCoreRegionalStepEndcapTrackingRegions',
0121 )
0122
0123 from RecoTracker.TkSeedGenerator.seedCreatorFromRegionConsecutiveHitsEDProducer_cff import seedCreatorFromRegionConsecutiveHitsEDProducer as _seedCreatorFromRegionConsecutiveHitsEDProducer
0124 jetCoreRegionalStepSeeds = _seedCreatorFromRegionConsecutiveHitsEDProducer.clone(
0125 seedingHitSets = 'jetCoreRegionalStepHitDoublets',
0126 forceKinematicWithRegionDirection = True
0127 )
0128 import RecoTracker.TkSeedGenerator.deepCoreSeedGenerator_cfi
0129 jetCoreRegionalStepSeedsBarrel = RecoTracker.TkSeedGenerator.deepCoreSeedGenerator_cfi.deepCoreSeedGenerator.clone(
0130 vertices = "firstStepPrimaryVertices",
0131 cores = "jetsForCoreTrackingBarrel"
0132 )
0133
0134 jetCoreRegionalStepSeedsEndcap = jetCoreRegionalStepSeeds.clone(
0135 seedingHitSets = 'jetCoreRegionalStepEndcapHitDoublets',
0136 )
0137
0138
0139 import TrackingTools.TrajectoryFiltering.TrajectoryFilter_cff
0140 jetCoreRegionalStepTrajectoryFilter = TrackingTools.TrajectoryFiltering.TrajectoryFilter_cff.CkfBaseTrajectoryFilter_block.clone(
0141 minimumNumberOfHits = 4,
0142 seedPairPenalty = 0,
0143 minPt = 0.1
0144 )
0145 jetCoreRegionalStepBarrelTrajectoryFilter = jetCoreRegionalStepTrajectoryFilter.clone(
0146 minimumNumberOfHits = 2,
0147 maxConsecLostHits = 2,
0148 maxLostHitsFraction = 1.1,
0149 seedPairPenalty = 0,
0150 minPt = 0.9
0151 )
0152 jetCoreRegionalStepEndcapTrajectoryFilter = jetCoreRegionalStepTrajectoryFilter.clone()
0153
0154
0155 from Configuration.Eras.Modifier_pp_on_XeXe_2017_cff import pp_on_XeXe_2017
0156 from Configuration.ProcessModifiers.pp_on_AA_cff import pp_on_AA
0157 (pp_on_XeXe_2017 | pp_on_AA).toModify(jetCoreRegionalStepTrajectoryFilter, minPt=5.0)
0158
0159 import TrackingTools.KalmanUpdators.Chi2MeasurementEstimator_cfi
0160 jetCoreRegionalStepChi2Est = TrackingTools.KalmanUpdators.Chi2MeasurementEstimator_cfi.Chi2MeasurementEstimator.clone(
0161 ComponentName = 'jetCoreRegionalStepChi2Est',
0162 nSigma = 3.0,
0163 MaxChi2 = 30.0
0164 )
0165
0166
0167 import RecoTracker.CkfPattern.GroupedCkfTrajectoryBuilder_cfi
0168
0169 CkfBaseTrajectoryFilter_block = RecoTracker.CkfPattern.GroupedCkfTrajectoryBuilder_cfi.CkfBaseTrajectoryFilter_block
0170 jetCoreRegionalStepTrajectoryBuilder = RecoTracker.CkfPattern.GroupedCkfTrajectoryBuilder_cfi.GroupedCkfTrajectoryBuilderIterativeDefault.clone(
0171 trajectoryFilter = dict(refToPSet_ = 'jetCoreRegionalStepTrajectoryFilter'),
0172 maxCand = 50,
0173 estimator = 'jetCoreRegionalStepChi2Est',
0174 maxDPhiForLooperReconstruction = 2.0,
0175 maxPtForLooperReconstruction = 0.7,
0176 )
0177 trackingNoLoopers.toModify(jetCoreRegionalStepTrajectoryBuilder,
0178 maxPtForLooperReconstruction = 0.0)
0179 jetCoreRegionalStepBarrelTrajectoryBuilder = RecoTracker.CkfPattern.GroupedCkfTrajectoryBuilder_cfi.GroupedCkfTrajectoryBuilderIterativeDefault.clone(
0180 trajectoryFilter = dict(refToPSet_ = 'jetCoreRegionalStepBarrelTrajectoryFilter'),
0181 maxCand = 30,
0182 estimator = 'jetCoreRegionalStepChi2Est',
0183 keepOriginalIfRebuildFails = True,
0184 lockHits = False,
0185 requireSeedHitsInRebuild = False,
0186 )
0187 trackingNoLoopers.toModify(jetCoreRegionalStepBarrelTrajectoryBuilder,
0188 maxPtForLooperReconstruction = cms.double(0.0))
0189 jetCoreRegionalStepEndcapTrajectoryBuilder = jetCoreRegionalStepTrajectoryBuilder.clone(
0190 trajectoryFilter = cms.PSet(refToPSet_ = cms.string('jetCoreRegionalStepEndcapTrajectoryFilter')),
0191 maxCand = 30,
0192
0193 )
0194 trackingNoLoopers.toModify(jetCoreRegionalStepEndcapTrajectoryBuilder,
0195 maxPtForLooperReconstruction = cms.double(0.0))
0196
0197 from TrackingTools.TrajectoryCleaning.TrajectoryCleanerBySharedHits_cfi import trajectoryCleanerBySharedHits
0198 jetCoreRegionalStepDeepCoreTrajectoryCleaner = trajectoryCleanerBySharedHits.clone(
0199 ComponentName = 'jetCoreRegionalStepDeepCoreTrajectoryCleaner',
0200 fractionShared = 0.45
0201 )
0202
0203
0204
0205
0206
0207
0208
0209
0210
0211
0212
0213
0214 import RecoTracker.CkfPattern.CkfTrackCandidates_cfi
0215 jetCoreRegionalStepTrackCandidates = RecoTracker.CkfPattern.CkfTrackCandidates_cfi.ckfTrackCandidatesIterativeDefault.clone(
0216 src = 'jetCoreRegionalStepSeeds',
0217 maxSeedsBeforeCleaning = 10000,
0218 TrajectoryBuilderPSet = dict(refToPSet_ = 'jetCoreRegionalStepTrajectoryBuilder'),
0219 NavigationSchool = 'SimpleNavigationSchool',
0220
0221
0222
0223 )
0224 jetCoreRegionalStepBarrelTrackCandidates = jetCoreRegionalStepTrackCandidates.clone(
0225 src = 'jetCoreRegionalStepSeedsBarrel',
0226 TrajectoryBuilderPSet = cms.PSet( refToPSet_ = cms.string('jetCoreRegionalStepBarrelTrajectoryBuilder')),
0227
0228
0229
0230 TrajectoryCleaner = 'jetCoreRegionalStepDeepCoreTrajectoryCleaner',
0231 doSeedingRegionRebuilding = True,
0232 )
0233 jetCoreRegionalStepEndcapTrackCandidates = jetCoreRegionalStepTrackCandidates.clone(
0234 src = 'jetCoreRegionalStepSeedsEndcap',
0235 TrajectoryBuilderPSet = cms.PSet( refToPSet_ = cms.string('jetCoreRegionalStepEndcapTrajectoryBuilder')),
0236
0237
0238
0239 )
0240
0241
0242 import RecoTracker.TrackProducer.TrackProducerIterativeDefault_cfi
0243 jetCoreRegionalStepTracks = RecoTracker.TrackProducer.TrackProducerIterativeDefault_cfi.TrackProducerIterativeDefault.clone(
0244 AlgorithmName = 'jetCoreRegionalStep',
0245 src = 'jetCoreRegionalStepTrackCandidates',
0246 Fitter = 'FlexibleKFFittingSmoother'
0247 )
0248 jetCoreRegionalStepBarrelTracks = jetCoreRegionalStepTracks.clone(
0249 src = 'jetCoreRegionalStepBarrelTrackCandidates',
0250 )
0251 jetCoreRegionalStepEndcapTracks = jetCoreRegionalStepTracks.clone(
0252 src = 'jetCoreRegionalStepEndcapTrackCandidates',
0253 )
0254
0255 from Configuration.Eras.Modifier_fastSim_cff import fastSim
0256 import RecoTracker.FinalTrackSelectors.trackListMerger_cfi
0257 _fastSim_jetCoreRegionalStepTracks = RecoTracker.FinalTrackSelectors.trackListMerger_cfi.trackListMerger.clone(
0258 TrackProducers = [],
0259 hasSelector = [],
0260 selectedTrackQuals = [],
0261 copyExtras = True
0262 )
0263 fastSim.toReplaceWith(jetCoreRegionalStepTracks,_fastSim_jetCoreRegionalStepTracks)
0264 from Configuration.ProcessModifiers.seedingDeepCore_cff import seedingDeepCore
0265 (seedingDeepCore & fastSim).toReplaceWith(jetCoreRegionalStepBarrelTracks,_fastSim_jetCoreRegionalStepTracks)
0266 (seedingDeepCore & fastSim).toReplaceWith(jetCoreRegionalStepEndcapTracks,_fastSim_jetCoreRegionalStepTracks)
0267
0268
0269 from RecoTracker.FinalTrackSelectors.TrackCutClassifier_cff import *
0270 jetCoreRegionalStep = TrackCutClassifier.clone(
0271 src = 'jetCoreRegionalStepTracks',
0272 mva = dict(
0273 minPixelHits = [1,1,1],
0274 maxChi2 = [9999.,9999.,9999.],
0275 maxChi2n = [1.6,1.0,0.7],
0276 minLayers = [3,5,5],
0277 min3DLayers = [1,2,3],
0278 maxLostLayers = [4,3,2],
0279 maxDz = [0.5,0.35,0.2],
0280 maxDr = [0.3,0.2,0.1]
0281 ),
0282 vertices = 'firstStepGoodPrimaryVertices'
0283 )
0284 jetCoreRegionalStepBarrel = jetCoreRegionalStep.clone(
0285 src = 'jetCoreRegionalStepBarrelTracks',
0286 mva = dict(
0287
0288 min3DLayers = [1,2,2],
0289 ),
0290 )
0291
0292 from RecoTracker.FinalTrackSelectors.TrackMVAClassifierPrompt_cfi import *
0293 trackingPhase1.toReplaceWith(jetCoreRegionalStep, TrackMVAClassifierPrompt.clone(
0294 mva = dict(GBRForestLabel = 'MVASelectorJetCoreRegionalStep_Phase1'),
0295 src = 'jetCoreRegionalStepTracks',
0296 qualityCuts = [-0.2,0.0,0.4]
0297 ))
0298
0299 trackingPhase1.toReplaceWith(jetCoreRegionalStepBarrel, jetCoreRegionalStep.clone(
0300 src = 'jetCoreRegionalStepBarrelTracks',
0301 ))
0302
0303 pp_on_AA.toModify(jetCoreRegionalStep, qualityCuts = [-0.2, 0.0, 0.8])
0304 pp_on_AA.toModify(jetCoreRegionalStepBarrel, qualityCuts = [-0.2, 0.0, 0.8])
0305
0306 from RecoTracker.FinalTrackSelectors.trackTfClassifier_cfi import *
0307 from RecoTracker.FinalTrackSelectors.trackSelectionTf_cfi import *
0308 from RecoTracker.FinalTrackSelectors.trackSelectionTf_CKF_cfi import *
0309 trackdnn.toReplaceWith(jetCoreRegionalStep, trackTfClassifier.clone(
0310 src = 'jetCoreRegionalStepTracks',
0311 qualityCuts = qualityCutDictionary.JetCoreRegionalStep.value()
0312 ))
0313 trackdnn.toReplaceWith(jetCoreRegionalStepBarrel, trackTfClassifier.clone(
0314 src = 'jetCoreRegionalStepBarrelTracks',
0315 qualityCuts = qualityCutDictionary.JetCoreRegionalStep.value()
0316 ))
0317
0318 fastSim.toModify(jetCoreRegionalStep,vertices = 'firstStepPrimaryVerticesBeforeMixing')
0319
0320
0321 jetCoreRegionalStepEndcap = jetCoreRegionalStep.clone(
0322 src = 'jetCoreRegionalStepEndcapTracks',
0323 )
0324
0325
0326 JetCoreRegionalStepTask = cms.Task(jetsForCoreTracking,
0327 firstStepGoodPrimaryVertices,
0328
0329 jetCoreRegionalStepSeedLayers,
0330 jetCoreRegionalStepTrackingRegions,
0331 jetCoreRegionalStepHitDoublets,
0332 jetCoreRegionalStepSeeds,
0333 jetCoreRegionalStepTrackCandidates,
0334 jetCoreRegionalStepTracks,
0335
0336 jetCoreRegionalStep)
0337 JetCoreRegionalStep = cms.Sequence(JetCoreRegionalStepTask)
0338
0339 JetCoreRegionalStepBarrelTask = cms.Task(jetsForCoreTrackingBarrel,
0340 firstStepGoodPrimaryVertices,
0341
0342 jetCoreRegionalStepSeedLayers,
0343 jetCoreRegionalStepSeedsBarrel,
0344 jetCoreRegionalStepBarrelTrackCandidates,
0345 jetCoreRegionalStepBarrelTracks,
0346 jetCoreRegionalStepBarrel)
0347
0348 JetCoreRegionalStepEndcapTask = cms.Task(jetsForCoreTrackingEndcap,
0349 firstStepGoodPrimaryVertices,
0350
0351 jetCoreRegionalStepSeedLayers,
0352 jetCoreRegionalStepEndcapTrackingRegions,
0353 jetCoreRegionalStepEndcapHitDoublets,
0354 jetCoreRegionalStepSeedsEndcap,
0355 jetCoreRegionalStepEndcapTrackCandidates,
0356 jetCoreRegionalStepEndcapTracks,
0357 jetCoreRegionalStepEndcap)
0358
0359
0360 import RecoTracker.FinalTrackSelectors.multiTrackSelector_cfi
0361 jetCoreRegionalStepSelector = RecoTracker.FinalTrackSelectors.multiTrackSelector_cfi.multiTrackSelector.clone(
0362 src = 'jetCoreRegionalStepTracks',
0363 trackSelectors= cms.VPSet(
0364 RecoTracker.FinalTrackSelectors.multiTrackSelector_cfi.looseMTS.clone(
0365 name = 'jetCoreRegionalStepLoose',
0366 chi2n_par = 2.0,
0367 res_par = ( 0.003, 0.002 ),
0368 minNumberLayers = 3,
0369 maxNumberLostLayers = 3,
0370 minNumber3DLayers = 3,
0371 d0_par1 = ( 0.8, 4.0 ),
0372 dz_par1 = ( 0.9, 4.0 ),
0373 d0_par2 = ( 0.6, 4.0 ),
0374 dz_par2 = ( 0.8, 4.0 )
0375 ),
0376 RecoTracker.FinalTrackSelectors.multiTrackSelector_cfi.tightMTS.clone(
0377 name = 'jetCoreRegionalStepTight',
0378 preFilterName = 'jetCoreRegionalStepLoose',
0379 chi2n_par = 1.4,
0380 res_par = ( 0.003, 0.002 ),
0381 minNumberLayers = 3,
0382 maxNumberLostLayers = 2,
0383 minNumber3DLayers = 3,
0384 d0_par1 = ( 0.7, 4.0 ),
0385 dz_par1 = ( 0.8, 4.0 ),
0386 d0_par2 = ( 0.5, 4.0 ),
0387 dz_par2 = ( 0.7, 4.0 )
0388 ),
0389 RecoTracker.FinalTrackSelectors.multiTrackSelector_cfi.highpurityMTS.clone(
0390 name = 'jetCoreRegionalStep',
0391 preFilterName = 'jetCoreRegionalStepTight',
0392 min_eta = -4.1,
0393 max_eta = 4.1,
0394 chi2n_par = 1.2,
0395 res_par = ( 0.003, 0.001 ),
0396 minNumberLayers = 3,
0397 maxNumberLostLayers = 2,
0398 minNumber3DLayers = 3,
0399 d0_par1 = ( 0.6, 4.0 ),
0400 dz_par1 = ( 0.7, 4.0 ),
0401 d0_par2 = ( 0.45, 4.0 ),
0402 dz_par2 = ( 0.55, 4.0 )
0403 ),
0404 ),
0405 )
0406 _JetCoreRegionalStepTask_trackingPhase2 = JetCoreRegionalStepTask.copy()
0407 _JetCoreRegionalStepTask_trackingPhase2.replace(jetCoreRegionalStep, jetCoreRegionalStepSelector)
0408 jetCoreInPhase2.toReplaceWith(JetCoreRegionalStepTask, _JetCoreRegionalStepTask_trackingPhase2)
0409 jetCoreInPhase2.toModify(jetCoreRegionalStepTracks, TrajectoryInEvent = True)
0410 jetCoreInPhase2.toModify(jetCoreRegionalStepTrajectoryBuilder, maxCand = 5)
0411
0412 from RecoTracker.FinalTrackSelectors.TrackCollectionMerger_cfi import *
0413 (seedingDeepCore & ~fastSim).toReplaceWith(jetCoreRegionalStepTracks, TrackCollectionMerger.clone(
0414 trackProducers = ["jetCoreRegionalStepBarrelTracks",
0415 "jetCoreRegionalStepEndcapTracks",],
0416 inputClassifiers = ["jetCoreRegionalStepBarrel",
0417 "jetCoreRegionalStepEndcap",],
0418 foundHitBonus = 100.0,
0419 lostHitPenalty = 1.0
0420 ))
0421
0422 seedingDeepCore.toReplaceWith(jetCoreRegionalStep, jetCoreRegionalStepTracks.clone())
0423
0424 seedingDeepCore.toReplaceWith(JetCoreRegionalStepTask, cms.Task(
0425 JetCoreRegionalStepBarrelTask,
0426 JetCoreRegionalStepEndcapTask,
0427 cms.Task(jetCoreRegionalStepTracks,jetCoreRegionalStep)
0428 ))
0429
0430
0431 fastSim.toReplaceWith(JetCoreRegionalStepTask,
0432 cms.Task(jetCoreRegionalStepTracks,
0433 jetCoreRegionalStep))
0434 (seedingDeepCore & fastSim).toReplaceWith(JetCoreRegionalStepTask,
0435 cms.Task(jetCoreRegionalStepBarrelTracks, jetCoreRegionalStepEndcapTracks,
0436 jetCoreRegionalStepTracks,
0437 jetCoreRegionalStepBarrel, jetCoreRegionalStepEndcap,
0438 jetCoreRegionalStep))