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
11 trackerClusterCheckPreSplitting = _trackerClusterCheck.clone(
12 PixelClusterCollectionLabel =
'siPixelClustersPreSplitting' 17 import RecoTracker.TkSeedingLayers.PixelLayerQuadruplets_cfi
18 initialStepSeedLayersPreSplitting = RecoTracker.TkSeedingLayers.PixelLayerTriplets_cfi.PixelLayerTriplets.clone(
19 FPix = dict(HitProducer =
'siPixelRecHitsPreSplitting'),
20 BPix = dict(HitProducer =
'siPixelRecHitsPreSplitting')
22 from Configuration.Eras.Modifier_trackingPhase1_cff
import trackingPhase1
23 trackingPhase1.toModify(initialStepSeedLayersPreSplitting,
24 layerList = RecoTracker.TkSeedingLayers.PixelLayerQuadruplets_cfi.PixelLayerQuadruplets.layerList.value()
28 from RecoTracker.TkTrackingRegions.globalTrackingRegionFromBeamSpot_cfi
import globalTrackingRegionFromBeamSpot
as _globalTrackingRegionFromBeamSpot
29 initialStepTrackingRegionsPreSplitting = _globalTrackingRegionFromBeamSpot.clone(RegionPSet = dict(
34 trackingPhase1.toModify(initialStepTrackingRegionsPreSplitting, RegionPSet = dict(ptMin = 0.5))
37 from RecoTracker.TkHitPairs.hitPairEDProducer_cfi
import hitPairEDProducer
as _hitPairEDProducer
38 initialStepHitDoubletsPreSplitting = _hitPairEDProducer.clone(
39 seedingLayers =
'initialStepSeedLayersPreSplitting',
40 trackingRegions =
'initialStepTrackingRegionsPreSplitting',
41 clusterCheck =
'trackerClusterCheckPreSplitting',
42 maxElement = 50000000,
43 produceIntermediateHitDoublets =
True,
45 from RecoTracker.PixelSeeding.pixelTripletHLTEDProducer_cfi
import pixelTripletHLTEDProducer
as _pixelTripletHLTEDProducer
47 import RecoTracker.PixelLowPtUtilities.LowPtClusterShapeSeedComparitor_cfi
48 initialStepHitTripletsPreSplitting = _pixelTripletHLTEDProducer.clone(
49 doublets =
'initialStepHitDoubletsPreSplitting',
50 produceSeedingHitSets =
True,
51 SeedComparitorPSet = RecoTracker.PixelLowPtUtilities.LowPtClusterShapeSeedComparitor_cfi.LowPtClusterShapeSeedComparitor.clone(
52 clusterShapeCacheSrc =
'siPixelClusterShapeCachePreSplitting' 55 from RecoTracker.PixelSeeding.caHitQuadrupletEDProducer_cfi
import caHitQuadrupletEDProducer
as _caHitQuadrupletEDProducer
56 trackingPhase1.toModify(initialStepHitDoubletsPreSplitting, layerPairs = [0,1,2])
57 initialStepHitQuadrupletsPreSplitting = _caHitQuadrupletEDProducer.clone(
58 doublets =
'initialStepHitDoubletsPreSplitting',
59 extraHitRPhitolerance = initialStepHitTripletsPreSplitting.extraHitRPhitolerance,
60 SeedComparitorPSet = initialStepHitTripletsPreSplitting.SeedComparitorPSet,
63 value1 = 200, value2 = 50,
65 useBendingCorrection =
True,
67 fitFastCircleChi2Cut =
True,
71 from RecoTracker.TkSeedGenerator.seedCreatorFromRegionConsecutiveHitsEDProducer_cff
import seedCreatorFromRegionConsecutiveHitsEDProducer
as _seedCreatorFromRegionConsecutiveHitsEDProducer
72 initialStepSeedsPreSplitting = _seedCreatorFromRegionConsecutiveHitsEDProducer.clone(
73 seedingHitSets =
'initialStepHitTripletsPreSplitting',
75 trackingPhase1.toModify(initialStepSeedsPreSplitting, seedingHitSets =
'initialStepHitQuadrupletsPreSplitting')
80 initialStepTrajectoryFilterBasePreSplitting = TrackingTools.TrajectoryFiltering.TrajectoryFilter_cff.CkfBaseTrajectoryFilter_block.clone(
81 minimumNumberOfHits = 4,
84 minGoodStripCharge = cms.PSet(refToPSet_ = cms.string(
'SiStripClusterChargeCutLoose'))
86 from Configuration.Eras.Modifier_tracker_apv_vfp30_2016_cff
import tracker_apv_vfp30_2016
87 _tracker_apv_vfp30_2016.toModify(initialStepTrajectoryFilterBasePreSplitting, maxCCCLostHits = 2)
88 import RecoTracker.PixelLowPtUtilities.StripSubClusterShapeTrajectoryFilter_cfi
89 initialStepTrajectoryFilterShapePreSplitting = RecoTracker.PixelLowPtUtilities.StripSubClusterShapeTrajectoryFilter_cfi.StripSubClusterShapeTrajectoryFilterTIX12.clone()
90 initialStepTrajectoryFilterPreSplitting = cms.PSet(
91 ComponentType = cms.string(
'CompositeTrajectoryFilter'),
93 cms.PSet( refToPSet_ = cms.string(
'initialStepTrajectoryFilterBasePreSplitting')),
94 cms.PSet( refToPSet_ = cms.string(
'initialStepTrajectoryFilterShapePreSplitting'))),
97 import RecoTracker.MeasurementDet.Chi2ChargeMeasurementEstimator_cfi
98 initialStepChi2EstPreSplitting = RecoTracker.MeasurementDet.Chi2ChargeMeasurementEstimator_cfi.Chi2ChargeMeasurementEstimator.clone(
99 ComponentName =
'initialStepChi2EstPreSplitting',
102 clusterChargeCut = cms.PSet(refToPSet_ = cms.string(
'SiStripClusterChargeCutLoose')),
104 _tracker_apv_vfp30_2016.toModify(initialStepChi2EstPreSplitting,
105 clusterChargeCut = dict(refToPSet_ =
'SiStripClusterChargeCutTiny')
109 initialStepTrajectoryBuilderPreSplitting = RecoTracker.CkfPattern.GroupedCkfTrajectoryBuilder_cfi.GroupedCkfTrajectoryBuilderIterativeDefault.clone(
110 trajectoryFilter = dict(refToPSet_ =
'initialStepTrajectoryFilterPreSplitting'),
111 alwaysUseInvalidHits =
True,
113 estimator =
'initialStepChi2Est',
117 initialStepTrackCandidatesPreSplitting = RecoTracker.CkfPattern.CkfTrackCandidates_cfi.ckfTrackCandidatesIterativeDefault.clone(
118 src =
'initialStepSeedsPreSplitting',
120 numHitsForSeedCleaner = 50,
121 onlyPixelHitsForSeedCleaner =
True,
122 TrajectoryBuilderPSet = dict(refToPSet_ =
'initialStepTrajectoryBuilderPreSplitting'),
123 doSeedingRegionRebuilding =
True,
124 useHitsSplitting =
True,
125 MeasurementTrackerEvent =
'MeasurementTrackerEventPreSplitting',
129 from RecoTracker.MkFit.mkFitGeometryESProducer_cfi
import mkFitGeometryESProducer
130 import RecoTracker.MkFit.mkFitSiPixelHitConverter_cfi
as _mkFitSiPixelHitConverter_cfi
131 import RecoTracker.MkFit.mkFitSiStripHitConverter_cfi
as _mkFitSiStripHitConverter_cfi
132 import RecoTracker.MkFit.mkFitEventOfHitsProducer_cfi
as _mkFitEventOfHitsProducer_cfi
133 import RecoTracker.MkFit.mkFitSeedConverter_cfi
as _mkFitSeedConverter_cfi
134 import RecoTracker.MkFit.mkFitIterationConfigESProducer_cfi
as _mkFitIterationConfigESProducer_cfi
135 import RecoTracker.MkFit.mkFitProducer_cfi
as _mkFitProducer_cfi
136 import RecoTracker.MkFit.mkFitOutputConverter_cfi
as _mkFitOutputConverter_cfi
137 mkFitSiPixelHitsPreSplitting = _mkFitSiPixelHitConverter_cfi.mkFitSiPixelHitConverter.clone(
138 hits =
'siPixelRecHitsPreSplitting' 140 mkFitSiStripHits = _mkFitSiStripHitConverter_cfi.mkFitSiStripHitConverter.clone()
141 mkFitEventOfHitsPreSplitting = _mkFitEventOfHitsProducer_cfi.mkFitEventOfHitsProducer.clone(
142 pixelHits =
'mkFitSiPixelHitsPreSplitting' 144 initialStepTrackCandidatesMkFitSeedsPreSplitting = _mkFitSeedConverter_cfi.mkFitSeedConverter.clone(
145 seeds =
'initialStepSeedsPreSplitting',
148 initialStepTrackCandidatesMkFitConfigPreSplitting = _mkFitIterationConfigESProducer_cfi.mkFitIterationConfigESProducer.clone(
149 ComponentName =
'initialStepTrackCandidatesMkFitConfigPreSplitting',
150 config =
'RecoTracker/MkFit/data/mkfit-phase1-initialStep.json',
152 initialStepTrackCandidatesMkFitPreSplitting = _mkFitProducer_cfi.mkFitProducer.clone(
153 pixelHits =
'mkFitSiPixelHitsPreSplitting',
154 eventOfHits =
'mkFitEventOfHitsPreSplitting',
155 seeds =
'initialStepTrackCandidatesMkFitSeedsPreSplitting',
156 config = (
'',
'initialStepTrackCandidatesMkFitConfigPreSplitting'),
158 trackingMkFitInitialStepPreSplitting.toReplaceWith(initialStepTrackCandidatesPreSplitting, _mkFitOutputConverter_cfi.mkFitOutputConverter.clone(
159 mkFitPixelHits =
'mkFitSiPixelHitsPreSplitting',
160 mkFitEventOfHits =
'mkFitEventOfHitsPreSplitting',
161 seeds =
'initialStepSeedsPreSplitting',
162 mkFitSeeds =
'initialStepTrackCandidatesMkFitSeedsPreSplitting',
163 tracks =
'initialStepTrackCandidatesMkFitPreSplitting',
167 import RecoTracker.TrackProducer.TrackProducerIterativeDefault_cfi
168 initialStepTracksPreSplitting = RecoTracker.TrackProducer.TrackProducerIterativeDefault_cfi.TrackProducerIterativeDefault.clone(
169 src =
'initialStepTrackCandidatesPreSplitting',
170 AlgorithmName =
'initialStep',
171 Fitter =
'FlexibleKFFittingSmoother',
172 NavigationSchool =
'',
173 MeasurementTrackerEvent =
'' 175 initialStepTracksPreSplitting.MeasurementTrackerEvent =
'MeasurementTrackerEventPreSplitting' 179 firstStepPrimaryVerticesPreSplitting = _offlinePrimaryVertices.clone(
180 TrackLabel =
'initialStepTracksPreSplitting',
181 vertexCollections = [_offlinePrimaryVertices.vertexCollections[0].
clone()]
183 from Configuration.Eras.Modifier_pp_on_XeXe_2017_cff
import pp_on_XeXe_2017
184 from Configuration.ProcessModifiers.pp_on_AA_cff
import pp_on_AA
185 (pp_on_XeXe_2017 | pp_on_AA).toModify(firstStepPrimaryVerticesPreSplitting,
186 TkFilterParameters = dict(
187 trackQuality =
'any',
188 maxNumTracksThreshold = 2**31-1
195 initialStepTrackRefsForJetsPreSplitting = initialStepTrackRefsForJets.clone(
196 src =
'initialStepTracksPreSplitting' 198 caloTowerForTrkPreSplitting = caloTowerForTrk.clone()
199 ak4CaloJetsForTrkPreSplitting = ak4CaloJetsForTrk.clone(
200 src =
'caloTowerForTrkPreSplitting',
201 srcPVs =
'firstStepPrimaryVerticesPreSplitting' 203 jetsForCoreTrackingPreSplitting = jetsForCoreTracking.clone(
204 src =
'ak4CaloJetsForTrkPreSplitting' 208 from RecoLocalTracker.SubCollectionProducers.jetCoreClusterSplitter_cfi
import jetCoreClusterSplitter
209 siPixelClusters = jetCoreClusterSplitter.clone(
210 pixelClusters =
'siPixelClustersPreSplitting',
211 vertices =
'firstStepPrimaryVerticesPreSplitting',
212 cores =
'jetsForCoreTrackingPreSplitting' 219 InitialStepPreSplittingTask = cms.Task(trackerClusterCheckPreSplitting,
220 initialStepSeedLayersPreSplitting,
221 initialStepTrackingRegionsPreSplitting,
222 initialStepHitDoubletsPreSplitting,
223 initialStepHitTripletsPreSplitting,
224 initialStepSeedsPreSplitting,
225 initialStepTrackCandidatesPreSplitting,
226 initialStepTracksPreSplitting,
227 firstStepPrimaryVerticesPreSplitting,
228 initialStepTrackRefsForJetsPreSplitting,
229 caloTowerForTrkPreSplitting,
230 ak4CaloJetsForTrkPreSplitting,
231 jetsForCoreTrackingPreSplitting,
234 MeasurementTrackerEvent,
235 siPixelClusterShapeCache)
236 InitialStepPreSplitting = cms.Sequence(InitialStepPreSplittingTask)
237 _InitialStepPreSplittingTask_trackingPhase1 = InitialStepPreSplittingTask.copy()
238 _InitialStepPreSplittingTask_trackingPhase1.replace(initialStepHitTripletsPreSplitting, cms.Task(initialStepHitTripletsPreSplitting,initialStepHitQuadrupletsPreSplitting))
239 trackingPhase1.toReplaceWith(InitialStepPreSplittingTask, _InitialStepPreSplittingTask_trackingPhase1.copyAndExclude([initialStepHitTripletsPreSplitting]))
242 _InitialStepPreSplittingTask_trackingMkFitCommon = InitialStepPreSplittingTask.copy()
243 _InitialStepPreSplittingTask_trackingMkFitCommon.add(mkFitSiStripHits, mkFitGeometryESProducer)
244 trackingMkFitCommon.toReplaceWith(InitialStepPreSplittingTask, _InitialStepPreSplittingTask_trackingMkFitCommon)
245 _InitialStepPreSplittingTask_trackingMkFit = InitialStepPreSplittingTask.copy()
246 _InitialStepPreSplittingTask_trackingMkFit.add(mkFitSiPixelHitsPreSplitting, mkFitEventOfHitsPreSplitting, initialStepTrackCandidatesMkFitSeedsPreSplitting, initialStepTrackCandidatesMkFitPreSplitting, initialStepTrackCandidatesMkFitConfigPreSplitting)
247 trackingMkFitInitialStepPreSplitting.toReplaceWith(InitialStepPreSplittingTask, _InitialStepPreSplittingTask_trackingMkFit)
265 from Configuration.Eras.Modifier_trackingLowPU_cff
import trackingLowPU
266 trackingLowPU.toReplaceWith(siPixelClusters, _siPixelClusters)
267 from Configuration.Eras.Modifier_trackingPhase2PU140_cff
import trackingPhase2PU140
268 trackingPhase2PU140.toReplaceWith(siPixelClusters, _siPixelClusters)
269 _InitialStepPreSplittingTask_LowPU_Phase2PU140 = cms.Task(
272 MeasurementTrackerEvent ,
273 siPixelClusterShapeCache
275 trackingLowPU.toReplaceWith(InitialStepPreSplittingTask, _InitialStepPreSplittingTask_LowPU_Phase2PU140)
276 trackingPhase2PU140.toReplaceWith(InitialStepPreSplittingTask, _InitialStepPreSplittingTask_LowPU_Phase2PU140)
TEveGeoShape * clone(const TEveElement *element, TEveElement *parent)