CMS 3D CMS Logo

DetachedTripletStep_cff.py
Go to the documentation of this file.
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
3 import RecoTracker.IterativeTracking.iterativeTkConfig as _cfg
4 
5 #for fastsim
6 from Configuration.Eras.Modifier_fastSim_cff import fastSim
7 from FastSimulation.Tracking.SeedingMigration import _hitSetProducerToFactoryPSet
8 
9 
10 ###############################################
11 # Low pT and detached tracks from pixel triplets
12 ###############################################
13 
14 # REMOVE HITS ASSIGNED TO GOOD TRACKS FROM PREVIOUS ITERATIONS
15 detachedTripletStepClusters = _cfg.clusterRemoverForIter("DetachedTripletStep")
16 for _eraName, _postfix, _era in _cfg.nonDefaultEras():
17  _era.toReplaceWith(detachedTripletStepClusters, _cfg.clusterRemoverForIter("DetachedTripletStep", _eraName, _postfix))
18 
19 # SEEDING LAYERS
21 detachedTripletStepSeedLayers = RecoTracker.TkSeedingLayers.PixelLayerTriplets_cfi.PixelLayerTriplets.clone()
22 detachedTripletStepSeedLayers.BPix.skipClusters = cms.InputTag('detachedTripletStepClusters')
23 detachedTripletStepSeedLayers.FPix.skipClusters = cms.InputTag('detachedTripletStepClusters')
24 _phase1LayerList = [
25  'BPix1+BPix2+BPix3',
26  'BPix2+BPix3+BPix4',
27 # 'BPix1+BPix3+BPix4', # has "hole", not tested
28 # 'BPix1+BPix2+BPix4', # has "hole", not tested
29  'BPix2+BPix3+FPix1_pos', 'BPix2+BPix3+FPix1_neg',
30 # 'BPix1+BPix2+FPix1_pos', 'BPix1+BPix2+FPix1_neg', # mostly fake tracks, lots of seeds
31 # 'BPix1+BPix3+FPix1_pos', 'BPix1+BPix3+FPix1_neg', # has "hole", not tested
32  'BPix2+FPix1_pos+FPix2_pos', 'BPix2+FPix1_neg+FPix2_neg',
33 # 'BPix1+FPix1_pos+FPix2_pos', 'BPix1+FPix1_neg+FPix2_neg', # mostly fake tracks, lots of seeds
34 # 'BPix1+BPix2+FPix2_pos', 'BPix1+BPix2+FPix2_neg', # has "hole", not tested
35  'FPix1_pos+FPix2_pos+FPix3_pos', 'FPix1_neg+FPix2_neg+FPix3_neg',
36 # 'BPix1+FPix2_pos+FPix3_pos', 'BPix1+FPix2_neg+FPix3_neg', # has "hole", not tested
37 # 'BPix1+FPix1_pos+FPix3_pos', 'BPix1+FPix1_neg+FPix3_neg' # has "hole", not tested
38  ]
39 from Configuration.Eras.Modifier_trackingPhase1_cff import trackingPhase1
40 trackingPhase1.toModify(detachedTripletStepSeedLayers, layerList=_phase1LayerList)
41 
42 # TrackingRegion
43 from RecoTracker.TkTrackingRegions.globalTrackingRegionFromBeamSpotFixedZ_cfi import globalTrackingRegionFromBeamSpotFixedZ as _globalTrackingRegionFromBeamSpotFixedZ
44 detachedTripletStepTrackingRegions = _globalTrackingRegionFromBeamSpotFixedZ.clone(RegionPSet = dict(
45  ptMin = 0.3,
46  originHalfLength = 15.0,
47  originRadius = 1.5
48 ))
49 trackingPhase1.toModify(detachedTripletStepTrackingRegions, RegionPSet = dict(ptMin = 0.25))
50 
51 from Configuration.Eras.Modifier_pp_on_XeXe_2017_cff import pp_on_XeXe_2017
52 from Configuration.Eras.Modifier_pp_on_AA_2018_cff import pp_on_AA_2018
53 from RecoTracker.TkTrackingRegions.globalTrackingRegionWithVertices_cff import globalTrackingRegionWithVertices as _globalTrackingRegionWithVertices
54 (pp_on_XeXe_2017 | pp_on_AA_2018).toReplaceWith(detachedTripletStepTrackingRegions,
55  _globalTrackingRegionWithVertices.clone(RegionPSet=dict(
56  fixedError = 2.5,
57  ptMin = 0.9,
58  originRadius = 1.5
59  )
60  )
61 )
62 from Configuration.Eras.Modifier_highBetaStar_2018_cff import highBetaStar_2018
63 highBetaStar_2018.toModify(detachedTripletStepTrackingRegions, RegionPSet = dict(ptMin = 0.05))
64 
65 
66 # seeding
67 from RecoTracker.TkHitPairs.hitPairEDProducer_cfi import hitPairEDProducer as _hitPairEDProducer
68 detachedTripletStepHitDoublets = _hitPairEDProducer.clone(
69  seedingLayers = "detachedTripletStepSeedLayers",
70  trackingRegions = "detachedTripletStepTrackingRegions",
71  maxElement = 0,
72  produceIntermediateHitDoublets = True,
73 )
74 from RecoPixelVertexing.PixelTriplets.pixelTripletLargeTipEDProducer_cfi import pixelTripletLargeTipEDProducer as _pixelTripletLargeTipEDProducer
76 detachedTripletStepHitTriplets = _pixelTripletLargeTipEDProducer.clone(
77  doublets = "detachedTripletStepHitDoublets",
78  produceSeedingHitSets = True,
79 )
80 from RecoTracker.TkSeedGenerator.seedCreatorFromRegionConsecutiveHitsTripletOnlyEDProducer_cff import seedCreatorFromRegionConsecutiveHitsTripletOnlyEDProducer as _seedCreatorFromRegionConsecutiveHitsTripletOnlyEDProducer
81 detachedTripletStepSeeds = _seedCreatorFromRegionConsecutiveHitsTripletOnlyEDProducer.clone(
82  seedingHitSets = "detachedTripletStepHitTriplets",
83  SeedComparitorPSet = dict(# FIXME: is this defined in any cfi that could be imported instead of copy-paste?
84  ComponentName = 'PixelClusterShapeSeedComparitor',
85  FilterAtHelixStage = cms.bool(False),
86  FilterPixelHits = cms.bool(True),
87  FilterStripHits = cms.bool(False),
88  ClusterShapeHitFilterName = cms.string('ClusterShapeHitFilter'),
89  ClusterShapeCacheSrc = cms.InputTag('siPixelClusterShapeCache')
90  ),
91 )
92 
93 from RecoPixelVertexing.PixelTriplets.caHitTripletEDProducer_cfi import caHitTripletEDProducer as _caHitTripletEDProducer
94 trackingPhase1.toModify(detachedTripletStepHitDoublets, layerPairs = [0,1]) # layer pairs (0,1), (1,2)
95 trackingPhase1.toReplaceWith(detachedTripletStepHitTriplets, _caHitTripletEDProducer.clone(
96  doublets = "detachedTripletStepHitDoublets",
97  extraHitRPhitolerance = detachedTripletStepHitTriplets.extraHitRPhitolerance,
98  maxChi2 = dict(
99  pt1 = 0.8, pt2 = 2,
100  value1 = 300 , value2 = 10,
101  ),
102  useBendingCorrection = True,
103  CAThetaCut = 0.001,
104  CAPhiCut = 0,
105  CAHardPtCut = 0.2,
106 ))
107 highBetaStar_2018.toModify(detachedTripletStepHitTriplets,CAThetaCut = 0.002,CAPhiCut = 0.1,CAHardPtCut = 0)
108 
109 import FastSimulation.Tracking.TrajectorySeedProducer_cfi
110 _fastSim_detachedTripletStepSeeds = FastSimulation.Tracking.TrajectorySeedProducer_cfi.trajectorySeedProducer.clone(
111  trackingRegions = "detachedTripletStepTrackingRegions",
112  hitMasks = cms.InputTag("detachedTripletStepMasks"),
113  seedFinderSelector = dict( pixelTripletGeneratorFactory = _hitSetProducerToFactoryPSet(detachedTripletStepHitTriplets),
114  layerList = detachedTripletStepSeedLayers.layerList.value())
115 )
116 #new for phase1
117 trackingPhase1.toModify(_fastSim_detachedTripletStepSeeds, seedFinderSelector = dict(
118  pixelTripletGeneratorFactory = None,
119  CAHitTripletGeneratorFactory = _hitSetProducerToFactoryPSet(detachedTripletStepHitTriplets),
120  #new parameters required for phase1 seeding
121  BPix = dict(
122  TTRHBuilder = 'WithoutRefit',
123  HitProducer = 'TrackingRecHitProducer',
124  ),
125  FPix = dict(
126  TTRHBuilder = 'WithoutRefit',
127  HitProducer = 'TrackingRecHitProducer',
128  ),
129  layerPairs = detachedTripletStepHitDoublets.layerPairs.value()
130  )
131 )
132 fastSim.toReplaceWith(detachedTripletStepSeeds,_fastSim_detachedTripletStepSeeds)
133 
134 # QUALITY CUTS DURING TRACK BUILDING
136 _detachedTripletStepTrajectoryFilterBase = TrackingTools.TrajectoryFiltering.TrajectoryFilter_cff.CkfBaseTrajectoryFilter_block.clone(
137 # maxLostHitsFraction = cms.double(1./10.),
138 # constantValueForLostHitsFractionFilter = cms.double(0.701),
139  minimumNumberOfHits = 3,
140  minPt = 0.075,
141 )
142 detachedTripletStepTrajectoryFilterBase = _detachedTripletStepTrajectoryFilterBase.clone(
143  maxCCCLostHits = 0,
144  minGoodStripCharge = cms.PSet(refToPSet_ = cms.string('SiStripClusterChargeCutLoose'))
145 )
146 from Configuration.Eras.Modifier_tracker_apv_vfp30_2016_cff import tracker_apv_vfp30_2016
147 _tracker_apv_vfp30_2016.toModify(detachedTripletStepTrajectoryFilterBase, maxCCCLostHits = 2)
148 from Configuration.Eras.Modifier_trackingLowPU_cff import trackingLowPU
149 trackingLowPU.toReplaceWith(detachedTripletStepTrajectoryFilterBase, _detachedTripletStepTrajectoryFilterBase.clone(
150  maxLostHitsFraction = 1./10.,
151  constantValueForLostHitsFractionFilter = 0.701,
152 ))
153 
154 for e in [pp_on_XeXe_2017, pp_on_AA_2018]:
155  e.toModify(detachedTripletStepTrajectoryFilterBase, minPt=0.9)
156 
157 import RecoPixelVertexing.PixelLowPtUtilities.StripSubClusterShapeTrajectoryFilter_cfi
158 detachedTripletStepTrajectoryFilterShape = RecoPixelVertexing.PixelLowPtUtilities.StripSubClusterShapeTrajectoryFilter_cfi.StripSubClusterShapeTrajectoryFilterTIX12.clone()
159 detachedTripletStepTrajectoryFilter = cms.PSet(
160  ComponentType = cms.string('CompositeTrajectoryFilter'),
161  filters = cms.VPSet(
162  cms.PSet( refToPSet_ = cms.string('detachedTripletStepTrajectoryFilterBase')),
163 # cms.PSet( refToPSet_ = cms.string('detachedTripletStepTrajectoryFilterShape'))
164  ),
165 )
166 
167 
168 import RecoTracker.MeasurementDet.Chi2ChargeMeasurementEstimator_cfi
169 detachedTripletStepChi2Est = RecoTracker.MeasurementDet.Chi2ChargeMeasurementEstimator_cfi.Chi2ChargeMeasurementEstimator.clone(
170  ComponentName = cms.string('detachedTripletStepChi2Est'),
171  nSigma = cms.double(3.0),
172  MaxChi2 = cms.double(9.0),
173  clusterChargeCut = cms.PSet(refToPSet_ = cms.string('SiStripClusterChargeCutTight')),
174 )
175 _tracker_apv_vfp30_2016.toModify(detachedTripletStepChi2Est,
176  clusterChargeCut = dict(refToPSet_ = "SiStripClusterChargeCutTiny")
177 )
178 
179 # TRACK BUILDING
181 detachedTripletStepTrajectoryBuilder = RecoTracker.CkfPattern.GroupedCkfTrajectoryBuilder_cfi.GroupedCkfTrajectoryBuilder.clone(
182  MeasurementTrackerName = '',
183  trajectoryFilter = cms.PSet(refToPSet_ = cms.string('detachedTripletStepTrajectoryFilter')),
184  maxCand = 3,
185  alwaysUseInvalidHits = True,
186  estimator = cms.string('detachedTripletStepChi2Est'),
187  maxDPhiForLooperReconstruction = cms.double(2.0),
188  maxPtForLooperReconstruction = cms.double(0.7)
189  )
190 trackingLowPU.toModify(detachedTripletStepTrajectoryBuilder,
191  maxCand = 2,
192  alwaysUseInvalidHits = False,
193 )
194 
195 # MAKING OF TRACK CANDIDATES
197 detachedTripletStepTrackCandidates = RecoTracker.CkfPattern.CkfTrackCandidates_cfi.ckfTrackCandidates.clone(
198  src = cms.InputTag('detachedTripletStepSeeds'),
199  clustersToSkip = cms.InputTag('detachedTripletStepClusters'),
200  ### these two parameters are relevant only for the CachingSeedCleanerBySharedInput
201  numHitsForSeedCleaner = cms.int32(50),
202  onlyPixelHitsForSeedCleaner = cms.bool(True),
203  TrajectoryBuilderPSet = cms.PSet(refToPSet_ = cms.string('detachedTripletStepTrajectoryBuilder')),
204  doSeedingRegionRebuilding = True,
205  useHitsSplitting = True
206  )
207 
208 from TrackingTools.TrajectoryCleaning.TrajectoryCleanerBySharedHits_cfi import trajectoryCleanerBySharedHits
209 detachedTripletStepTrajectoryCleanerBySharedHits = trajectoryCleanerBySharedHits.clone(
210  ComponentName = cms.string('detachedTripletStepTrajectoryCleanerBySharedHits'),
211  fractionShared = cms.double(0.13),
212  allowSharedFirstHit = cms.bool(True)
213  )
214 detachedTripletStepTrackCandidates.TrajectoryCleaner = 'detachedTripletStepTrajectoryCleanerBySharedHits'
215 trackingLowPU.toModify(detachedTripletStepTrajectoryCleanerBySharedHits, fractionShared = 0.19)
216 
217 import FastSimulation.Tracking.TrackCandidateProducer_cfi
218 _fastSim_detachedTripletStepTrackCandidates = FastSimulation.Tracking.TrackCandidateProducer_cfi.trackCandidateProducer.clone(
219  src = cms.InputTag("detachedTripletStepSeeds"),
220  MinNumberOfCrossedLayers = 3,
221  hitMasks = cms.InputTag("detachedTripletStepMasks")
222  )
223 fastSim.toReplaceWith(detachedTripletStepTrackCandidates,_fastSim_detachedTripletStepTrackCandidates)
224 
225 
226 # TRACK FITTING
228 detachedTripletStepTracks = RecoTracker.TrackProducer.TrackProducer_cfi.TrackProducer.clone(
229  AlgorithmName = cms.string('detachedTripletStep'),
230  src = 'detachedTripletStepTrackCandidates',
231  Fitter = cms.string('FlexibleKFFittingSmoother')
232  )
233 fastSim.toModify(detachedTripletStepTracks,TTRHBuilder = 'WithoutRefit')
234 
235 # TRACK SELECTION AND QUALITY FLAG SETTING.
236 
237 
240 detachedTripletStepClassifier1 = TrackMVAClassifierDetached.clone()
241 detachedTripletStepClassifier1.src = 'detachedTripletStepTracks'
242 detachedTripletStepClassifier1.mva.GBRForestLabel = 'MVASelectorIter3_13TeV'
243 detachedTripletStepClassifier1.qualityCuts = [-0.5,0.0,0.5]
244 fastSim.toModify(detachedTripletStepClassifier1,vertices = "firstStepPrimaryVerticesBeforeMixing")
245 
246 detachedTripletStepClassifier2 = TrackMVAClassifierPrompt.clone()
247 detachedTripletStepClassifier2.src = 'detachedTripletStepTracks'
248 detachedTripletStepClassifier2.mva.GBRForestLabel = 'MVASelectorIter0_13TeV'
249 detachedTripletStepClassifier2.qualityCuts = [-0.2,0.0,0.4]
250 fastSim.toModify(detachedTripletStepClassifier2,vertices = "firstStepPrimaryVerticesBeforeMixing")
251 
253 detachedTripletStep = ClassifierMerger.clone()
254 detachedTripletStep.inputClassifiers=['detachedTripletStepClassifier1','detachedTripletStepClassifier2']
255 
256 trackingPhase1.toReplaceWith(detachedTripletStep, detachedTripletStepClassifier1.clone(
257  mva = dict(GBRForestLabel = 'MVASelectorDetachedTripletStep_Phase1'),
258  qualityCuts = [-0.2,0.3,0.8],
259 ))
260 highBetaStar_2018.toModify(detachedTripletStep,qualityCuts = [-0.5,0.0,0.5])
261 
262 # For LowPU
263 import RecoTracker.FinalTrackSelectors.multiTrackSelector_cfi
264 detachedTripletStepSelector = RecoTracker.FinalTrackSelectors.multiTrackSelector_cfi.multiTrackSelector.clone(
265  src = 'detachedTripletStepTracks',
266  useAnyMVA = cms.bool(False),
267  GBRForestLabel = cms.string('MVASelectorIter3'),
268  trackSelectors = [
269  RecoTracker.FinalTrackSelectors.multiTrackSelector_cfi.looseMTS.clone(
270  name = 'detachedTripletStepVtxLoose',
271  chi2n_par = 1.6,
272  res_par = ( 0.003, 0.001 ),
273  minNumberLayers = 3,
274  d0_par1 = ( 1.2, 3.0 ),
275  dz_par1 = ( 1.2, 3.0 ),
276  d0_par2 = ( 1.3, 3.0 ),
277  dz_par2 = ( 1.3, 3.0 )
278  ),
279  RecoTracker.FinalTrackSelectors.multiTrackSelector_cfi.looseMTS.clone(
280  name = 'detachedTripletStepTrkLoose',
281  chi2n_par = 0.7,
282  res_par = ( 0.003, 0.001 ),
283  minNumberLayers = 3,
284  d0_par1 = ( 1.6, 4.0 ),
285  dz_par1 = ( 1.6, 4.0 ),
286  d0_par2 = ( 1.6, 4.0 ),
287  dz_par2 = ( 1.6, 4.0 )
288  ),
289  RecoTracker.FinalTrackSelectors.multiTrackSelector_cfi.tightMTS.clone(
290  name = 'detachedTripletStepVtxTight',
291  preFilterName = 'detachedTripletStepVtxLoose',
292  chi2n_par = 0.7,
293  res_par = ( 0.003, 0.001 ),
294  minNumberLayers = 3,
295  maxNumberLostLayers = 1,
296  minNumber3DLayers = 3,
297  d0_par1 = ( 0.95, 3.0 ),
298  dz_par1 = ( 0.9, 3.0 ),
299  d0_par2 = ( 1.0, 3.0 ),
300  dz_par2 = ( 1.0, 3.0 )
301  ),
302  RecoTracker.FinalTrackSelectors.multiTrackSelector_cfi.tightMTS.clone(
303  name = 'detachedTripletStepTrkTight',
304  preFilterName = 'detachedTripletStepTrkLoose',
305  chi2n_par = 0.5,
306  res_par = ( 0.003, 0.001 ),
307  minNumberLayers = 5,
308  maxNumberLostLayers = 1,
309  minNumber3DLayers = 3,
310  d0_par1 = ( 1.1, 4.0 ),
311  dz_par1 = ( 1.1, 4.0 ),
312  d0_par2 = ( 1.1, 4.0 ),
313  dz_par2 = ( 1.1, 4.0 )
314  ),
315  RecoTracker.FinalTrackSelectors.multiTrackSelector_cfi.highpurityMTS.clone(
316  name = 'detachedTripletStepVtx',
317  preFilterName = 'detachedTripletStepVtxTight',
318  chi2n_par = 0.7,
319  res_par = ( 0.003, 0.001 ),
320  minNumberLayers = 3,
321  maxNumberLostLayers = 1,
322  minNumber3DLayers = 3,
323  d0_par1 = ( 0.85, 3.0 ),
324  dz_par1 = ( 0.8, 3.0 ),
325  d0_par2 = ( 0.9, 3.0 ),
326  dz_par2 = ( 0.9, 3.0 )
327  ),
328  RecoTracker.FinalTrackSelectors.multiTrackSelector_cfi.highpurityMTS.clone(
329  name = 'detachedTripletStepTrk',
330  preFilterName = 'detachedTripletStepTrkTight',
331  chi2n_par = 0.4,
332  res_par = ( 0.003, 0.001 ),
333  minNumberLayers = 5,
334  maxNumberLostLayers = 1,
335  minNumber3DLayers = 4,
336  d0_par1 = ( 1.0, 4.0 ),
337  dz_par1 = ( 1.0, 4.0 ),
338  d0_par2 = ( 1.0, 4.0 ),
339  dz_par2 = ( 1.0, 4.0 )
340  )
341  ] #end of vpset
342 ) #end of clone
343 
344 from RecoTracker.FinalTrackSelectors.trackAlgoPriorityOrder_cfi import trackAlgoPriorityOrder
346 trackingLowPU.toReplaceWith(detachedTripletStep, RecoTracker.FinalTrackSelectors.trackListMerger_cfi.trackListMerger.clone(
347  TrackProducers = [
348  'detachedTripletStepTracks',
349  'detachedTripletStepTracks',
350  ],
351  hasSelector = [1,1],
352  selectedTrackQuals = [
353  cms.InputTag("detachedTripletStepSelector","detachedTripletStepVtx"),
354  cms.InputTag("detachedTripletStepSelector","detachedTripletStepTrk")
355  ],
356  setsToMerge = [cms.PSet( tLists=cms.vint32(0,1), pQual=cms.bool(True) )],
357  writeOnlyTrkQuals =True
358 ))
359 
360 DetachedTripletStepTask = cms.Task(detachedTripletStepClusters,
361  detachedTripletStepSeedLayers,
362  detachedTripletStepTrackingRegions,
363  detachedTripletStepHitDoublets,
364  detachedTripletStepHitTriplets,
365  detachedTripletStepSeeds,
366  detachedTripletStepTrackCandidates,
367  detachedTripletStepTracks,
368  detachedTripletStepClassifier1,detachedTripletStepClassifier2,
369  detachedTripletStep)
370 DetachedTripletStep = cms.Sequence(DetachedTripletStepTask)
371 _DetachedTripletStepTask_LowPU = DetachedTripletStepTask.copyAndExclude([detachedTripletStepClassifier2])
372 _DetachedTripletStepTask_LowPU.replace(detachedTripletStepClassifier1, detachedTripletStepSelector)
373 trackingLowPU.toReplaceWith(DetachedTripletStepTask, _DetachedTripletStepTask_LowPU)
374 
375 # fast tracking mask producer
376 from FastSimulation.Tracking.FastTrackerRecHitMaskProducer_cfi import maskProducerFromClusterRemover
377 detachedTripletStepMasks = maskProducerFromClusterRemover(detachedTripletStepClusters)
378 fastSim.toReplaceWith(DetachedTripletStepTask,
379  cms.Task(detachedTripletStepMasks
380  ,detachedTripletStepTrackingRegions
381  ,detachedTripletStepSeeds
382  ,detachedTripletStepTrackCandidates
383  ,detachedTripletStepTracks
384  ,detachedTripletStepClassifier1,detachedTripletStepClassifier2
385  ,detachedTripletStep
386  ) )
COSMIC RUN SETTING.
def _hitSetProducerToFactoryPSet(producer)