CMS 3D CMS Logo

hiLowPtQuadStep_cff.py
Go to the documentation of this file.
1 from __future__ import absolute_import
3 from .HIPixelTripletSeeds_cff import *
4 from .HIPixel3PrimTracks_cfi import *
5 
6 hiLowPtQuadStepClusters = cms.EDProducer("HITrackClusterRemover",
7  clusterLessSolution = cms.bool(True),
8  trajectories = cms.InputTag("hiGlobalPrimTracks"),
9  overrideTrkQuals = cms.InputTag('hiInitialStepSelector','hiInitialStep'),
10  TrackQuality = cms.string('highPurity'),
11  minNumberOfLayersWithMeasBeforeFiltering = cms.int32(0),
12  pixelClusters = cms.InputTag("siPixelClusters"),
13  stripClusters = cms.InputTag("siStripClusters"),
14  Common = cms.PSet(
15  maxChi2 = cms.double(9.0),
16  ),
17  Strip = cms.PSet(
18  #Yen-Jie's mod to preserve merged clusters
19  maxSize = cms.uint32(2),
20  maxChi2 = cms.double(9.0)
21  )
22 )
23 
24 # SEEDING LAYERS
25 # Using 4 layers layerlist
26 hiLowPtQuadStepSeedLayers = hiPixelLayerQuadruplets.clone()
27 hiLowPtQuadStepSeedLayers.BPix.skipClusters = cms.InputTag('hiLowPtQuadStepClusters')
28 hiLowPtQuadStepSeedLayers.FPix.skipClusters = cms.InputTag('hiLowPtQuadStepClusters')
29 
30 # SEEDS
31 from RecoTracker.TkTrackingRegions.globalTrackingRegionWithVertices_cfi import globalTrackingRegionWithVertices as _globalTrackingRegionWithVertices
32 from RecoTracker.TkHitPairs.hitPairEDProducer_cfi import hitPairEDProducer as _hitPairEDProducer
33 from RecoPixelVertexing.PixelTriplets.pixelTripletHLTEDProducer_cfi import pixelTripletHLTEDProducer as _pixelTripletHLTEDProducer
39 
40 hiLowPtQuadStepTrackingRegions = _globalTrackingRegionWithVertices.clone(RegionPSet=dict(
41  precise = True,
42  useMultipleScattering = False,
43  useFakeVertices = False,
44  beamSpot = "offlineBeamSpot",
45  useFixedError = True,
46  nSigmaZ = 4.0,
47  sigmaZVertex = 4.0,
48  fixedError = 0.5,
49  VertexCollection = "hiSelectedPixelVertex",
50  ptMin = 0.3,#0.2 for pp
51  useFoundVertices = True,
52  originRadius = 0.02 #0.02 for pp
53 ))
54 hiLowPtQuadStepTracksHitDoubletsCA = _hitPairEDProducer.clone(
55  clusterCheck = "",
56  seedingLayers = "hiLowPtQuadStepSeedLayers",
57  trackingRegions = "hiLowPtQuadStepTrackingRegions",
58  maxElement = 50000000,
59  produceIntermediateHitDoublets = True,
60  layerPairs = [0,1,2]
61 )
62 
63 import RecoPixelVertexing.PixelLowPtUtilities.LowPtClusterShapeSeedComparitor_cfi
64 from RecoPixelVertexing.PixelTriplets.caHitQuadrupletEDProducer_cfi import caHitQuadrupletEDProducer as _caHitQuadrupletEDProducer
65 hiLowPtQuadStepTracksHitQuadrupletsCA = _caHitQuadrupletEDProducer.clone(
66  doublets = "hiLowPtQuadStepTracksHitDoubletsCA",
67  extraHitRPhitolerance = 0.0,
68  SeedComparitorPSet = RecoPixelVertexing.PixelLowPtUtilities.LowPtClusterShapeSeedComparitor_cfi.LowPtClusterShapeSeedComparitor.clone(),
69  maxChi2 = dict(
70  pt1 = 0.7, pt2 = 2,
71  value1 = 1000, value2 = 150,
72  ),
73  useBendingCorrection = True,
74  fitFastCircle = True,
75  fitFastCircleChi2Cut = True,
76  CAThetaCut = 0.0017,
77  CAPhiCut = 0.3,
78 )
79 
80 
81 hiLowPtQuadStepPixelTracksFilter = hiFilter.clone(
82  nSigmaTipMaxTolerance = 0,
83  lipMax = 1.0,
84  tipMax = 1.0,
85  ptMin = 0.4, #seeding region is 0.3
86 )
87 hiLowPtQuadStepPixelTracks = cms.EDProducer("PixelTrackProducer",
88 
89  passLabel = cms.string('Pixel detached tracks with vertex constraint'),
90 
91  # Ordered Hits
92  SeedingHitSets = cms.InputTag("hiLowPtQuadStepTracksHitQuadrupletsCA"),
93 
94  # Fitter
95  Fitter = cms.InputTag("pixelFitterByHelixProjections"),
96 
97  # Filter
98  Filter = cms.InputTag("hiLowPtQuadStepPixelTracksFilter"),
99 
100  # Cleaner
101  Cleaner = cms.string("trackCleaner")
102 )
103 
104 
105 import RecoPixelVertexing.PixelLowPtUtilities.TrackSeeds_cfi
106 hiLowPtQuadStepSeeds = RecoPixelVertexing.PixelLowPtUtilities.TrackSeeds_cfi.pixelTrackSeeds.clone(
107  InputCollection = 'hiLowPtQuadStepPixelTracks'
108  )
109 
110 # QUALITY CUTS DURING TRACK BUILDING
112 hiLowPtQuadStepTrajectoryFilter = TrackingTools.TrajectoryFiltering.TrajectoryFilter_cff.CkfBaseTrajectoryFilter_block.clone(
113  #maxLostHits = 1,
114  minimumNumberOfHits = 3,#3 for pp
115  minPt = cms.double(0.075),# 0.075 for pp
116  #constantValueForLostHitsFractionFilter = cms.double(0.701)
117  )
118 
120 hiLowPtQuadStepChi2Est = TrackingTools.KalmanUpdators.Chi2MeasurementEstimator_cfi.Chi2MeasurementEstimator.clone(
121  ComponentName = cms.string('hiLowPtQuadStepChi2Est'),
122  nSigma = cms.double(3.0),
123  MaxChi2 = cms.double(9.0)
124  )
125 
126 
127 # TRACK BUILDING
129 hiLowPtQuadStepTrajectoryBuilder = RecoTracker.CkfPattern.GroupedCkfTrajectoryBuilder_cfi.GroupedCkfTrajectoryBuilder.clone(
130  MeasurementTrackerName = '',
131  trajectoryFilter = cms.PSet(refToPSet_ = cms.string('hiLowPtQuadStepTrajectoryFilter')),
132  maxCand = 4,#4 for pp
133  estimator = cms.string('hiLowPtQuadStepChi2Est'),
134  maxDPhiForLooperReconstruction = cms.double(2.0),#2.0 for pp
135  # 0.63 GeV is the maximum pT for a charged particle to loop within the 1.1m radius
136  # of the outermost Tracker barrel layer (B=3.8T)
137  maxPtForLooperReconstruction = cms.double(0.7),# 0.7 for pp
138  alwaysUseInvalidHits = cms.bool(False)
139  )
140 
141 # MAKING OF TRACK CANDIDATES
142 
143 # Trajectory cleaner in default
144 
146 hiLowPtQuadStepTrackCandidates = RecoTracker.CkfPattern.CkfTrackCandidates_cfi.ckfTrackCandidates.clone(
147  src = cms.InputTag('hiLowPtQuadStepSeeds'),
148 
149  numHitsForSeedCleaner = cms.int32(50),
150  onlyPixelHitsForSeedCleaner = cms.bool(True),
151  TrajectoryBuilderPSet = cms.PSet(refToPSet_ = cms.string('hiLowPtQuadStepTrajectoryBuilder')),
152  TrajectoryBuilder = cms.string('hiLowPtQuadStepTrajectoryBuilder'),
153  clustersToSkip = cms.InputTag('hiLowPtQuadStepClusters'),
154  doSeedingRegionRebuilding = True,
155  useHitsSplitting = True
156  )
157 
158 
159 # TRACK FITTING
161 hiLowPtQuadStepTracks = RecoTracker.TrackProducer.TrackProducer_cfi.TrackProducer.clone(
162  src = 'hiLowPtQuadStepTrackCandidates',
163  AlgorithmName = cms.string('lowPtQuadStep'),
164  Fitter=cms.string('FlexibleKFFittingSmoother')
165  )
166 
167 # Final selection
169 hiLowPtQuadStepSelector = RecoHI.HiTracking.hiMultiTrackSelector_cfi.hiMultiTrackSelector.clone(
170  src='hiLowPtQuadStepTracks',
171  useAnyMVA = cms.bool(True),
172  GBRForestLabel = cms.string('HIMVASelectorIter8'),#FIXME MVA for new iteration
173  GBRForestVars = cms.vstring(['chi2perdofperlayer', 'nhits', 'nlayers', 'eta']),
174  trackSelectors= cms.VPSet(
175  RecoHI.HiTracking.hiMultiTrackSelector_cfi.hiLooseMTS.clone(
176  name = 'hiLowPtQuadStepLoose',
177  applyAdaptedPVCuts = cms.bool(False),
178  useMVA = cms.bool(False),
179  ), #end of pset
180  RecoHI.HiTracking.hiMultiTrackSelector_cfi.hiTightMTS.clone(
181  name = 'hiLowPtQuadStepTight',
182  preFilterName = 'hiLowPtQuadStepLoose',
183  applyAdaptedPVCuts = cms.bool(False),
184  useMVA = cms.bool(True),
185  minMVA = cms.double(-0.2)
186  ),
187  RecoHI.HiTracking.hiMultiTrackSelector_cfi.hiHighpurityMTS.clone(
188  name = 'hiLowPtQuadStep',
189  preFilterName = 'hiLowPtQuadStepTight',
190  applyAdaptedPVCuts = cms.bool(False),
191  useMVA = cms.bool(True),
192  minMVA = cms.double(-0.09)
193  ),
194  ) #end of vpset
195  ) #end of clone
196 from Configuration.Eras.Modifier_trackingPhase1_cff import trackingPhase1
197 trackingPhase1.toModify(hiLowPtQuadStepSelector, useAnyMVA = cms.bool(False))
198 trackingPhase1.toModify(hiLowPtQuadStepSelector, trackSelectors= cms.VPSet(
199  RecoHI.HiTracking.hiMultiTrackSelector_cfi.hiLooseMTS.clone(
200  name = 'hiLowPtQuadStepLoose',
201  applyAdaptedPVCuts = cms.bool(False),
202  useMVA = cms.bool(False),
203  ), #end of pset
204  RecoHI.HiTracking.hiMultiTrackSelector_cfi.hiTightMTS.clone(
205  name = 'hiLowPtQuadStepTight',
206  preFilterName = 'hiLowPtQuadStepLoose',
207  applyAdaptedPVCuts = cms.bool(False),
208  useMVA = cms.bool(False),
209  minMVA = cms.double(-0.2)
210  ),
211  RecoHI.HiTracking.hiMultiTrackSelector_cfi.hiHighpurityMTS.clone(
212  name = 'hiLowPtQuadStep',
213  preFilterName = 'hiLowPtQuadStepTight',
214  applyAdaptedPVCuts = cms.bool(False),
215  useMVA = cms.bool(False),
216  minMVA = cms.double(-0.09)
217  ),
218  ) #end of vpset
219 )
220 
222 hiLowPtQuadStepQual = RecoTracker.FinalTrackSelectors.trackListMerger_cfi.trackListMerger.clone(
223  TrackProducers=cms.VInputTag(cms.InputTag('hiLowPtQuadStepTracks')),
224  hasSelector=cms.vint32(1),
225  selectedTrackQuals = cms.VInputTag(cms.InputTag("hiLowPtQuadStepSelector","hiLowPtQuadStep")),
226  copyExtras = True,
227  makeReKeyedSeeds = cms.untracked.bool(False),
228  )
229 
230 
231 hiLowPtQuadStepTask = cms.Task(hiLowPtQuadStepClusters,
232  hiLowPtQuadStepSeedLayers,
233  hiLowPtQuadStepTrackingRegions,
234  hiLowPtQuadStepTracksHitDoubletsCA,
235  hiLowPtQuadStepTracksHitQuadrupletsCA,
236  pixelFitterByHelixProjections,
237  hiLowPtQuadStepPixelTracksFilter,
238  hiLowPtQuadStepPixelTracks,
239  hiLowPtQuadStepSeeds,
240  hiLowPtQuadStepTrackCandidates,
241  hiLowPtQuadStepTracks,
242  hiLowPtQuadStepSelector,
243  hiLowPtQuadStepQual)
244 hiLowPtQuadStep = cms.Sequence(hiLowPtQuadStepTask)
245 
TrajectoryFilter_cff
GroupedCkfTrajectoryBuilder_cfi
hiMultiTrackSelector_cfi
ClusterShapeHitFilterESProducer_cfi
HITrackingRegionProducer_cfi
pp iterative tracking modified for hiOffline reco (the vertex is the one reconstructed in HI) 3rd ste...
trackListMerger_cfi
TrackProducer_cfi
globalTrackingRegionWithVertices_cfi
pixelFitterByHelixProjections_cfi
LowPtQuadStep_cff
trackCleaner_cfi
HIPixelTrackFilter_cff
CkfTrackCandidates_cfi
Chi2MeasurementEstimator_cfi