CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
hiLowPtTripletStep_cff.py
Go to the documentation of this file.
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
36 hiLowPtTripletStepPixelTracks = cms.EDProducer("PixelTrackProducer",
37 
38  passLabel = cms.string('Pixel primary tracks with vertex constraint'),
39 
40  # Region
41  RegionFactoryPSet = cms.PSet(
42  ComponentName = cms.string("GlobalTrackingRegionWithVerticesProducer"),
43  RegionPSet = cms.PSet(
44  precise = cms.bool(True),
45  beamSpot = cms.InputTag("offlineBeamSpot"),
46  useFixedError = cms.bool(False),
47  nSigmaZ = cms.double(4.0),
48  sigmaZVertex = cms.double(4.0),
49  fixedError = cms.double(0.2),
50  VertexCollection = cms.InputTag("hiSelectedVertex"),
51  ptMin = cms.double(0.4),
52  useFoundVertices = cms.bool(True),
53  originRadius = cms.double(0.02)
54  )
55  ),
56 
57  # Ordered Hits
58  OrderedHitsFactoryPSet = cms.PSet(
59  ComponentName = cms.string( "StandardHitTripletGenerator" ),
60  SeedingLayers = cms.InputTag( "PixelLayerTriplets" ),
61  GeneratorPSet = cms.PSet(
62  PixelTripletHLTGenerator
63  )
64  ),
65 
66  # Fitter
67  FitterPSet = cms.PSet(
68  ComponentName = cms.string('PixelFitterByHelixProjections'),
69  TTRHBuilder = cms.string('TTRHBuilderWithoutAngle4PixelTriplets')
70  ),
71 
72  # Filter
73  useFilterWithES = cms.bool( True ),
74  FilterPSet = cms.PSet(
75  nSigmaLipMaxTolerance = cms.double(4.0),
76  chi2 = cms.double(1000.0),
77  ComponentName = cms.string('HIPixelTrackFilter'),
78  nSigmaTipMaxTolerance = cms.double(4.0),
79  clusterShapeCacheSrc = cms.InputTag("siPixelClusterShapeCache"),
80  VertexCollection = cms.InputTag("hiSelectedVertex"),
81  useClusterShape = cms.bool(False),
82  lipMax = cms.double(0),
83  tipMax = cms.double(0),
84  ptMin = cms.double(0.4)
85  ),
86 
87  # Cleaner
88  CleanerPSet = cms.PSet(
89  ComponentName = cms.string( "TrackCleaner" )
90  )
91 )
92 
93 hiLowPtTripletStepPixelTracks.OrderedHitsFactoryPSet.GeneratorPSet.maxElement = cms.uint32(5000000)
94 hiLowPtTripletStepPixelTracks.OrderedHitsFactoryPSet.SeedingLayers = cms.InputTag('hiLowPtTripletStepSeedLayers')
95 
96 hiLowPtTripletStepPixelTracks.OrderedHitsFactoryPSet.GeneratorPSet.SeedComparitorPSet = RecoPixelVertexing.PixelLowPtUtilities.LowPtClusterShapeSeedComparitor_cfi.LowPtClusterShapeSeedComparitor
97 
98 
99 
100 import RecoPixelVertexing.PixelLowPtUtilities.TrackSeeds_cfi
101 hiLowPtTripletStepSeeds = RecoPixelVertexing.PixelLowPtUtilities.TrackSeeds_cfi.pixelTrackSeeds.clone(
102  InputCollection = 'hiLowPtTripletStepPixelTracks'
103  )
104 
105 
106 # QUALITY CUTS DURING TRACK BUILDING
108 hiLowPtTripletStepTrajectoryFilter = TrackingTools.TrajectoryFiltering.TrajectoryFilter_cff.CkfBaseTrajectoryFilter_block.clone(
109  maxLostHits = 1,
110  minimumNumberOfHits = 6,
111  minPt = 0.4
112  )
113 
115 hiLowPtTripletStepChi2Est = TrackingTools.KalmanUpdators.Chi2MeasurementEstimatorESProducer_cfi.Chi2MeasurementEstimator.clone(
116  ComponentName = cms.string('hiLowPtTripletStepChi2Est'),
117  nSigma = cms.double(3.0),
118  MaxChi2 = cms.double(9.0)
119  )
120 
121 # TRACK BUILDING
123 hiLowPtTripletStepTrajectoryBuilder = RecoTracker.CkfPattern.GroupedCkfTrajectoryBuilder_cfi.GroupedCkfTrajectoryBuilder.clone(
124  MeasurementTrackerName = '',
125  trajectoryFilter = cms.PSet(refToPSet_ = cms.string('hiLowPtTripletStepTrajectoryFilter')),
126  maxCand = 3,
127  estimator = cms.string('hiLowPtTripletStepChi2Est'),
128  maxDPhiForLooperReconstruction = cms.double(2.0),
129  # 0.63 GeV is the maximum pT for a charged particle to loop within the 1.1m radius
130  # of the outermost Tracker barrel layer (with B=3.8T)
131  maxPtForLooperReconstruction = cms.double(0.7)
132  )
133 
134 # MAKING OF TRACK CANDIDATES
136 hiLowPtTripletStepTrackCandidates = RecoTracker.CkfPattern.CkfTrackCandidates_cfi.ckfTrackCandidates.clone(
137  src = cms.InputTag('hiLowPtTripletStepSeeds'),
138  ### these two parameters are relevant only for the CachingSeedCleanerBySharedInput
139  numHitsForSeedCleaner = cms.int32(50),
140  onlyPixelHitsForSeedCleaner = cms.bool(True),
141  TrajectoryBuilderPSet = cms.PSet(refToPSet_ = cms.string('hiLowPtTripletStepTrajectoryBuilder')),
142  clustersToSkip = cms.InputTag('hiLowPtTripletStepClusters'),
143  doSeedingRegionRebuilding = True,
144  useHitsSplitting = True
145  )
146 
147 # TRACK FITTING
149 hiLowPtTripletStepTracks = RecoTracker.TrackProducer.TrackProducer_cfi.TrackProducer.clone(
150  src = 'hiLowPtTripletStepTrackCandidates',
151  AlgorithmName = cms.string('lowPtTripletStep'),
152  Fitter=cms.string('FlexibleKFFittingSmoother')
153  )
154 
155 
156 
157 # Final selection
158 import RecoHI.HiTracking.hiMultiTrackSelector_cfi
159 hiLowPtTripletStepSelector = RecoHI.HiTracking.hiMultiTrackSelector_cfi.hiMultiTrackSelector.clone(
160  src='hiLowPtTripletStepTracks',
161  trackSelectors= cms.VPSet(
162  RecoHI.HiTracking.hiMultiTrackSelector_cfi.hiLooseMTS.clone(
163  name = 'hiLowPtTripletStepLoose',
164  ), #end of pset
165  RecoHI.HiTracking.hiMultiTrackSelector_cfi.hiTightMTS.clone(
166  name = 'hiLowPtTripletStepTight',
167  preFilterName = 'hiLowPtTripletStepLoose',
168  ),
169  RecoHI.HiTracking.hiMultiTrackSelector_cfi.hiHighpurityMTS.clone(
170  name = 'hiLowPtTripletStep',
171  preFilterName = 'hiLowPtTripletStepTight',
172  # min_nhits = 14
173  ),
174  ) #end of vpset
175  ) #end of clone
176 
177 
178 import RecoTracker.FinalTrackSelectors.trackListMerger_cfi
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  hiLowPtTripletStepPixelTracks*hiLowPtTripletStepSeeds*
193  hiLowPtTripletStepTrackCandidates*
194  hiLowPtTripletStepTracks*
195  hiLowPtTripletStepSelector*
196  hiLowPtTripletStepQual
197  )
pp iterative tracking modified for hiOffline reco (the vertex is the one reconstructed in HI) 3rd ste...