File indexing completed on 2024-04-06 12:09:13
0001 import FWCore.ParameterSet.Config as cms
0002 import FWCore.Utilities.FileUtils as FileUtils
0003 from FWCore.ParameterSet.VarParsing import VarParsing
0004
0005 options = VarParsing('analysis')
0006 options.register('inputTag',
0007 'LayerRot_9p43e-6',
0008 VarParsing.multiplicity.singleton,
0009 VarParsing.varType.string,
0010 "input tag")
0011 options.register('inputFile',
0012 '/store/relval/CMSSW_14_0_0_pre1/RelValZMM_14/GEN-SIM-RECO/133X_mcRun3_2023_realistic_v3-v1/2590000/586487a4-71be-4b23-b5a4-5662fab803c9.root',
0013 VarParsing.multiplicity.singleton,
0014 VarParsing.varType.string,
0015 "input file")
0016 options.register('isAlCaReco',
0017 False,
0018 VarParsing.multiplicity.singleton,
0019 VarParsing.varType.bool,
0020 "is alcareco input file?")
0021 options.register('isUnitTest',
0022 False,
0023 VarParsing.multiplicity.singleton,
0024 VarParsing.varType.bool,
0025 "is this configuration run in unit test?")
0026 options.parseArguments()
0027
0028 from Configuration.Eras.Era_Run3_cff import Run3
0029 process = cms.Process("TrackingResolution", Run3)
0030
0031
0032
0033
0034 process.load('Configuration.StandardSequences.Services_cff')
0035 process.load('SimGeneral.HepPDTESSource.pythiapdt_cfi')
0036 process.load('FWCore.MessageService.MessageLogger_cfi')
0037 process.MessageLogger.cerr.FwkReport.reportEvery = (100 if options.isUnitTest else 100000)
0038 process.load('Configuration.EventContent.EventContent_cff')
0039 process.load('Configuration.StandardSequences.GeometryRecoDB_cff')
0040 process.load('Configuration.StandardSequences.MagneticField_cff')
0041 process.load('DQMOffline.Configuration.DQMOffline_cff')
0042 process.load('Configuration.StandardSequences.EndOfProcess_cff')
0043 process.load('Configuration.StandardSequences.FrontierConditions_GlobalTag_cff')
0044
0045
0046
0047
0048 process.load("RecoVertex.BeamSpotProducer.BeamSpot_cfi")
0049
0050
0051
0052
0053 process.load("RecoTracker.MeasurementDet.MeasurementTrackerEventProducer_cfi")
0054 if(options.isAlCaReco):
0055
0056 process.MeasurementTrackerEvent.pixelClusterProducer = "ALCARECOTkAlDiMuon"
0057 process.MeasurementTrackerEvent.stripClusterProducer = "ALCARECOTkAlDiMuon"
0058 process.MeasurementTrackerEvent.inactivePixelDetectorLabels = cms.VInputTag()
0059 process.MeasurementTrackerEvent.inactiveStripDetectorLabels = cms.VInputTag()
0060
0061 process.maxEvents = cms.untracked.PSet(
0062 input = cms.untracked.int32(10 if options.isUnitTest else -1)
0063 )
0064
0065
0066
0067
0068
0069
0070
0071
0072 readFiles = cms.untracked.vstring(options.inputFile)
0073 process.source = cms.Source("PoolSource",fileNames = readFiles)
0074
0075 process.options = cms.untracked.PSet()
0076
0077
0078
0079
0080 process.DQMoutput = cms.OutputModule("DQMRootOutputModule",
0081 dataset = cms.untracked.PSet(
0082 dataTier = cms.untracked.string('DQMIO'),
0083 filterName = cms.untracked.string('')
0084 ),
0085 fileName = cms.untracked.string('file:step1_DQM_'+options.inputTag+'_'+('fromALCA' if options.isAlCaReco else 'fromRECO' )+'.root'),
0086 outputCommands = process.DQMEventContent.outputCommands,
0087 splitLevel = cms.untracked.int32(0)
0088 )
0089
0090
0091
0092
0093 from Configuration.AlCa.GlobalTag import GlobalTag
0094
0095 process.GlobalTag = GlobalTag(process.GlobalTag, "125X_mcRun3_2022_design_v6", '')
0096 process.GlobalTag.toGet = cms.VPSet(cms.PSet(connect = cms.string("frontier://FrontierProd/CMS_CONDITIONS"),
0097 record = cms.string('TrackerAlignmentRcd'),
0098 tag = cms.string(options.inputTag)))
0099
0100
0101
0102
0103 process.load("DQM.TrackingMonitorSource.shortTrackResolution_cff")
0104 process.load("RecoTracker.TrackProducer.TrackRefitters_cff")
0105 import RecoTracker.TrackProducer.TrackRefitters_cff
0106 process.LongTracksRefit = process.TrackRefitter.clone(
0107 src = 'SingleLongTrackProducer',
0108 TrajectoryInEvent = True,
0109 TTRHBuilder = "WithAngleAndTemplate",
0110 NavigationSchool = ''
0111 )
0112
0113 process.ShortTrackCandidates3.src = cms.InputTag("LongTracksRefit")
0114 process.ShortTrackCandidates4.src = cms.InputTag("LongTracksRefit")
0115 process.ShortTrackCandidates5.src = cms.InputTag("LongTracksRefit")
0116 process.ShortTrackCandidates6.src = cms.InputTag("LongTracksRefit")
0117 process.ShortTrackCandidates7.src = cms.InputTag("LongTracksRefit")
0118 process.ShortTrackCandidates8.src = cms.InputTag("LongTracksRefit")
0119
0120 process.SingleLongTrackProducer.matchMuons = cms.InputTag("muons")
0121 if(options.isAlCaReco):
0122 process.SingleLongTrackProducer.requiredDr = cms.double(-9999.)
0123 process.SingleLongTrackProducer.allTracks = cms.InputTag("ALCARECOTkAlDiMuon")
0124
0125
0126
0127
0128 process.analysis_step = cms.Path(process.offlineBeamSpot *
0129 process.MeasurementTrackerEvent *
0130 process.SingleLongTrackProducer *
0131 process.LongTracksRefit *
0132 process.ShortTrackCandidates3 *
0133 process.ShortTrackCandidates4 *
0134 process.ShortTrackCandidates5 *
0135 process.ShortTrackCandidates6 *
0136 process.ShortTrackCandidates7 *
0137 process.ShortTrackCandidates8 *
0138 process.RefittedShortTracks3 *
0139 process.RefittedShortTracks4 *
0140 process.RefittedShortTracks5 *
0141 process.RefittedShortTracks6 *
0142 process.RefittedShortTracks7 *
0143 process.RefittedShortTracks8 *
0144 process.trackingResolution)
0145
0146
0147
0148
0149 process.endjob_step = cms.EndPath(process.endOfProcess)
0150 process.DQMoutput_step = cms.EndPath(process.DQMoutput)
0151
0152 process.schedule = cms.Schedule(process.analysis_step, process.endjob_step, process.DQMoutput_step)
0153
0154
0155
0156
0157 process.options.numberOfThreads = 8