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