1 import FWCore.ParameterSet.Config
as cms
5 hiLowPtTripletStepClusters = cms.EDProducer(
"HITrackClusterRemover",
6 clusterLessSolution= cms.bool(
True),
7 oldClusterRemovalInfo = cms.InputTag(
"hiDetachedTripletStepClusters"),
8 trajectories = cms.InputTag(
"hiDetachedTripletStepTracks"),
9 overrideTrkQuals = cms.InputTag(
"hiDetachedTripletStepSelector",
"hiDetachedTripletStep"),
10 TrackQuality = cms.string(
'highPurity'),
11 minNumberOfLayersWithMeasBeforeFiltering = cms.int32(0),
12 pixelClusters = cms.InputTag(
"siPixelClusters"),
13 stripClusters = cms.InputTag(
"siStripClusters"),
15 maxChi2 = cms.double(9.0)
19 maxSize = cms.uint32(2),
20 maxChi2 = cms.double(9.0)
27 hiLowPtTripletStepSeedLayers = RecoTracker.TkSeedingLayers.PixelLayerTriplets_cfi.PixelLayerTriplets.clone(
28 BPix = dict(skipClusters = cms.InputTag(
'hiLowPtTripletStepClusters')),
29 FPix = dict(skipClusters = cms.InputTag(
'hiLowPtTripletStepClusters'))
34 from RecoTracker.TkHitPairs.hitPairEDProducer_cfi
import hitPairEDProducer
as _hitPairEDProducer
35 from RecoTracker.PixelSeeding.pixelTripletHLTEDProducer_cfi
import pixelTripletHLTEDProducer
as _pixelTripletHLTEDProducer
37 import RecoTracker.PixelLowPtUtilities.LowPtClusterShapeSeedComparitor_cfi
43 hiLowPtTripletStepTrackingRegions = _globalTrackingRegionWithVertices.clone(RegionPSet=dict(
45 useMultipleScattering =
False,
46 useFakeVertices =
False,
47 beamSpot =
"offlineBeamSpot",
48 useFixedError =
False,
52 VertexCollection =
"hiSelectedPixelVertex",
54 useFoundVertices =
True,
57 hiLowPtTripletStepTracksHitDoublets = _hitPairEDProducer.clone(
59 seedingLayers =
"hiLowPtTripletStepSeedLayers",
60 trackingRegions =
"hiLowPtTripletStepTrackingRegions",
61 maxElement = 50000000,
62 produceIntermediateHitDoublets =
True,
64 import RecoTracker.PixelLowPtUtilities.LowPtClusterShapeSeedComparitor_cfi
65 hiLowPtTripletStepTracksHitTriplets = _pixelTripletHLTEDProducer.clone(
66 doublets =
"hiLowPtTripletStepTracksHitDoublets",
68 SeedComparitorPSet = RecoTracker.PixelLowPtUtilities.LowPtClusterShapeSeedComparitor_cfi.LowPtClusterShapeSeedComparitor.clone(),
69 produceSeedingHitSets =
True,
72 from RecoTracker.PixelSeeding.caHitTripletEDProducer_cfi
import caHitTripletEDProducer
as _caHitTripletEDProducer
73 hiLowPtTripletStepTracksHitDoubletsCA = hiLowPtTripletStepTracksHitDoublets.clone(
76 hiLowPtTripletStepTracksHitTripletsCA = _caHitTripletEDProducer.clone(
77 doublets =
"hiLowPtTripletStepTracksHitDoubletsCA",
78 extraHitRPhitolerance = hiLowPtTripletStepTracksHitTriplets.extraHitRPhitolerance,
79 SeedComparitorPSet = hiLowPtTripletStepTracksHitTriplets.SeedComparitorPSet,
82 value1 = 70 , value2 = 8,
84 useBendingCorrection =
True,
89 hiLowPtTripletStepPixelTracksFilter = hiFilter.clone(
90 nSigmaLipMaxTolerance = 4.0,
91 nSigmaTipMaxTolerance = 4.0,
96 import RecoTracker.PixelTrackFitting.pixelTracks_cfi
as _mod
98 hiLowPtTripletStepPixelTracks = _mod.pixelTracks.clone(
99 passLabel =
'Pixel primary tracks with vertex constraint',
101 SeedingHitSets =
"hiLowPtTripletStepTracksHitTriplets",
103 Fitter =
"pixelFitterByHelixProjections",
105 Filter =
"hiLowPtTripletStepPixelTracksFilter",
107 Cleaner =
"trackCleaner" 109 from Configuration.Eras.Modifier_trackingPhase1_cff
import trackingPhase1
110 trackingPhase1.toModify(hiLowPtTripletStepPixelTracks,
111 SeedingHitSets =
"hiLowPtTripletStepTracksHitTripletsCA" 115 import RecoTracker.PixelLowPtUtilities.TrackSeeds_cfi
116 hiLowPtTripletStepSeeds = RecoTracker.PixelLowPtUtilities.TrackSeeds_cfi.pixelTrackSeeds.clone(
117 InputCollection =
'hiLowPtTripletStepPixelTracks' 123 hiLowPtTripletStepTrajectoryFilter = TrackingTools.TrajectoryFiltering.TrajectoryFilter_cff.CkfBaseTrajectoryFilter_block.clone(
125 minimumNumberOfHits = 6,
130 hiLowPtTripletStepChi2Est = TrackingTools.KalmanUpdators.Chi2MeasurementEstimator_cfi.Chi2MeasurementEstimator.clone(
131 ComponentName =
'hiLowPtTripletStepChi2Est',
138 hiLowPtTripletStepTrajectoryBuilder = RecoTracker.CkfPattern.GroupedCkfTrajectoryBuilder_cfi.GroupedCkfTrajectoryBuilder.clone(
139 trajectoryFilter = dict(refToPSet_ =
'hiLowPtTripletStepTrajectoryFilter'),
141 estimator =
'hiLowPtTripletStepChi2Est',
142 maxDPhiForLooperReconstruction = 2.0,
145 maxPtForLooperReconstruction = 0.7,
150 hiLowPtTripletStepTrackCandidates = RecoTracker.CkfPattern.CkfTrackCandidates_cfi.ckfTrackCandidates.clone(
151 src =
'hiLowPtTripletStepSeeds',
153 numHitsForSeedCleaner = 50,
154 onlyPixelHitsForSeedCleaner =
True,
155 TrajectoryBuilderPSet = dict(refToPSet_ =
'hiLowPtTripletStepTrajectoryBuilder'),
156 clustersToSkip =
'hiLowPtTripletStepClusters',
157 doSeedingRegionRebuilding =
True,
158 useHitsSplitting =
True 163 hiLowPtTripletStepTracks = RecoTracker.TrackProducer.TrackProducer_cfi.TrackProducer.clone(
164 src =
'hiLowPtTripletStepTrackCandidates',
165 AlgorithmName =
'lowPtTripletStep',
166 Fitter=
'FlexibleKFFittingSmoother' 173 hiLowPtTripletStepSelector = RecoHI.HiTracking.hiMultiTrackSelector_cfi.hiMultiTrackSelector.clone(
174 src =
'hiLowPtTripletStepTracks',
176 GBRForestLabel =
'HIMVASelectorIter5',
177 GBRForestVars = [
'chi2perdofperlayer',
'dxyperdxyerror',
'dzperdzerror',
'relpterr',
'nhits',
'nlayers',
'eta'],
178 trackSelectors= cms.VPSet(
179 RecoHI.HiTracking.hiMultiTrackSelector_cfi.hiLooseMTS.clone(
180 name =
'hiLowPtTripletStepLoose',
183 RecoHI.HiTracking.hiMultiTrackSelector_cfi.hiTightMTS.clone(
184 name =
'hiLowPtTripletStepTight',
185 preFilterName =
'hiLowPtTripletStepLoose',
189 RecoHI.HiTracking.hiMultiTrackSelector_cfi.hiHighpurityMTS.clone(
190 name =
'hiLowPtTripletStep',
191 preFilterName =
'hiLowPtTripletStepTight',
197 from Configuration.Eras.Modifier_trackingPhase1_cff
import trackingPhase1
198 trackingPhase1.toModify(hiLowPtTripletStepSelector, useAnyMVA =
False)
199 trackingPhase1.toModify(hiLowPtTripletStepSelector, trackSelectors= cms.VPSet(
200 RecoHI.HiTracking.hiMultiTrackSelector_cfi.hiLooseMTS.clone(
201 name =
'hiLowPtTripletStepLoose',
204 RecoHI.HiTracking.hiMultiTrackSelector_cfi.hiTightMTS.clone(
205 name =
'hiLowPtTripletStepTight',
206 preFilterName =
'hiLowPtTripletStepLoose',
210 RecoHI.HiTracking.hiMultiTrackSelector_cfi.hiHighpurityMTS.clone(
211 name =
'hiLowPtTripletStep',
212 preFilterName =
'hiLowPtTripletStepTight',
221 hiLowPtTripletStepQual = RecoTracker.FinalTrackSelectors.trackListMerger_cfi.trackListMerger.clone(
222 TrackProducers = [
'hiLowPtTripletStepTracks'],
224 selectedTrackQuals = [
"hiLowPtTripletStepSelector:hiLowPtTripletStep"],
226 makeReKeyedSeeds = cms.untracked.bool(
False),
232 hiLowPtTripletStepTask = cms.Task(hiLowPtTripletStepClusters,
233 hiLowPtTripletStepSeedLayers,
234 hiLowPtTripletStepTrackingRegions,
235 hiLowPtTripletStepTracksHitDoublets,
236 hiLowPtTripletStepTracksHitTriplets,
237 pixelFitterByHelixProjections,
238 hiLowPtTripletStepPixelTracksFilter,
239 hiLowPtTripletStepPixelTracks,hiLowPtTripletStepSeeds,
240 hiLowPtTripletStepTrackCandidates,
241 hiLowPtTripletStepTracks,
242 hiLowPtTripletStepSelector,
243 hiLowPtTripletStepQual
245 hiLowPtTripletStep = cms.Sequence(hiLowPtTripletStepTask)
246 hiLowPtTripletStepTask_Phase1 = hiLowPtTripletStepTask.copy()
247 hiLowPtTripletStepTask_Phase1.replace(hiLowPtTripletStepTracksHitDoublets, hiLowPtTripletStepTracksHitDoubletsCA)
248 hiLowPtTripletStepTask_Phase1.replace(hiLowPtTripletStepTracksHitTriplets, hiLowPtTripletStepTracksHitTripletsCA)
249 trackingPhase1.toReplaceWith(hiLowPtTripletStepTask, hiLowPtTripletStepTask_Phase1)
pp iterative tracking modified for hiOffline reco (the vertex is the one reconstructed in HI) 3rd ste...