CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
PixelLessStep_cff.py
Go to the documentation of this file.
2 
3 ##########################################################################
4 # Large impact parameter tracking using TIB/TID/TEC stereo layer seeding #
5 ##########################################################################
6 
8 pixelLessStepClusters = trackClusterRemover.clone(
9  maxChi2 = cms.double(9.0),
10  trajectories = cms.InputTag("mixedTripletStepTracks"),
11  pixelClusters = cms.InputTag("siPixelClusters"),
12  stripClusters = cms.InputTag("siStripClusters"),
13  oldClusterRemovalInfo = cms.InputTag("mixedTripletStepClusters"),
14  overrideTrkQuals = cms.InputTag('mixedTripletStep'),
15  TrackQuality = cms.string('highPurity'),
16  minNumberOfLayersWithMeasBeforeFiltering = cms.int32(0),
17 )
18 
19 # SEEDING LAYERS
21 pixelLessStepSeedLayers = cms.EDProducer("SeedingLayersEDProducer",
22  layerList = cms.vstring(
23  #TIB
24  'TIB1+TIB2+MTIB3',
25  #TIB+TID
26  'TIB1+TIB2+MTID1_pos','TIB1+TIB2+MTID1_neg',
27  #TID
28  'TID1_pos+TID2_pos+TID3_pos','TID1_neg+TID2_neg+TID3_neg',#ring 1-2 (matched)
29  'TID1_pos+TID2_pos+MTID3_pos','TID1_neg+TID2_neg+MTID3_neg',#ring 3 (mono)
30  'TID1_pos+TID2_pos+MTEC1_pos','TID1_neg+TID2_neg+MTEC1_neg',
31  #TID+TEC RING 1-3
32  'TID2_pos+TID3_pos+TEC1_pos','TID2_neg+TID3_neg+TEC1_neg',#ring 1-2 (matched)
33  'TID2_pos+TID3_pos+MTEC1_pos','TID2_neg+TID3_neg+MTEC1_neg',#ring 3 (mono)
34  #TEC RING 1-3
35  'TEC1_pos+TEC2_pos+TEC3_pos', 'TEC1_neg+TEC2_neg+TEC3_neg',
36  'TEC1_pos+TEC2_pos+MTEC3_pos','TEC1_neg+TEC2_neg+MTEC3_neg',
37  'TEC1_pos+TEC2_pos+TEC4_pos', 'TEC1_neg+TEC2_neg+TEC4_neg',
38  'TEC1_pos+TEC2_pos+MTEC4_pos','TEC1_neg+TEC2_neg+MTEC4_neg',
39  'TEC2_pos+TEC3_pos+TEC4_pos', 'TEC2_neg+TEC3_neg+TEC4_neg',
40  'TEC2_pos+TEC3_pos+MTEC4_pos','TEC2_neg+TEC3_neg+MTEC4_neg',
41  'TEC2_pos+TEC3_pos+TEC5_pos', 'TEC2_neg+TEC3_neg+TEC5_neg',
42  'TEC2_pos+TEC3_pos+TEC6_pos', 'TEC2_neg+TEC3_neg+TEC6_neg',
43  'TEC3_pos+TEC4_pos+TEC5_pos', 'TEC3_neg+TEC4_neg+TEC5_neg',
44  'TEC3_pos+TEC4_pos+MTEC5_pos','TEC3_neg+TEC4_neg+MTEC5_neg',
45  'TEC3_pos+TEC5_pos+TEC6_pos', 'TEC3_neg+TEC5_neg+TEC6_neg',
46  'TEC4_pos+TEC5_pos+TEC6_pos', 'TEC4_neg+TEC5_neg+TEC6_neg'
47  ),
48  TIB = cms.PSet(
49  TTRHBuilder = cms.string('WithTrackAngle'), clusterChargeCut = cms.PSet(refToPSet_ = cms.string('SiStripClusterChargeCutTight')),
50  matchedRecHits = cms.InputTag("siStripMatchedRecHits","matchedRecHit"),
51  skipClusters = cms.InputTag('pixelLessStepClusters')
52  ),
53  MTIB = cms.PSet(
54  TTRHBuilder = cms.string('WithTrackAngle'), clusterChargeCut = cms.PSet(refToPSet_ = cms.string('SiStripClusterChargeCutTight')),
55  skipClusters = cms.InputTag('pixelLessStepClusters'),
56  rphiRecHits = cms.InputTag("siStripMatchedRecHits","rphiRecHit")
57  ),
58  TID = cms.PSet(
59  matchedRecHits = cms.InputTag("siStripMatchedRecHits","matchedRecHit"),
60  skipClusters = cms.InputTag('pixelLessStepClusters'),
61  useRingSlector = cms.bool(True),
62  TTRHBuilder = cms.string('WithTrackAngle'), clusterChargeCut = cms.PSet(refToPSet_ = cms.string('SiStripClusterChargeCutTight')),
63  minRing = cms.int32(1),
64  maxRing = cms.int32(2)
65  ),
66  MTID = cms.PSet(
67  rphiRecHits = cms.InputTag("siStripMatchedRecHits","rphiRecHit"),
68  skipClusters = cms.InputTag('pixelLessStepClusters'),
69  useRingSlector = cms.bool(True),
70  TTRHBuilder = cms.string('WithTrackAngle'), clusterChargeCut = cms.PSet(refToPSet_ = cms.string('SiStripClusterChargeCutTight')),
71  minRing = cms.int32(3),
72  maxRing = cms.int32(3)
73  ),
74  TEC = cms.PSet(
75  matchedRecHits = cms.InputTag("siStripMatchedRecHits","matchedRecHit"),
76  skipClusters = cms.InputTag('pixelLessStepClusters'),
77  useRingSlector = cms.bool(True),
78  TTRHBuilder = cms.string('WithTrackAngle'), clusterChargeCut = cms.PSet(refToPSet_ = cms.string('SiStripClusterChargeCutTight')),
79  minRing = cms.int32(1),
80  maxRing = cms.int32(2)
81  ),
82  MTEC = cms.PSet(
83  rphiRecHits = cms.InputTag("siStripMatchedRecHits","rphiRecHit"),
84  skipClusters = cms.InputTag('pixelLessStepClusters'),
85  useRingSlector = cms.bool(True),
86  TTRHBuilder = cms.string('WithTrackAngle'), clusterChargeCut = cms.PSet(refToPSet_ = cms.string('SiStripClusterChargeCutTight')),
87  minRing = cms.int32(3),
88  maxRing = cms.int32(3)
89  )
90 )
91 
92 # SEEDS
94 pixelLessStepSeeds = RecoTracker.TkSeedGenerator.GlobalSeedsFromTriplets_cff.globalSeedsFromTriplets.clone()
95 #OrderedHitsFactory
96 pixelLessStepSeeds.OrderedHitsFactoryPSet.SeedingLayers = 'pixelLessStepSeedLayers'
97 pixelLessStepSeeds.OrderedHitsFactoryPSet.ComponentName = 'StandardMultiHitGenerator'
98 import RecoTracker.TkSeedGenerator.MultiHitGeneratorFromChi2_cfi
99 pixelLessStepSeeds.OrderedHitsFactoryPSet.GeneratorPSet = RecoTracker.TkSeedGenerator.MultiHitGeneratorFromChi2_cfi.MultiHitGeneratorFromChi2.clone()
100 #SeedCreator
101 pixelLessStepSeeds.SeedCreatorPSet.ComponentName = 'SeedFromConsecutiveHitsTripletOnlyCreator'
102 #RegionFactory
103 pixelLessStepSeeds.RegionFactoryPSet.RegionPSet.ptMin = 0.4
104 pixelLessStepSeeds.RegionFactoryPSet.RegionPSet.originHalfLength = 12.0
105 pixelLessStepSeeds.RegionFactoryPSet.RegionPSet.originRadius = 1.0
106 #SeedComparitor
108 pixelLessStepClusterShapeHitFilter = RecoPixelVertexing.PixelLowPtUtilities.ClusterShapeHitFilterESProducer_cfi.ClusterShapeHitFilterESProducer.clone(
109  ComponentName = cms.string('pixelLessStepClusterShapeHitFilter'),
110  PixelShapeFile= cms.string('RecoPixelVertexing/PixelLowPtUtilities/data/pixelShape.par'),
111  doStripShapeCut = cms.bool(False),
112  clusterChargeCut = cms.PSet(refToPSet_ = cms.string('SiStripClusterChargeCutTight'))
113  )
114 import RecoPixelVertexing.PixelLowPtUtilities.StripSubClusterShapeSeedFilter_cfi
115 pixelLessStepSeeds.SeedComparitorPSet = cms.PSet(
116  ComponentName = cms.string('CombinedSeedComparitor'),
117  mode = cms.string("and"),
118  comparitors = cms.VPSet(
119  cms.PSet(
120  ComponentName = cms.string('PixelClusterShapeSeedComparitor'),
121  FilterAtHelixStage = cms.bool(True),
122  FilterPixelHits = cms.bool(False),
123  FilterStripHits = cms.bool(True),
124  ClusterShapeHitFilterName = cms.string('pixelLessStepClusterShapeHitFilter'),
125  ClusterShapeCacheSrc = cms.InputTag("siPixelClusterShapeCache") # not really needed here since FilterPixelHits=False
126  ),
127  RecoPixelVertexing.PixelLowPtUtilities.StripSubClusterShapeSeedFilter_cfi.StripSubClusterShapeSeedFilter.clone()
128  )
129  )
130 
131 # QUALITY CUTS DURING TRACK BUILDING
133 pixelLessStepTrajectoryFilter = TrackingTools.TrajectoryFiltering.TrajectoryFilter_cff.CkfBaseTrajectoryFilter_block.clone(
134  maxLostHits = 0,
135  minimumNumberOfHits = 4,
136  minPt = 0.1
137  )
138 
139 import RecoTracker.MeasurementDet.Chi2ChargeMeasurementEstimatorESProducer_cfi
140 pixelLessStepChi2Est = RecoTracker.MeasurementDet.Chi2ChargeMeasurementEstimatorESProducer_cfi.Chi2ChargeMeasurementEstimator.clone(
141  ComponentName = cms.string('pixelLessStepChi2Est'),
142  nSigma = cms.double(3.0),
143  MaxChi2 = cms.double(9.0),
144  clusterChargeCut = cms.PSet(refToPSet_ = cms.string('SiStripClusterChargeCutTight'))
145 )
146 
147 # TRACK BUILDING
149 pixelLessStepTrajectoryBuilder = RecoTracker.CkfPattern.GroupedCkfTrajectoryBuilder_cfi.GroupedCkfTrajectoryBuilder.clone(
150  MeasurementTrackerName = '',
151  trajectoryFilter = cms.PSet(refToPSet_ = cms.string('pixelLessStepTrajectoryFilter')),
152  minNrOfHitsForRebuild = 4,
153  maxCand = 2,
154  alwaysUseInvalidHits = False,
155  estimator = cms.string('pixelLessStepChi2Est'),
156  maxDPhiForLooperReconstruction = cms.double(2.0),
157  maxPtForLooperReconstruction = cms.double(0.7)
158  )
159 
160 # MAKING OF TRACK CANDIDATES
162 pixelLessStepTrackCandidates = RecoTracker.CkfPattern.CkfTrackCandidates_cfi.ckfTrackCandidates.clone(
163  src = cms.InputTag('pixelLessStepSeeds'),
164  clustersToSkip = cms.InputTag('pixelLessStepClusters'),
165  ### these two parameters are relevant only for the CachingSeedCleanerBySharedInput
166  numHitsForSeedCleaner = cms.int32(50),
167  #onlyPixelHitsForSeedCleaner = cms.bool(True),
168  TrajectoryBuilderPSet = cms.PSet(refToPSet_ = cms.string('pixelLessStepTrajectoryBuilder'))
169 )
170 
171 from TrackingTools.TrajectoryCleaning.TrajectoryCleanerBySharedHits_cfi import trajectoryCleanerBySharedHits
172 pixelLessStepTrajectoryCleanerBySharedHits = trajectoryCleanerBySharedHits.clone(
173  ComponentName = cms.string('pixelLessStepTrajectoryCleanerBySharedHits'),
174  fractionShared = cms.double(0.11),
175  allowSharedFirstHit = cms.bool(True)
176  )
177 pixelLessStepTrackCandidates.TrajectoryCleaner = 'pixelLessStepTrajectoryCleanerBySharedHits'
178 
179 
180 # TRACK FITTING
182 pixelLessStepTracks = RecoTracker.TrackProducer.TrackProducer_cfi.TrackProducer.clone(
183  src = 'pixelLessStepTrackCandidates',
184  AlgorithmName = cms.string('pixelLessStep'),
185  Fitter = cms.string('FlexibleKFFittingSmoother')
186  )
187 
188 import RecoTracker.FinalTrackSelectors.multiTrackSelector_cfi
189 from RecoTracker.IterativeTracking.MixedTripletStep_cff import mixedTripletStepSelector
190 pixelLessStepSelector = RecoTracker.FinalTrackSelectors.multiTrackSelector_cfi.multiTrackSelector.clone(
191  src='pixelLessStepTracks',
192  useAnyMVA = cms.bool(True),
193  GBRForestLabel = cms.string('MVASelectorIter5_13TeV_v0'),
194  trackSelectors= cms.VPSet(
195  RecoTracker.FinalTrackSelectors.multiTrackSelector_cfi.looseMTS.clone(
196  name = 'pixelLessStepLoose',
197  chi2n_par = 9999,
198  useMVA = cms.bool(True),
199  minMVA = cms.double(0.4),
200  #chi2n_par = 0.4,
201  #res_par = ( 0.003, 0.001 ),
202  #minNumberLayers = 4,
203  #maxNumberLostLayers = 1,
204  #minNumber3DLayers = 3,
205  d0_par1 = ( 1.2, 4.0 ),
206  dz_par1 = ( 1.2, 4.0 ),
207  d0_par2 = ( 1.2, 4.0 ),
208  dz_par2 = ( 1.2, 4.0 )
209  ),
210  RecoTracker.FinalTrackSelectors.multiTrackSelector_cfi.tightMTS.clone(
211  name = 'pixelLessStepTight',
212  preFilterName = 'pixelLessStepLoose',
213  chi2n_par = 0.3,
214  res_par = ( 0.003, 0.001 ),
215  minNumberLayers = 4,
216  maxNumberLostLayers = 0,
217  minNumber3DLayers = 3,
218  d0_par1 = ( 0.9, 4.0 ),
219  dz_par1 = ( 0.9, 4.0 ),
220  d0_par2 = ( 0.9, 4.0 ),
221  dz_par2 = ( 0.9, 4.0 )
222  ),
223  RecoTracker.FinalTrackSelectors.multiTrackSelector_cfi.looseMTS.clone(
224  name = 'pixelLessStep',
225  preFilterName = 'pixelLessStepLoose',
226  chi2n_par = cms.double(9999),
227  useMVA = cms.bool(True),
228  minMVA = cms.double(0.4),
229  qualityBit = cms.string('highPurity'),
230  keepAllTracks = cms.bool(True),
231  #chi2n_par = 0.2,
232  #res_par = ( 0.003, 0.001 ),
233  #minNumberLayers = 4,
234  #maxNumberLostLayers = 0,
235  #minNumber3DLayers = 3,
236  #max_minMissHitOutOrIn = 2,
237  #max_lostHitFraction = 1.0,
238  d0_par1 = ( 0.7, 4.0 ),
239  dz_par1 = ( 0.7, 4.0 ),
240  d0_par2 = ( 0.7, 4.0 ),
241  dz_par2 = ( 0.7, 4.0 )
242  ),
243  mixedTripletStepSelector.trackSelectors[4].clone(
244  name = 'pixelLessStepVtx',
245  preFilterName=cms.string(''),
246  keepAllTracks = cms.bool(False)
247  ),
248  mixedTripletStepSelector.trackSelectors[5].clone(
249  name = 'pixelLessStepTrk',
250  preFilterName=cms.string(''),
251  keepAllTracks = cms.bool(False)
252  )
253  ) #end of vpset
254  ) #end of clone
255 
256 # need to merge the three sets
257 import RecoTracker.FinalTrackSelectors.trackListMerger_cfi
258 pixelLessStep = RecoTracker.FinalTrackSelectors.trackListMerger_cfi.trackListMerger.clone(
259  TrackProducers = cms.VInputTag(cms.InputTag("pixelLessStepTracks"),
260  cms.InputTag("pixelLessStepTracks"),
261  cms.InputTag("pixelLessStepTracks")),
262  hasSelector=cms.vint32(1,1,1),
263  shareFrac=cms.double(0.11),
264  indivShareFrac=cms.vdouble(0.11,0.11,0.11),
265  selectedTrackQuals = cms.VInputTag(cms.InputTag("pixelLessStepSelector","pixelLessStep"),
266  cms.InputTag("pixelLessStepSelector","pixelLessStepVtx"),
267  cms.InputTag("pixelLessStepSelector","pixelLessStepTrk")),
268  setsToMerge = cms.VPSet( cms.PSet( tLists=cms.vint32(0,1,2), pQual=cms.bool(True) )),
269  writeOnlyTrkQuals=cms.bool(True)
270 )
271 
272 PixelLessStep = cms.Sequence(pixelLessStepClusters*
273  pixelLessStepSeedLayers*
274  pixelLessStepSeeds*
275  pixelLessStepTrackCandidates*
276  pixelLessStepTracks*
277  pixelLessStepSelector*
278  pixelLessStep)
TEveGeoShape * clone(const TEveElement *element, TEveElement *parent)
Definition: eve_macros.cc:135