Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2022-01-11 23:32:41

0001 import FWCore.ParameterSet.Config as cms
0002 
0003 # This step runs over all clusters
0004 
0005 # run only if there are high pT jets
0006 from RecoJets.JetProducers.TracksForJets_cff import trackRefsForJets
0007 hiInitialStepTrackRefsForJets = trackRefsForJets.clone(src = 'hiGlobalPrimTracks')
0008 
0009 #change this to import Bkg substracted Heavy Ion jets:
0010 from RecoHI.HiJetAlgos.hiCaloJetsForTrk_cff import *
0011 
0012 hiJetsForCoreTracking = cms.EDFilter("CandPtrSelector", src = cms.InputTag("akPu4CaloJetsSelected"), cut = cms.string("pt > 30 && abs(eta) < 2.5"))
0013 
0014 # care only at tracks from main PV
0015 hiFirstStepGoodPrimaryVertices = cms.EDFilter("PrimaryVertexObjectFilter",
0016      filterParams = cms.PSet(
0017              minNdof = cms.double(25.0),
0018              maxZ = cms.double(15.0),
0019              maxRho = cms.double(2.0)
0020      ),
0021      src=cms.InputTag('hiSelectedPixelVertex')
0022 )
0023 
0024 import RecoTracker.TkSeedingLayers.seedingLayersEDProducer_cfi as _mod
0025 
0026 # SEEDING LAYERS
0027 hiJetCoreRegionalStepSeedLayers = _mod.seedingLayersEDProducer.clone(
0028     layerList = ['BPix1+BPix2+BPix3', 
0029                  'BPix1+BPix2+FPix1_pos', 
0030                  'BPix1+BPix2+FPix1_neg', 
0031                  'BPix1+FPix1_pos+FPix2_pos', 
0032                  'BPix1+FPix1_neg+FPix2_neg',
0033                  'BPix1+BPix2+TIB1', 
0034                  'BPix1+BPix3+TIB1', 
0035                  'BPix2+BPix3+TIB1', 
0036                 ],
0037     TIB = dict(
0038         matchedRecHits = cms.InputTag("siStripMatchedRecHits","matchedRecHit"),
0039         TTRHBuilder = cms.string('WithTrackAngle'), 
0040         clusterChargeCut = cms.PSet(
0041             refToPSet_ = cms.string('SiStripClusterChargeCutNone')
0042         )
0043     ),
0044     BPix = dict(
0045         useErrorsFromParam = cms.bool(True),
0046         hitErrorRPhi = cms.double(0.0027),
0047         hitErrorRZ = cms.double(0.006),
0048         TTRHBuilder = cms.string('WithTrackAngle'),
0049         HitProducer = cms.string('siPixelRecHits'),
0050     ),
0051     FPix = dict(
0052         useErrorsFromParam = cms.bool(True),
0053         hitErrorRPhi = cms.double(0.0051),
0054         hitErrorRZ = cms.double(0.0036),
0055         TTRHBuilder = cms.string('WithTrackAngle'),
0056         HitProducer = cms.string('siPixelRecHits'),
0057     )
0058 )
0059 from Configuration.Eras.Modifier_trackingPhase1_cff import trackingPhase1
0060 trackingPhase1.toModify(hiJetCoreRegionalStepSeedLayers, layerList = ['BPix1+BPix2+BPix3',
0061     'BPix2+BPix3+BPix4',
0062     'BPix1+BPix3+BPix4',
0063     'BPix1+BPix2+BPix4',
0064     'BPix2+BPix3+FPix1_pos',
0065     'BPix2+BPix3+FPix1_neg',
0066     'BPix1+BPix2+FPix1_pos',
0067     'BPix1+BPix2+FPix1_neg',
0068     'BPix2+FPix1_pos+FPix2_pos',
0069     'BPix2+FPix1_neg+FPix2_neg',
0070     'BPix1+FPix1_pos+FPix2_pos',
0071     'BPix1+FPix1_neg+FPix2_neg',
0072     'FPix1_pos+FPix2_pos+FPix3_pos',
0073     'FPix1_neg+FPix2_neg+FPix3_neg',#up to here, same as what is in RecoTracker/TkSeedingLayers/python/PixelLayerTriplets_cfi.py for phase 1
0074     'BPix1+BPix2+TIB1',#use TIB1 to try to recover tracks w/ 2 hits missing in pix barrel
0075     'BPix1+BPix3+TIB1',
0076     'BPix1+BPix4+TIB1',
0077     'BPix2+BPix3+TIB1',
0078     'BPix2+BPix4+TIB1',
0079     'BPix3+BPix4+TIB1',
0080     ]
0081 )
0082 
0083 # SEEDS
0084 import RecoTracker.TkSeedGenerator.GlobalSeedsFromTriplets_cff
0085 hiJetCoreRegionalStepSeeds = RecoTracker.TkSeedGenerator.GlobalSeedsFromTriplets_cff.globalSeedsFromTriplets.clone(
0086     RegionFactoryPSet = cms.PSet(
0087     ComponentName = cms.string( "TauRegionalPixelSeedGenerator" ),#not so nice to depend on RecoTau...
0088     RegionPSet = cms.PSet(
0089        precise = cms.bool( True ),
0090        useMultipleScattering = cms.bool(False),
0091        useFakeVertices       = cms.bool(False),
0092        originRadius = cms.double( 0.2 ),
0093        ptMin = cms.double( 15. ),
0094        originHalfLength = cms.double( 0.2 ),
0095        deltaPhiRegion = cms.double( 0.30 ), 
0096        deltaEtaRegion = cms.double( 0.30 ), 
0097        JetSrc = cms.InputTag( "hiJetsForCoreTracking" ),
0098        vertexSrc = cms.InputTag( "hiFirstStepGoodPrimaryVertices" ),
0099        measurementTrackerName = cms.InputTag( "MeasurementTrackerEvent" ),
0100        howToUseMeasurementTracker = cms.string( "Never" )
0101     )
0102     ),
0103     OrderedHitsFactoryPSet = dict(SeedingLayers = 'hiJetCoreRegionalStepSeedLayers'),
0104     SeedComparitorPSet     = dict(ComponentName = 'none'),
0105     SeedCreatorPSet        = dict(forceKinematicWithRegionDirection = True),
0106     ClusterCheckPSet       = dict(doClusterCheck = False)
0107 )
0108 # QUALITY CUTS DURING TRACK BUILDING
0109 import TrackingTools.TrajectoryFiltering.TrajectoryFilter_cff
0110 hiJetCoreRegionalStepTrajectoryFilter = TrackingTools.TrajectoryFiltering.TrajectoryFilter_cff.CkfBaseTrajectoryFilter_block.clone(
0111     minimumNumberOfHits = 6,
0112     minPt = 10.0
0113 )
0114 
0115 import TrackingTools.KalmanUpdators.Chi2MeasurementEstimator_cfi
0116 hiJetCoreRegionalStepChi2Est = TrackingTools.KalmanUpdators.Chi2MeasurementEstimator_cfi.Chi2MeasurementEstimator.clone(
0117     ComponentName = 'hiJetCoreRegionalStepChi2Est',
0118     nSigma  = 3.0,
0119     MaxChi2 = 30.0
0120 )
0121 
0122 # TRACK BUILDING
0123 import RecoTracker.CkfPattern.GroupedCkfTrajectoryBuilder_cfi
0124 hiJetCoreRegionalStepTrajectoryBuilder = RecoTracker.CkfPattern.GroupedCkfTrajectoryBuilder_cfi.GroupedCkfTrajectoryBuilder.clone(
0125     trajectoryFilter = dict(refToPSet_ = 'hiJetCoreRegionalStepTrajectoryFilter'),
0126     maxCand = 50,
0127     estimator = 'hiJetCoreRegionalStepChi2Est',
0128     maxDPhiForLooperReconstruction = 2.0,
0129     maxPtForLooperReconstruction = 0.7,
0130 )
0131 
0132 # MAKING OF TRACK CANDIDATES
0133 import RecoTracker.CkfPattern.CkfTrackCandidates_cfi
0134 hiJetCoreRegionalStepTrackCandidates = RecoTracker.CkfPattern.CkfTrackCandidates_cfi.ckfTrackCandidates.clone(
0135     src = 'hiJetCoreRegionalStepSeeds',
0136     maxSeedsBeforeCleaning = 10000,
0137     TrajectoryBuilderPSet = dict(refToPSet_ = 'hiJetCoreRegionalStepTrajectoryBuilder'),
0138     NavigationSchool = 'SimpleNavigationSchool',
0139 )
0140 
0141 
0142 # TRACK FITTING
0143 import RecoTracker.TrackProducer.TrackProducer_cfi
0144 hiJetCoreRegionalStepTracks = RecoTracker.TrackProducer.TrackProducer_cfi.TrackProducer.clone(
0145     AlgorithmName = 'jetCoreRegionalStep',
0146     src = 'hiJetCoreRegionalStepTrackCandidates',
0147     Fitter = 'FlexibleKFFittingSmoother'
0148 )
0149 
0150 # Final selection
0151 import RecoHI.HiTracking.hiMultiTrackSelector_cfi
0152 hiJetCoreRegionalStepSelector = RecoHI.HiTracking.hiMultiTrackSelector_cfi.hiMultiTrackSelector.clone(
0153     src = 'hiJetCoreRegionalStepTracks',
0154     vertices = "hiFirstStepGoodPrimaryVertices",
0155     trackSelectors= cms.VPSet(
0156         RecoHI.HiTracking.hiMultiTrackSelector_cfi.hiLooseMTS.clone(
0157             name = 'hiJetCoreRegionalStepLoose',
0158             ), #end of pset
0159         RecoHI.HiTracking.hiMultiTrackSelector_cfi.hiTightMTS.clone(
0160             name = 'hiJetCoreRegionalStepTight',
0161             preFilterName = 'hiJetCoreRegionalStepLoose',
0162             ),
0163         RecoHI.HiTracking.hiMultiTrackSelector_cfi.hiHighpurityMTS.clone(
0164             name = 'hiJetCoreRegionalStep',
0165             preFilterName = 'hiJetCoreRegionalStepTight',
0166             min_nhits = 14
0167             ),
0168         ) #end of vpset
0169     ) #end of clone
0170 
0171 # Final sequence
0172 hiJetCoreRegionalStepTask = cms.Task(
0173                                    hiCaloJetsForTrkTask,hiJetsForCoreTracking,
0174                                    hiFirstStepGoodPrimaryVertices,
0175                                    hiJetCoreRegionalStepSeedLayers,
0176                                    hiJetCoreRegionalStepSeeds,
0177                                    hiJetCoreRegionalStepTrackCandidates,
0178                                    hiJetCoreRegionalStepTracks,
0179                                    hiJetCoreRegionalStepSelector)
0180 hiJetCoreRegionalStep = cms.Sequence(hiJetCoreRegionalStepTask)