CMS 3D CMS Logo

hiLowPtTripletStep_cff.py
Go to the documentation of this file.
1 import FWCore.ParameterSet.Config as cms
2 
3 
4 # NEW CLUSTERS (remove previously used clusters)
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"),
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 
25 # SEEDING LAYERS
27 hiLowPtTripletStepSeedLayers = RecoTracker.TkSeedingLayers.PixelLayerTriplets_cfi.PixelLayerTriplets.clone()
28 hiLowPtTripletStepSeedLayers.BPix.skipClusters = cms.InputTag('hiLowPtTripletStepClusters')
29 hiLowPtTripletStepSeedLayers.FPix.skipClusters = cms.InputTag('hiLowPtTripletStepClusters')
30 
31 # SEEDS
32 from RecoTracker.TkTrackingRegions.globalTrackingRegionWithVertices_cfi import globalTrackingRegionWithVertices as _globalTrackingRegionWithVertices
33 from RecoTracker.TkHitPairs.hitPairEDProducer_cfi import hitPairEDProducer as _hitPairEDProducer
34 from RecoPixelVertexing.PixelTriplets.pixelTripletHLTEDProducer_cfi import pixelTripletHLTEDProducer as _pixelTripletHLTEDProducer
36 import RecoPixelVertexing.PixelLowPtUtilities.LowPtClusterShapeSeedComparitor_cfi
41 
42 hiLowPtTripletStepTrackingRegions = _globalTrackingRegionWithVertices.clone(RegionPSet=dict(
43  precise = True,
44  useMultipleScattering = False,
45  useFakeVertices = False,
46  beamSpot = "offlineBeamSpot",
47  useFixedError = False,
48  nSigmaZ = 4.0,
49  sigmaZVertex = 4.0,
50  fixedError = 0.2,
51  VertexCollection = "hiSelectedVertex",
52  ptMin = 0.4,
53  useFoundVertices = True,
54  originRadius = 0.02
55 ))
56 hiLowPtTripletStepTracksHitDoublets = _hitPairEDProducer.clone(
57  clusterCheck = "",
58  seedingLayers = "hiLowPtTripletStepSeedLayers",
59  trackingRegions = "hiLowPtTripletStepTrackingRegions",
60  maxElement = 0,
61  produceIntermediateHitDoublets = True,
62 )
63 hiLowPtTripletStepTracksHitTriplets = _pixelTripletHLTEDProducer.clone(
64  doublets = "hiLowPtTripletStepTracksHitDoublets",
65  maxElement = 5000000,
66  SeedComparitorPSet = RecoPixelVertexing.PixelLowPtUtilities.LowPtClusterShapeSeedComparitor_cfi.LowPtClusterShapeSeedComparitor.clone(),
67  produceSeedingHitSets = True,
68 )
69 hiLowPtTripletStepPixelTracksFilter = hiFilter.clone(
70  nSigmaLipMaxTolerance = 4.0,
71  nSigmaTipMaxTolerance = 4.0,
72  lipMax = 0,
73  ptMin = 0.4,
74 )
75 hiLowPtTripletStepPixelTracks = cms.EDProducer("PixelTrackProducer",
76 
77  passLabel = cms.string('Pixel primary tracks with vertex constraint'),
78 
79  # Ordered Hits
80  SeedingHitSets = cms.InputTag("hiLowPtTripletStepTracksHitTriplets"),
81 
82  # Fitter
83  Fitter = cms.InputTag("pixelFitterByHelixProjections"),
84 
85  # Filter
86  Filter = cms.InputTag("hiLowPtTripletStepPixelTracksFilter"),
87 
88  # Cleaner
89  Cleaner = cms.string("trackCleaner")
90 )
91 
92 
93 import RecoPixelVertexing.PixelLowPtUtilities.TrackSeeds_cfi
94 hiLowPtTripletStepSeeds = RecoPixelVertexing.PixelLowPtUtilities.TrackSeeds_cfi.pixelTrackSeeds.clone(
95  InputCollection = 'hiLowPtTripletStepPixelTracks'
96  )
97 
98 
99 # QUALITY CUTS DURING TRACK BUILDING
101 hiLowPtTripletStepTrajectoryFilter = TrackingTools.TrajectoryFiltering.TrajectoryFilter_cff.CkfBaseTrajectoryFilter_block.clone(
102  maxLostHits = 1,
103  minimumNumberOfHits = 6,
104  minPt = 0.4
105  )
106 
108 hiLowPtTripletStepChi2Est = TrackingTools.KalmanUpdators.Chi2MeasurementEstimator_cfi.Chi2MeasurementEstimator.clone(
109  ComponentName = cms.string('hiLowPtTripletStepChi2Est'),
110  nSigma = cms.double(3.0),
111  MaxChi2 = cms.double(9.0)
112  )
113 
114 # TRACK BUILDING
116 hiLowPtTripletStepTrajectoryBuilder = RecoTracker.CkfPattern.GroupedCkfTrajectoryBuilder_cfi.GroupedCkfTrajectoryBuilder.clone(
117  MeasurementTrackerName = '',
118  trajectoryFilter = cms.PSet(refToPSet_ = cms.string('hiLowPtTripletStepTrajectoryFilter')),
119  maxCand = 3,
120  estimator = cms.string('hiLowPtTripletStepChi2Est'),
121  maxDPhiForLooperReconstruction = cms.double(2.0),
122  # 0.63 GeV is the maximum pT for a charged particle to loop within the 1.1m radius
123  # of the outermost Tracker barrel layer (with B=3.8T)
124  maxPtForLooperReconstruction = cms.double(0.7)
125  )
126 
127 # MAKING OF TRACK CANDIDATES
129 hiLowPtTripletStepTrackCandidates = RecoTracker.CkfPattern.CkfTrackCandidates_cfi.ckfTrackCandidates.clone(
130  src = cms.InputTag('hiLowPtTripletStepSeeds'),
131  ### these two parameters are relevant only for the CachingSeedCleanerBySharedInput
132  numHitsForSeedCleaner = cms.int32(50),
133  onlyPixelHitsForSeedCleaner = cms.bool(True),
134  TrajectoryBuilderPSet = cms.PSet(refToPSet_ = cms.string('hiLowPtTripletStepTrajectoryBuilder')),
135  clustersToSkip = cms.InputTag('hiLowPtTripletStepClusters'),
136  doSeedingRegionRebuilding = True,
137  useHitsSplitting = True
138  )
139 
140 # TRACK FITTING
142 hiLowPtTripletStepTracks = RecoTracker.TrackProducer.TrackProducer_cfi.TrackProducer.clone(
143  src = 'hiLowPtTripletStepTrackCandidates',
144  AlgorithmName = cms.string('lowPtTripletStep'),
145  Fitter=cms.string('FlexibleKFFittingSmoother')
146  )
147 
148 
149 
150 # Final selection
152 hiLowPtTripletStepSelector = RecoHI.HiTracking.hiMultiTrackSelector_cfi.hiMultiTrackSelector.clone(
153  src='hiLowPtTripletStepTracks',
154  useAnyMVA = cms.bool(True),
155  GBRForestLabel = cms.string('HIMVASelectorIter5'),
156  GBRForestVars = cms.vstring(['chi2perdofperlayer', 'dxyperdxyerror', 'dzperdzerror', 'relpterr', 'nhits', 'nlayers', 'eta']),
157  trackSelectors= cms.VPSet(
158  RecoHI.HiTracking.hiMultiTrackSelector_cfi.hiLooseMTS.clone(
159  name = 'hiLowPtTripletStepLoose',
160  useMVA = cms.bool(False)
161  ), #end of pset
162  RecoHI.HiTracking.hiMultiTrackSelector_cfi.hiTightMTS.clone(
163  name = 'hiLowPtTripletStepTight',
164  preFilterName = 'hiLowPtTripletStepLoose',
165  useMVA = cms.bool(True),
166  minMVA = cms.double(-0.58)
167  ),
168  RecoHI.HiTracking.hiMultiTrackSelector_cfi.hiHighpurityMTS.clone(
169  name = 'hiLowPtTripletStep',
170  preFilterName = 'hiLowPtTripletStepTight',
171  useMVA = cms.bool(True),
172  minMVA = cms.double(0.35)
173  ),
174  ) #end of vpset
175  ) #end of clone
176 
177 
179 hiLowPtTripletStepQual = RecoTracker.FinalTrackSelectors.trackListMerger_cfi.trackListMerger.clone(
180  TrackProducers = cms.VInputTag(cms.InputTag('hiLowPtTripletStepTracks')),
181  hasSelector=cms.vint32(1),
182  selectedTrackQuals = cms.VInputTag(cms.InputTag("hiLowPtTripletStepSelector","hiLowPtTripletStep")),
183  copyExtras = True,
184  makeReKeyedSeeds = cms.untracked.bool(False),
185  #writeOnlyTrkQuals = True
186  )
187 
188 # Final sequence
189 
190 hiLowPtTripletStep = cms.Sequence(hiLowPtTripletStepClusters*
191  hiLowPtTripletStepSeedLayers*
192  hiLowPtTripletStepTrackingRegions*
193  hiLowPtTripletStepTracksHitDoublets*
194  hiLowPtTripletStepTracksHitTriplets*
195  pixelFitterByHelixProjections*
196  hiLowPtTripletStepPixelTracksFilter*
197  hiLowPtTripletStepPixelTracks*hiLowPtTripletStepSeeds*
198  hiLowPtTripletStepTrackCandidates*
199  hiLowPtTripletStepTracks*
200  hiLowPtTripletStepSelector*
201  hiLowPtTripletStepQual
202  )
pp iterative tracking modified for hiOffline reco (the vertex is the one reconstructed in HI) 3rd ste...