CMS 3D CMS Logo

Classes | Functions | Variables
jetTools Namespace Reference

Classes

class  AddJetCollection
 
class  AddJetID
 
class  SetTagInfos
 
class  SwitchJetCollection
 
class  UpdateJetCollection
 

Functions

def checkJetCorrectionsFormat (jetCorrections)
 
def deprecatedOptionOutputModule (obj)
 
def jetAlgo (algo)
 
def rerunningIVF ()
 
def rerunningIVFMiniAOD ()
 
def setupBTagging (process, jetSource, pfCandidates, explicitJTA, pvSource, svSource, elSource, muSource, runIVF, tightBTagNTkHits, loadStdRecoBTag, svClustering, fatJets, groomedFatJets, algo, rParam, btagDiscriminators, btagInfos, patJets, labelName, btagPrefix, postfix)
 
def setupJetCorrections (process, knownModules, jetCorrections, jetSource, pvSource, patJets, labelName, postfix)
 
def setupSVClustering (btagInfo, svClustering, algo, rParam, fatJets=cms.InputTag(''), groomedFatJets=cms.InputTag(''))
 
def undefinedLabelName (obj)
 
def unsupportedJetAlgorithm (obj)
 

Variables

 supportedJetAlgos
 dictionary with supported jet clustering algorithms More...
 

Function Documentation

def jetTools.checkJetCorrectionsFormat (   jetCorrections)

Definition at line 17 of file jetTools.py.

Referenced by jetTools.AddJetCollection.toolCode(), and jetTools.UpdateJetCollection.toolCode().

17 def checkJetCorrectionsFormat(jetCorrections):
18  ## check for the correct format
19  if type(jetCorrections) != type(('PAYLOAD-LABEL',['CORRECTION-LEVEL-A','CORRECTION-LEVEL-B'], 'MET-LABEL')):
20  raise ValueError("In addJetCollection: 'jetCorrections' must be 'None' (as a python value w/o quotation marks), or of type ('PAYLOAD-LABEL', ['CORRECTION-LEVEL-A', \
21  'CORRECTION-LEVEL-B', ...], 'MET-LABEL'). Note that 'MET-LABEL' can be set to 'None' (as a string in quotation marks) in case you do not want to apply MET(Type1) \
22  corrections.")
23 
24 
def checkJetCorrectionsFormat(jetCorrections)
Definition: jetTools.py:17
def jetTools.deprecatedOptionOutputModule (   obj)

Definition at line 1840 of file jetTools.py.

1841  sys.stderr.write("-------------------------------------------------------\n")
1842  sys.stderr.write(" Error: the option 'outputModule' is not supported\n")
1843  sys.stderr.write(" anymore by:\n")
1844  sys.stderr.write(" " + obj._label + "\n")
1845  sys.stderr.write(" please use 'outputModules' now and specify the\n")
1846  sys.stderr.write(" names of all needed OutModules in there\n")
1847  sys.stderr.write(" (default: ['out'])\n")
1848  sys.stderr.write("-------------------------------------------------------\n")
1849  raise KeyError("Unsupported option 'outputModule' used in '"+obj._label+"'")
1850 
def deprecatedOptionOutputModule(obj)
Definition: jetTools.py:1840
def jetTools.jetAlgo (   algo)

Definition at line 3 of file jetTools.py.

3 def jetAlgo( algo ):
4 
5  # print 'PF2PAT: selecting jet algorithm ', algo
6 
7  if algo == 'IC5':
8  import RecoJets.JetProducers.ic5PFJets_cfi as jetConfig
9  jetAlgo = jetConfig.iterativeCone5PFJets.clone()
10  elif algo == 'AK4':
11  import RecoJets.JetProducers.ak4PFJets_cfi as jetConfig
12  jetAlgo = jetConfig.ak4PFJets.clone()
13  elif algo == 'AK7':
14  import RecoJets.JetProducers.ak4PFJets_cfi as jetConfig
15  jetAlgo = jetConfig.ak4PFJets.clone()
16  jetAlgo.rParam = cms.double(0.7)
17  jetAlgo.doAreaFastjet = cms.bool(False)
18 
19  jetAlgo.src = 'pfNoElectronJME'
20  jetAlgo.doPVCorrection = True
21  jetAlgo.srcPVs = cms.InputTag("goodOfflinePrimaryVertices")
22  return jetAlgo
23 
def jetAlgo(algo)
Definition: jetTools.py:3
def jetTools.rerunningIVF ( )

Definition at line 1867 of file jetTools.py.

Referenced by setupBTagging().

1868  sys.stderr.write("-------------------------------------------------------------------\n")
1869  sys.stderr.write(" Warning: You are attempting to remake the IVF secondary vertices\n")
1870  sys.stderr.write(" already produced by the standard reconstruction. This\n")
1871  sys.stderr.write(" option is not enabled by default so please use it only if\n")
1872  sys.stderr.write(" you know what you are doing.\n")
1873  sys.stderr.write("-------------------------------------------------------------------\n")
1874 
def rerunningIVF()
Definition: jetTools.py:1867
def jetTools.rerunningIVFMiniAOD ( )

Definition at line 1875 of file jetTools.py.

Referenced by setupBTagging().

1876  sys.stderr.write("-------------------------------------------------------------------\n")
1877  sys.stderr.write(" Warning: You are attempting to remake IVF secondary vertices from\n")
1878  sys.stderr.write(" MiniAOD. If that was your intention, note that secondary\n")
1879  sys.stderr.write(" vertices remade from MiniAOD will have somewhat degraded\n")
1880  sys.stderr.write(" performance compared to those remade from RECO/AOD.\n")
1881  sys.stderr.write("-------------------------------------------------------------------\n")
1882 
def rerunningIVFMiniAOD()
Definition: jetTools.py:1875
def jetTools.setupBTagging (   process,
  jetSource,
  pfCandidates,
  explicitJTA,
  pvSource,
  svSource,
  elSource,
  muSource,
  runIVF,
  tightBTagNTkHits,
  loadStdRecoBTag,
  svClustering,
  fatJets,
  groomedFatJets,
  algo,
  rParam,
  btagDiscriminators,
  btagInfos,
  patJets,
  labelName,
  btagPrefix,
  postfix 
)

Definition at line 236 of file jetTools.py.

References helpers.addToProcessAndTask(), any(), clone(), helpers.getPatAlgosToolsTask(), list(), helpers.loadWithPrefix(), rerunningIVF(), rerunningIVFMiniAOD(), and setupSVClustering().

Referenced by setupSVClustering(), jetTools.AddJetCollection.toolCode(), and jetTools.UpdateJetCollection.toolCode().

236  algo, rParam, btagDiscriminators, btagInfos, patJets, labelName, btagPrefix, postfix):
237 
238  task = getPatAlgosToolsTask(process)
239 
240  ## expand the btagDiscriminators to remove the meta taggers and substitute the equivalent sources
241  discriminators = set(btagDiscriminators)
242  present_meta = discriminators.intersection(set(supportedMetaDiscr.keys()))
243  discriminators -= present_meta
244  for meta_tagger in present_meta:
245  for src in supportedMetaDiscr[meta_tagger]:
246  discriminators.add(src)
247  btagDiscriminators = list(discriminators)
248 
249  ## expand tagInfos to what is explicitly required by user + implicit
250  ## requirements that come in from one or the other discriminator
251  requiredTagInfos = list(btagInfos)
252  for btagDiscr in btagDiscriminators :
253  for tagInfoList in supportedBtagDiscr[btagDiscr] :
254  for requiredTagInfo in tagInfoList :
255  tagInfoCovered = False
256  for tagInfo in requiredTagInfos :
257  if requiredTagInfo == tagInfo :
258  tagInfoCovered = True
259  break
260  if not tagInfoCovered :
261  requiredTagInfos.append(requiredTagInfo)
262  ## load sequences and setups needed for btagging
263  if hasattr( process, 'candidateJetProbabilityComputer' ) == False :
264  if loadStdRecoBTag: # also loading modules already run in the standard reconstruction
265  process.load("RecoBTag.ImpactParameter.impactParameter_cff")
266  task.add(process.impactParameterTask)
267  process.load("RecoBTag.SecondaryVertex.secondaryVertex_cff")
268  task.add(process.secondaryVertexTask)
269  process.load("RecoBTag.SoftLepton.softLepton_cff")
270  task.add(process.softLeptonTask)
271  process.load("RecoBTag.Combined.combinedMVA_cff")
272  task.add(process.combinedMVATask)
273  process.load("RecoBTag.CTagging.cTagging_cff")
274  task.add(process.cTaggingTask)
275  else: # to prevent loading of modules already run in the standard reconstruction
276  process.load("RecoBTag.ImpactParameter.impactParameter_EventSetup_cff")
277  process.load("RecoBTag.SecondaryVertex.secondaryVertex_EventSetup_cff")
278  process.load("RecoBTag.SoftLepton.softLepton_EventSetup_cff")
279  process.load("RecoBTag.Combined.combinedMVA_EventSetup_cff")
280  process.load("RecoBTag.CTagging.cTagging_EventSetup_cff")
282  import RecoJets.JetProducers.caTopTaggers_cff as toptag
283 
284  if tightBTagNTkHits:
285  if not runIVF:
286  sys.stderr.write("-------------------------------------------------------------------\n")
287  sys.stderr.write(" Warning: For a complete switch to the legacy tight b-tag track\n")
288  sys.stderr.write(" selection, please also enable the \'runIVF\' switch.\n")
289  sys.stderr.write("-------------------------------------------------------------------\n")
290  if btagPrefix == '':
291  sys.stderr.write("-------------------------------------------------------------------\n")
292  sys.stderr.write(" Warning: With the tight b-tag track selection enabled, it is\n")
293  sys.stderr.write(" advisable to set \'btagPrefix\' to a non-empty string to\n")
294  sys.stderr.write(" avoid unintentional modifications to the default\n")
295  sys.stderr.write(" b tagging setup that might be loaded in the same job.\n")
296  sys.stderr.write("-------------------------------------------------------------------\n")
297 
298  ## define c tagging CvsL SV source (for now tied to the default SV source
299  ## in the first part of the module label, product instance label and process name)
300  svSourceCvsL = copy.deepcopy(svSource)
301  svSourceCvsL.setModuleLabel(svSource.getModuleLabel()+'CvsL')
302 
303  ## check if and under what conditions to re-run IVF
304  runIVFforCTagOnly = False
305  ivfcTagInfos = ['pfInclusiveSecondaryVertexFinderCvsLTagInfos', 'pfInclusiveSecondaryVertexFinderNegativeCvsLTagInfos']
306  ## if MiniAOD and running c tagging
307  if pvSource.getModuleLabel() == 'offlineSlimmedPrimaryVertices' and any(i in requiredTagInfos for i in ivfcTagInfos) and not runIVF:
308  runIVFforCTagOnly = True
309  runIVF = True
310  sys.stderr.write("-------------------------------------------------------------------\n")
311  sys.stderr.write(" Info: To run c tagging on MiniAOD, c-tag-specific IVF secondary\n")
312  sys.stderr.write(" vertices will be remade.\n")
313  sys.stderr.write("-------------------------------------------------------------------\n")
314  ## adjust svSources
315  if runIVF and btagPrefix != '':
316  if runIVFforCTagOnly:
317  svSourceCvsL.setModuleLabel(btagPrefix+svSourceCvsL.getModuleLabel())
318  else:
319  svSource.setModuleLabel(btagPrefix+svSource.getModuleLabel())
320  svSourceCvsL.setModuleLabel(btagPrefix+svSourceCvsL.getModuleLabel())
321 
322  ## setup all required btagInfos : we give a dedicated treatment for different
323  ## types of tagInfos here. A common treatment is possible but might require a more
324  ## general approach anyway in coordination with the btagging POG.
325 
326  runNegativeVertexing = False
327  runNegativeCvsLVertexing = False
328  for btagInfo in requiredTagInfos:
329  if btagInfo == 'pfInclusiveSecondaryVertexFinderNegativeTagInfos' or btagInfo == 'pfNegativeDeepFlavourTagInfos':
330  runNegativeVertexing = True
331  if btagInfo == 'pfInclusiveSecondaryVertexFinderNegativeCvsLTagInfos':
332  runNegativeCvsLVertexing = True
333 
334  if runNegativeVertexing or runNegativeCvsLVertexing:
335  import RecoVertex.AdaptiveVertexFinder.inclusiveNegativeVertexing_cff as NegVertex
336 
337  if runNegativeVertexing:
338  addToProcessAndTask(btagPrefix+'inclusiveCandidateNegativeVertexFinder'+labelName+postfix,
339  NegVertex.inclusiveCandidateNegativeVertexFinder.clone(primaryVertices = pvSource,tracks=pfCandidates),
340  process, task)
341  addToProcessAndTask(btagPrefix+'candidateNegativeVertexMerger'+labelName+postfix,
342  NegVertex.candidateNegativeVertexMerger.clone(secondaryVertices = cms.InputTag(btagPrefix+'inclusiveCandidateNegativeVertexFinder'+labelName+postfix)),
343  process, task)
344  addToProcessAndTask(btagPrefix+'candidateNegativeVertexArbitrator'+labelName+postfix,
345  NegVertex.candidateNegativeVertexArbitrator.clone( secondaryVertices = cms.InputTag(btagPrefix+'candidateNegativeVertexMerger'+labelName+postfix)
346  ,primaryVertices = pvSource
347  ,tracks=pfCandidates),
348  process, task)
349  addToProcessAndTask(btagPrefix+'inclusiveCandidateNegativeSecondaryVertices'+labelName+postfix,
350  NegVertex.inclusiveCandidateNegativeSecondaryVertices.clone(secondaryVertices = cms.InputTag(btagPrefix+'candidateNegativeVertexArbitrator'+labelName+postfix)),
351  process, task)
352 
353  if runNegativeCvsLVertexing:
354  addToProcessAndTask(btagPrefix+'inclusiveCandidateNegativeVertexFinderCvsL'+labelName+postfix,
355  NegVertex.inclusiveCandidateNegativeVertexFinderCvsL.clone(primaryVertices = pvSource,tracks=pfCandidates),
356  process, task)
357  addToProcessAndTask(btagPrefix+'candidateNegativeVertexMergerCvsL'+labelName+postfix,
358  NegVertex.candidateNegativeVertexMergerCvsL.clone(secondaryVertices = cms.InputTag(btagPrefix+'inclusiveCandidateNegativeVertexFinderCvsL'+labelName+postfix)),
359  process, task)
360  addToProcessAndTask(btagPrefix+'candidateNegativeVertexArbitratorCvsL'+labelName+postfix,
361  NegVertex.candidateNegativeVertexArbitratorCvsL.clone( secondaryVertices = cms.InputTag(btagPrefix+'candidateNegativeVertexMergerCvsL'+labelName+postfix)
362  ,primaryVertices = pvSource
363  ,tracks=pfCandidates),
364  process, task)
365  addToProcessAndTask(btagPrefix+'inclusiveCandidateNegativeSecondaryVerticesCvsL'+labelName+postfix,
366  NegVertex.inclusiveCandidateNegativeSecondaryVerticesCvsL.clone(secondaryVertices = cms.InputTag(btagPrefix+'candidateNegativeVertexArbitratorCvsL'+labelName+postfix)),
367  process, task)
368 
369 
370  acceptedTagInfos = list()
371  for btagInfo in requiredTagInfos:
372  if hasattr(btag,btagInfo):
373  if btagInfo == 'pfImpactParameterTagInfos':
374  addToProcessAndTask(btagPrefix+btagInfo+labelName+postfix,
375  btag.pfImpactParameterTagInfos.clone(jets = jetSource,primaryVertex=pvSource,candidates=pfCandidates),
376  process, task)
377  if explicitJTA:
378  _btagInfo = getattr(process, btagPrefix+btagInfo+labelName+postfix)
379  _btagInfo.explicitJTA = cms.bool(explicitJTA)
380  if tightBTagNTkHits:
381  _btagInfo = getattr(process, btagPrefix+btagInfo+labelName+postfix)
382  _btagInfo.minimumNumberOfPixelHits = cms.int32(2)
383  _btagInfo.minimumNumberOfHits = cms.int32(8)
384  if btagInfo == 'pfImpactParameterAK8TagInfos':
385  addToProcessAndTask(btagPrefix+btagInfo+labelName+postfix,
386  btag.pfImpactParameterAK8TagInfos.clone(jets = jetSource,primaryVertex=pvSource,candidates=pfCandidates),
387  process, task)
388  if explicitJTA:
389  _btagInfo = getattr(process, btagPrefix+btagInfo+labelName+postfix)
390  _btagInfo.explicitJTA = cms.bool(explicitJTA)
391  if tightBTagNTkHits:
392  _btagInfo = getattr(process, btagPrefix+btagInfo+labelName+postfix)
393  _btagInfo.minimumNumberOfPixelHits = cms.int32(2)
394  _btagInfo.minimumNumberOfHits = cms.int32(8)
395  if btagInfo == 'pfImpactParameterCA15TagInfos':
396  addToProcessAndTask(btagPrefix+btagInfo+labelName+postfix,
397  btag.pfImpactParameterCA15TagInfos.clone(jets = jetSource,primaryVertex=pvSource,candidates=pfCandidates),
398  process, task)
399  if explicitJTA:
400  _btagInfo = getattr(process, btagPrefix+btagInfo+labelName+postfix)
401  _btagInfo.explicitJTA = cms.bool(explicitJTA)
402  if tightBTagNTkHits:
403  _btagInfo = getattr(process, btagPrefix+btagInfo+labelName+postfix)
404  _btagInfo.minimumNumberOfPixelHits = cms.int32(2)
405  _btagInfo.minimumNumberOfHits = cms.int32(8)
406  if btagInfo == 'pfSecondaryVertexTagInfos':
407  addToProcessAndTask(btagPrefix+btagInfo+labelName+postfix,
408  btag.pfSecondaryVertexTagInfos.clone(
409  trackIPTagInfos = cms.InputTag(btagPrefix+'pfImpactParameterTagInfos'+labelName+postfix)),
410  process, task)
411  if tightBTagNTkHits:
412  _btagInfo = getattr(process, btagPrefix+btagInfo+labelName+postfix)
413  _btagInfo.trackSelection.pixelHitsMin = cms.uint32(2)
414  _btagInfo.trackSelection.totalHitsMin = cms.uint32(8)
415  if btagInfo == 'pfDeepCSVTagInfos':
416  addToProcessAndTask(btagPrefix+btagInfo+labelName+postfix,
417  btag.pfDeepCSVTagInfos.clone(
418  svTagInfos = cms.InputTag(btagPrefix+'pfInclusiveSecondaryVertexFinderTagInfos'+labelName+postfix)),
419  process, task)
420  if svClustering or fatJets != cms.InputTag(''):
421  setupSVClustering(getattr(process, btagPrefix+btagInfo+labelName+postfix), svClustering, algo, rParam, fatJets, groomedFatJets)
422  if btagInfo == 'pfDeepCSVNegativeTagInfos':
423  addToProcessAndTask(btagPrefix+btagInfo+labelName+postfix,
424  btag.pfDeepCSVNegativeTagInfos.clone(
425  svTagInfos = cms.InputTag(btagPrefix+'pfInclusiveSecondaryVertexFinderNegativeTagInfos'+labelName+postfix)),
426  process, task)
427  if svClustering or fatJets != cms.InputTag(''):
428  setupSVClustering(getattr(process, btagPrefix+btagInfo+labelName+postfix), svClustering, algo, rParam, fatJets, groomedFatJets)
429  if btagInfo == 'pfDeepCSVPositiveTagInfos':
430  addToProcessAndTask(btagPrefix+btagInfo+labelName+postfix,
431  btag.pfDeepCSVPositiveTagInfos.clone(
432  svTagInfos = cms.InputTag(btagPrefix+'pfInclusiveSecondaryVertexFinderTagInfos'+labelName+postfix)),
433  process, task)
434  if svClustering or fatJets != cms.InputTag(''):
435  setupSVClustering(getattr(process, btagPrefix+btagInfo+labelName+postfix), svClustering, algo, rParam, fatJets, groomedFatJets)
436  if btagInfo == 'pfDeepCMVATagInfos':
437  addToProcessAndTask(btagPrefix+btagInfo+labelName+postfix,
438  btag.pfDeepCMVATagInfos.clone(
439  deepNNTagInfos = cms.InputTag(btagPrefix+'pfDeepCSVTagInfos'+labelName+postfix),
440  ipInfoSrc = cms.InputTag(btagPrefix+"pfImpactParameterTagInfos"+labelName+postfix),
441  muInfoSrc = cms.InputTag(btagPrefix+"softPFMuonsTagInfos"+labelName+postfix),
442  elInfoSrc = cms.InputTag(btagPrefix+"softPFElectronsTagInfos"+labelName+postfix)),
443  process, task)
444  if svClustering or fatJets != cms.InputTag(''):
445  setupSVClustering(getattr(process, btagPrefix+btagInfo+labelName+postfix), svClustering, algo, rParam, fatJets, groomedFatJets)
446  if btagInfo == 'pfDeepCMVANegativeTagInfos':
447  addToProcessAndTask(btagPrefix+btagInfo+labelName+postfix,
448  btag.pfDeepCMVATagInfos.clone(
449  deepNNTagInfos = cms.InputTag(btagPrefix+'pfDeepCSVTagInfos'+labelName+postfix),
450  ipInfoSrc = cms.InputTag(btagPrefix+"pfImpactParameterTagInfos"+labelName+postfix),
451  muInfoSrc = cms.InputTag(btagPrefix+"softPFMuonsTagInfos"+labelName+postfix),
452  elInfoSrc = cms.InputTag(btagPrefix+"softPFElectronsTagInfos"+labelName+postfix)),
453  process, task)
454  if svClustering or fatJets != cms.InputTag(''):
455  setupSVClustering(getattr(process, btagPrefix+btagInfo+labelName+postfix), svClustering, algo, rParam, fatJets, groomedFatJets)
456  if btagInfo == 'pfDeepCMVAPositiveTagInfos':
457  addToProcessAndTask(btagPrefix+btagInfo+labelName+postfix,
458  btag.pfDeepCMVATagInfos.clone(
459  deepNNTagInfos = cms.InputTag(btagPrefix+'pfDeepCSVTagInfos'+labelName+postfix),
460  ipInfoSrc = cms.InputTag(btagPrefix+"pfImpactParameterTagInfos"+labelName+postfix),
461  muInfoSrc = cms.InputTag(btagPrefix+"softPFMuonsTagInfos"+labelName+postfix),
462  elInfoSrc = cms.InputTag(btagPrefix+"softPFElectronsTagInfos"+labelName+postfix)),
463  process, task)
464  if svClustering or fatJets != cms.InputTag(''):
465  setupSVClustering(getattr(process, btagPrefix+btagInfo+labelName+postfix), svClustering, algo, rParam, fatJets, groomedFatJets)
466  if btagInfo == 'pfInclusiveSecondaryVertexFinderTagInfos':
467  addToProcessAndTask(btagPrefix+btagInfo+labelName+postfix,
468  btag.pfInclusiveSecondaryVertexFinderTagInfos.clone(
469  trackIPTagInfos = cms.InputTag(btagPrefix+'pfImpactParameterTagInfos'+labelName+postfix),
470  extSVCollection=svSource),
471  process, task)
472  if svClustering or fatJets != cms.InputTag(''):
473  setupSVClustering(getattr(process, btagPrefix+btagInfo+labelName+postfix), svClustering, algo, rParam, fatJets, groomedFatJets)
474  if btagInfo == 'pfInclusiveSecondaryVertexFinderAK8TagInfos':
475  addToProcessAndTask(btagPrefix+btagInfo+labelName+postfix,
476  btag.pfInclusiveSecondaryVertexFinderAK8TagInfos.clone(
477  trackIPTagInfos = cms.InputTag(btagPrefix+'pfImpactParameterAK8TagInfos'+labelName+postfix),
478  extSVCollection=svSource),
479  process, task)
480  if svClustering or fatJets != cms.InputTag(''):
481  setupSVClustering(getattr(process, btagPrefix+btagInfo+labelName+postfix), svClustering, algo, rParam, fatJets, groomedFatJets)
482  if btagInfo == 'pfBoostedDoubleSVAK8TagInfos':
483  addToProcessAndTask(btagPrefix+btagInfo+labelName+postfix,
484  btag.pfBoostedDoubleSVAK8TagInfos.clone(
485  svTagInfos = cms.InputTag(btagPrefix+'pfInclusiveSecondaryVertexFinderAK8TagInfos'+labelName+postfix)),
486  process, task)
487  if btagInfo == 'pfInclusiveSecondaryVertexFinderCA15TagInfos':
488  addToProcessAndTask(btagPrefix+btagInfo+labelName+postfix,
489  btag.pfInclusiveSecondaryVertexFinderCA15TagInfos.clone(
490  trackIPTagInfos = cms.InputTag(btagPrefix+'pfImpactParameterCA15TagInfos'+labelName+postfix),
491  extSVCollection=svSource),
492  process, task)
493  if svClustering or fatJets != cms.InputTag(''):
494  setupSVClustering(getattr(process, btagPrefix+btagInfo+labelName+postfix), svClustering, algo, rParam, fatJets, groomedFatJets)
495  if btagInfo == 'pfBoostedDoubleSVCA15TagInfos':
496  addToProcessAndTask(btagPrefix+btagInfo+labelName+postfix,
497  btag.pfBoostedDoubleSVCA15TagInfos.clone(
498  svTagInfos = cms.InputTag(btagPrefix+'pfInclusiveSecondaryVertexFinderCA15TagInfos'+labelName+postfix)),
499  process, task)
500  if btagInfo == 'pfInclusiveSecondaryVertexFinderCvsLTagInfos':
501  addToProcessAndTask(btagPrefix+btagInfo+labelName+postfix,
502  btag.pfInclusiveSecondaryVertexFinderCvsLTagInfos.clone(
503  trackIPTagInfos = cms.InputTag(btagPrefix+'pfImpactParameterTagInfos'+labelName+postfix),
504  extSVCollection=svSourceCvsL),
505  process, task)
506  if svClustering or fatJets != cms.InputTag(''):
507  setupSVClustering(getattr(process, btagPrefix+btagInfo+labelName+postfix), svClustering, algo, rParam, fatJets, groomedFatJets)
508  if btagInfo == 'pfInclusiveSecondaryVertexFinderNegativeCvsLTagInfos':
509  addToProcessAndTask(btagPrefix+btagInfo+labelName+postfix,
510  btag.pfInclusiveSecondaryVertexFinderNegativeCvsLTagInfos.clone(
511  trackIPTagInfos = cms.InputTag(btagPrefix+'pfImpactParameterTagInfos'+labelName+postfix),
512  extSVCollection = btagPrefix+'inclusiveCandidateNegativeSecondaryVerticesCvsL'+labelName+postfix),
513  process, task)
514  if svClustering or fatJets != cms.InputTag(''):
515  setupSVClustering(getattr(process, btagPrefix+btagInfo+labelName+postfix), svClustering, algo, rParam, fatJets, groomedFatJets)
516  if btagInfo == 'pfGhostTrackVertexTagInfos':
517  addToProcessAndTask(btagPrefix+btagInfo+labelName+postfix,
518  btag.pfGhostTrackVertexTagInfos.clone(
519  trackIPTagInfos = cms.InputTag(btagPrefix+'pfImpactParameterTagInfos'+labelName+postfix)),
520  process, task)
521  if btagInfo == 'pfSecondaryVertexNegativeTagInfos':
522  addToProcessAndTask(btagPrefix+btagInfo+labelName+postfix,
523  btag.pfSecondaryVertexNegativeTagInfos.clone(
524  trackIPTagInfos = cms.InputTag(btagPrefix+'pfImpactParameterTagInfos'+labelName+postfix)),
525  process, task)
526  if tightBTagNTkHits:
527  _btagInfo = getattr(process, btagPrefix+btagInfo+labelName+postfix)
528  _btagInfo.trackSelection.pixelHitsMin = cms.uint32(2)
529  _btagInfo.trackSelection.totalHitsMin = cms.uint32(8)
530  if btagInfo == 'pfInclusiveSecondaryVertexFinderNegativeTagInfos':
531  addToProcessAndTask(btagPrefix+btagInfo+labelName+postfix,
532  btag.pfInclusiveSecondaryVertexFinderNegativeTagInfos.clone(
533  trackIPTagInfos = cms.InputTag(btagPrefix+'pfImpactParameterTagInfos'+labelName+postfix),
534  extSVCollection=cms.InputTag(btagPrefix+'inclusiveCandidateNegativeSecondaryVertices'+labelName+postfix)),
535  process, task)
536  if svClustering or fatJets != cms.InputTag(''):
537  setupSVClustering(getattr(process, btagPrefix+btagInfo+labelName+postfix), svClustering, algo, rParam, fatJets, groomedFatJets)
538  if btagInfo == 'impactParameterTagInfos':
539  addToProcessAndTask(btagPrefix+btagInfo+labelName+postfix,
540  btag.impactParameterTagInfos.clone(
541  jetTracks = cms.InputTag('jetTracksAssociatorAtVertex'+labelName+postfix),
542  primaryVertex=pvSource),
543  process, task)
544  if btagInfo == 'secondaryVertexTagInfos':
545  addToProcessAndTask(btagPrefix+btagInfo+labelName+postfix,
546  btag.secondaryVertexTagInfos.clone(
547  trackIPTagInfos = cms.InputTag(btagPrefix+'impactParameterTagInfos'+labelName+postfix)),
548  process, task)
549  if btagInfo == 'inclusiveSecondaryVertexFinderTagInfos':
550  addToProcessAndTask(btagPrefix+btagInfo+labelName+postfix,
551  btag.inclusiveSecondaryVertexFinderTagInfos.clone(
552  trackIPTagInfos = cms.InputTag(btagPrefix+'impactParameterTagInfos'+labelName+postfix)),
553  process, task)
554  if svClustering or fatJets != cms.InputTag(''):
555  setupSVClustering(getattr(process, btagPrefix+btagInfo+labelName+postfix), svClustering, algo, rParam, fatJets, groomedFatJets)
556  if btagInfo == 'inclusiveSecondaryVertexFinderFilteredTagInfos':
557  addToProcessAndTask(btagPrefix+btagInfo+labelName+postfix,
558  btag.inclusiveSecondaryVertexFinderFilteredTagInfos.clone(
559  trackIPTagInfos = cms.InputTag(btagPrefix+'impactParameterTagInfos'+labelName+postfix)),
560  process, task)
561  if svClustering or fatJets != cms.InputTag(''):
562  setupSVClustering(getattr(process, btagPrefix+btagInfo+labelName+postfix), svClustering, algo, rParam, fatJets, groomedFatJets)
563  if btagInfo == 'secondaryVertexNegativeTagInfos':
564  addToProcessAndTask(btagPrefix+btagInfo+labelName+postfix,
565  btag.secondaryVertexNegativeTagInfos.clone(
566  trackIPTagInfos = cms.InputTag(btagPrefix+'impactParameterTagInfos'+labelName+postfix)),
567  process, task)
568  if btagInfo == 'inclusiveSecondaryVertexFinderNegativeTagInfos':
569  addToProcessAndTask(btagPrefix+btagInfo+labelName+postfix,
570  btag.inclusiveSecondaryVertexFinderNegativeTagInfos.clone(
571  trackIPTagInfos = cms.InputTag(btagPrefix+'impactParameterTagInfos'+labelName+postfix)),
572  process, task)
573  if svClustering or fatJets != cms.InputTag(''):
574  setupSVClustering(getattr(process, btagPrefix+btagInfo+labelName+postfix), svClustering, algo, rParam, fatJets, groomedFatJets)
575  if btagInfo == 'inclusiveSecondaryVertexFinderFilteredNegativeTagInfos':
576  addToProcessAndTask(btagPrefix+btagInfo+labelName+postfix,
577  btag.inclusiveSecondaryVertexFinderFilteredNegativeTagInfos.clone(
578  trackIPTagInfos = cms.InputTag(btagPrefix+'impactParameterTagInfos'+labelName+postfix)),
579  process, task)
580  if svClustering or fatJets != cms.InputTag(''):
581  setupSVClustering(getattr(process, btagPrefix+btagInfo+labelName+postfix), svClustering, algo, rParam, fatJets, groomedFatJets)
582  if btagInfo == 'softMuonTagInfos':
583  addToProcessAndTask(btagPrefix+btagInfo+labelName+postfix,
584  btag.softMuonTagInfos.clone(jets = jetSource, primaryVertex=pvSource),
585  process, task)
586  if btagInfo == 'softPFMuonsTagInfos':
587  addToProcessAndTask(btagPrefix+btagInfo+labelName+postfix,
588  btag.softPFMuonsTagInfos.clone(jets = jetSource, primaryVertex=pvSource, muons=muSource),
589  process, task)
590  if btagInfo == 'softPFElectronsTagInfos':
591  addToProcessAndTask(btagPrefix+btagInfo+labelName+postfix,
592  btag.softPFElectronsTagInfos.clone(jets = jetSource, primaryVertex=pvSource, electrons=elSource),
593  process, task)
594 
595  if 'DeepFlavourTagInfos' in btagInfo:
596  svUsed = svSource
597  if btagInfo == 'pfNegativeDeepFlavourTagInfos':
598  deep_csv_tag_infos = 'pfDeepCSVNegativeTagInfos'
599  svUsed = cms.InputTag(btagPrefix+'inclusiveCandidateNegativeSecondaryVertices'+labelName+postfix)
600  flip = True
601  else:
602  deep_csv_tag_infos = 'pfDeepCSVTagInfos'
603  flip = False
604  # use right input tags when running with RECO PF candidates, which actually
605  # depens of wether jets were slimmed or not (check for s/S-limmed in name)
606  if not ('limmed' in jetSource.value()):
607  puppi_value_map = cms.InputTag("puppi")
608  vertex_associator = cms.InputTag("primaryVertexAssociation","original")
609  else:
610  puppi_value_map = cms.InputTag("")
611  vertex_associator = cms.InputTag("")
612  addToProcessAndTask(btagPrefix+btagInfo+labelName+postfix,
613  btag.pfDeepFlavourTagInfos.clone(
614  jets = jetSource,
615  vertices=pvSource,
616  secondary_vertices=svUsed,
617  shallow_tag_infos = cms.InputTag(btagPrefix+deep_csv_tag_infos+labelName+postfix),
618  puppi_value_map = puppi_value_map,
619  vertex_associator = vertex_associator,
620  flip = flip),
621  process, task)
622 
623  if btagInfo == 'pfDeepDoubleXTagInfos':
624  addToProcessAndTask(btagPrefix+btagInfo+labelName+postfix,
625  btag.pfDeepDoubleXTagInfos.clone(
626  jets = jetSource,
627  vertices=pvSource,
628  secondary_vertices=svSource,
629  shallow_tag_infos = cms.InputTag(btagPrefix+'pfBoostedDoubleSVAK8TagInfos'+labelName+postfix),
630  ),
631  process, task)
632 
633  if btagInfo == 'pfDeepDoubleXTagInfosNopt':
634  addToProcessAndTask(btagPrefix+btagInfo+labelName+postfix,
635  btag.pfDeepDoubleXTagInfosNopt.clone(
636  jets = jetSource,
637  vertices=pvSource,
638  secondary_vertices=svSource,
639  shallow_tag_infos = cms.InputTag(btagPrefix+'pfBoostedDoubleSVAK8TagInfos'+labelName+postfix),
640  ),
641  process, task)
642 
643  if btagInfo == 'pfDeepBoostedJetTagInfos':
644  if pfCandidates.value() == 'packedPFCandidates':
645  # case 1: running over jets whose daughters are PackedCandidates (only via updateJetCollection for now)
646  jetSrcName = jetSource.value().lower()
647  if 'updated' in jetSrcName:
648  puppi_value_map = ""
649  vertex_associator = ""
650  fix_daughters = False
651  if 'withpuppidaughter' in jetSrcName:
652  # special case for Puppi jets reclustered from MiniAOD by analyzers
653  # need to specify 'WithPuppiDaughters' in the postfix when calling updateJetCollection
654  # daughters of these jets are already scaled by their puppi weights
655  has_puppi_weighted_daughters = True
656  else:
657  # default case for updating jet collection stored in MiniAOD, e.g., slimmedJetsAK8
658  # daughters are links to the original PackedCandidates, so NOT scaled by their puppi weights yet
659  has_puppi_weighted_daughters = False
660  if 'slimmed' in jetSrcName:
661  # when adding DeepBoostedJetTag to slimmedJetsAK8 in the MiniAOD step
662  fix_daughters = True
663  else:
664  raise ValueError("Invalid jet collection: %s. pfDeepBoostedJetTagInfos only supports running via updateJetCollection." % jetSource.value())
665  elif pfCandidates.value() == 'particleFlow':
666  raise ValueError("Running pfDeepBoostedJetTagInfos with reco::PFCandidates is currently not supported.")
667  # case 2: running on new jet collection whose daughters are PFCandidates (e.g., cluster jets in RECO/AOD)
668  # daughters are the particles used in jet clustering, so already scaled by their puppi weights
669  # Uncomment the lines below after running pfDeepBoostedJetTagInfos with reco::PFCandidates becomes supported
670 # has_puppi_weighted_daughters = True
671 # fix_daughters = False
672 # puppi_value_map = "puppi"
673 # vertex_associator = "primaryVertexAssociation:original"
674  else:
675  raise ValueError("Invalid pfCandidates collection: %s." % pfCandidates.value())
676  addToProcessAndTask(btagPrefix+btagInfo+labelName+postfix,
677  btag.pfDeepBoostedJetTagInfos.clone(
678  jets = jetSource,
679  vertices = pvSource,
680  secondary_vertices = svSource,
681  has_puppi_weighted_daughters = has_puppi_weighted_daughters,
682  fix_daughters = fix_daughters,
683  puppi_value_map = puppi_value_map,
684  vertex_associator = vertex_associator,
685  ),
686  process, task)
687 
688  acceptedTagInfos.append(btagInfo)
689  elif hasattr(toptag, btagInfo) :
690  acceptedTagInfos.append(btagInfo)
691  else:
692  print ' --> %s ignored, since not available via RecoBTag.Configuration.RecoBTag_cff!'%(btagInfo)
693  ## setup all required btagDiscriminators
694  acceptedBtagDiscriminators = list()
695  for discriminator_name in btagDiscriminators :
696  btagDiscr = discriminator_name.split(':')[0] #split input tag to get the producer label
697  #print discriminator_name, '-->', btagDiscr
698  if hasattr(btag,btagDiscr):
699  newDiscr = btagPrefix+btagDiscr+labelName+postfix #new discriminator name
700  if hasattr(process, newDiscr):
701  pass
702  elif hasattr(getattr(btag, btagDiscr), 'tagInfos'):
704  newDiscr,
705  getattr(btag, btagDiscr).clone(
706  tagInfos = cms.VInputTag(
707  *[ cms.InputTag(btagPrefix+x+labelName+postfix) \
708  for x in supportedBtagDiscr[discriminator_name][0] ]
709  )
710  ),
711  process,
712  task
713  )
714  elif hasattr(getattr(btag, btagDiscr), 'src'):
716  newDiscr,
717  getattr(btag, btagDiscr).clone(
718  src = cms.InputTag(btagPrefix+supportedBtagDiscr[discriminator_name][0][0]+labelName+postfix)
719  ),
720  process,
721  task
722  )
723  else:
724  raise ValueError('I do not know how to update %s it does not have neither "tagInfos" nor "src" attributes' % btagDiscr)
725  acceptedBtagDiscriminators.append(discriminator_name)
726  else:
727  print ' --> %s ignored, since not available via RecoBTag.Configuration.RecoBTag_cff!'%(btagDiscr)
728  #update meta-taggers, if any
729  for meta_tagger in present_meta:
730  btagDiscr = meta_tagger.split(':')[0] #split input tag to get the producer label
731  #print discriminator_name, '-->', btagDiscr
732  if hasattr(btag,btagDiscr):
733  newDiscr = btagPrefix+btagDiscr+labelName+postfix #new discriminator name
734  if hasattr(process, newDiscr):
735  pass
736  else:
738  newDiscr,
739  getattr(btag, btagDiscr).clone(),
740  process,
741  task
742  )
743  for dependency in supportedMetaDiscr[meta_tagger]:
744  if ':' in dependency:
745  new_dep = btagPrefix+dependency.split(':')[0]+labelName+postfix+':'+dependency.split(':')[1]
746  else:
747  new_dep = btagPrefix+dependency+labelName+postfix
748  replace = MassSearchReplaceAnyInputTagVisitor(dependency, new_dep)
749  replace.doIt(getattr(process, newDiscr), newDiscr)
750  acceptedBtagDiscriminators.append(meta_tagger)
751  else:
752  print ' --> %s ignored, since not available via RecoBTag.Configuration.RecoBTag_cff!'%(btagDiscr)
753 
754  ## replace corresponding tags for pat jet production
755  patJets.tagInfoSources = cms.VInputTag( *[ cms.InputTag(btagPrefix+x+labelName+postfix) for x in acceptedTagInfos ] )
756  patJets.discriminatorSources = cms.VInputTag(*[
757  cms.InputTag(btagPrefix+x+labelName+postfix) \
758  if ':' not in x else \
759  cms.InputTag(btagPrefix+x.split(':')[0]+labelName+postfix+':'+x.split(':')[1]) \
760  for x in acceptedBtagDiscriminators
761  ])
762  if len(acceptedBtagDiscriminators) > 0 :
763  patJets.addBTagInfo = True
764  ## if re-running IVF
765  if runIVF:
766  if not tightBTagNTkHits:
767  if pvSource.getModuleLabel() == 'offlineSlimmedPrimaryVertices': ## MiniAOD case
768  if not runIVFforCTagOnly: rerunningIVFMiniAOD()
769  else:
770  rerunningIVF()
771  from PhysicsTools.PatAlgos.tools.helpers import loadWithPrefix
772  ivfbTagInfos = ['pfInclusiveSecondaryVertexFinderTagInfos', 'pfInclusiveSecondaryVertexFinderAK8TagInfos', 'pfInclusiveSecondaryVertexFinderCA15TagInfos']
773  if any(i in acceptedTagInfos for i in ivfbTagInfos) and not runIVFforCTagOnly:
774  if not hasattr( process, btagPrefix+'inclusiveCandidateVertexFinder' ):
775  loadWithPrefix(process, 'RecoVertex.AdaptiveVertexFinder.inclusiveVertexing_cff', btagPrefix, task.label())
776  if tightBTagNTkHits:
777  if hasattr( process, btagPrefix+'inclusiveCandidateVertexFinder' ):
778  _temp = getattr(process, btagPrefix+'inclusiveCandidateVertexFinder')
779  _temp.minHits = cms.uint32(8)
780  ## MiniAOD case
781  if pvSource.getModuleLabel() == 'offlineSlimmedPrimaryVertices':
782  if hasattr( process, btagPrefix+'inclusiveCandidateVertexFinder' ):
783  _temp = getattr(process, btagPrefix+'inclusiveCandidateVertexFinder')
784  _temp.primaryVertices = pvSource
785  _temp.tracks = pfCandidates
786  if hasattr( process, btagPrefix+'candidateVertexArbitrator' ):
787  _temp = getattr(process, btagPrefix+'candidateVertexArbitrator')
788  _temp.primaryVertices = pvSource
789  _temp.tracks = pfCandidates
790  if hasattr( process, btagPrefix+'inclusiveCandidateSecondaryVertices' ) and not hasattr( process, svSource.getModuleLabel() ):
791  addToProcessAndTask(svSource.getModuleLabel(),
792  getattr(process, btagPrefix+'inclusiveCandidateSecondaryVertices').clone(),
793  process, task)
794  if any(i in acceptedTagInfos for i in ivfcTagInfos):
795  if not hasattr( process, btagPrefix+'inclusiveCandidateVertexFinderCvsL' ):
796  loadWithPrefix(process, 'RecoVertex.AdaptiveVertexFinder.inclusiveVertexing_cff', btagPrefix, task.label())
797  if tightBTagNTkHits:
798  if hasattr( process, btagPrefix+'inclusiveCandidateVertexFinderCvsL' ):
799  _temp = getattr(process, btagPrefix+'inclusiveCandidateVertexFinderCvsL')
800  _temp.minHits = cms.uint32(8)
801  ## MiniAOD case
802  if pvSource.getModuleLabel() == 'offlineSlimmedPrimaryVertices':
803  if hasattr( process, btagPrefix+'inclusiveCandidateVertexFinderCvsL' ):
804  _temp = getattr(process, btagPrefix+'inclusiveCandidateVertexFinderCvsL')
805  _temp.primaryVertices = pvSource
806  _temp.tracks = pfCandidates
807  if hasattr( process, btagPrefix+'candidateVertexArbitratorCvsL' ):
808  _temp = getattr(process, btagPrefix+'candidateVertexArbitratorCvsL')
809  _temp.primaryVertices = pvSource
810  _temp.tracks = pfCandidates
811  if hasattr( process, btagPrefix+'inclusiveCandidateSecondaryVerticesCvsL' ) and not hasattr( process, svSourceCvsL.getModuleLabel() ):
812  addToProcessAndTask(svSourceCvsL.getModuleLabel(),
813  getattr(process, btagPrefix+'inclusiveCandidateSecondaryVerticesCvsL').clone(),
814  process, task)
815  if 'inclusiveSecondaryVertexFinderTagInfos' in acceptedTagInfos:
816  if not hasattr( process, 'inclusiveVertexing' ):
817  process.load( 'RecoVertex.AdaptiveVertexFinder.inclusiveVertexing_cff' )
818  task.add(process.inclusiveVertexingTask)
819  task.add(process.inclusiveCandidateVertexingTask)
820  task.add(process.inclusiveCandidateVertexingCvsLTask)
821  if 'inclusiveSecondaryVertexFinderFilteredTagInfos' in acceptedTagInfos:
822  if not hasattr( process, 'inclusiveVertexing' ):
823  process.load( 'RecoVertex.AdaptiveVertexFinder.inclusiveVertexing_cff' )
824  task.add(process.inclusiveVertexingTask)
825  task.add(process.inclusiveCandidateVertexingTask)
826  task.add(process.inclusiveCandidateVertexingCvsLTask)
827  if 'inclusiveSecondaryVertexFinderFilteredTagInfos' in acceptedTagInfos:
828  if not hasattr( process, 'inclusiveSecondaryVerticesFiltered' ):
829  process.load( 'RecoBTag.SecondaryVertex.inclusiveSecondaryVerticesFiltered_cfi' )
830  task.add(process.inclusiveSecondaryVerticesFiltered)
831  task.add(process.bVertexFilter)
832  if not hasattr( process, 'bToCharmDecayVertexMerged' ):
833  process.load( 'RecoBTag.SecondaryVertex.bToCharmDecayVertexMerger_cfi' )
834  task.add(process.bToCharmDecayVertexMerged)
835  if 'caTopTagInfos' in acceptedTagInfos :
836  patJets.addTagInfos = True
837  if not hasattr( process, 'caTopTagInfos' ) and not hasattr( process, 'caTopTagInfosAK8' ):
838  process.load( 'RecoJets.JetProducers.caTopTaggers_cff' )
839  task.add(process.caTopTaggersTask)
840 
bool any(const std::vector< T > &v, const T &what)
Definition: ECalSD.cc:37
def addToProcessAndTask(label, module, process, task)
Definition: helpers.py:27
def setupSVClustering(btagInfo, svClustering, algo, rParam, fatJets=cms.InputTag(''), groomedFatJets=cms.InputTag(''))
Definition: jetTools.py:224
def loadWithPrefix(process, moduleName, prefix='', loadedProducersAndFilters=None)
Definition: helpers.py:42
def rerunningIVF()
Definition: jetTools.py:1867
TEveGeoShape * clone(const TEveElement *element, TEveElement *parent)
Definition: eve_macros.cc:135
def rerunningIVFMiniAOD()
Definition: jetTools.py:1875
def getPatAlgosToolsTask(process)
Definition: helpers.py:12
How EventSelector::AcceptEvent() decides whether to accept an event for output otherwise it is excluding the probing of A single or multiple positive and the trigger will pass if any such matching triggers are PASS or EXCEPTION[A criterion thatmatches no triggers at all is detected and causes a throw.] A single negative with an expectation of appropriate bit checking in the decision and the trigger will pass if any such matching triggers are FAIL or EXCEPTION A wildcarded negative criterion that matches more than one trigger in the trigger list("!*","!HLTx*"if it matches 2 triggers or more) will accept the event if all the matching triggers are FAIL.It will reject the event if any of the triggers are PASS or EXCEPTION(this matches the behavior of"!*"before the partial wildcard feature was incorporated).Triggers which are in the READY state are completely ignored.(READY should never be returned since the trigger paths have been run
def jetTools.setupJetCorrections (   process,
  knownModules,
  jetCorrections,
  jetSource,
  pvSource,
  patJets,
  labelName,
  postfix 
)

Definition at line 25 of file jetTools.py.

References helpers.addToProcessAndTask(), clone(), KineDebug3.count(), helpers.getPatAlgosToolsTask(), and split.

Referenced by jetTools.AddJetCollection.toolCode(), and jetTools.UpdateJetCollection.toolCode().

25 def setupJetCorrections(process, knownModules, jetCorrections, jetSource, pvSource, patJets, labelName, postfix):
26 
27  task = getPatAlgosToolsTask(process)
28 
29  ## determine type of jet constituents from jetSource; supported
30  ## jet constituent types are calo, pf, jpt, for pf also particleflow
31  ## is aloowed as part of the jetSource label, which might be used
32  ## in CommonTools.ParticleFlow
33  _type="NONE"
34  if jetCorrections[0].count('PF')>0:
35  _type='PF'
36  elif jetCorrections[0].count('Calo')>0:
37  _type='Calo'
38  elif jetCorrections[0].count('JPT')>0:
39  _type='JPT'
40  else:
41  raise TypeError("In addJetCollection: Jet energy corrections are only supported for PF, JPT and Calo jets.")
42  from PhysicsTools.PatAlgos.recoLayer0.jetCorrFactors_cfi import patJetCorrFactors
43  if 'patJetCorrFactors'+labelName+postfix in knownModules :
44  _newPatJetCorrFactors=getattr(process, 'patJetCorrFactors'+labelName+postfix)
45  _newPatJetCorrFactors.src=jetSource
46  _newPatJetCorrFactors.primaryVertices=pvSource
47  else:
48  addToProcessAndTask('patJetCorrFactors'+labelName+postfix,
49  patJetCorrFactors.clone(src=jetSource, primaryVertices=pvSource),
50  process, task)
51  _newPatJetCorrFactors=getattr(process, "patJetCorrFactors"+labelName+postfix)
52  _newPatJetCorrFactors.payload=jetCorrections[0]
53  _newPatJetCorrFactors.levels=jetCorrections[1]
54  ## check whether L1Offset or L1FastJet is part of levels
55  error=False
56  for x in jetCorrections[1]:
57  if x == 'L1Offset' :
58  if not error :
59  _newPatJetCorrFactors.useNPV=True
60  _newPatJetCorrFactors.primaryVertices='offlinePrimaryVertices'
61  _newPatJetCorrFactors.useRho=False
62  ## we set this to True now as a L1 correction type should appear only once
63  ## otherwise levels is miss configured
64  error=True
65  else:
66  raise ValueError, "In addJetCollection: Correction levels for jet energy corrections are miss configured. An L1 correction type should appear not more than \
67  once. Check the list of correction levels you requested to be applied: ", jetCorrections[1]
68  if x == 'L1FastJet' :
69  if not error :
70  if _type == "JPT" :
71  raise TypeError("In addJetCollection: L1FastJet corrections are only supported for PF and Calo jets.")
72  ## configure module
73  _newPatJetCorrFactors.useRho=True
74  if "PF" in _type :
75  _newPatJetCorrFactors.rho=cms.InputTag('fixedGridRhoFastjetAll')
76  else :
77  _newPatJetCorrFactors.rho=cms.InputTag('fixedGridRhoFastjetAllCalo')
78  ## we set this to True now as a L1 correction type should appear only once
79  ## otherwise levels is miss configured
80  error=True
81  else:
82  raise ValueError, "In addJetCollection: Correction levels for jet energy corrections are miss configured. An L1 correction type should appear not more than \
83  once. Check the list of correction levels you requested to be applied: ", jetCorrections[1]
84  patJets.jetCorrFactorsSource=cms.VInputTag(cms.InputTag('patJetCorrFactors'+labelName+postfix))
85  ## configure MET(Type1) corrections
86  if jetCorrections[2].lower() != 'none' and jetCorrections[2] != '':
87  if not jetCorrections[2].lower() == 'type-1' and not jetCorrections[2].lower() == 'type-2':
88  raise valueError, "In addJetCollection: Wrong choice of MET corrections for new jet collection. Possible choices are None (or empty string), Type-1, Type-2 (i.e.\
89  Type-1 and Type-2 corrections applied). This choice is not case sensitive. Your choice was: ", jetCorrections[2]
90  if _type == "JPT":
91  raise ValueError("In addJecCollection: MET(type1) corrections are not supported for JPTJets. Please set the MET-LABEL to \"None\" (as string in quatiation \
92  marks) and use raw tcMET together with JPTJets.")
93  ## set up jet correctors for MET corrections
94  process.load( "JetMETCorrections.Configuration.JetCorrectorsAllAlgos_cff") # FIXME: This adds a lot of garbage
95  # I second the FIXME comment on the last line. When I counted it, this brought in 344 EDProducers
96  # to be available to run unscheduled. All jet correctors, probably some small fraction of which
97  # are actually used.
98  task.add(process.jetCorrectorsAllAlgosTask)
99  _payloadType = jetCorrections[0].split(_type)[0].lower()+_type
100  if "PF" in _type :
101  addToProcessAndTask(jetCorrections[0]+'L1FastJet',
102  getattr(process, _payloadType+'L1FastjetCorrector').clone(srcRho=cms.InputTag('fixedGridRhoFastjetAll')),
103  process, task)
104  else :
105  addToProcessAndTask(jetCorrections[0]+'L1FastJet',
106  getattr(process, _payloadType+'L1FastjetCorrector').clone(srcRho=cms.InputTag('fixedGridRhoFastjetAllCalo')),
107  process, task)
108  addToProcessAndTask(jetCorrections[0]+'L1Offset', getattr(process, _payloadType+'L1OffsetCorrector').clone(), process, task)
109  addToProcessAndTask(jetCorrections[0]+'L2Relative', getattr(process, _payloadType+'L2RelativeCorrector').clone(), process, task)
110  addToProcessAndTask(jetCorrections[0]+'L3Absolute', getattr(process, _payloadType+'L3AbsoluteCorrector').clone(), process, task)
111  addToProcessAndTask(jetCorrections[0]+'L2L3Residual', getattr(process, _payloadType+'ResidualCorrector').clone(), process, task)
112  addToProcessAndTask(jetCorrections[0]+'CombinedCorrector',
113  cms.EDProducer( 'ChainedJetCorrectorProducer', correctors = cms.VInputTag()),
114  process, task)
115  for x in jetCorrections[1]:
116  if x != 'L1FastJet' and x != 'L1Offset' and x != 'L2Relative' and x != 'L3Absolute' and x != 'L2L3Residual':
117  raise ValueError('In addJetCollection: Unsupported JEC for MET(Type1). Currently supported jet correction levels are L1FastJet, L1Offset, L2Relative, L3Asolute, L2L3Residual. Requested was: %s'%(x))
118  else:
119  _corrector = _payloadType
120  if x == 'L1FastJet':
121  _corrector += 'L1Fastjet'
122  elif x == 'L2L3Residual':
123  _corrector += 'Residual'
124  else:
125  _corrector += x
126  _corrector += 'Corrector'
127  getattr(process, jetCorrections[0]+'CombinedCorrector').correctors.append(cms.InputTag(_corrector))
128 
129  ## set up MET(Type1) correction modules
130  _labelCorrName = labelName
131  if labelName != '':
132  _labelCorrName = 'For' + labelName
133  if _type == 'Calo':
134  from JetMETCorrections.Type1MET.correctionTermsCaloMet_cff import corrCaloMetType1
135  from JetMETCorrections.Type1MET.correctionTermsCaloMet_cff import corrCaloMetType2
136  from JetMETCorrections.Type1MET.correctedMet_cff import caloMetT1
137  from JetMETCorrections.Type1MET.correctedMet_cff import caloMetT1T2
139  jetCorrections[0]+_labelCorrName+'JetMETcorr'+postfix,
140  corrCaloMetType1.clone(src=jetSource,srcMET = "caloMetM",jetCorrLabel = cms.InputTag(jetCorrections[0]+'CombinedCorrector')),
141  process, task)
143  jetCorrections[0]+_labelCorrName+'JetMETcorr2'+postfix,
144  corrCaloMetType2.clone(srcUnclEnergySums = cms.VInputTag(
145  cms.InputTag(jetCorrections[0]+_labelCorrName+'JetMETcorr'+postfix, 'type2'),
146  cms.InputTag(jetCorrections[0]+_labelCorrName+'JetMETcorr'+postfix, 'offset'),
147  cms.InputTag('muCaloMetCorr'))),
148  process, task)
150  jetCorrections[0]+_labelCorrName+'Type1CorMet'+postfix,
151  caloMetT1.clone(src = "caloMetM", srcCorrections = cms.VInputTag(
152  cms.InputTag(jetCorrections[0]+_labelCorrName+'JetMETcorr'+postfix, 'type1'))),
153  process, task)
154  addToProcessAndTask(jetCorrections[0]+_labelCorrName+'Type1p2CorMet'+postfix,
155  caloMetT1T2.clone(src = "caloMetM", srcCorrections = cms.VInputTag(
156  cms.InputTag(jetCorrections[0]+_labelCorrName+'JetMETcorr'+postfix, 'type1'),
157  cms.InputTag(jetCorrections[0]+_labelCorrName+'JetMETcorr2'+postfix))),
158  process, task)
159  elif _type == 'PF':
161  from JetMETCorrections.Type1MET.correctionTermsPfMetType1Type2_cff import pfCandsNotInJetsPtrForMetCorr
162  from JetMETCorrections.Type1MET.correctionTermsPfMetType1Type2_cff import pfCandsNotInJetsForMetCorr
166  from JetMETCorrections.Type1MET.correctedMet_cff import pfMetT1
167  from JetMETCorrections.Type1MET.correctedMet_cff import pfMetT1T2
168  addToProcessAndTask(jetCorrections[0]+_labelCorrName+'pfJetsPtrForMetCorr'+postfix,
169  pfJetsPtrForMetCorr.clone(src = jetSource), process, task)
171  jetCorrections[0]+_labelCorrName+'pfCandsNotInJetsPtrForMetCorr'+postfix,
172  pfCandsNotInJetsPtrForMetCorr.clone(topCollection = jetCorrections[0]+_labelCorrName+'pfJetsPtrForMetCorr'+postfix),
173  process, task)
175  jetCorrections[0]+_labelCorrName+'pfCandsNotInJetsForMetCorr'+postfix,
176  pfCandsNotInJetsForMetCorr.clone(src = jetCorrections[0]+_labelCorrName+'pfCandsNotInJetsPtrForMetCorr'+postfix),
177  process, task)
179  jetCorrections[0]+_labelCorrName+'CandMETcorr'+postfix,
180  pfCandMETcorr.clone(src = cms.InputTag(jetCorrections[0]+_labelCorrName+'pfCandsNotInJetsForMetCorr'+postfix)),
181  process, task)
183  jetCorrections[0]+_labelCorrName+'JetMETcorr'+postfix,
184  corrPfMetType1.clone(src = jetSource, jetCorrLabel = cms.InputTag(jetCorrections[0]+'CombinedCorrector')),
185  process, task) # FIXME: Originally w/o jet corrections?
186  addToProcessAndTask(jetCorrections[0]+_labelCorrName+'corrPfMetType2'+postfix,
187  corrPfMetType2.clone(srcUnclEnergySums = cms.VInputTag(
188  cms.InputTag(jetCorrections[0]+_labelCorrName+'JetMETcorr'+postfix, 'type2'),
189  cms.InputTag(jetCorrections[0]+_labelCorrName+'JetMETcorr'+postfix, 'offset'),
190  cms.InputTag(jetCorrections[0]+_labelCorrName+'CandMETcorr'+postfix))),
191  process, task)
192  addToProcessAndTask(jetCorrections[0]+_labelCorrName+'Type1CorMet'+postfix,
193  pfMetT1.clone(srcCorrections = cms.VInputTag(
194  cms.InputTag(jetCorrections[0]+_labelCorrName+'JetMETcorr'+postfix, 'type1'))),
195  process, task)
196  addToProcessAndTask(jetCorrections[0]+_labelCorrName+'Type1p2CorMet'+postfix,
197  pfMetT1T2.clone(srcCorrections = cms.VInputTag(
198  cms.InputTag(jetCorrections[0]+_labelCorrName+'JetMETcorr'+postfix, 'type1'),
199  jetCorrections[0]+_labelCorrName+'corrPfMetType2'+postfix)),
200  process, task)
201 
202  ## common configuration for Calo and PF
203  if ('L1FastJet' in jetCorrections[1] or 'L1Fastjet' in jetCorrections[1]):
204  getattr(process,jetCorrections[0]+_labelCorrName+'JetMETcorr'+postfix).offsetCorrLabel = cms.InputTag(jetCorrections[0]+'L1FastJet')
205  #FIXME: What is wrong here?
206  #elif ('L1Offset' in jetCorrections[1]):
207  #getattr(process,jetCorrections[0]+_labelCorrName+'JetMETcorr'+postfix).offsetCorrLabel = cms.InputTag(jetCorrections[0]+'L1Offset')
208  else:
209  getattr(process,jetCorrections[0]+_labelCorrName+'JetMETcorr'+postfix).offsetCorrLabel = cms.InputTag('')
210 
212  if jetCorrections[2].lower() == 'type-1':
213  addToProcessAndTask('patMETs'+labelName+postfix,
214  patMETs.clone(metSource = cms.InputTag(jetCorrections[0]+_labelCorrName+'Type1CorMet'+postfix),
215  addMuonCorrections = False),
216  process, task)
217  elif jetCorrections[2].lower() == 'type-2':
218  addToProcessAndTask('patMETs'+labelName+postfix,
219  patMETs.clone(metSource = cms.InputTag(jetCorrections[0]+_labelCorrName+'Type1p2CorMet'+postfix),
220  addMuonCorrections = False),
221  process, task)
222 
223 
def addToProcessAndTask(label, module, process, task)
Definition: helpers.py:27
def setupJetCorrections(process, knownModules, jetCorrections, jetSource, pvSource, patJets, labelName, postfix)
Definition: jetTools.py:25
TEveGeoShape * clone(const TEveElement *element, TEveElement *parent)
Definition: eve_macros.cc:135
def getPatAlgosToolsTask(process)
Definition: helpers.py:12
double split
Definition: MVATrainer.cc:139
def jetTools.setupSVClustering (   btagInfo,
  svClustering,
  algo,
  rParam,
  fatJets = cms.InputTag(''),
  groomedFatJets = cms.InputTag('') 
)

Definition at line 224 of file jetTools.py.

References setupBTagging().

Referenced by setupBTagging().

224 def setupSVClustering(btagInfo, svClustering, algo, rParam, fatJets=cms.InputTag(''), groomedFatJets=cms.InputTag('')):
225  btagInfo.useSVClustering = cms.bool(svClustering)
226  btagInfo.jetAlgorithm = cms.string(algo)
227  btagInfo.rParam = cms.double(rParam)
228  ## if the jet is actually a subjet
229  if fatJets != cms.InputTag(''):
230  btagInfo.fatJets = fatJets
231  if groomedFatJets != cms.InputTag(''):
232  btagInfo.groomedFatJets = groomedFatJets
233 
234 
def setupSVClustering(btagInfo, svClustering, algo, rParam, fatJets=cms.InputTag(''), groomedFatJets=cms.InputTag(''))
Definition: jetTools.py:224
def jetTools.undefinedLabelName (   obj)

Definition at line 1851 of file jetTools.py.

Referenced by jetTools.AddJetCollection.toolCode().

1852  sys.stderr.write("-------------------------------------------------------\n")
1853  sys.stderr.write(" Error: the jet 'labelName' is not defined.\n")
1854  sys.stderr.write(" All added jets must have 'labelName' defined.\n")
1855  sys.stderr.write("-------------------------------------------------------\n")
1856  raise KeyError("Undefined jet 'labelName' used in '"+obj._label+"'")
1857 
def undefinedLabelName(obj)
Definition: jetTools.py:1851
def jetTools.unsupportedJetAlgorithm (   obj)

Definition at line 1858 of file jetTools.py.

Referenced by jetTools.AddJetCollection.toolCode(), and jetTools.UpdateJetCollection.toolCode().

1859  sys.stderr.write("-------------------------------------------------------\n")
1860  sys.stderr.write(" Error: Unsupported jet algorithm detected.\n")
1861  sys.stderr.write(" The supported algorithms are:\n")
1862  for key in supportedJetAlgos.keys():
1863  sys.stderr.write(" " + key.upper() + ", " + key.lower() + ": " + supportedJetAlgos[key] + "\n")
1864  sys.stderr.write("-------------------------------------------------------\n")
1865  raise KeyError("Unsupported jet algorithm used in '"+obj._label+"'")
1866 
def unsupportedJetAlgorithm(obj)
Definition: jetTools.py:1858

Variable Documentation

jetTools.supportedJetAlgos

dictionary with supported jet clustering algorithms

Definition at line 11 of file jetTools.py.