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