CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
hiSecondPixelTripletStep_cff.py
Go to the documentation of this file.
2 
3 
4 # NEW CLUSTERS (remove previously used clusters)
5 hiSecondPixelTripletClusters = cms.EDProducer("HITrackClusterRemover",
6  clusterLessSolution= cms.bool(True),
7  trajectories = cms.InputTag("hiGlobalPrimTracks"),
8  overrideTrkQuals = cms.InputTag("hiInitialStepSelector","hiInitialStep"),
9  TrackQuality = cms.string('highPurity'),
10  minNumberOfLayersWithMeasBeforeFiltering = cms.int32(0),
11  pixelClusters = cms.InputTag("siPixelClusters"),
12  stripClusters = cms.InputTag("siStripClusters"),
13  Common = cms.PSet(
14  maxChi2 = cms.double(9.0)
15  ),
16  Strip = cms.PSet(
17  #Yen-Jie's mod to preserve merged clusters
18  maxSize = cms.uint32(2),
19  maxChi2 = cms.double(9.0)
20  )
21  )
22 
23 
24 # SEEDING LAYERS
26 hiSecondPixelTripletSeedLayers = RecoTracker.TkSeedingLayers.PixelLayerTriplets_cfi.PixelLayerTriplets.clone()
27 hiSecondPixelTripletSeedLayers.BPix.skipClusters = cms.InputTag('hiSecondPixelTripletClusters')
28 hiSecondPixelTripletSeedLayers.FPix.skipClusters = cms.InputTag('hiSecondPixelTripletClusters')
29 
30 # SEEDS
35 hiSecondPixelTracks = cms.EDProducer("PixelTrackProducer",
36 
37  passLabel = cms.string('Pixel primary tracks with vertex constraint'),
38 
39  # Region
40  RegionFactoryPSet = cms.PSet(
41  ComponentName = cms.string("GlobalTrackingRegionWithVerticesProducer"),
42  RegionPSet = cms.PSet(
43  precise = cms.bool(True),
44  beamSpot = cms.InputTag("offlineBeamSpot"),
45  useFixedError = cms.bool(False),
46  nSigmaZ = cms.double(4.0),
47  sigmaZVertex = cms.double(4.0),
48  fixedError = cms.double(0.2),
49  VertexCollection = cms.InputTag("hiSelectedVertex"),
50  ptMin = cms.double(0.4),
51  useFoundVertices = cms.bool(True),
52  originRadius = cms.double(0.02)
53  )
54  ),
55 
56  # Ordered Hits
57  OrderedHitsFactoryPSet = cms.PSet(
58  ComponentName = cms.string( "StandardHitTripletGenerator" ),
59  SeedingLayers = cms.InputTag( "PixelLayerTriplets" ),
60  GeneratorPSet = cms.PSet(
61  PixelTripletHLTGenerator
62  )
63  ),
64 
65  # Fitter
66  FitterPSet = cms.PSet(
67  ComponentName = cms.string('PixelFitterByHelixProjections'),
68  TTRHBuilder = cms.string('TTRHBuilderWithoutAngle4PixelTriplets')
69  ),
70 
71  # Filter
72  useFilterWithES = cms.bool( True ),
73  FilterPSet = cms.PSet(
74  nSigmaLipMaxTolerance = cms.double(4.0),
75  chi2 = cms.double(1000.0),
76  ComponentName = cms.string('HIPixelTrackFilter'),
77  nSigmaTipMaxTolerance = cms.double(4.0),
78  clusterShapeCacheSrc = cms.InputTag("siPixelClusterShapeCache"),
79  VertexCollection = cms.InputTag("hiSelectedVertex"),
80  useClusterShape = cms.bool(False),
81  lipMax = cms.double(0),
82  tipMax = cms.double(0),
83  ptMin = cms.double(0.4)
84  ),
85 
86  # Cleaner
87  CleanerPSet = cms.PSet(
88  ComponentName = cms.string( "TrackCleaner" )
89  )
90 )
91 
92 hiSecondPixelTracks.OrderedHitsFactoryPSet.GeneratorPSet.maxElement = cms.uint32(5000000)
93 hiSecondPixelTracks.OrderedHitsFactoryPSet.SeedingLayers = cms.InputTag('hiSecondPixelTripletSeedLayers')
94 
95 hiSecondPixelTracks.OrderedHitsFactoryPSet.GeneratorPSet.SeedComparitorPSet = RecoPixelVertexing.PixelLowPtUtilities.LowPtClusterShapeSeedComparitor_cfi.LowPtClusterShapeSeedComparitor
96 
97 
98 
99 import RecoPixelVertexing.PixelLowPtUtilities.TrackSeeds_cfi
100 hiSecondPixelTrackSeeds = RecoPixelVertexing.PixelLowPtUtilities.TrackSeeds_cfi.pixelTrackSeeds.clone(
101  InputCollection = 'hiSecondPixelTracks'
102  )
103 
104 
105 
107 from RecoTracker.TkTrackingRegions.GlobalTrackingRegionFromBeamSpot_cfi import RegionPsetFomBeamSpotBlock
108 hiSecondPixelTripletSeeds = RecoTracker.TkSeedGenerator.GlobalSeedsFromTriplets_cff.globalSeedsFromTriplets.clone(
109  RegionFactoryPSet = RegionPsetFomBeamSpotBlock.clone(
110  ComponentName = cms.string('GlobalTrackingRegionWithVerticesProducer'),
111  RegionPSet = cms.PSet(
112  precise = cms.bool(True),
113  beamSpot = cms.InputTag("offlineBeamSpot"),
114  useFixedError = cms.bool(False), #def value True
115  nSigmaZ = cms.double(4.0),
116  sigmaZVertex = cms.double(4.0), #def value 3
117  fixedError = cms.double(0.2),
118  VertexCollection = cms.InputTag("hiSelectedVertex"),
119  ptMin = cms.double(0.4),
120  useFoundVertices = cms.bool(True),
121  originRadius = cms.double(0.02)
122  )
123  )
124 )
125 
126 hiSecondPixelTripletSeeds.OrderedHitsFactoryPSet.SeedingLayers = 'hiSecondPixelTripletSeedLayers'
127 hiSecondPixelTripletSeeds.OrderedHitsFactoryPSet.GeneratorPSet.maxElement = 5000000
128 hiSecondPixelTripletSeeds.ClusterCheckPSet.MaxNumberOfPixelClusters = 5000000
129 hiSecondPixelTripletSeeds.ClusterCheckPSet.MaxNumberOfCosmicClusters = 50000000
130 del hiSecondPixelTripletSeeds.ClusterCheckPSet.cut
131 
133 import RecoPixelVertexing.PixelLowPtUtilities.LowPtClusterShapeSeedComparitor_cfi
134 hiSecondPixelTripletSeeds.OrderedHitsFactoryPSet.GeneratorPSet.SeedComparitorPSet = RecoPixelVertexing.PixelLowPtUtilities.LowPtClusterShapeSeedComparitor_cfi.LowPtClusterShapeSeedComparitor
135 
136 # QUALITY CUTS DURING TRACK BUILDING
138 hiSecondPixelTripletTrajectoryFilter = TrackingTools.TrajectoryFiltering.TrajectoryFilter_cff.CkfBaseTrajectoryFilter_block.clone(
139  maxLostHits = 1,
140  minimumNumberOfHits = 6,
141  minPt = 0.4
142  )
143 
145 hiSecondPixelTripletChi2Est = TrackingTools.KalmanUpdators.Chi2MeasurementEstimatorESProducer_cfi.Chi2MeasurementEstimator.clone(
146  ComponentName = cms.string('hiSecondPixelTripletChi2Est'),
147  nSigma = cms.double(3.0),
148  MaxChi2 = cms.double(9.0)
149  )
150 
151 # TRACK BUILDING
153 hiSecondPixelTripletTrajectoryBuilder = RecoTracker.CkfPattern.GroupedCkfTrajectoryBuilder_cfi.GroupedCkfTrajectoryBuilder.clone(
154  MeasurementTrackerName = '',
155  trajectoryFilter = cms.PSet(refToPSet_ = cms.string('hiSecondPixelTripletTrajectoryFilter')),
156  maxCand = 3,
157  estimator = cms.string('hiSecondPixelTripletChi2Est'),
158  maxDPhiForLooperReconstruction = cms.double(2.0),
159  # 0.63 GeV is the maximum pT for a charged particle to loop within the 1.1m radius
160  # of the outermost Tracker barrel layer (with B=3.8T)
161  maxPtForLooperReconstruction = cms.double(0.7)
162  )
163 
164 # MAKING OF TRACK CANDIDATES
166 hiSecondPixelTripletTrackCandidates = RecoTracker.CkfPattern.CkfTrackCandidates_cfi.ckfTrackCandidates.clone(
167  src = cms.InputTag('hiSecondPixelTrackSeeds'),
168  ### these two parameters are relevant only for the CachingSeedCleanerBySharedInput
169  numHitsForSeedCleaner = cms.int32(50),
170  onlyPixelHitsForSeedCleaner = cms.bool(True),
171  TrajectoryBuilderPSet = cms.PSet(refToPSet_ = cms.string('hiSecondPixelTripletTrajectoryBuilder')),
172  clustersToSkip = cms.InputTag('hiSecondPixelTripletClusters'),
173  doSeedingRegionRebuilding = True,
174  useHitsSplitting = True
175  )
176 
177 # TRACK FITTING
179 hiSecondPixelTripletGlobalPrimTracks = RecoTracker.TrackProducer.TrackProducer_cfi.TrackProducer.clone(
180  src = 'hiSecondPixelTripletTrackCandidates',
181  AlgorithmName = cms.string('lowPtTripletStep'),
182  Fitter=cms.string('FlexibleKFFittingSmoother')
183  )
184 
185 
186 
187 # Final selection
188 import RecoHI.HiTracking.hiMultiTrackSelector_cfi
189 hiSecondPixelTripletStepSelector = RecoHI.HiTracking.hiMultiTrackSelector_cfi.hiMultiTrackSelector.clone(
190  src='hiSecondPixelTripletGlobalPrimTracks',
191  trackSelectors= cms.VPSet(
192  RecoHI.HiTracking.hiMultiTrackSelector_cfi.hiLooseMTS.clone(
193  name = 'hiSecondPixelTripletStepLoose',
194  ), #end of pset
195  RecoHI.HiTracking.hiMultiTrackSelector_cfi.hiTightMTS.clone(
196  name = 'hiSecondPixelTripletStepTight',
197  preFilterName = 'hiSecondPixelTripletStepLoose',
198  ),
199  RecoHI.HiTracking.hiMultiTrackSelector_cfi.hiHighpurityMTS.clone(
200  name = 'hiSecondPixelTripletStep',
201  preFilterName = 'hiSecondPixelTripletStepTight',
202  # min_nhits = 14
203  ),
204  ) #end of vpset
205  ) #end of clone
206 
207 
208 import RecoTracker.FinalTrackSelectors.trackListMerger_cfi
209 hiSecondQual = RecoTracker.FinalTrackSelectors.trackListMerger_cfi.trackListMerger.clone(
210  TrackProducers = cms.VInputTag(cms.InputTag('hiSecondPixelTripletGlobalPrimTracks')),
211  hasSelector=cms.vint32(1),
212  selectedTrackQuals = cms.VInputTag(cms.InputTag("hiSecondPixelTripletStepSelector","hiSecondPixelTripletStep")),
213  copyExtras = True,
214  makeReKeyedSeeds = cms.untracked.bool(False),
215  #writeOnlyTrkQuals = True
216  )
217 
218 # Final sequence
219 
220 hiSecondPixelTripletStep = cms.Sequence(hiSecondPixelTripletClusters*
221  hiSecondPixelTripletSeedLayers*
222  hiSecondPixelTracks*hiSecondPixelTrackSeeds*
223  hiSecondPixelTripletTrackCandidates*
224  hiSecondPixelTripletGlobalPrimTracks*
225  hiSecondPixelTripletStepSelector
226  *hiSecondQual
227  )
pp iterative tracking modified for hiOffline reco (the vertex is the one reconstructed in HI) 3rd ste...