File indexing completed on 2024-11-27 03:18:01
0001 import FWCore.ParameterSet.Config as cms
0002
0003 from PhysicsTools.PatAlgos.tools.ConfigToolBase import *
0004 import PhysicsTools.PatAlgos.tools.helpers as configtools
0005 from PhysicsTools.PatAlgos.tools.helpers import getPatAlgosToolsTask, addToProcessAndTask, addTaskToProcess
0006 from PhysicsTools.PatAlgos.tools.jetTools import switchJetCollection
0007 import CommonTools.CandAlgos.candPtrProjector_cfi as _mod
0008 from PhysicsTools.PatUtils.tools.pfforTrkMET_cff import *
0009 import JetMETCorrections.Type1MET.BadPFCandidateJetsEEnoiseProducer_cfi as _modbad
0010 import JetMETCorrections.Type1MET.UnclusteredBlobProducer_cfi as _modunc
0011
0012
0013 def isValidInputTag(input):
0014 input_str = input
0015 if isinstance(input, cms.InputTag):
0016 input_str = input.value()
0017 if input is None or input_str == '""':
0018 return False
0019 else:
0020 return True
0021
0022
0023 class RunMETCorrectionsAndUncertainties(ConfigToolBase):
0024
0025 _label='RunMETCorrectionsAndUncertainties'
0026 _defaultParameters=dicttypes.SortedKeysDict()
0027
0028 def __init__(self):
0029 ConfigToolBase.__init__(self)
0030
0031 self.addParameter(self._defaultParameters, 'metType', "PF",
0032 "Type of considered MET (only PF and Puppi supported so far)", Type=str, allowedValues = ["PF","Puppi"])
0033 self.addParameter(self._defaultParameters, 'correctionLevel', [""],
0034 "level of correction : available corrections for pfMet are T0, T1, T2, Txy and Smear)",
0035 allowedValues=["T0","T1","T2","Txy","Smear",""])
0036 self.addParameter(self._defaultParameters, 'computeUncertainties', True,
0037 "enable/disable the uncertainty computation", Type=bool)
0038 self.addParameter(self._defaultParameters, 'produceIntermediateCorrections', False,
0039 "enable/disable the production of all correction schemes (only for the most common)", Type=bool)
0040
0041
0042 self.addParameter(self._defaultParameters, 'electronCollection', cms.InputTag('selectedPatElectrons'),
0043 "Input electron collection", Type=cms.InputTag, acceptNoneValue=True)
0044 self.addParameter(self._defaultParameters, 'photonCollection', cms.InputTag('selectedPatPhotons'),
0045 "Input photon collection", Type=cms.InputTag, acceptNoneValue=True)
0046 self.addParameter(self._defaultParameters, 'muonCollection', cms.InputTag('selectedPatMuons'),
0047 "Input muon collection", Type=cms.InputTag, acceptNoneValue=True)
0048 self.addParameter(self._defaultParameters, 'tauCollection', cms.InputTag('selectedPatTaus'),
0049 "Input tau collection", Type=cms.InputTag, acceptNoneValue=True)
0050 self.addParameter(self._defaultParameters, 'jetCollectionUnskimmed', cms.InputTag('patJets'),
0051 "Input unskimmed jet collection for T1 MET computation", Type=cms.InputTag, acceptNoneValue=True)
0052
0053
0054 self.addParameter(self._defaultParameters, 'pfCandCollection', cms.InputTag('particleFlow'),
0055 "pf Candidate collection", Type=cms.InputTag, acceptNoneValue=True)
0056
0057
0058 self.addParameter(self._defaultParameters, 'autoJetCleaning', 'LepClean',
0059 "Enable the jet cleaning for the uncertainty computation: Full for tau/photons/jet cleaning, Partial for jet cleaning, LepClean for jet cleaning with muon and electrons only, None or Manual for no cleaning",
0060 allowedValues = ["Full","Partial","LepClean","None"])
0061 self.addParameter(self._defaultParameters, 'jetFlavor', 'AK4PFchs',
0062 "Use AK4PF/AK4PFchs for PFJets,AK4Calo for CaloJets", Type=str, allowedValues = ["AK4PF","AK4PFchs","AK4PFPuppi","CaloJets"])
0063 self.addParameter(self._defaultParameters, 'jetCorrectionType', 'L1L2L3-L1',
0064 "Use L1L2L3-L1 for the standard L1 removal / L1L2L3-RC for the random-cone correction", Type=str, allowedValues = ["L1L2L3-L1","L1L2L3-RC"])
0065
0066
0067 self.addParameter(self._defaultParameters, 'jetCorLabelUpToL3', "ak4PFCHSL1FastL2L3Corrector", "Use ak4PFL1FastL2L3Corrector (ak4PFCHSL1FastL2L3Corrector) for PFJets with (without) charged hadron subtraction, ak4CaloL1FastL2L3Corrector for CaloJets", Type=str)
0068 self.addParameter(self._defaultParameters, 'jetCorLabelL3Res', "ak4PFCHSL1FastL2L3ResidualCorrector", "Use ak4PFL1FastL2L3ResidualCorrector (ak4PFCHSL1FastL2L3ResidualCorrector) for PFJets with (without) charged hadron subtraction, ak4CaloL1FastL2L3ResidualCorrector for CaloJets", Type=str)
0069
0070 self.addParameter(self._defaultParameters, 'jecUncertaintyFile', '',
0071 "Extra JES uncertainty file", Type=str)
0072 self.addParameter(self._defaultParameters, 'jecUncertaintyTag', None,
0073 "JES uncertainty Tag", acceptNoneValue=True)
0074
0075
0076 self.addParameter(self._defaultParameters, 'manualJetConfig', False,
0077 "Enable jet configuration options", Type=bool)
0078 self.addParameter(self._defaultParameters, 'jetSelection', 'pt>15 && abs(eta)<9.9',
0079 "Advanced jet kinematic selection", Type=str)
0080
0081
0082 self.addParameter(self._defaultParameters, 'recoMetFromPFCs', False,
0083 "Recompute the MET from scratch using the pfCandidate collection", Type=bool)
0084 self.addParameter(self._defaultParameters, 'reapplyJEC', True,
0085 "Flag to enable/disable JEC update", Type=bool)
0086 self.addParameter(self._defaultParameters, 'reclusterJets', False,
0087 "Flag to enable/disable the jet reclustering", Type=bool)
0088 self.addParameter(self._defaultParameters, 'computeMETSignificance', True,
0089 "Flag to enable/disable the MET significance computation", Type=bool)
0090 self.addParameter(self._defaultParameters, 'CHS', False,
0091 "Flag to enable/disable the CHS jets", Type=bool)
0092
0093
0094 self.addParameter(self._defaultParameters, 'runOnData', False,
0095 "Switch for data/MC processing", Type=bool)
0096 self.addParameter(self._defaultParameters, 'onMiniAOD', False,
0097 "Switch on miniAOD configuration", Type=bool)
0098
0099
0100 self.addParameter(self._defaultParameters,'fixEE2017', False,
0101 "Exclude jets and PF candidates with EE noise characteristics (fix for 2017 run)", Type=bool)
0102 self.addParameter(self._defaultParameters,'fixEE2017Params', {'userawPt': True, 'ptThreshold': 50.0, 'minEtaThreshold': 2.65, 'maxEtaThreshold': 3.139},
0103 "Parameters dict for fixEE2017: userawPt, ptThreshold, minEtaThreshold, maxEtaThreshold", Type=dict)
0104 self.addParameter(self._defaultParameters, 'extractDeepMETs', False,
0105 "Extract DeepMETs from miniAOD, instead of recomputing them.", Type=bool)
0106
0107
0108 self.addParameter(self._defaultParameters, 'Puppi', False,
0109 "Puppi algorithm (private)", Type=bool)
0110 self.addParameter(self._defaultParameters, 'puppiProducerLabel', 'puppi',
0111 "PuppiProducer module for jet clustering label name", Type=str)
0112 self.addParameter(self._defaultParameters, 'puppiProducerForMETLabel', 'puppiNoLep',
0113 "PuppiProducer module for MET clustering label name", Type=str)
0114 self.addParameter(self._defaultParameters, 'addToPatDefaultSequence', False,
0115 "Flag to enable/disable that metUncertaintySequence is inserted into patDefaultSequence", Type=bool)
0116 self.addParameter(self._defaultParameters, 'postfix', '',
0117 "Technical parameter to identify the resulting sequences/tasks and its corresponding modules (allows multiple calls in a job)", Type=str)
0118
0119 self.addParameter(self._defaultParameters, 'campaign', '', 'Production campaign', Type=str)
0120 self.addParameter(self._defaultParameters, 'era', '', 'Era e.g. 2018, 2017B, ...', Type=str)
0121
0122
0123 self._parameters = copy.deepcopy(self._defaultParameters)
0124 self._comment = ""
0125
0126
0127 def getDefaultParameters(self):
0128 return self._defaultParameters
0129
0130
0131
0132
0133
0134 def __call__(self, process,
0135 metType =None,
0136 correctionLevel =None,
0137 computeUncertainties =None,
0138 produceIntermediateCorrections = None,
0139 electronCollection =None,
0140 photonCollection =None,
0141 muonCollection =None,
0142 tauCollection =None,
0143 jetCollectionUnskimmed =None,
0144 pfCandCollection =None,
0145 autoJetCleaning =None,
0146 jetFlavor =None,
0147 jetCorr =None,
0148 jetCorLabelUpToL3 =None,
0149 jetCorLabelL3Res =None,
0150 jecUncertaintyFile =None,
0151 jecUncertaintyTag =None,
0152 addToPatDefaultSequence =None,
0153 manualJetConfig =None,
0154 jetSelection =None,
0155 recoMetFromPFCs =None,
0156 reapplyJEC =None,
0157 reclusterJets =None,
0158 computeMETSignificance =None,
0159 CHS =None,
0160 puppiProducerLabel =None,
0161 puppiProducerForMETLabel = None,
0162 runOnData =None,
0163 onMiniAOD =None,
0164 fixEE2017 =None,
0165 fixEE2017Params =None,
0166 extractDeepMETs =None,
0167 campaign =None,
0168 era =None,
0169 postfix =None):
0170 electronCollection = self.initializeInputTag(electronCollection, 'electronCollection')
0171 photonCollection = self.initializeInputTag(photonCollection, 'photonCollection')
0172 muonCollection = self.initializeInputTag(muonCollection, 'muonCollection')
0173 tauCollection = self.initializeInputTag(tauCollection, 'tauCollection')
0174 jetCollectionUnskimmed = self.initializeInputTag(jetCollectionUnskimmed, 'jetCollectionUnskimmed')
0175 pfCandCollection = self.initializeInputTag(pfCandCollection, 'pfCandCollection')
0176 if metType is None :
0177 metType = self._defaultParameters['metType'].value
0178 if correctionLevel is None :
0179 correctionLevel = self._defaultParameters['correctionLevel'].value
0180 if computeUncertainties is None :
0181 computeUncertainties = self._defaultParameters['computeUncertainties'].value
0182 if produceIntermediateCorrections is None :
0183 produceIntermediateCorrections = self._defaultParameters['produceIntermediateCorrections'].value
0184 if electronCollection is None :
0185 electronCollection = self._defaultParameters['electronCollection'].value
0186 if photonCollection is None :
0187 photonCollection = self._defaultParameters['photonCollection'].value
0188 if muonCollection is None :
0189 muonCollection = self._defaultParameters['muonCollection'].value
0190 if tauCollection is None :
0191 tauCollection = self._defaultParameters['tauCollection'].value
0192 if jetCollectionUnskimmed is None :
0193 jetCollectionUnskimmed = self._defaultParameters['jetCollectionUnskimmed'].value
0194 if pfCandCollection is None :
0195 pfCandCollection = self._defaultParameters['pfCandCollection'].value
0196 if autoJetCleaning is None :
0197 autoJetCleaning = self._defaultParameters['autoJetCleaning'].value
0198 if jetFlavor is None :
0199 jetFlavor = self._defaultParameters['jetFlavor'].value
0200 if jetCorr is None :
0201 jetCorr = self._defaultParameters['jetCorrectionType'].value
0202 if jetCorLabelUpToL3 is None:
0203 jetCorLabelUpToL3 = self._defaultParameters['jetCorLabelUpToL3'].value
0204 if jetCorLabelL3Res is None:
0205 jetCorLabelL3Res = self._defaultParameters['jetCorLabelL3Res'].value
0206 if jecUncertaintyFile is None:
0207 jecUncertaintyFile = self._defaultParameters['jecUncertaintyFile'].value
0208 if jecUncertaintyTag is None:
0209 jecUncertaintyTag = self._defaultParameters['jecUncertaintyTag'].value
0210 if addToPatDefaultSequence is None :
0211 addToPatDefaultSequence = self._defaultParameters['addToPatDefaultSequence'].value
0212 if manualJetConfig is None :
0213 manualJetConfig = self._defaultParameters['manualJetConfig'].value
0214 if jetSelection is None :
0215 jetSelection = self._defaultParameters['jetSelection'].value
0216 recoMetFromPFCsIsNone = (recoMetFromPFCs is None)
0217 if recoMetFromPFCs is None :
0218 recoMetFromPFCs = self._defaultParameters['recoMetFromPFCs'].value
0219 if reapplyJEC is None :
0220 reapplyJEC = self._defaultParameters['reapplyJEC'].value
0221 reclusterJetsIsNone = (reclusterJets is None)
0222 if reclusterJets is None :
0223 reclusterJets = self._defaultParameters['reclusterJets'].value
0224 if computeMETSignificance is None :
0225 computeMETSignificance = self._defaultParameters['computeMETSignificance'].value
0226 if CHS is None :
0227 CHS = self._defaultParameters['CHS'].value
0228 if puppiProducerLabel is None:
0229 puppiProducerLabel = self._defaultParameters['puppiProducerLabel'].value
0230 if puppiProducerForMETLabel is None:
0231 puppiProducerForMETLabel = self._defaultParameters['puppiProducerForMETLabel'].value
0232 if runOnData is None :
0233 runOnData = self._defaultParameters['runOnData'].value
0234 if onMiniAOD is None :
0235 onMiniAOD = self._defaultParameters['onMiniAOD'].value
0236 if postfix is None :
0237 postfix = self._defaultParameters['postfix'].value
0238 if fixEE2017 is None :
0239 fixEE2017 = self._defaultParameters['fixEE2017'].value
0240 if fixEE2017Params is None :
0241 fixEE2017Params = self._defaultParameters['fixEE2017Params'].value
0242 if extractDeepMETs is None :
0243 extractDeepMETs = self._defaultParameters['extractDeepMETs'].value
0244 if campaign is None :
0245 campaign = self._defaultParameters['campaign'].value
0246 if era is None :
0247 era = self._defaultParameters['era'].value
0248
0249 self.setParameter('metType',metType),
0250 self.setParameter('correctionLevel',correctionLevel),
0251 self.setParameter('computeUncertainties',computeUncertainties),
0252 self.setParameter('produceIntermediateCorrections',produceIntermediateCorrections),
0253 self.setParameter('electronCollection',electronCollection),
0254 self.setParameter('photonCollection',photonCollection),
0255 self.setParameter('muonCollection',muonCollection),
0256 self.setParameter('tauCollection',tauCollection),
0257 self.setParameter('jetCollectionUnskimmed',jetCollectionUnskimmed),
0258 self.setParameter('pfCandCollection',pfCandCollection),
0259
0260 self.setParameter('autoJetCleaning',autoJetCleaning),
0261 self.setParameter('jetFlavor',jetFlavor),
0262
0263
0264 self.setParameter('jecUncertaintyFile',jecUncertaintyFile),
0265 self.setParameter('jecUncertaintyTag',jecUncertaintyTag),
0266
0267 self.setParameter('addToPatDefaultSequence',addToPatDefaultSequence),
0268 self.setParameter('jetSelection',jetSelection),
0269 self.setParameter('recoMetFromPFCs',recoMetFromPFCs),
0270 self.setParameter('reclusterJets',reclusterJets),
0271 self.setParameter('computeMETSignificance',computeMETSignificance),
0272 self.setParameter('reapplyJEC',reapplyJEC),
0273 self.setParameter('CHS',CHS),
0274 self.setParameter('puppiProducerLabel',puppiProducerLabel),
0275 self.setParameter('puppiProducerForMETLabel',puppiProducerForMETLabel),
0276 self.setParameter('runOnData',runOnData),
0277 self.setParameter('onMiniAOD',onMiniAOD),
0278 self.setParameter('postfix',postfix),
0279 self.setParameter('fixEE2017',fixEE2017),
0280 self.setParameter('fixEE2017Params',fixEE2017Params),
0281 self.setParameter('extractDeepMETs',extractDeepMETs),
0282 self.setParameter('campaign',campaign),
0283 self.setParameter('era',era),
0284
0285
0286 if metType == "Puppi":
0287 self.setParameter('CHS',False),
0288
0289
0290
0291 self.setParameter('Puppi',self._defaultParameters['Puppi'].value)
0292 if metType == "Puppi":
0293 self.setParameter('metType',"PF")
0294 self.setParameter('Puppi',True)
0295
0296
0297 if manualJetConfig:
0298 self.setParameter('CHS',CHS)
0299 self.setParameter('jetCorLabelUpToL3',jetCorLabelUpToL3)
0300 self.setParameter('jetCorLabelL3Res',jetCorLabelL3Res)
0301 self.setParameter('reclusterJets',reclusterJets)
0302 else:
0303
0304 self.jetConfiguration()
0305
0306
0307
0308 if fixEE2017:
0309 if recoMetFromPFCsIsNone: self.setParameter('recoMetFromPFCs',True)
0310 if reclusterJetsIsNone: self.setParameter('reclusterJets',False)
0311
0312
0313 if recoMetFromPFCs and reclusterJetsIsNone and not fixEE2017:
0314 self.setParameter('reclusterJets',True)
0315
0316 self.apply(process)
0317
0318
0319 def toolCode(self, process):
0320
0321
0322
0323
0324
0325 metType = self._parameters['metType'].value
0326 correctionLevel = self._parameters['correctionLevel'].value
0327 computeUncertainties = self._parameters['computeUncertainties'].value
0328 produceIntermediateCorrections = self._parameters['produceIntermediateCorrections'].value
0329
0330 electronCollection = self._parameters['electronCollection'].value
0331 photonCollection = self._parameters['photonCollection'].value
0332 muonCollection = self._parameters['muonCollection'].value
0333 tauCollection = self._parameters['tauCollection'].value
0334 jetCollectionUnskimmed = self._parameters['jetCollectionUnskimmed'].value
0335 pfCandCollection = self._parameters['pfCandCollection'].value
0336
0337 jetSelection = self._parameters['jetSelection'].value
0338 autoJetCleaning = self._parameters['autoJetCleaning'].value
0339 jetFlavor = self._parameters['jetFlavor'].value
0340 jetCorLabelUpToL3 = self._parameters['jetCorLabelUpToL3'].value
0341 jetCorLabelL3Res = self._parameters['jetCorLabelL3Res'].value
0342 jecUncertaintyFile = self._parameters['jecUncertaintyFile'].value
0343 jecUncertaintyTag = self._parameters['jecUncertaintyTag'].value
0344
0345 recoMetFromPFCs = self._parameters['recoMetFromPFCs'].value
0346 reapplyJEC = self._parameters['reapplyJEC'].value
0347 reclusterJets = self._parameters['reclusterJets'].value
0348 computeMETSignificance = self._parameters['computeMETSignificance'].value
0349 extractDeepMETs = self._parameters['extractDeepMETs'].value
0350
0351 fixEE2017 = self._parameters['fixEE2017'].value
0352 fixEE2017Params = self._parameters['fixEE2017Params'].value
0353 campaign = self._parameters['campaign'].value
0354 era = self._parameters['era'].value
0355
0356 onMiniAOD = self._parameters['onMiniAOD'].value
0357 addToPatDefaultSequence = self._parameters['addToPatDefaultSequence'].value
0358 postfix = self._parameters['postfix'].value
0359
0360
0361 jetUncInfos = {
0362 "jCorrPayload":jetFlavor,
0363 "jCorLabelUpToL3":jetCorLabelUpToL3,
0364 "jCorLabelL3Res":jetCorLabelL3Res,
0365 "jecUncFile":jecUncertaintyFile,
0366 "jecUncTag":"Uncertainty"
0367 }
0368
0369
0370 if (jecUncertaintyFile!="" and jecUncertaintyTag==None):
0371 jetUncInfos[ "jecUncTag" ] = ""
0372
0373 elif (jecUncertaintyTag!=None):
0374 jetUncInfos[ "jecUncTag" ] = jecUncertaintyTag
0375
0376
0377
0378
0379
0380
0381 patMetModuleTask = cms.Task()
0382
0383
0384 if fixEE2017:
0385 pfCandCollection, jetCollectionUnskimmed = self.runFixEE2017(process,
0386 fixEE2017Params,
0387 jetCollectionUnskimmed,
0388 pfCandCollection,
0389 [electronCollection,muonCollection,tauCollection,photonCollection],
0390 patMetModuleTask,
0391 postfix,
0392 )
0393
0394
0395 if recoMetFromPFCs:
0396 self.recomputeRawMetFromPfcs(process,
0397 pfCandCollection,
0398 onMiniAOD,
0399 patMetModuleTask,
0400 postfix)
0401
0402 elif onMiniAOD:
0403 self.extractMET(process, "raw", patMetModuleTask, postfix)
0404
0405
0406 if reclusterJets:
0407 jetCollectionUnskimmed = self.ak4JetReclustering(process, pfCandCollection,
0408 patMetModuleTask, postfix)
0409
0410
0411 if onMiniAOD:
0412 if not reclusterJets and reapplyJEC:
0413 jetCollectionUnskimmed = self.updateJECs(process, jetCollectionUnskimmed, patMetModuleTask, postfix)
0414
0415
0416
0417
0418 jetCollection = self.getJetCollectionForCorsAndUncs(process,
0419 jetCollectionUnskimmed,
0420 jetSelection,
0421 autoJetCleaning,
0422 patMetModuleTask,
0423 postfix)
0424
0425 if onMiniAOD:
0426
0427 self.miniAODConfigurationPre(process, patMetModuleTask, pfCandCollection, postfix)
0428 else:
0429 from PhysicsTools.PatUtils.pfeGammaToCandidate_cfi import pfeGammaToCandidate
0430 addToProcessAndTask("pfeGammaToCandidate", pfeGammaToCandidate.clone(
0431 electrons = copy.copy(electronCollection),
0432 photons = copy.copy(photonCollection)),
0433 process, patMetModuleTask)
0434 if hasattr(process,"patElectrons") and process.patElectrons.electronSource == cms.InputTag("reducedEgamma","reducedGedGsfElectrons"):
0435 process.pfeGammaToCandidate.electron2pf = "reducedEgamma:reducedGsfElectronPfCandMap"
0436 if hasattr(process,"patPhotons") and process.patPhotons.photonSource == cms.InputTag("reducedEgamma","reducedGedPhotons"):
0437 process.pfeGammaToCandidate.photon2pf = "reducedEgamma:reducedPhotonPfCandMap"
0438
0439
0440 self.produceMET(process, metType, patMetModuleTask, postfix)
0441
0442
0443
0444
0445 if onMiniAOD:
0446 self.miniAODConfiguration(process,
0447 pfCandCollection,
0448 jetCollection,
0449 patMetModuleTask,
0450 postfix
0451 )
0452
0453
0454
0455
0456
0457 patMetCorrectionTask = cms.Task()
0458 metModName = self.getCorrectedMET(process, metType, correctionLevel,
0459 produceIntermediateCorrections,
0460 jetCollection,
0461 patMetCorrectionTask, postfix )
0462
0463
0464
0465 if "T1" in metModName:
0466 getattr(process,"patPFMetT1T2Corr"+postfix).src = jetCollection
0467 getattr(process,"patPFMetT2Corr"+postfix).src = jetCollection
0468
0469 if self._parameters["Puppi"].value:
0470 getattr(process,"patPFMetT1T2Corr"+postfix).offsetCorrLabel = cms.InputTag("")
0471 getattr(process,"patPFMetT2Corr"+postfix).offsetCorrLabel = cms.InputTag("")
0472 if "Smear" in metModName:
0473 getattr(process,"patSmearedJets"+postfix).src = jetCollection
0474 if self._parameters["Puppi"].value:
0475 getattr(process,"patPFMetT1T2SmearCorr"+postfix).offsetCorrLabel = cms.InputTag("")
0476
0477
0478
0479
0480
0481
0482 patMetUncertaintyTask = cms.Task()
0483 if not hasattr(process, "patMetUncertaintyTask"+postfix):
0484 if self._parameters["Puppi"].value:
0485 patMetUncertaintyTask.add(cms.Task(getattr(process, "ak4PFPuppiL1FastL2L3CorrectorTask"),getattr(process, "ak4PFPuppiL1FastL2L3ResidualCorrectorTask")))
0486 else:
0487 patMetUncertaintyTask.add(cms.Task(getattr(process, "ak4PFCHSL1FastL2L3CorrectorTask"),getattr(process, "ak4PFCHSL1FastL2L3ResidualCorrectorTask")))
0488 patShiftedModuleTask = cms.Task()
0489 if computeUncertainties:
0490 self.getMETUncertainties(process, metType, metModName,
0491 electronCollection,
0492 photonCollection,
0493 muonCollection,
0494 tauCollection,
0495 pfCandCollection,
0496 jetCollection,
0497 jetUncInfos,
0498 patMetUncertaintyTask,
0499 postfix)
0500
0501
0502
0503
0504
0505
0506 addTaskToProcess(process, "patMetCorrectionTask"+postfix, patMetCorrectionTask)
0507 addTaskToProcess(process, "patMetUncertaintyTask"+postfix, patMetUncertaintyTask)
0508 addTaskToProcess(process, "patShiftedModuleTask"+postfix, patShiftedModuleTask)
0509 addTaskToProcess(process, "patMetModuleTask"+postfix, patMetModuleTask)
0510
0511
0512 fullPatMetTask = cms.Task()
0513 fullPatMetTask.add(getattr(process, "patMetModuleTask"+postfix))
0514 fullPatMetTask.add(getattr(process, "patMetCorrectionTask"+postfix))
0515 fullPatMetTask.add(getattr(process, "patMetUncertaintyTask"+postfix))
0516 fullPatMetTask.add(getattr(process, "patShiftedModuleTask"+postfix))
0517
0518
0519 if hasattr(process, "patCaloMet"):
0520 fullPatMetTask.add(getattr(process, "patCaloMet"))
0521
0522 if hasattr(process, "deepMETsResolutionTune"):
0523 fullPatMetTask.add(getattr(process, "deepMETsResolutionTune"))
0524 if hasattr(process, "deepMETsResponseTune"):
0525 fullPatMetTask.add(getattr(process, "deepMETsResponseTune"))
0526
0527 if hasattr(process, "slimmedMETs"+postfix):
0528 fullPatMetTask.add(getattr(process, "slimmedMETs"+postfix))
0529
0530
0531 addTaskToProcess(process, "fullPatMetTask"+postfix, fullPatMetTask)
0532
0533
0534 task = getPatAlgosToolsTask(process)
0535 task.add(getattr(process,"fullPatMetTask"+postfix))
0536
0537
0538
0539
0540
0541 self.miniAODConfigurationPost(process, postfix)
0542
0543
0544 if addToPatDefaultSequence:
0545 if not hasattr(process, "patDefaultSequence"):
0546 raise ValueError("PAT default sequence is not defined !!")
0547 process.patDefaultSequence += getattr(process, "fullPatMetSequence"+postfix)
0548
0549
0550 def produceMET(self, process, metType, patMetModuleTask, postfix):
0551
0552 produceMET_task, produceMET_label = cms.Task(), "produceMET_task{}".format(postfix)
0553
0554
0555 if metType == "PF" and not hasattr(process, 'pat'+metType+'Met'):
0556 process.load("PhysicsTools.PatUtils.patPFMETCorrections_cff")
0557 produceMET_task.add(process.producePatPFMETCorrectionsTask)
0558 produceMET_task.add(process.patPFMetT2SmearCorrTask)
0559 produceMET_task.add(process.patPFMetTxyCorrTask)
0560 produceMET_task.add(process.jetCorrectorsTask)
0561
0562
0563 _myPatMet = 'pat'+metType+'Met'+postfix
0564
0565 if postfix != "" and metType == "PF" and not hasattr(process, _myPatMet):
0566 noClonesTmp = [ "particleFlowDisplacedVertex", "pfCandidateToVertexAssociation" ]
0567
0568
0569 configtools.cloneProcessingSnippetTask(process, getattr(process,"producePatPFMETCorrectionsTask"), postfix, noClones = noClonesTmp)
0570 produceMET_task.add(getattr(process,"producePatPFMETCorrectionsTask"+postfix))
0571
0572 addToProcessAndTask(_myPatMet, getattr(process,'patPFMet').clone(), process, produceMET_task)
0573
0574 getattr(process, _myPatMet).metSource = cms.InputTag("pfMet"+postfix)
0575 getattr(process, _myPatMet).srcPFCands = copy.copy(self.getvalue("pfCandCollection"))
0576
0577 if self.getvalue("Puppi"):
0578 getattr(process, _myPatMet).srcWeights = self._parameters['puppiProducerForMETLabel'].value
0579
0580 if metType == "PF":
0581 getattr(process, _myPatMet).srcLeptons = \
0582 cms.VInputTag(copy.copy(self.getvalue("electronCollection")) if self.getvalue("onMiniAOD") else
0583 cms.InputTag("pfeGammaToCandidate","electrons"),
0584 copy.copy(self.getvalue("muonCollection")),
0585 copy.copy(self.getvalue("photonCollection")) if self.getvalue("onMiniAOD") else
0586 cms.InputTag("pfeGammaToCandidate","photons"))
0587
0588 if self.getvalue("runOnData"):
0589 getattr(process, _myPatMet).addGenMET = False
0590
0591
0592 produceMET_task.add(getattr(process, _myPatMet ))
0593
0594
0595 addTaskToProcess(process, produceMET_label, produceMET_task)
0596
0597
0598 patMetModuleTask.add(getattr(process, produceMET_label))
0599
0600
0601 def getCorrectedMET(self, process, metType, correctionLevel, produceIntermediateCorrections,
0602 jetCollection, metModuleTask, postfix):
0603
0604
0605 getCorrectedMET_task, getCorrectedMET_label = cms.Task(), "getCorrectedMET_task{}".format(postfix)
0606
0607 metModName = "pat"+metType+"Met"+postfix
0608
0609
0610
0611 corTypeNames = {
0612 "T0":"T0pc",
0613 "T1":"T1",
0614 "T2":"T2",
0615 "Txy":"Txy",
0616 "Smear":"Smear",
0617 }
0618
0619
0620
0621 for corType in correctionLevel:
0622 if corType not in corTypeNames.keys():
0623 if corType != "":
0624 raise ValueError(corType+" is not a proper MET correction name! Aborting the MET correction production")
0625 else:
0626 return metModName
0627
0628
0629 corTypeTaskNames = {
0630 "T0": "patPFMetT0CorrTask"+postfix,
0631 "T1": "patPFMetT1T2CorrTask"+postfix,
0632 "T2": "patPFMetT2CorrTask"+postfix,
0633 "Txy": "patPFMetTxyCorrTask"+postfix,
0634 "Smear": "patPFMetSmearCorrTask"+postfix,
0635 "T2Smear": "patPFMetT2SmearCorrTask"+postfix
0636 }
0637
0638
0639
0640 if postfix != "":
0641 noClonesTmp = [ "particleFlowDisplacedVertex", "pfCandidateToVertexAssociation" ]
0642 if not hasattr(process, "patPFMetT0CorrTask"+postfix):
0643 configtools.cloneProcessingSnippetTask(process, getattr(process,"patPFMetT0CorrTask"), postfix, noClones = noClonesTmp)
0644 getCorrectedMET_task.add(getattr(process,"patPFMetT0CorrTask"+postfix))
0645 if not hasattr(process, "patPFMetT1T2CorrTask"+postfix):
0646 configtools.cloneProcessingSnippetTask(process, getattr(process,"patPFMetT1T2CorrTask"), postfix)
0647 getCorrectedMET_task.add(getattr(process,"patPFMetT1T2CorrTask"+postfix))
0648 if not hasattr(process, "patPFMetT2CorrTask"+postfix):
0649 configtools.cloneProcessingSnippetTask(process, getattr(process,"patPFMetT2CorrTask"), postfix)
0650 getCorrectedMET_task.add(getattr(process,"patPFMetT2CorrTask"+postfix))
0651 if not hasattr(process, "patPFMetTxyCorrTask"+postfix):
0652 configtools.cloneProcessingSnippetTask(process, getattr(process,"patPFMetTxyCorrTask"), postfix)
0653 getCorrectedMET_task.add(getattr(process,"patPFMetTxyCorrTask"+postfix))
0654 if not hasattr(process, "patPFMetSmearCorrTask"+postfix):
0655 configtools.cloneProcessingSnippetTask(process, getattr(process,"patPFMetSmearCorrTask"), postfix)
0656 getCorrectedMET_task.add(getattr(process,"patPFMetSmearCorrTask"+postfix))
0657 if not hasattr(process, "patPFMetT2SmearCorrTask"+postfix):
0658 configtools.cloneProcessingSnippetTask(process, getattr(process,"patPFMetT2SmearCorrTask"), postfix)
0659 getCorrectedMET_task.add(getattr(process,"patPFMetT2SmearCorrTask"+postfix))
0660
0661
0662 corTypeTasks = {}
0663 for corType in corTypeTaskNames.keys():
0664 corTypeTasks[corType] = getattr(process, corTypeTaskNames[corType] )
0665
0666
0667 corTypeTags = {
0668 "T0":['patPFMetT0Corr'+postfix,''],
0669 "T1":['patPFMetT1T2Corr'+postfix, 'type1'],
0670 "T2":['patPFMetT2Corr'+postfix, 'type2'],
0671 "Txy": ['patPFMetTxyCorr'+postfix,''],
0672 "Smear":['patPFMetT1T2SmearCorr'+postfix, 'type1'],
0673 "T2Smear":['patPFMetT2SmearCorr'+postfix, 'type2']
0674 }
0675
0676
0677 correctionScheme=""
0678 correctionProducts = []
0679 correctionTasks = []
0680 for corType in correctionLevel:
0681 correctionScheme += corTypeNames[corType]
0682 correctionProducts.append(cms.InputTag(corTypeTags[corType][0],corTypeTags[corType][1]))
0683 correctionTasks.append(corTypeTasks[corType])
0684
0685
0686 if "T2" in correctionLevel and "Smear" in correctionLevel:
0687 correctionProducts.append(cms.InputTag(corTypeTags["T2Smear"][0],corTypeTags["T2Smear"][1]))
0688 correctionTasks.append(corTypeTasks["T2Smear"])
0689
0690
0691 if "T1" in correctionLevel and "Smear" in correctionLevel:
0692 correctionProducts.remove(cms.InputTag(corTypeTags["T1"][0],corTypeTags["T1"][1]))
0693
0694
0695 if "Txy" in correctionLevel:
0696 datamc = "DATA" if self.getvalue("runOnData") else "MC"
0697 self.tuneTxyParameters(process, correctionScheme, postfix, datamc, self.getvalue("campaign"), self.getvalue("era"))
0698 getattr(process, "patPFMetTxyCorr"+postfix).srcPFlow = self._parameters["pfCandCollection"].value
0699 if self.getvalue("Puppi"):
0700 getattr(process, "patPFMetTxyCorr"+postfix).srcWeights = self._parameters['puppiProducerForMETLabel'].value
0701
0702
0703
0704 if "T1" in correctionLevel:
0705 _myPatMet = "pat"+metType+"Met"+postfix
0706 getattr(process, _myPatMet).computeMETSignificance = cms.bool(self.getvalue("computeMETSignificance"))
0707 getattr(process, _myPatMet).srcPFCands = copy.copy(self.getvalue("pfCandCollection"))
0708 getattr(process, _myPatMet).srcLeptons = \
0709 cms.VInputTag(copy.copy(self.getvalue("electronCollection")) if self.getvalue("onMiniAOD") else
0710 cms.InputTag("pfeGammaToCandidate","electrons"),
0711 copy.copy(self.getvalue("muonCollection")),
0712 copy.copy(self.getvalue("photonCollection")) if self.getvalue("onMiniAOD") else
0713 cms.InputTag("pfeGammaToCandidate","photons"))
0714 if postfix=="NoHF":
0715 getattr(process, _myPatMet).computeMETSignificance = cms.bool(False)
0716 if self.getvalue("runOnData"):
0717 from RecoMET.METProducers.METSignificanceParams_cfi import METSignificanceParams_Data
0718 getattr(process, _myPatMet).parameters = METSignificanceParams_Data
0719 if self.getvalue("Puppi"):
0720 getattr(process, _myPatMet).srcWeights = self._parameters['puppiProducerForMETLabel'].value
0721 getattr(process, _myPatMet).srcJets = cms.InputTag('cleanedPatJets'+postfix)
0722 getattr(process, _myPatMet).srcJetSF = 'AK4PFPuppi'
0723 getattr(process, _myPatMet).srcJetResPt = 'AK4PFPuppi_pt'
0724 getattr(process, _myPatMet).srcJetResPhi = 'AK4PFPuppi_phi'
0725
0726
0727 if not self._parameters["onMiniAOD"].value and not postfix=="NoHF":
0728 _myPatMet = "patMETs"+postfix
0729 getattr(process, _myPatMet).computeMETSignificance = cms.bool(self.getvalue("computeMETSignificance"))
0730 getattr(process, _myPatMet).srcPFCands=copy.copy(self.getvalue("pfCandCollection"))
0731 getattr(process, _myPatMet).srcLeptons = \
0732 cms.VInputTag(copy.copy(self.getvalue("electronCollection")) if self.getvalue("onMiniAOD") else
0733 cms.InputTag("pfeGammaToCandidate","electrons"),
0734 copy.copy(self.getvalue("muonCollection")),
0735 copy.copy(self.getvalue("photonCollection")) if self.getvalue("onMiniAOD") else
0736 cms.InputTag("pfeGammaToCandidate","photons"))
0737 if self.getvalue("Puppi"):
0738 getattr(process, _myPatMet).srcWeights = self._parameters['puppiProducerForMETLabel'].value
0739
0740 if hasattr(process, "patCaloMet"):
0741 getattr(process, "patCaloMet").computeMETSignificance = cms.bool(False)
0742
0743
0744 if "T1" in correctionLevel and not self._parameters["CHS"].value and not self._parameters["Puppi"].value:
0745 addToProcessAndTask("corrPfMetType1"+postfix, getattr(process, "corrPfMetType1" ).clone(), process, getCorrectedMET_task)
0746 getattr(process, "corrPfMetType1"+postfix).src = cms.InputTag("ak4PFJets"+postfix)
0747 getattr(process, "corrPfMetType1"+postfix).jetCorrLabel = cms.InputTag("ak4PFL1FastL2L3Corrector")
0748 getattr(process, "corrPfMetType1"+postfix).jetCorrLabelRes = cms.InputTag("ak4PFL1FastL2L3ResidualCorrector")
0749 getattr(process, "corrPfMetType1"+postfix).offsetCorrLabel = cms.InputTag("ak4PFL1FastjetCorrector")
0750 getattr(process, "basicJetsForMet"+postfix).offsetCorrLabel = cms.InputTag("ak4PFL1FastjetCorrector")
0751
0752 if "T1" in correctionLevel and self._parameters["Puppi"].value:
0753 addToProcessAndTask("corrPfMetType1"+postfix, getattr(process, "corrPfMetType1" ).clone(), process, getCorrectedMET_task)
0754 getattr(process, "corrPfMetType1"+postfix).src = cms.InputTag("ak4PFJets"+postfix)
0755 getattr(process, "corrPfMetType1"+postfix).jetCorrLabel = cms.InputTag("ak4PFPuppiL1FastL2L3Corrector")
0756 getattr(process, "corrPfMetType1"+postfix).jetCorrLabelRes = cms.InputTag("ak4PFPuppiL1FastL2L3ResidualCorrector")
0757 getattr(process, "corrPfMetType1"+postfix).offsetCorrLabel = cms.InputTag("ak4PFPuppiL1FastjetCorrector")
0758 getattr(process, "basicJetsForMet"+postfix).offsetCorrLabel = cms.InputTag("L1FastJet")
0759
0760 if "T1" in correctionLevel and self._parameters["CHS"].value and self._parameters["reclusterJets"].value:
0761 getattr(process, "corrPfMetType1"+postfix).src = cms.InputTag("ak4PFJetsCHS"+postfix)
0762
0763
0764 metModName = "pat"+metType+"Met"+correctionScheme+postfix
0765
0766 taskName=""
0767 corMetProducer=None
0768
0769
0770 if metType == "PF":
0771 corMetProducer = cms.EDProducer("CorrectedPATMETProducer",
0772 src = cms.InputTag('pat'+metType+'Met' + postfix),
0773 srcCorrections = cms.VInputTag(correctionProducts)
0774 )
0775 taskName="getCorrectedMET_task"
0776
0777 addToProcessAndTask(metModName, corMetProducer, process, getCorrectedMET_task)
0778
0779
0780 if not hasattr(process, getCorrectedMET_label):
0781 for corTask in correctionTasks:
0782 getCorrectedMET_task.add(corTask)
0783
0784 else:
0785 for corType in corTypeTaskNames.keys():
0786 if not configtools.contains(getCorrectedMET_task, corTypeTags[corType][0]) and corType in correctionLevel:
0787 getCorrectedMET_task.add(corTypeTasks[corType])
0788
0789
0790 getCorrectedMET_task.add(getattr(process, metModName))
0791
0792 addTaskToProcess(process, getCorrectedMET_label, getCorrectedMET_task)
0793 metModuleTask.add(getattr(process, getCorrectedMET_label))
0794
0795
0796
0797 if produceIntermediateCorrections:
0798 interMets = self.addIntermediateMETs(process, metType, correctionLevel, correctionScheme, corTypeTags,corTypeNames, postfix)
0799 for met in interMets.keys():
0800 addToProcessAndTask(met, interMets[met], process, getCorrectedMET_task)
0801
0802 return metModName
0803
0804
0805
0806 def addIntermediateMETs(self, process, metType, correctionLevel, corScheme, corTags, corNames, postfix):
0807 interMets = {}
0808
0809
0810 if len(correctionLevel) == 1:
0811 return interMets
0812
0813
0814 nCor=len(correctionLevel)+1
0815 ids = [0]*nCor
0816 for i in range(nCor**nCor):
0817 tmp=i
0818 exists=False
0819 corName=""
0820 corrections = []
0821 for j in range(nCor):
0822 ids[j] = tmp%nCor
0823 tmp = tmp//nCor
0824
0825 if j != 0 and ids[j-1] < ids[j]:
0826 exists=True
0827 for k in range(0,j):
0828 if ids[k] == ids[j] and ids[k]!=0:
0829 exists=True
0830
0831 if exists or sum(ids[j] for j in range(nCor))==0:
0832 continue
0833
0834 for cor in range(nCor):
0835 cid = ids[nCor-cor-1]
0836 cKey = correctionLevel[cid-1]
0837 if cid ==0:
0838 continue
0839 else :
0840 corName += corNames[cKey]
0841 corrections.append( cms.InputTag(corTags[ cKey ][0], corTags[ cKey ][1]) )
0842
0843 if corName == corScheme:
0844 continue
0845
0846 corName='pat'+metType+'Met' + corName + postfix
0847 if configtools.contains(getattr(process,"getCorrectedMET_task"+postfix), corName ) and hasattr(process, corName):
0848 continue
0849
0850 interMets[corName] = cms.EDProducer("CorrectedPATMETProducer",
0851 src = cms.InputTag('pat'+metType+'Met' + postfix),
0852 srcCorrections = cms.VInputTag(corrections)
0853 )
0854
0855
0856 return interMets
0857
0858
0859
0860 def getMETUncertainties(self, process, metType, metModName, electronCollection,
0861 photonCollection, muonCollection, tauCollection,
0862 pfCandCollection, jetCollection, jetUncInfos,
0863 patMetUncertaintyTask,
0864 postfix):
0865
0866 getMETUncertainties_task, getMETUncertainties_label = cms.Task(), "getMETUncertainties_task{}".format(postfix)
0867
0868
0869
0870
0871 if not isValidInputTag(jetCollection):
0872 print("INFO : jet collection %s does not exists, no energy resolution shifting will be performed in MET uncertainty tools" % jetCollection)
0873 else:
0874 preId=""
0875 if "Smear" in metModName:
0876 preId="Smeared"
0877
0878 metJERUncModules = self.getVariations(process, metModName, "Jet",preId, jetCollection, "Res", patMetUncertaintyTask, postfix=postfix )
0879
0880 for mod in metJERUncModules.keys():
0881 addToProcessAndTask(mod, metJERUncModules[mod], process, getMETUncertainties_task)
0882
0883
0884
0885
0886 if not hasattr(process, "pfCandsForUnclusteredUnc"+postfix):
0887
0888
0889 pfCandsNoJets = _mod.candPtrProjector.clone(
0890 src = pfCandCollection,
0891 veto = jetCollection,
0892 )
0893 addToProcessAndTask("pfCandsNoJets"+postfix, pfCandsNoJets, process, getMETUncertainties_task)
0894
0895
0896 pfCandsNoJetsNoEle = _mod.candPtrProjector.clone(
0897 src = "pfCandsNoJets"+postfix,
0898 veto = electronCollection,
0899 )
0900 if not self.getvalue("onMiniAOD"):
0901 pfCandsNoJetsNoEle.veto = "pfeGammaToCandidate:electrons"
0902 addToProcessAndTask("pfCandsNoJetsNoEle"+postfix, pfCandsNoJetsNoEle, process, getMETUncertainties_task)
0903
0904
0905 pfCandsNoJetsNoEleNoMu = _mod.candPtrProjector.clone(
0906 src = "pfCandsNoJetsNoEle"+postfix,
0907 veto = muonCollection,
0908 )
0909 addToProcessAndTask("pfCandsNoJetsNoEleNoMu"+postfix, pfCandsNoJetsNoEleNoMu, process, getMETUncertainties_task)
0910
0911
0912 pfCandsNoJetsNoEleNoMuNoTau = _mod.candPtrProjector.clone(
0913 src = "pfCandsNoJetsNoEleNoMu"+postfix,
0914 veto = tauCollection,
0915 )
0916 addToProcessAndTask("pfCandsNoJetsNoEleNoMuNoTau"+postfix, pfCandsNoJetsNoEleNoMuNoTau, process, getMETUncertainties_task)
0917
0918
0919 pfCandsForUnclusteredUnc = _mod.candPtrProjector.clone(
0920 src = "pfCandsNoJetsNoEleNoMuNoTau"+postfix,
0921 veto = photonCollection,
0922 )
0923 if not self.getvalue("onMiniAOD"):
0924 pfCandsForUnclusteredUnc.veto = "pfeGammaToCandidate:photons"
0925 addToProcessAndTask("pfCandsForUnclusteredUnc"+postfix, pfCandsForUnclusteredUnc, process, getMETUncertainties_task)
0926
0927
0928
0929
0930
0931
0932
0933
0934
0935
0936 pfElectrons = cms.EDFilter("CandPtrSelector",
0937 src = electronCollection,
0938 cut = cms.string("pt > 5 && isPF && gsfTrack.isAvailable() && gsfTrack.hitPattern().numberOfLostHits(\'MISSING_INNER_HITS\') < 2")
0939 )
0940 addToProcessAndTask("pfElectrons"+postfix, pfElectrons, process, getMETUncertainties_task)
0941
0942
0943
0944 pfTaus = cms.EDFilter("PATTauRefSelector",
0945 src = tauCollection,
0946 cut = cms.string('pt > 18.0 & abs(eta) < 2.6 & tauID("decayModeFinding") > 0.5 & isPFTau')
0947 )
0948 addToProcessAndTask("pfTaus"+postfix, pfTaus, process, getMETUncertainties_task)
0949
0950
0951
0952 pfMuons = cms.EDFilter("CandPtrSelector",
0953 src = muonCollection,
0954 cut = cms.string("pt > 5.0 && isPFMuon && abs(eta) < 2.4")
0955 )
0956 addToProcessAndTask("pfMuons"+postfix, pfMuons, process, getMETUncertainties_task)
0957
0958
0959
0960 if self._parameters["Puppi"].value or not self._parameters["onMiniAOD"].value:
0961 cutforpfNoPileUp = cms.string("")
0962 else:
0963 cutforpfNoPileUp = cms.string("fromPV > 1")
0964
0965 pfNoPileUp = cms.EDFilter("CandPtrSelector",
0966 src = pfCandCollection,
0967 cut = cutforpfNoPileUp
0968 )
0969 addToProcessAndTask("pfNoPileUp"+postfix, pfNoPileUp, process, getMETUncertainties_task)
0970
0971 pfPhotons = cms.EDFilter("CandPtrSelector",
0972 src = pfCandCollection if self._parameters["Puppi"].value or not self._parameters["onMiniAOD"].value else cms.InputTag("pfCHS"),
0973 cut = cms.string("abs(pdgId) = 22")
0974 )
0975 addToProcessAndTask("pfPhotons"+postfix, pfPhotons, process, getMETUncertainties_task)
0976
0977
0978
0979 electronCollection = cms.InputTag("pfElectrons"+postfix)
0980 muonCollection = cms.InputTag("pfMuons"+postfix)
0981 tauCollection = cms.InputTag("pfTaus"+postfix)
0982 photonCollection = cms.InputTag("pfPhotons"+postfix)
0983
0984
0985 objectCollections = { "Jet":jetCollection,
0986 "Electron":electronCollection,
0987 "Photon":photonCollection,
0988 "Muon":muonCollection,
0989 "Unclustered":cms.InputTag("pfCandsForUnclusteredUnc"+postfix),
0990 "Tau":tauCollection,
0991 }
0992
0993 for obj in objectCollections.keys():
0994 if not isValidInputTag(objectCollections[obj]):
0995 print("INFO : %s collection %s does not exists, no energy scale shifting will be performed in MET uncertainty tools" %(obj, objectCollections[obj]))
0996 else:
0997 metObjUncModules = self.getVariations(process, metModName, obj,"", objectCollections[obj], "En", patMetUncertaintyTask, jetUncInfos, postfix )
0998
0999
1000 for mod in metObjUncModules.keys():
1001 addToProcessAndTask(mod, metObjUncModules[mod], process, getMETUncertainties_task)
1002
1003
1004 addTaskToProcess(process, getMETUncertainties_label, getMETUncertainties_task)
1005 patMetUncertaintyTask.add(getattr(process, getMETUncertainties_label))
1006
1007
1008 def createEnergyScaleShiftedUpModule(self, process,identifier, objectCollection,
1009 varyByNsigmas, jetUncInfos=None, postfix=""):
1010
1011 shiftedModuleUp = None
1012
1013 if identifier == "Electron":
1014 shiftedModuleUp = cms.EDProducer("ShiftedParticleProducer",
1015 src = objectCollection,
1016 uncertainty = cms.string('((abs(y)<1.479)?(0.006+0*x):(0.015+0*x))'),
1017 shiftBy = cms.double(+1.*varyByNsigmas),
1018 srcWeights = cms.InputTag("")
1019 )
1020
1021 if identifier == "Photon":
1022 shiftedModuleUp = cms.EDProducer("ShiftedParticleProducer",
1023 src = objectCollection,
1024 uncertainty = cms.string('((abs(y)<1.479)?(0.01+0*x):(0.025+0*x))'),
1025 shiftBy = cms.double(+1.*varyByNsigmas),
1026 srcWeights = cms.InputTag("")
1027 )
1028
1029 if identifier == "Muon":
1030 shiftedModuleUp = cms.EDProducer("ShiftedParticleProducer",
1031 src = objectCollection,
1032 uncertainty = cms.string('((x<100)?(0.002+0*y):(0.05+0*y))'),
1033 shiftBy = cms.double(+1.*varyByNsigmas),
1034 srcWeights = cms.InputTag("")
1035 )
1036
1037 if identifier == "Tau":
1038 shiftedModuleUp = cms.EDProducer("ShiftedParticleProducer",
1039 src = objectCollection,
1040 uncertainty = cms.string('0.03+0*x*y'),
1041 shiftBy = cms.double(+1.*varyByNsigmas),
1042 srcWeights = cms.InputTag("")
1043 )
1044
1045 if identifier == "Unclustered":
1046 shiftedModuleUp = cms.EDProducer("ShiftedParticleProducer",
1047 src = objectCollection,
1048 binning = cms.VPSet(
1049
1050 cms.PSet(
1051 binSelection = cms.string('charge!=0'),
1052 binUncertainty = cms.string('sqrt(pow(0.00009*x,2)+pow(0.0085/sqrt(sin(2*atan(exp(-y)))),2))')
1053 ),
1054
1055 cms.PSet(
1056 binSelection = cms.string('pdgId==130'),
1057 energyDependency = cms.bool(True),
1058 binUncertainty = cms.string('((abs(y)<1.3)?(min(0.25,sqrt(0.64/x+0.0025))):(min(0.30,sqrt(1.0/x+0.0016))))')
1059 ),
1060
1061 cms.PSet(
1062 binSelection = cms.string('pdgId==22'),
1063 energyDependency = cms.bool(True),
1064 binUncertainty = cms.string('sqrt(0.0009/x+0.000001)+0*y')
1065 ),
1066
1067 cms.PSet(
1068 binSelection = cms.string('pdgId==1 || pdgId==2'),
1069 energyDependency = cms.bool(True),
1070 binUncertainty = cms.string('sqrt(1./x+0.0025)+0*y')
1071 ),
1072 ),
1073 shiftBy = cms.double(+1.*varyByNsigmas),
1074 srcWeights = cms.InputTag("")
1075 )
1076
1077 if identifier == "Jet":
1078 moduleType="ShiftedPATJetProducer"
1079
1080 if jetUncInfos["jecUncFile"] == "":
1081 shiftedModuleUp = cms.EDProducer(moduleType,
1082 src = objectCollection,
1083 jetCorrUncertaintyTag = cms.string(jetUncInfos["jecUncTag"] ),
1084 addResidualJES = cms.bool(True),
1085 jetCorrLabelUpToL3 = cms.InputTag(jetUncInfos["jCorLabelUpToL3"] ),
1086 jetCorrLabelUpToL3Res = cms.InputTag(jetUncInfos["jCorLabelL3Res"] ),
1087 jetCorrPayloadName = cms.string(jetUncInfos["jCorrPayload"] ),
1088 shiftBy = cms.double(+1.*varyByNsigmas),
1089 )
1090 else:
1091 shiftedModuleUp = cms.EDProducer(moduleType,
1092 src = objectCollection,
1093 jetCorrInputFileName = cms.FileInPath(jetUncInfos["jecUncFile"] ),
1094 jetCorrUncertaintyTag = cms.string(jetUncInfos["jecUncTag"] ),
1095 addResidualJES = cms.bool(True),
1096 jetCorrLabelUpToL3 = cms.InputTag(jetUncInfos["jCorLabelUpToL3"] ),
1097 jetCorrLabelUpToL3Res = cms.InputTag(jetUncInfos["jCorLabelL3Res"] ),
1098 jetCorrPayloadName = cms.string(jetUncInfos["jCorrPayload"] ),
1099 shiftBy = cms.double(+1.*varyByNsigmas),
1100 )
1101
1102
1103 return shiftedModuleUp
1104
1105
1106
1107
1108
1109
1110 def removePostfix(self, name, postfix):
1111
1112 if postfix=="":
1113 return name
1114
1115 baseName = name
1116 if baseName[-len(postfix):] == postfix:
1117 baseName = baseName[0:-len(postfix)]
1118 else:
1119 raise Exception("Tried to remove postfix %s from %s, but it wasn't there" % (postfix, baseName))
1120
1121 return baseName
1122
1123
1124 def tuneTxyParameters(self, process, corScheme, postfix, datamc="", campaign="", era="" ):
1125 corSchemes = ["Txy", "T1Txy", "T0pcTxy", "T0pcT1Txy", "T1T2Txy", "T0pcT1T2Txy", "T1SmearTxy", "T1T2SmearTxy", "T0pcT1SmearTxy", "T0pcT1T2SmearTxy"]
1126 import PhysicsTools.PatUtils.patPFMETCorrections_cff as metCors
1127 xyTags = {}
1128 for corScheme_ in corSchemes:
1129 xyTags["{}_{}".format(corScheme_,"50ns")]=getattr(metCors,"{}_{}_{}".format("patMultPhiCorrParams",corScheme_,"50ns"))
1130 xyTags["{}_{}".format(corScheme_,"25ns")]=getattr(metCors,"{}_{}_{}".format("patMultPhiCorrParams",corScheme_,"25ns"))
1131 if datamc!="" and campaign!="" and era!="":
1132 if not self.getvalue("Puppi"):
1133 xyTags["{}_{}_{}_{}".format(corScheme_,campaign,datamc,era)]=getattr(metCors,"{}_{}{}{}".format("patMultPhiCorrParams",campaign,datamc,era))
1134 else:
1135 xyTags["{}_{}_{}_{}".format(corScheme_,campaign,datamc,era)]=getattr(metCors,"{}_{}{}{}".format("patMultPhiCorrParams_Puppi",campaign,datamc,era))
1136
1137 if datamc!="" and campaign!="" and era!="":
1138 getattr(process, "patPFMetTxyCorr"+postfix).parameters = xyTags["{}_{}_{}_{}".format(corScheme,campaign,datamc,era)]
1139 else:
1140 getattr(process, "patPFMetTxyCorr"+postfix).parameters = xyTags[corScheme+"_25ns"]
1141
1142
1143
1144 def getVariations(self, process, metModName, identifier,preId, objectCollection, varType,
1145 patMetUncertaintyTask, jetUncInfos=None, postfix="" ):
1146
1147
1148 varyByNsigmas=1
1149
1150
1151 baseName = self.removePostfix(metModName, postfix)
1152
1153
1154 shiftedMetProducers = {preId+identifier+varType+'Up':None, preId+identifier+varType+'Down':None}
1155
1156
1157 shiftedCollModules = {'Up':None, 'Down':None}
1158
1159 if identifier=="Jet" and varType=="Res":
1160 smear=False
1161 if "Smear" in metModName:
1162 smear=True
1163 else:
1164 smear=True
1165 varyByNsigmas=101
1166
1167 shiftedCollModules['Up'] = self.createShiftedJetResModule(process, smear, objectCollection, +1.*varyByNsigmas,
1168 "Up", postfix)
1169 shiftedCollModules['Down'] = self.createShiftedJetResModule(process, smear, objectCollection, -1.*varyByNsigmas,
1170 "Down", postfix)
1171 else:
1172 shiftedCollModules['Up'] = self.createEnergyScaleShiftedUpModule(process, identifier, objectCollection, varyByNsigmas, jetUncInfos, postfix)
1173 shiftedCollModules['Down'] = shiftedCollModules['Up'].clone( shiftBy = cms.double(-1.*varyByNsigmas) )
1174
1175 if identifier=="Jet" and varType=="Res":
1176 smear=False
1177 if "Smear" in metModName:
1178 objectCollection=cms.InputTag("selectedPatJetsForMetT1T2SmearCorr"+postfix)
1179
1180
1181
1182
1183 if identifier=="Jet" and varType=="Res" and self._parameters["runOnData"].value:
1184 shiftedMetProducers = self.copyCentralMETProducer(process, shiftedCollModules, identifier, metModName, varType, postfix)
1185 else:
1186 shiftedMetProducers = self.createShiftedModules(process, shiftedCollModules, identifier, preId, objectCollection,
1187 metModName, varType, patMetUncertaintyTask, postfix)
1188
1189 return shiftedMetProducers
1190
1191
1192
1193 def copyCentralMETProducer(self, process, shiftedCollModules, identifier, metModName, varType, postfix):
1194
1195
1196 shiftedMetProducers = {}
1197 baseName = self.removePostfix(metModName, postfix)
1198 for mod in shiftedCollModules.keys():
1199 modName = baseName+identifier+varType+mod+postfix
1200 shiftedMETModule = getattr(process, metModName).clone()
1201 shiftedMetProducers[ modName ] = shiftedMETModule
1202
1203 return shiftedMetProducers
1204
1205
1206
1207 def createShiftedJetResModule(self, process, smear, objectCollection, varyByNsigmas, varDir, postfix ):
1208
1209 smearedJetModule = self.createSmearedJetModule(process, objectCollection, smear, varyByNsigmas, varDir, postfix)
1210
1211 return smearedJetModule
1212
1213
1214
1215 def createShiftedModules(self, process, shiftedCollModules, identifier, preId, objectCollection,
1216 metModName, varType, patMetUncertaintyTask, postfix):
1217
1218 createShiftedModules_task, createShiftedModules_label = cms.Task(), "createShiftedModules_task{}".format(postfix)
1219
1220 shiftedMetProducers = {}
1221
1222
1223 baseName = self.removePostfix(metModName, postfix)
1224
1225
1226 for mod in shiftedCollModules.keys():
1227 modName = "shiftedPat"+preId+identifier+varType+mod+postfix
1228 if (identifier=="Photon" or identifier=="Unclustered") and self.getvalue("Puppi"):
1229 shiftedCollModules[mod].srcWeights = self._parameters['puppiProducerForMETLabel'].value
1230 if not hasattr(process, modName):
1231 addToProcessAndTask(modName, shiftedCollModules[mod], process, createShiftedModules_task)
1232
1233
1234 modName = "shiftedPat"+preId+identifier+varType+mod+postfix
1235
1236
1237 if "PF" in metModName:
1238
1239 shiftedMETCorrModule = self.createShiftedMETModule(process, objectCollection, modName)
1240 if (identifier=="Photon" or identifier=="Unclustered") and self.getvalue("Puppi"):
1241 shiftedMETCorrModule.srcWeights = self._parameters['puppiProducerForMETLabel'].value
1242 modMETShiftName = "shiftedPatMETCorr"+preId+identifier+varType+mod+postfix
1243 if not hasattr(process, modMETShiftName):
1244 addToProcessAndTask(modMETShiftName, shiftedMETCorrModule, process, createShiftedModules_task)
1245
1246
1247 modName = baseName+identifier+varType+mod+postfix
1248 shiftedMETModule = getattr(process, metModName).clone(
1249 src = cms.InputTag( metModName ),
1250 srcCorrections = cms.VInputTag( cms.InputTag(modMETShiftName) )
1251 )
1252 shiftedMetProducers[ modName ] = shiftedMETModule
1253
1254
1255
1256 addTaskToProcess(process, createShiftedModules_label, createShiftedModules_task)
1257 patMetUncertaintyTask.add(getattr(process, createShiftedModules_label))
1258
1259 return shiftedMetProducers
1260
1261
1262
1263 def createShiftedMETModule(self, process, originCollection, shiftedCollection):
1264
1265 shiftedModule = cms.EDProducer("ShiftedParticleMETcorrInputProducer",
1266 srcOriginal = originCollection,
1267 srcShifted = cms.InputTag(shiftedCollection),
1268 srcWeights = cms.InputTag("")
1269 )
1270
1271 return shiftedModule
1272
1273
1274 def createSmearedJetModule(self, process, jetCollection, smear, varyByNsigmas, varDir, postfix):
1275
1276 smearedJetModule = None
1277
1278 modName = "pat"
1279 selJetModName= "selectedPatJetsForMetT1T2"
1280 if smear:
1281 modName += "SmearedJets"
1282 selJetModName += "SmearCorr"
1283 else:
1284 modName += "Jets"
1285
1286
1287 if varDir != "":
1288 modName += "Res"+varDir
1289 selJetModName += "Res"+varDir
1290
1291 modName += postfix
1292 selJetModName += postfix
1293
1294 genJetsCollection=cms.InputTag('ak4GenJetsNoNu')
1295 if self._parameters["onMiniAOD"].value:
1296 genJetsCollection=cms.InputTag("slimmedGenJets")
1297
1298 if self._parameters["Puppi"].value:
1299 getattr(process, "patSmearedJets"+postfix).algo = 'AK4PFPuppi'
1300 getattr(process, "patSmearedJets"+postfix).algopt = 'AK4PFPuppi_pt'
1301
1302 if "PF" == self._parameters["metType"].value:
1303 smearedJetModule = getattr(process, "patSmearedJets"+postfix).clone(
1304 src = jetCollection,
1305 enabled = cms.bool(smear),
1306 variation = cms.int32( int(varyByNsigmas) ),
1307 genJets = genJetsCollection,
1308 )
1309
1310 return smearedJetModule
1311
1312
1313
1314 def labelsInSequence(process, sequenceLabel, postfix=""):
1315 result = [ m.label()[:-len(postfix)] for m in listModules( getattr(process,sequenceLabel+postfix))]
1316 result.extend([ m.label()[:-len(postfix)] for m in listSequences( getattr(process,sequenceLabel+postfix))] )
1317 if postfix == "":
1318 result = [ m.label() for m in listModules( getattr(process,sequenceLabel+postfix))]
1319 result.extend([ m.label() for m in listSequences( getattr(process,sequenceLabel+postfix))] )
1320 return result
1321
1322 def initializeInputTag(self, input, default):
1323 retVal = None
1324 if input is None:
1325 retVal = self._defaultParameters[default].value
1326 elif isinstance(input, str):
1327 retVal = cms.InputTag(input)
1328 else:
1329 retVal = input
1330 return retVal
1331
1332
1333 def recomputeRawMetFromPfcs(self, process, pfCandCollection, onMiniAOD, patMetModuleTask, postfix):
1334
1335 recomputeRawMetFromPfcs_task, recomputeRawMetFromPfcs_label = cms.Task(), "recomputeRawMetFromPfcs_task{}".format(postfix)
1336
1337
1338 if not hasattr(process, "pfMet"+postfix) and self._parameters["metType"].value == "PF":
1339
1340
1341
1342 from RecoMET.METProducers.pfMet_cfi import pfMet
1343 addToProcessAndTask("pfMet"+postfix, pfMet.clone(), process, recomputeRawMetFromPfcs_task)
1344 getattr(process, "pfMet"+postfix).src = pfCandCollection
1345 getattr(process, "pfMet"+postfix).calculateSignificance = False
1346
1347 if self.getvalue("Puppi"):
1348 getattr(process, "pfMet"+postfix).applyWeight = True
1349 getattr(process, "pfMet"+postfix).srcWeights = self._parameters['puppiProducerForMETLabel'].value
1350
1351
1352 if not hasattr(process, "patMETs"+postfix) and self._parameters["metType"].value == "PF":
1353 process.load("PhysicsTools.PatAlgos.producersLayer1.metProducer_cff")
1354 recomputeRawMetFromPfcs_task.add(process.makePatMETsTask)
1355 configtools.cloneProcessingSnippetTask(process, getattr(process,"patMETCorrectionsTask"), postfix)
1356 recomputeRawMetFromPfcs_task.add(getattr(process,"patMETCorrectionsTask"+postfix))
1357
1358
1359 if not onMiniAOD:
1360
1361 getattr(process, "pfMetT1"+postfix).src = cms.InputTag("pfMet"+postfix)
1362 recomputeRawMetFromPfcs_task.add(getattr(process, "pfMetT1"+postfix))
1363 _myPatMet = 'patMETs'+postfix
1364 addToProcessAndTask(_myPatMet, getattr(process,'patMETs' ).clone(), process, recomputeRawMetFromPfcs_task)
1365 getattr(process, _myPatMet).metSource = cms.InputTag("pfMetT1"+postfix)
1366 getattr(process, _myPatMet).computeMETSignificance = cms.bool(self.getvalue("computeMETSignificance"))
1367 getattr(process, _myPatMet).srcLeptons = \
1368 cms.VInputTag(copy.copy(self.getvalue("electronCollection")) if self.getvalue("onMiniAOD") else
1369 cms.InputTag("pfeGammaToCandidate","electrons"),
1370 copy.copy(self.getvalue("muonCollection")),
1371 copy.copy(self.getvalue("photonCollection")) if self.getvalue("onMiniAOD") else
1372 cms.InputTag("pfeGammaToCandidate","photons"))
1373 if postfix=="NoHF":
1374 getattr(process, _myPatMet).computeMETSignificance = cms.bool(False)
1375
1376 if self.getvalue("Puppi"):
1377 getattr(process, _myPatMet).srcWeights = self._parameters['puppiProducerForMETLabel'].value
1378 getattr(process, _myPatMet).srcJets = cms.InputTag('cleanedPatJets'+postfix)
1379 getattr(process, _myPatMet).srcJetSF = cms.string('AK4PFPuppi')
1380 getattr(process, _myPatMet).srcJetResPt = cms.string('AK4PFPuppi_pt')
1381 getattr(process, _myPatMet).srcJetResPhi = cms.string('AK4PFPuppi_phi')
1382
1383
1384 addTaskToProcess(process, recomputeRawMetFromPfcs_label, recomputeRawMetFromPfcs_task)
1385
1386
1387 patMetModuleTask.add(getattr(process, recomputeRawMetFromPfcs_label))
1388
1389
1390 def extractMET(self, process, correctionLevel, patMetModuleTask, postfix):
1391
1392 extractMET_task, extractMET_label = cms.Task(), "extractMET_task{}".format(postfix)
1393
1394
1395 pfMet = cms.EDProducer("RecoMETExtractor",
1396 metSource= cms.InputTag("slimmedMETs" if not self._parameters["Puppi"].value else "slimmedMETsPuppi",processName=cms.InputTag.skipCurrentProcess()),
1397 correctionLevel = cms.string(correctionLevel)
1398 )
1399 if(correctionLevel=="raw"):
1400 addToProcessAndTask("pfMet"+postfix, pfMet, process, extractMET_task)
1401 else:
1402 addToProcessAndTask("met"+correctionLevel+postfix, pfMet, process, extractMET_task)
1403
1404 if not hasattr(process, "genMetExtractor"+postfix) and not self._parameters["runOnData"].value:
1405 genMetExtractor = cms.EDProducer("GenMETExtractor",
1406 metSource= cms.InputTag("slimmedMETs",processName=cms.InputTag.skipCurrentProcess())
1407 )
1408 addToProcessAndTask("genMetExtractor"+postfix, genMetExtractor, process, extractMET_task)
1409
1410
1411 addTaskToProcess(process, extractMET_label, extractMET_task)
1412
1413
1414 patMetModuleTask.add(getattr(process, extractMET_label))
1415
1416
1417 def updateJECs(self,process,jetCollection, patMetModuleTask, postfix):
1418 from PhysicsTools.PatAlgos.producersLayer1.jetUpdater_cff import updatedPatJetCorrFactors
1419
1420 patJetCorrFactorsReapplyJEC = updatedPatJetCorrFactors.clone(
1421 src = jetCollection if not self._parameters["Puppi"].value else cms.InputTag("slimmedJetsPuppi"),
1422 levels = ['L1FastJet',
1423 'L2Relative',
1424 'L3Absolute'],
1425 payload = 'AK4PFchs' if not self._parameters["Puppi"].value else 'AK4PFPuppi' )
1426
1427 if self._parameters["runOnData"].value:
1428 patJetCorrFactorsReapplyJEC.levels.append("L2L3Residual")
1429
1430 from PhysicsTools.PatAlgos.producersLayer1.jetUpdater_cff import updatedPatJets
1431 patJetsReapplyJEC = updatedPatJets.clone(
1432 jetSource = jetCollection if not self._parameters["Puppi"].value else cms.InputTag("slimmedJetsPuppi"),
1433 jetCorrFactorsSource = cms.VInputTag(cms.InputTag("patJetCorrFactorsReapplyJEC"+postfix))
1434 )
1435
1436
1437 updateJECs_task, updateJECs_label = cms.Task(), "updateJECs_task{}".format(postfix)
1438 addToProcessAndTask("patJetCorrFactorsReapplyJEC"+postfix, patJetCorrFactorsReapplyJEC, process, updateJECs_task)
1439 addToProcessAndTask("patJetsReapplyJEC"+postfix, patJetsReapplyJEC.clone(), process, updateJECs_task)
1440
1441
1442 addTaskToProcess(process, updateJECs_label, updateJECs_task)
1443
1444
1445 patMetModuleTask.add(getattr(process, updateJECs_label))
1446
1447 return cms.InputTag("patJetsReapplyJEC"+postfix)
1448
1449
1450 def getJetCollectionForCorsAndUncs(self, process, jetCollectionUnskimmed,
1451 jetSelection, autoJetCleaning,patMetModuleTask, postfix):
1452
1453 getJetCollectionForCorsAndUncs_task, getJetCollectionForCorsAndUncs_label = cms.Task(), "getJetCollectionForCorsAndUncs_task{}".format(postfix)
1454
1455 basicJetsForMet = cms.EDProducer("PATJetCleanerForType1MET",
1456 src = jetCollectionUnskimmed,
1457 jetCorrEtaMax = cms.double(9.9),
1458 jetCorrLabel = cms.InputTag("L3Absolute"),
1459 jetCorrLabelRes = cms.InputTag("L2L3Residual"),
1460 offsetCorrLabel = cms.InputTag("L1FastJet"),
1461 skipEM = cms.bool(True),
1462 skipEMfractionThreshold = cms.double(0.9),
1463 skipMuonSelection = cms.string('isGlobalMuon | isStandAloneMuon'),
1464 skipMuons = cms.bool(True),
1465 type1JetPtThreshold = cms.double(15.0)
1466 )
1467 addToProcessAndTask("basicJetsForMet"+postfix, basicJetsForMet, process, getJetCollectionForCorsAndUncs_task)
1468
1469 from PhysicsTools.PatAlgos.selectionLayer1.jetSelector_cfi import selectedPatJets
1470 jetSelector = selectedPatJets.clone(
1471 src = cms.InputTag("basicJetsForMet"+postfix),
1472 cut = cms.string(jetSelection)
1473 )
1474 addToProcessAndTask("jetSelectorForMet"+postfix, jetSelector, process, getJetCollectionForCorsAndUncs_task)
1475
1476 jetCollection = self.jetCleaning(process, "jetSelectorForMet"+postfix, autoJetCleaning, patMetModuleTask, postfix)
1477 addTaskToProcess(process, getJetCollectionForCorsAndUncs_label, getJetCollectionForCorsAndUncs_task)
1478 patMetModuleTask.add(getattr(process, getJetCollectionForCorsAndUncs_label))
1479
1480 return jetCollection
1481
1482
1483 def ak4JetReclustering(self,process, pfCandCollection, patMetModuleTask, postfix):
1484
1485 ak4JetReclustering_task, ak4JetReclustering_label = cms.Task(), "ak4JetReclustering_task{}".format(postfix)
1486
1487 chs = self._parameters["CHS"].value
1488 jetColName="ak4PFJets"
1489 CHSname=""
1490 pfCandColl=pfCandCollection
1491 if chs:
1492 CHSname="chs"
1493 jetColName="ak4PFJetsCHS"
1494 if self._parameters["onMiniAOD"].value:
1495 pfCandColl = cms.InputTag("pfCHS")
1496 else:
1497 addToProcessAndTask("tmpPFCandCollPtr"+postfix,
1498 cms.EDProducer("PFCandidateFwdPtrProducer",
1499 src = pfCandCollection ),
1500 process, ak4JetReclustering_task)
1501 process.load("CommonTools.ParticleFlow.pfNoPileUpJME_cff")
1502 ak4JetReclustering_task.add(process.pfNoPileUpJMETask)
1503 configtools.cloneProcessingSnippetTask(process, getattr(process,"pfNoPileUpJMETask"), postfix)
1504 ak4JetReclustering_task.add(getattr(process,"pfNoPileUpJMETask"+postfix))
1505 getattr(process, "pfPileUpJME"+postfix).PFCandidates = "tmpPFCandCollPtr"+postfix
1506 getattr(process, "pfNoPileUpJME"+postfix).bottomCollection = "tmpPFCandCollPtr"+postfix
1507 pfCandColl = "pfNoPileUpJME"+postfix
1508
1509 jetColName+=postfix
1510 if not hasattr(process, jetColName):
1511 from RecoJets.JetProducers.ak4PFJets_cfi import ak4PFJets
1512
1513 addToProcessAndTask(jetColName, ak4PFJets.clone(), process, ak4JetReclustering_task)
1514 getattr(process, jetColName).src = pfCandColl
1515 getattr(process, jetColName).doAreaFastjet = True
1516
1517
1518 if self._parameters["Puppi"].value:
1519 getattr(process, jetColName).srcWeights = cms.InputTag(self._parameters['puppiProducerLabel'].value)
1520 getattr(process, jetColName).applyWeight = True
1521
1522 corLevels=['L1FastJet', 'L2Relative', 'L3Absolute']
1523 if self._parameters["runOnData"].value:
1524 corLevels.append("L2L3Residual")
1525
1526 switchJetCollection(process,
1527 jetSource = cms.InputTag(jetColName),
1528 jetCorrections = ('AK4PF'+CHSname, corLevels , ''),
1529 postfix=postfix
1530 )
1531
1532 getattr(process,"patJets"+postfix).addGenJetMatch = False
1533 getattr(process,"patJets"+postfix).addGenPartonMatch = False
1534 getattr(process,"patJets"+postfix).addPartonJetMatch = False
1535 getattr(process,"patJets"+postfix).embedGenPartonMatch = False
1536 getattr(process,"patJets"+postfix).embedGenJetMatch = False
1537 if self._parameters['onMiniAOD'].value:
1538 del getattr(process,"patJets"+postfix).JetFlavourInfoSource
1539 del getattr(process,"patJets"+postfix).JetPartonMapSource
1540 del getattr(process,"patJets"+postfix).genPartonMatch
1541 del getattr(process,"patJets"+postfix).genJetMatch
1542 getattr(process,"patJets"+postfix).getJetMCFlavour = False
1543
1544 getattr(process,"patJetCorrFactors"+postfix).src=cms.InputTag(jetColName)
1545 getattr(process,"patJetCorrFactors"+postfix).primaryVertices= cms.InputTag("offlineSlimmedPrimaryVertices")
1546 if self._parameters["Puppi"].value:
1547 getattr(process,"patJetCorrFactors"+postfix).payload=cms.string('AK4PFPuppi')
1548
1549
1550 addTaskToProcess(process, ak4JetReclustering_label, ak4JetReclustering_task)
1551
1552
1553 patMetModuleTask.add(getattr(process, ak4JetReclustering_label))
1554
1555 return cms.InputTag("patJets"+postfix)
1556
1557
1558
1559 def miniAODConfigurationPre(self, process, patMetModuleTask, pfCandCollection, postfix):
1560
1561 miniAODConfigurationPre_task, miniAODConfigurationPre_label = cms.Task(), "miniAODConfigurationPre_task{}".format(postfix)
1562
1563
1564
1565 self.extractMET(process, "rawCalo", miniAODConfigurationPre_task, postfix)
1566 caloMetName="metrawCalo" if hasattr(process,"metrawCalo") else "metrawCalo"+postfix
1567 from PhysicsTools.PatAlgos.tools.metTools import addMETCollection
1568 addMETCollection(process,
1569 labelName = "patCaloMet",
1570 metSource = caloMetName
1571 )
1572 getattr(process,"patCaloMet").addGenMET = False
1573 miniAODConfigurationPre_task.add(getattr(process,"patCaloMet"))
1574
1575
1576
1577 if self._parameters["extractDeepMETs"].value:
1578 self.extractMET(process, "rawDeepResponseTune", miniAODConfigurationPre_task, postfix)
1579 deepMetResponseTuneName = "metrawDeepResponseTune" if hasattr(process, "metrawDeepResponseTune") else "metrawDeepResponseTune"+postfix
1580 addMETCollection(process,
1581 labelName = "deepMETsResponseTune",
1582 metSource = deepMetResponseTuneName
1583 )
1584 getattr(process, "deepMETsResponseTune").addGenMET = False
1585 getattr(process, "deepMETsResponseTune").computeMETSignificance = cms.bool(False)
1586 miniAODConfigurationPre_task.add(getattr(process, "deepMETsResponseTune"))
1587
1588 self.extractMET(process, "rawDeepResolutionTune", miniAODConfigurationPre_task, postfix)
1589 deepMetResolutionTuneName = "metrawDeepResolutionTune" if hasattr(process, "metrawDeepResolutionTune") else "metrawDeepResolutionTune"+postfix
1590 addMETCollection(process,
1591 labelName = "deepMETsResolutionTune",
1592 metSource = deepMetResolutionTuneName
1593 )
1594 getattr(process, "deepMETsResolutionTune").addGenMET = False
1595 getattr(process, "deepMETsResolutionTune").computeMETSignificance = cms.bool(False)
1596 miniAODConfigurationPre_task.add(getattr(process, "deepMETsResolutionTune"))
1597
1598
1599
1600
1601
1602 from CommonTools.ParticleFlow.pfCHS_cff import pfCHS
1603 addToProcessAndTask("pfCHS", pfCHS.clone(), process, miniAODConfigurationPre_task)
1604 from RecoMET.METProducers.pfMet_cfi import pfMet
1605 pfMetCHS = pfMet.clone(src = "pfCHS")
1606 addToProcessAndTask("pfMetCHS", pfMetCHS, process, miniAODConfigurationPre_task)
1607
1608 addMETCollection(process,
1609 labelName = "patCHSMet",
1610 metSource = "pfMetCHS"
1611 )
1612
1613 process.patCHSMet.computeMETSignificant = cms.bool(False)
1614 process.patCHSMet.addGenMET = cms.bool(False)
1615 miniAODConfigurationPre_task.add(process.patCHSMet)
1616
1617
1618 pfTrk = chargedPackedCandsForTkMet.clone()
1619 addToProcessAndTask("pfTrk", pfTrk, process, miniAODConfigurationPre_task)
1620 pfMetTrk = pfMet.clone(src = 'pfTrk')
1621 addToProcessAndTask("pfMetTrk", pfMetTrk, process, miniAODConfigurationPre_task)
1622
1623 addMETCollection(process,
1624 labelName = "patTrkMet",
1625 metSource = "pfMetTrk"
1626 )
1627
1628 process.patTrkMet.computeMETSignificant = cms.bool(False)
1629 process.patTrkMet.addGenMET = cms.bool(False)
1630 miniAODConfigurationPre_task.add(process.patTrkMet)
1631
1632 addTaskToProcess(process, miniAODConfigurationPre_label, miniAODConfigurationPre_task)
1633
1634 patMetModuleTask.add(getattr(process, miniAODConfigurationPre_label))
1635
1636 def miniAODConfigurationPost(self, process, postfix):
1637
1638 if self._parameters["metType"].value == "PF":
1639 if hasattr(process, "patPFMetTxyCorr"+postfix):
1640 getattr(process, "patPFMetTxyCorr"+postfix).vertexCollection = cms.InputTag("offlineSlimmedPrimaryVertices")
1641
1642 if self._parameters['computeUncertainties'].value and not self._parameters["runOnData"]:
1643 getattr(process, "shiftedPatJetResDown"+postfix).genJets = cms.InputTag("slimmedGenJets")
1644 getattr(process, "shiftedPatJetResUp"+postfix).genJets = cms.InputTag("slimmedGenJets")
1645
1646
1647 def miniAODConfiguration(self, process, pfCandCollection, jetCollection,
1648 patMetModuleTask, postfix ):
1649
1650 if self._parameters["metType"].value == "PF":
1651 if "T1" in self._parameters['correctionLevel'].value:
1652 getattr(process, "patPFMet"+postfix).srcJets = jetCollection
1653 getattr(process, "patPFMet"+postfix).addGenMET = False
1654 if not self._parameters["runOnData"].value:
1655 getattr(process, "patPFMet"+postfix).addGenMET = True
1656 getattr(process, "patPFMet"+postfix).genMETSource = cms.InputTag("genMetExtractor"+postfix)
1657
1658
1659 if "Smear" in self._parameters['correctionLevel'].value:
1660 getattr(process, "patSmearedJets"+postfix).genJets = cms.InputTag("slimmedGenJets")
1661
1662
1663 if not hasattr(process, "slimmedMETs"+postfix) and self._parameters["metType"].value == "PF":
1664
1665
1666 miniAODConfiguration_task, miniAODConfiguration_label = cms.Task(), "miniAODConfiguration_task{}".format(postfix)
1667
1668 from PhysicsTools.PatAlgos.slimming.slimmedMETs_cfi import slimmedMETs
1669 addToProcessAndTask("slimmedMETs"+postfix, slimmedMETs.clone(), process, miniAODConfiguration_task)
1670 getattr(process,"slimmedMETs"+postfix).src = cms.InputTag("patPFMetT1"+postfix)
1671 getattr(process,"slimmedMETs"+postfix).rawVariation = cms.InputTag("patPFMet"+postfix)
1672 getattr(process,"slimmedMETs"+postfix).t1Uncertainties = cms.InputTag("patPFMetT1%s"+postfix)
1673 getattr(process,"slimmedMETs"+postfix).t01Variation = cms.InputTag("patPFMetT0pcT1"+postfix)
1674 getattr(process,"slimmedMETs"+postfix).t1SmearedVarsAndUncs = cms.InputTag("patPFMetT1Smear%s"+postfix)
1675
1676 getattr(process,"slimmedMETs"+postfix).tXYUncForRaw = cms.InputTag("patPFMetTxy"+postfix)
1677 getattr(process,"slimmedMETs"+postfix).tXYUncForT1 = cms.InputTag("patPFMetT1Txy"+postfix)
1678 getattr(process,"slimmedMETs"+postfix).tXYUncForT01 = cms.InputTag("patPFMetT0pcT1Txy"+postfix)
1679 getattr(process,"slimmedMETs"+postfix).tXYUncForT1Smear = cms.InputTag("patPFMetT1SmearTxy"+postfix)
1680 getattr(process,"slimmedMETs"+postfix).tXYUncForT01Smear = cms.InputTag("patPFMetT0pcT1SmearTxy"+postfix)
1681
1682 getattr(process,"slimmedMETs"+postfix).runningOnMiniAOD = True
1683 getattr(process,"slimmedMETs"+postfix).t01Variation = cms.InputTag("slimmedMETs" if not self._parameters["Puppi"].value else "slimmedMETsPuppi",processName=cms.InputTag.skipCurrentProcess())
1684
1685 if hasattr(process, "deepMETsResolutionTune") and hasattr(process, "deepMETsResponseTune"):
1686
1687
1688 getattr(process,"slimmedMETs"+postfix).addDeepMETs = True
1689
1690
1691
1692 del getattr(process,"slimmedMETs"+postfix).tXYUncForRaw
1693 del getattr(process,"slimmedMETs"+postfix).tXYUncForT01
1694 del getattr(process,"slimmedMETs"+postfix).tXYUncForT1Smear
1695 del getattr(process,"slimmedMETs"+postfix).tXYUncForT01Smear
1696
1697
1698
1699 addTaskToProcess(process, miniAODConfiguration_label, miniAODConfiguration_task)
1700
1701
1702 patMetModuleTask.add(getattr(process, miniAODConfiguration_label))
1703
1704 def jetConfiguration(self):
1705
1706 jetFlavor = self._parameters["jetFlavor"].value
1707 jetCorr = self._parameters["jetCorrectionType"].value
1708
1709 jetCorLabelUpToL3Name="ak4PF"
1710 jetCorLabelL3ResName="ak4PF"
1711
1712
1713 if "chs" in jetFlavor:
1714 self.setParameter("CHS",True)
1715 jetCorLabelUpToL3Name += "CHS"
1716 jetCorLabelL3ResName += "CHS"
1717 elif "Puppi" in jetFlavor:
1718 self.setParameter("CHS",False)
1719 jetCorLabelUpToL3Name += "Puppi"
1720 jetCorLabelL3ResName += "Puppi"
1721
1722 else:
1723 self.setParameter("CHS",False)
1724
1725
1726 if jetCorr == "L1L2L3-L1":
1727 jetCorLabelUpToL3Name += "L1FastL2L3Corrector"
1728 jetCorLabelL3ResName += "L1FastL2L3ResidualCorrector"
1729 elif jetCorr == "L1L2L3-RC":
1730 jetCorLabelUpToL3Name += "L1FastL2L3Corrector"
1731 jetCorLabelL3ResName += "L1FastL2L3ResidualCorrector"
1732
1733 self.setParameter("jetCorLabelUpToL3",jetCorLabelUpToL3Name )
1734 self.setParameter("jetCorLabelL3Res",jetCorLabelL3ResName )
1735
1736
1737 def jetCleaning(self, process, jetCollectionName, autoJetCleaning, jetProductionTask, postfix ):
1738
1739 if autoJetCleaning == "None" or autoJetCleaning == "Manual" :
1740 return cms.InputTag(jetCollectionName)
1741
1742
1743 electronCollection = self._parameters["electronCollection"].value
1744 muonCollection = self._parameters["muonCollection"].value
1745 photonCollection = self._parameters["photonCollection"].value
1746 tauCollection = self._parameters["tauCollection"].value
1747
1748 jetCleaning_task, jetCleaning_label = cms.Task(), "jetCleaning_task{}".format(postfix)
1749
1750 if autoJetCleaning == "Full" :
1751 if isValidInputTag(tauCollection):
1752 process.load("PhysicsTools.PatAlgos.cleaningLayer1.tauCleaner_cfi")
1753 jetCleaning_task.add(process.cleanPatTaus)
1754 cleanPatTauProducer = getattr(process, "cleanPatTaus").clone(
1755 src = tauCollection
1756 )
1757 cleanPatTauProducer.checkOverlaps.electrons.src = electronCollection
1758 cleanPatTauProducer.checkOverlaps.muons.src = muonCollection
1759 addToProcessAndTask("cleanedPatTaus"+postfix, cleanPatTauProducer, process, jetCleaning_task)
1760 tauCollection = cms.InputTag("cleanedPatTaus"+postfix)
1761
1762 if isValidInputTag(photonCollection):
1763 process.load("PhysicsTools.PatAlgos.cleaningLayer1.photonCleaner_cfi")
1764 jetCleaning_task.add(process.cleanPatPhotons)
1765 cleanPatPhotonProducer = getattr(process, "cleanPatPhotons").clone(
1766 src = photonCollection
1767 )
1768 cleanPatPhotonProducer.checkOverlaps.electrons.src = electronCollection
1769 addToProcessAndTask("cleanedPatPhotons"+postfix, cleanPatPhotonProducer, process, jetCleaning_task)
1770 photonCollection = cms.InputTag("cleanedPatPhotons"+postfix)
1771
1772
1773 have_cleanPatJets = hasattr(process, "cleanPatJets")
1774 process.load("PhysicsTools.PatAlgos.cleaningLayer1.jetCleaner_cfi")
1775 cleanPatJetProducer = getattr(process, "cleanPatJets").clone(
1776 src = cms.InputTag(jetCollectionName)
1777 )
1778
1779 if not have_cleanPatJets:
1780 del process.cleanPatJets
1781 cleanPatJetProducer.checkOverlaps.muons.src = muonCollection
1782 cleanPatJetProducer.checkOverlaps.electrons.src = electronCollection
1783 if isValidInputTag(photonCollection) and autoJetCleaning != "LepClean":
1784 cleanPatJetProducer.checkOverlaps.photons.src = photonCollection
1785 else:
1786 del cleanPatJetProducer.checkOverlaps.photons
1787
1788 if isValidInputTag(tauCollection) and autoJetCleaning != "LepClean":
1789 cleanPatJetProducer.checkOverlaps.taus.src = tauCollection
1790 else:
1791 del cleanPatJetProducer.checkOverlaps.taus
1792
1793
1794 del cleanPatJetProducer.checkOverlaps.tkIsoElectrons
1795
1796 addToProcessAndTask("cleanedPatJets"+postfix, cleanPatJetProducer, process, jetCleaning_task)
1797
1798 addTaskToProcess(process, jetCleaning_label, jetCleaning_task)
1799 jetProductionTask.add(getattr(process, jetCleaning_label))
1800 return cms.InputTag("cleanedPatJets"+postfix)
1801
1802
1803 def runFixEE2017(self, process, params, jets, cands, goodcolls, patMetModuleTask, postfix):
1804
1805 runFixEE2017_task, runFixEE2017_label = cms.Task(), "runFixEE2017_task{}".format(postfix)
1806
1807 pfCandidateJetsWithEEnoise = _modbad.BadPFCandidateJetsEEnoiseProducer.clone(
1808 jetsrc = jets,
1809 userawPt = params["userawPt"],
1810 ptThreshold = params["ptThreshold"],
1811 minEtaThreshold = params["minEtaThreshold"],
1812 maxEtaThreshold = params["maxEtaThreshold"],
1813 )
1814 addToProcessAndTask("pfCandidateJetsWithEEnoise"+postfix, pfCandidateJetsWithEEnoise, process, runFixEE2017_task)
1815
1816 pfcandidateClustered = cms.EDProducer("CandViewMerger",
1817 src = cms.VInputTag(goodcolls+[jets])
1818 )
1819 addToProcessAndTask("pfcandidateClustered"+postfix, pfcandidateClustered, process, runFixEE2017_task)
1820
1821 pfcandidateForUnclusteredUnc = _mod.candPtrProjector.clone(
1822 src = cands,
1823 veto = "pfcandidateClustered"+postfix,
1824 )
1825 addToProcessAndTask("pfcandidateForUnclusteredUnc"+postfix, pfcandidateForUnclusteredUnc, process, runFixEE2017_task)
1826
1827 badUnclustered = cms.EDFilter("CandPtrSelector",
1828 src = cms.InputTag("pfcandidateForUnclusteredUnc"+postfix),
1829 cut = cms.string("abs(eta) > "+str(params["minEtaThreshold"])+" && abs(eta) < "+str(params["maxEtaThreshold"])),
1830 )
1831 addToProcessAndTask("badUnclustered"+postfix, badUnclustered, process, runFixEE2017_task)
1832
1833 blobUnclustered = cms.EDProducer("UnclusteredBlobProducer",
1834 candsrc = cms.InputTag("badUnclustered"+postfix),
1835 )
1836 addToProcessAndTask("blobUnclustered"+postfix, blobUnclustered, process, runFixEE2017_task)
1837
1838 superbad = cms.EDProducer("CandViewMerger",
1839 src = cms.VInputTag(
1840 cms.InputTag("blobUnclustered"+postfix),
1841 cms.InputTag("pfCandidateJetsWithEEnoise"+postfix,"bad"),
1842 )
1843 )
1844 addToProcessAndTask("superbad"+postfix, superbad, process, runFixEE2017_task)
1845
1846 pfCandidatesGoodEE2017 = _mod.candPtrProjector.clone(
1847 src = cands,
1848 veto = "superbad"+postfix,
1849 )
1850 addToProcessAndTask("pfCandidatesGoodEE2017"+postfix, pfCandidatesGoodEE2017, process, runFixEE2017_task)
1851
1852
1853 addTaskToProcess(process, runFixEE2017_label, runFixEE2017_task)
1854
1855
1856 patMetModuleTask.add(getattr(process, runFixEE2017_label))
1857
1858
1859 return (cms.InputTag("pfCandidatesGoodEE2017"+postfix), cms.InputTag("pfCandidateJetsWithEEnoise"+postfix,"good"))
1860
1861
1862 runMETCorrectionsAndUncertainties = RunMETCorrectionsAndUncertainties()
1863
1864
1865
1866
1867
1868 def runMetCorAndUncForMiniAODProduction(process, metType="PF",
1869 jetCollUnskimmed="patJets",
1870 photonColl="selectedPatPhotons",
1871 electronColl="selectedPatElectrons",
1872 muonColl="selectedPatMuons",
1873 tauColl="selectedPatTaus",
1874 pfCandColl = "particleFlow",
1875 jetCleaning="LepClean",
1876 jetSelection="pt>15 && abs(eta)<9.9",
1877 jecUnFile="",
1878 jetFlavor="AK4PFchs",
1879 recoMetFromPFCs=False,
1880 postfix=""):
1881
1882 runMETCorrectionsAndUncertainties = RunMETCorrectionsAndUncertainties()
1883
1884
1885 runMETCorrectionsAndUncertainties(process, metType=metType,
1886 correctionLevel=["T0","T1","T2","Smear","Txy"],
1887 computeUncertainties=False,
1888 produceIntermediateCorrections=True,
1889 addToPatDefaultSequence=False,
1890 jetCollectionUnskimmed=jetCollUnskimmed,
1891 photonCollection=photonColl,
1892 electronCollection=electronColl,
1893 muonCollection=muonColl,
1894 tauCollection=tauColl,
1895 pfCandCollection =pfCandColl,
1896 autoJetCleaning=jetCleaning,
1897 jecUncertaintyFile=jecUnFile,
1898 jetSelection=jetSelection,
1899 jetFlavor=jetFlavor,
1900 recoMetFromPFCs=recoMetFromPFCs,
1901 postfix=postfix
1902 )
1903
1904
1905 runMETCorrectionsAndUncertainties(process, metType=metType,
1906 correctionLevel=["T1"],
1907 computeUncertainties=True,
1908 produceIntermediateCorrections=False,
1909 addToPatDefaultSequence=False,
1910 jetCollectionUnskimmed=jetCollUnskimmed,
1911 photonCollection=photonColl,
1912 electronCollection=electronColl,
1913 muonCollection=muonColl,
1914 tauCollection=tauColl,
1915 pfCandCollection =pfCandColl,
1916 autoJetCleaning=jetCleaning,
1917 jecUncertaintyFile=jecUnFile,
1918 jetSelection=jetSelection,
1919 jetFlavor=jetFlavor,
1920 recoMetFromPFCs=recoMetFromPFCs,
1921 postfix=postfix
1922 )
1923
1924
1925 runMETCorrectionsAndUncertainties(process, metType=metType,
1926 correctionLevel=["T1","Smear"],
1927 computeUncertainties=True,
1928 produceIntermediateCorrections=False,
1929 addToPatDefaultSequence=False,
1930 jetCollectionUnskimmed=jetCollUnskimmed,
1931 photonCollection=photonColl,
1932 electronCollection=electronColl,
1933 muonCollection=muonColl,
1934 tauCollection=tauColl,
1935 pfCandCollection =pfCandColl,
1936 autoJetCleaning=jetCleaning,
1937 jecUncertaintyFile=jecUnFile,
1938 jetSelection=jetSelection,
1939 jetFlavor=jetFlavor,
1940 recoMetFromPFCs=recoMetFromPFCs,
1941 postfix=postfix,
1942 )
1943
1944
1945
1946
1947
1948 def runMetCorAndUncFromMiniAOD(process, metType="PF",
1949 jetCollUnskimmed="slimmedJets",
1950 photonColl="slimmedPhotons",
1951 electronColl="slimmedElectrons",
1952 muonColl="slimmedMuons",
1953 tauColl="slimmedTaus",
1954 pfCandColl = "packedPFCandidates",
1955 jetFlavor="AK4PFchs",
1956 jetCleaning="LepClean",
1957 isData=False,
1958 manualJetConfig=False,
1959 reclusterJets=None,
1960 jetSelection="pt>15 && abs(eta)<9.9",
1961 recoMetFromPFCs=None,
1962 jetCorLabelL3="ak4PFCHSL1FastL2L3Corrector",
1963 jetCorLabelRes="ak4PFCHSL1FastL2L3ResidualCorrector",
1964
1965 CHS=False,
1966 puppiProducerLabel="puppi",
1967 puppiProducerForMETLabel="puppiNoLep",
1968 reapplyJEC=True,
1969 jecUncFile="",
1970 computeMETSignificance=True,
1971 fixEE2017=False,
1972 fixEE2017Params=None,
1973 extractDeepMETs=False,
1974 campaign="",
1975 era="",
1976 postfix=""):
1977
1978 runMETCorrectionsAndUncertainties = RunMETCorrectionsAndUncertainties()
1979
1980
1981 runMETCorrectionsAndUncertainties(process, metType=metType,
1982 correctionLevel=["T1"],
1983 computeUncertainties=True,
1984 produceIntermediateCorrections=False,
1985 addToPatDefaultSequence=False,
1986 jetCollectionUnskimmed=jetCollUnskimmed,
1987 electronCollection=electronColl,
1988 muonCollection=muonColl,
1989 tauCollection=tauColl,
1990 photonCollection=photonColl,
1991 pfCandCollection =pfCandColl,
1992 runOnData=isData,
1993 onMiniAOD=True,
1994 reapplyJEC=reapplyJEC,
1995 reclusterJets=reclusterJets,
1996 jetSelection=jetSelection,
1997 recoMetFromPFCs=recoMetFromPFCs,
1998 autoJetCleaning=jetCleaning,
1999 manualJetConfig=manualJetConfig,
2000 jetFlavor=jetFlavor,
2001 jetCorLabelUpToL3=jetCorLabelL3,
2002 jetCorLabelL3Res=jetCorLabelRes,
2003 jecUncertaintyFile=jecUncFile,
2004 computeMETSignificance=computeMETSignificance,
2005 CHS=CHS,
2006 puppiProducerLabel=puppiProducerLabel,
2007 puppiProducerForMETLabel=puppiProducerForMETLabel,
2008 postfix=postfix,
2009 fixEE2017=fixEE2017,
2010 fixEE2017Params=fixEE2017Params,
2011 extractDeepMETs=extractDeepMETs,
2012 campaign=campaign,
2013 era=era,
2014 )
2015
2016
2017 runMETCorrectionsAndUncertainties(process, metType=metType,
2018 correctionLevel=["T1","Txy"],
2019 computeUncertainties=False,
2020 produceIntermediateCorrections=True,
2021 addToPatDefaultSequence=False,
2022 jetCollectionUnskimmed=jetCollUnskimmed,
2023 electronCollection=electronColl,
2024 muonCollection=muonColl,
2025 tauCollection=tauColl,
2026 photonCollection=photonColl,
2027 pfCandCollection =pfCandColl,
2028 runOnData=isData,
2029 onMiniAOD=True,
2030 reapplyJEC=reapplyJEC,
2031 reclusterJets=reclusterJets,
2032 jetSelection=jetSelection,
2033 recoMetFromPFCs=recoMetFromPFCs,
2034 autoJetCleaning=jetCleaning,
2035 manualJetConfig=manualJetConfig,
2036 jetFlavor=jetFlavor,
2037 jetCorLabelUpToL3=jetCorLabelL3,
2038 jetCorLabelL3Res=jetCorLabelRes,
2039 jecUncertaintyFile=jecUncFile,
2040 computeMETSignificance=computeMETSignificance,
2041 CHS=CHS,
2042 puppiProducerLabel=puppiProducerLabel,
2043 puppiProducerForMETLabel=puppiProducerForMETLabel,
2044 postfix=postfix,
2045 fixEE2017=fixEE2017,
2046 fixEE2017Params=fixEE2017Params,
2047 extractDeepMETs=extractDeepMETs,
2048 campaign=campaign,
2049 era=era,
2050 )
2051
2052 runMETCorrectionsAndUncertainties(process, metType=metType,
2053 correctionLevel=["T1","Smear"],
2054 computeUncertainties=True,
2055 produceIntermediateCorrections=False,
2056 addToPatDefaultSequence=False,
2057 jetCollectionUnskimmed=jetCollUnskimmed,
2058 electronCollection=electronColl,
2059 muonCollection=muonColl,
2060 tauCollection=tauColl,
2061 photonCollection=photonColl,
2062 pfCandCollection =pfCandColl,
2063 runOnData=isData,
2064 onMiniAOD=True,
2065 reapplyJEC=reapplyJEC,
2066 reclusterJets=reclusterJets,
2067 jetSelection=jetSelection,
2068 recoMetFromPFCs=recoMetFromPFCs,
2069 autoJetCleaning=jetCleaning,
2070 manualJetConfig=manualJetConfig,
2071 jetFlavor=jetFlavor,
2072 jetCorLabelUpToL3=jetCorLabelL3,
2073 jetCorLabelL3Res=jetCorLabelRes,
2074 jecUncertaintyFile=jecUncFile,
2075 computeMETSignificance=computeMETSignificance,
2076 CHS=CHS,
2077 puppiProducerLabel=puppiProducerLabel,
2078 puppiProducerForMETLabel=puppiProducerForMETLabel,
2079 postfix=postfix,
2080 fixEE2017=fixEE2017,
2081 fixEE2017Params=fixEE2017Params,
2082 extractDeepMETs=extractDeepMETs,
2083 campaign=campaign,
2084 era=era,
2085 )