1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
|
import json
import os
import FWCore.ParameterSet.Config as cms
import FWCore.PythonUtilities.LumiList as LumiList
from Alignment.OfflineValidation.TkAlAllInOneTool.defaultInputFiles_cff import filesDefaultMC_NoPU
from FWCore.ParameterSet.VarParsing import VarParsing
from Alignment.OfflineValidation.TkAlAllInOneTool.utils import _byteify
###################################################################
# Define process
###################################################################
process = cms.Process("PrimaryVertexResolution")
###################################################################
# Argument parsing
###################################################################
options = VarParsing()
options.register("config", "", VarParsing.multiplicity.singleton, VarParsing.varType.string , "AllInOne config")
options.parseArguments()
###################################################################
# Read in AllInOne config in JSON format
###################################################################
if options.config == "":
config = {"validation": {},
"alignment": {}}
else:
with open(options.config, "r") as configFile:
config = json.load(configFile)
###################################################################
# Read filenames from given TXT file and define input source
###################################################################
readFiles = []
if "dataset" in config["validation"]:
with open(config["validation"]["dataset"], "r") as datafiles:
for fileName in datafiles.readlines():
readFiles.append(fileName.replace("\n", ""))
process.source = cms.Source("PoolSource",
fileNames = cms.untracked.vstring(readFiles),
skipEvents = cms.untracked.uint32(0)
)
else:
print(">>>>>>>>>> SplitV_cfg.py: msg%-i: config not specified! Loading default MC simulation -> filesDefaultMC_NoPU!")
process.source = cms.Source("PoolSource",
fileNames = filesDefaultMC_NoPU,
skipEvents = cms.untracked.uint32(0)
)
###################################################################
# Get good lumi section and load data or handle MC
###################################################################
if "goodlumi" in config["validation"]:
if os.path.isfile(config["validation"]["goodlumi"]):
goodLumiSecs = cms.untracked.VLuminosityBlockRange(LumiList.LumiList(filename = config["validation"]["goodlumi"]).getCMSSWString().split(','))
else:
print("Does not exist: {}. Continue without good lumi section file.")
goodLumiSecs = cms.untracked.VLuminosityBlockRange()
else:
goodLumiSecs = cms.untracked.VLuminosityBlockRange()
###################################################################
# Runs and events
###################################################################
runboundary = config["validation"].get("runboundary", 1)
isMultipleRuns=False
if(isinstance(runboundary, (list, tuple))):
isMultipleRuns=True
print("Multiple Runs are selected")
if(isMultipleRuns):
process.source.firstRun = cms.untracked.uint32(runboundary[0])
else:
process.source.firstRun = cms.untracked.uint32(runboundary)
###################################################################
# Default set to 1 for unit tests
###################################################################
process.maxEvents = cms.untracked.PSet(input = cms.untracked.int32(config["validation"].get("maxevents", 1)))
###################################################################
# Bookeeping
###################################################################
process.options = cms.untracked.PSet(
wantSummary = cms.untracked.bool(False),
Rethrow = cms.untracked.vstring("ProductNotFound"), # make this exception fatal
fileMode = cms.untracked.string('NOMERGE'), # no ordering needed, but calls endRun/beginRun etc. at file boundaries
)
###################################################################
# Messages
###################################################################
process.load("FWCore.MessageLogger.MessageLogger_cfi")
process.MessageLogger.cerr.FwkReport.reportEvery = 1000
process.MessageLogger.cout.enableStatistics = cms.untracked.bool(True)
###################################################################
# Basic modules
###################################################################
process.load("RecoVertex.BeamSpotProducer.BeamSpot_cff")
process.load("Configuration.Geometry.GeometryDB_cff")
process.load('Configuration.StandardSequences.Services_cff')
process.load("Configuration.StandardSequences.MagneticField_cff")
####################################################################
# Load and Configure Track refitter
####################################################################
process.load("RecoTracker.TrackProducer.TrackRefitters_cff")
process.TrackRefitter.src = config["validation"].get("trackcollection", "generalTracks")
process.TrackRefitter.TTRHBuilder = config["validation"].get("tthrbuilder", "WithAngleAndTemplate")
process.TrackRefitter.NavigationSchool = ""
####################################################################
# Global tag
####################################################################
process.load("Configuration.StandardSequences.FrontierConditions_GlobalTag_cff")
from Configuration.AlCa.GlobalTag import GlobalTag
process.GlobalTag = GlobalTag(process.GlobalTag, config["alignment"].get("globaltag", "105X_upgrade2018_realistic_v4"))
####################################################################
# Load conditions if wished
####################################################################
if "conditions" in config["alignment"]:
from CalibTracker.Configuration.Common.PoolDBESSource_cfi import poolDBESSource
for condition in config["alignment"]["conditions"]:
setattr(process, "conditionsIn{}".format(condition), poolDBESSource.clone(
connect = cms.string(str(config["alignment"]["conditions"][condition]["connect"])),
toGet = cms.VPSet(
cms.PSet(
record = cms.string(str(condition)),
tag = cms.string(str(config["alignment"]["conditions"][condition]["tag"]))
)
)
)
)
setattr(process, "prefer_conditionsIn{}".format(condition), cms.ESPrefer("PoolDBESSource", "conditionsIn{}".format(condition)))
###################################################################
# The trigger filter module
###################################################################
from HLTrigger.HLTfilters.triggerResultsFilter_cfi import *
process.theHLTFilter = triggerResultsFilter.clone(
triggerConditions = cms.vstring(config["validation"].get("triggerBits", "*")),
hltResults = cms.InputTag( "TriggerResults", "", "HLT" ),
l1tResults = cms.InputTag( "" ),
throw = cms.bool(False)
)
###################################################################
# PV refit
###################################################################
process.load("TrackingTools.TransientTrack.TransientTrackBuilder_cfi")
from RecoVertex.PrimaryVertexProducer.OfflinePrimaryVertices_cfi import offlinePrimaryVertices
process.offlinePrimaryVerticesFromRefittedTrks = offlinePrimaryVertices.clone()
process.offlinePrimaryVerticesFromRefittedTrks.TrackLabel = cms.InputTag("TrackRefitter")
process.offlinePrimaryVerticesFromRefittedTrks.vertexCollections.maxDistanceToBeam = 1
process.offlinePrimaryVerticesFromRefittedTrks.TkFilterParameters.maxNormalizedChi2 = 20
process.offlinePrimaryVerticesFromRefittedTrks.TkFilterParameters.minSiliconLayersWithHits = 5
process.offlinePrimaryVerticesFromRefittedTrks.TkFilterParameters.maxD0Significance = 5.0
# as it was prior to https://github.com/cms-sw/cmssw/commit/c8462ae4313b6be3bbce36e45373aa6e87253c59
process.offlinePrimaryVerticesFromRefittedTrks.TkFilterParameters.maxD0Error = 1.0
process.offlinePrimaryVerticesFromRefittedTrks.TkFilterParameters.maxDzError = 1.0
process.offlinePrimaryVerticesFromRefittedTrks.TkFilterParameters.minPixelLayersWithHits = 2
# Use compressions settings of TFile
# see https://root.cern.ch/root/html534/TFile.html#TFile:SetCompressionSettings
# settings = 100 * algorithm + level
# level is from 1 (small) to 9 (large compression)
# algo: 1 (ZLIB), 2 (LMZA)
# see more about compression & performance: https://root.cern.ch/root/html534/guides/users-guide/InputOutput.html#compression-and-performance
compressionSettings = 207
###################################################################
# The PV resolution module
###################################################################
from Alignment.OfflineValidation.splitVertexResolution_cfi import splitVertexResolution as _splitVertexResolution
process.PrimaryVertexResolution = _splitVertexResolution.clone(
compressionSettings = compressionSettings,
storeNtuple = False,
vtxCollection = "offlinePrimaryVerticesFromRefittedTrks",
trackCollection = "TrackRefitter",
minVertexNdf = 10.,
minVertexMeanWeight = 0.5,
runControl = config["validation"].get("runControl", False),
runControlNumber = [runboundary]
)
process.TFileService = cms.Service("TFileService",
fileName = cms.string("{}/SplitV.root".format(config.get("output", os.getcwd()))),
closeFileFast = cms.untracked.bool(True)
)
print("Saving the output at %s" % process.TFileService.fileName.value())
###################################################################
# Beamspot compatibility check
###################################################################
from RecoVertex.BeamSpotProducer.beamSpotCompatibilityChecker_cfi import beamSpotCompatibilityChecker
process.BeamSpotChecker = beamSpotCompatibilityChecker.clone(
bsFromEvent = "offlineBeamSpot::RECO", # source of the event beamspot (in the ALCARECO files)
bsFromDB = "offlineBeamSpot", # source of the DB beamspot (from Global Tag) NOTE: only if dbFromEvent is True!
warningThr = config["validation"].get("bsIncompatibleWarnThresh", 3), # significance threshold to emit a warning message
errorThr = config["validation"].get("bsIncompatibleErrThresh", 5), # significance threshold to abort the job
)
process.theValidSequence = cms.Sequence(process.offlineBeamSpot +
process.BeamSpotChecker +
process.TrackRefitter +
process.offlinePrimaryVerticesFromRefittedTrks +
process.PrimaryVertexResolution)
HLTSel = config["validation"].get("HLTselection", False)
if (HLTSel):
process.p = cms.Path(process.theHLTFilter + process.theValidSequence)
else:
process.p = cms.Path(process.theValidSequence)
print("Done")
|