File indexing completed on 2024-04-12 23:01:04
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('scenario',
0007 '0',
0008 VarParsing.multiplicity.singleton,
0009 VarParsing.varType.string,
0010 "Name of input misalignment scenario")
0011 options.parseArguments()
0012
0013 valid_scenarios = ['-10e-6','-8e-6','-6e-6','-4e-6','-2e-6','0','2e-6','4e-6','6e-6','8e-6','10e-6']
0014
0015 if options.scenario not in valid_scenarios:
0016 print("Error: Invalid scenario specified. Please choose from the following list: ")
0017 print(valid_scenarios)
0018 exit(1)
0019
0020 process = cms.Process("TrackingResolution")
0021
0022
0023
0024
0025 process.load('Configuration.StandardSequences.Services_cff')
0026 process.load('SimGeneral.HepPDTESSource.pythiapdt_cfi')
0027 process.load('FWCore.MessageService.MessageLogger_cfi')
0028 process.MessageLogger.cerr.FwkReport.reportEvery = 100000
0029 process.load('Configuration.EventContent.EventContent_cff')
0030 process.load('Configuration.StandardSequences.GeometryRecoDB_cff')
0031 process.load('Configuration.StandardSequences.MagneticField_cff')
0032 process.load('Configuration.StandardSequences.FrontierConditions_GlobalTag_cff')
0033
0034
0035
0036
0037 process.load("RecoVertex.BeamSpotProducer.BeamSpot_cfi")
0038
0039
0040
0041
0042 process.load("RecoTracker.MeasurementDet.MeasurementTrackerEventProducer_cfi")
0043 process.MeasurementTrackerEvent.pixelClusterProducer = "ALCARECOTkAlDiMuon"
0044 process.MeasurementTrackerEvent.stripClusterProducer = "ALCARECOTkAlDiMuon"
0045 process.MeasurementTrackerEvent.inactivePixelDetectorLabels = cms.VInputTag()
0046 process.MeasurementTrackerEvent.inactiveStripDetectorLabels = cms.VInputTag()
0047
0048 process.maxEvents = cms.untracked.PSet(
0049 input = cms.untracked.int32(1000000)
0050 )
0051
0052
0053
0054
0055
0056
0057
0058 readFiles = cms.untracked.vstring('/store/mc/Run3Winter23Reco/DYJetsToMuMu_M-50_TuneCP5_13p6TeV-madgraphMLM-pythia8/ALCARECO/TkAlDiMuonAndVertex-TRKDesignNoPU_AlcaRecoTRKMu_designGaussSigmaZ4cm_125X_mcRun3_2022_design_v6-v1/60000/d3af17a5-2409-4551-9c3d-00deb2f3f64f.root')
0059 process.source = cms.Source("PoolSource",fileNames = readFiles)
0060
0061 process.options = cms.untracked.PSet()
0062
0063
0064
0065
0066 process.TFileService = cms.Service("TFileService",
0067 fileName = cms.string("shortenedTrackResolution_LayerRotation_"+options.scenario+".root"))
0068
0069
0070
0071
0072 from Configuration.AlCa.GlobalTag import GlobalTag
0073 process.GlobalTag = GlobalTag(process.GlobalTag, "125X_mcRun3_2022_design_v6", '')
0074 if (options.scenario=='null'):
0075 print("null scenario, do nothing")
0076 pass
0077 else:
0078 process.GlobalTag.toGet = cms.VPSet(cms.PSet(connect = cms.string("frontier://FrontierPrep/CMS_CONDITIONS"),
0079 record = cms.string('TrackerAlignmentRcd'),
0080 tag = cms.string("LayerRotation_"+options.scenario)))
0081
0082
0083
0084
0085 process.load("DQM.TrackingMonitorSource.shortTrackResolution_cff")
0086
0087
0088
0089
0090 process.load("RecoTracker.TrackProducer.TrackRefitters_cff")
0091 import RecoTracker.TrackProducer.TrackRefitters_cff
0092 process.LongTracksRefit = process.TrackRefitter.clone(
0093 src = 'SingleLongTrackProducer',
0094 TrajectoryInEvent = True,
0095 TTRHBuilder = "WithAngleAndTemplate",
0096 NavigationSchool = ''
0097 )
0098
0099 process.ShortTrackCandidates3.src = cms.InputTag("LongTracksRefit")
0100 process.ShortTrackCandidates4.src = cms.InputTag("LongTracksRefit")
0101 process.ShortTrackCandidates5.src = cms.InputTag("LongTracksRefit")
0102 process.ShortTrackCandidates6.src = cms.InputTag("LongTracksRefit")
0103 process.ShortTrackCandidates7.src = cms.InputTag("LongTracksRefit")
0104 process.ShortTrackCandidates8.src = cms.InputTag("LongTracksRefit")
0105
0106 process.SingleLongTrackProducer.requiredDr = cms.double(-9999.)
0107 process.SingleLongTrackProducer.matchMuons = cms.InputTag("muons")
0108 process.SingleLongTrackProducer.allTracks = cms.InputTag("ALCARECOTkAlDiMuon")
0109
0110
0111
0112
0113 from Alignment.OfflineValidation.shortenedTrackValidation_cfi import shortenedTrackValidation as _shortenedTrackValidation
0114 process.ShortenedTrackValidation = _shortenedTrackValidation.clone(folderName = "ShortTrackResolution",
0115 hitsRemainInput = ["3","4","5","6","7","8"],
0116 minTracksEtaInput = 0.0,
0117 maxTracksEtaInput = 2.2,
0118 minTracksPtInput = 15.0,
0119 maxTracksPtInput = 99999.9,
0120 maxDrInput = 0.01,
0121 tracksInputTag = "LongTracksRefit",
0122 tracksRerecoInputTag = ["RefittedShortTracks3",
0123 "RefittedShortTracks4",
0124 "RefittedShortTracks5",
0125 "RefittedShortTracks6",
0126 "RefittedShortTracks7",
0127 "RefittedShortTracks8"])
0128
0129
0130
0131
0132 process.analysis_step = cms.Path(process.offlineBeamSpot *
0133 process.MeasurementTrackerEvent *
0134 process.SingleLongTrackProducer *
0135 process.LongTracksRefit *
0136 process.ShortTrackCandidates3 *
0137 process.ShortTrackCandidates4 *
0138 process.ShortTrackCandidates5 *
0139 process.ShortTrackCandidates6 *
0140 process.ShortTrackCandidates7 *
0141 process.ShortTrackCandidates8 *
0142 process.RefittedShortTracks3 *
0143 process.RefittedShortTracks4 *
0144 process.RefittedShortTracks5 *
0145 process.RefittedShortTracks6 *
0146 process.RefittedShortTracks7 *
0147 process.RefittedShortTracks8 *
0148 process.ShortenedTrackValidation)
0149
0150
0151
0152
0153 process.options.numberOfThreads = 8