CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
Classes | Functions | Variables
heepElectronID_tools Namespace Reference

Classes

class  HEEP_WorkingPoint_V1
 
class  HEEP_WorkingPoint_V2
 

Functions

def addHEEPProducersToSeq
 
def addHEEPVMProducerToSeq
 
def configureHEEPElectronID_V51
 
def configureHEEPElectronID_V60
 
def configureHEEPElectronID_V60_80XAOD
 
def configureHEEPElectronID_V61
 
def configureHEEPElectronID_V70
 
def configureHEEPElectronID_V71
 
def psetGsfEleDEtaInSeedCut
 
def psetGsfEleDPhiInCut
 
def psetGsfEleDxyCut
 
def psetGsfEleEcalDrivenCut
 
def psetGsfEleEmHadD1IsoRhoCut
 
def psetGsfEleFull5x5E2x5OverE5x5Cut
 
def psetGsfEleFull5x5E2x5OverE5x5WithSatCut
 
def psetGsfEleFull5x5SigmaIEtaIEtaCut
 
def psetGsfEleFull5x5SigmaIEtaIEtaWithSatCut
 
def psetGsfEleHadronicOverEMLinearCut
 
def psetGsfEleMissingHitsCut
 
def psetGsfEleSCEtaMultiRangeCut
 
def psetGsfEleTrkPtFall16IsoCut
 
def psetGsfEleTrkPtIsoCut
 
def psetGsfEleTrkPtIsoRhoCut
 
def psetGsfEleTrkPtNoJetCoreIsoCut
 
def psetMinEcalEtCut
 
def psetMinPtCut
 

Variables

float ebCutOff = 1.479
 
float ebMax = 1.4442
 
float eeMin = 1.566
 

Function Documentation

def heepElectronID_tools.addHEEPProducersToSeq (   process,
  seq,
  useMiniAOD,
  task = None 
)
currently no additional products are needed to run the HEEP ID so this 
function simply remains as a placeholder

Definition at line 576 of file heepElectronID_tools.py.

Referenced by vid_id_tools.setupVIDElectronSelection().

577 def addHEEPProducersToSeq(process,seq,useMiniAOD, task=None):
578  '''currently no additional products are needed to run the HEEP ID so this
579  function simply remains as a placeholder
580  '''
581  return
def heepElectronID_tools.addHEEPVMProducerToSeq (   process,
  seq,
  useMiniAOD,
  task = None 
)

Definition at line 548 of file heepElectronID_tools.py.

549 def addHEEPVMProducerToSeq(process,seq,useMiniAOD, task=None):
550  newTask = cms.Task()
551  seq.associate(newTask)
552  if task is not None:
553  task.add(newTask)
554 
555  process.load("RecoEgamma.ElectronIdentification.heepIdVarValueMapProducer_cfi")
556  newTask.add(process.heepIDVarValueMaps)
557 
558  if useMiniAOD==False:
559  process.load("TrackingTools.TransientTrack.TransientTrackBuilder_cfi")
560  process.load("PhysicsTools.PatAlgos.slimming.primaryVertexAssociation_cfi")
561  process.load("PhysicsTools.PatAlgos.slimming.offlineSlimmedPrimaryVertices_cfi")
562  process.load("PhysicsTools.PatAlgos.slimming.packedPFCandidates_cfi")
563  from PhysicsTools.PatAlgos.slimming.packedPFCandidates_cfi import packedPFCandidates
564  process.packedCandsForTkIso = packedPFCandidates.clone()
565  process.packedCandsForTkIso.PuppiSrc=cms.InputTag("")
566  process.packedCandsForTkIso.PuppiNoLepSrc=cms.InputTag("")
567 
568  process.load("PhysicsTools.PatAlgos.slimming.lostTracks_cfi")
569  from PhysicsTools.PatAlgos.slimming.lostTracks_cfi import lostTracks
570  process.lostTracksForTkIso = lostTracks.clone()
571  process.lostTracksForTkIso.packedPFCandidates =cms.InputTag("packedCandsForTkIso")
572  newTask.add(process.primaryVertexAssociation,
573  process.offlineSlimmedPrimaryVertices,
574  process.packedCandsForTkIso,
575  process.lostTracksForTkIso)
def heepElectronID_tools.configureHEEPElectronID_V51 (   wpEB,
  wpEE 
)
This function configures the full cms.PSet for a VID ID and returns it.
The inputs: two objects of the type HEEP_WorkingPoint_V1, one
containing the cuts for the Barrel (EB) and the other one for the Endcap (EE).

Definition at line 396 of file heepElectronID_tools.py.

References psetGsfEleDEtaInSeedCut(), psetGsfEleDPhiInCut(), psetGsfEleDxyCut(), psetGsfEleEcalDrivenCut(), psetGsfEleEmHadD1IsoRhoCut(), psetGsfEleFull5x5E2x5OverE5x5Cut(), psetGsfEleFull5x5SigmaIEtaIEtaCut(), psetGsfEleHadronicOverEMLinearCut(), psetGsfEleMissingHitsCut(), psetGsfEleSCEtaMultiRangeCut(), psetGsfEleTrkPtIsoCut(), and psetMinPtCut().

397 def configureHEEPElectronID_V51(wpEB, wpEE):
398  """
399  This function configures the full cms.PSet for a VID ID and returns it.
400  The inputs: two objects of the type HEEP_WorkingPoint_V1, one
401  containing the cuts for the Barrel (EB) and the other one for the Endcap (EE).
402  """
403  parameterSet = cms.PSet(
404  idName = cms.string("heepElectronID-HEEPV51"),
405  cutFlow = cms.VPSet(
406  psetMinPtCut(cutValue=35.), #0
408  psetGsfEleDEtaInSeedCut(wpEB,wpEE), #2
409  psetGsfEleDPhiInCut(wpEB,wpEE), #3
410  psetGsfEleFull5x5SigmaIEtaIEtaCut(wpEB,wpEE), #4
411  psetGsfEleFull5x5E2x5OverE5x5Cut(wpEB,wpEE), #5
412  psetGsfEleHadronicOverEMLinearCut(wpEB,wpEE), #6
413  psetGsfEleTrkPtIsoCut(wpEB,wpEE), #7
414  psetGsfEleEmHadD1IsoRhoCut(wpEB,wpEE), #8
415  psetGsfEleDxyCut(wpEB,wpEE), #9
416  psetGsfEleMissingHitsCut(wpEB,wpEE), #10,
417  psetGsfEleEcalDrivenCut(wpEB,wpEE) #11
418  )
419  )
420  return parameterSet
def heepElectronID_tools.configureHEEPElectronID_V60 (   wpEB,
  wpEE 
)
This function configures the full cms.PSet for a VID ID and returns it.
The inputs: two objects of the type HEEP_WorkingPoint_V1, one
containing the cuts for the Barrel (EB) and the other one for the Endcap (EE).

Definition at line 421 of file heepElectronID_tools.py.

References psetGsfEleDEtaInSeedCut(), psetGsfEleDPhiInCut(), psetGsfEleDxyCut(), psetGsfEleEcalDrivenCut(), psetGsfEleEmHadD1IsoRhoCut(), psetGsfEleFull5x5E2x5OverE5x5Cut(), psetGsfEleFull5x5SigmaIEtaIEtaCut(), psetGsfEleHadronicOverEMLinearCut(), psetGsfEleMissingHitsCut(), psetGsfEleSCEtaMultiRangeCut(), psetGsfEleTrkPtIsoCut(), and psetMinPtCut().

422 def configureHEEPElectronID_V60(wpEB, wpEE):
423  """
424  This function configures the full cms.PSet for a VID ID and returns it.
425  The inputs: two objects of the type HEEP_WorkingPoint_V1, one
426  containing the cuts for the Barrel (EB) and the other one for the Endcap (EE).
427  """
428  parameterSet = cms.PSet(
429  idName = cms.string("heepElectronID-HEEPV60"),
430  cutFlow = cms.VPSet(
431  psetMinPtCut(cutValue=35.), #0
433  psetGsfEleDEtaInSeedCut(wpEB,wpEE), #2
434  psetGsfEleDPhiInCut(wpEB,wpEE), #3
435  psetGsfEleFull5x5SigmaIEtaIEtaCut(wpEB,wpEE), #4
436  psetGsfEleFull5x5E2x5OverE5x5Cut(wpEB,wpEE), #5
437  psetGsfEleHadronicOverEMLinearCut(wpEB,wpEE), #6
438  psetGsfEleTrkPtIsoCut(wpEB,wpEE), #7
439  psetGsfEleEmHadD1IsoRhoCut(wpEB,wpEE), #8
440  psetGsfEleDxyCut(wpEB,wpEE), #9
441  psetGsfEleMissingHitsCut(wpEB,wpEE), #10,
442  psetGsfEleEcalDrivenCut(wpEB,wpEE) #11
443  )
444  )
445  return parameterSet
446 
def heepElectronID_tools.configureHEEPElectronID_V60_80XAOD (   idName,
  wpEB,
  wpEE 
)
This function configures the full cms.PSet for a VID ID and returns it.
The inputs: two objects of the type HEEP_WorkingPoint_V1, one
containing the cuts for the Barrel (EB) and the other one for the Endcap (EE).

Definition at line 498 of file heepElectronID_tools.py.

References psetGsfEleDEtaInSeedCut(), psetGsfEleDPhiInCut(), psetGsfEleDxyCut(), psetGsfEleEcalDrivenCut(), psetGsfEleEmHadD1IsoRhoCut(), psetGsfEleFull5x5E2x5OverE5x5WithSatCut(), psetGsfEleFull5x5SigmaIEtaIEtaWithSatCut(), psetGsfEleHadronicOverEMLinearCut(), psetGsfEleMissingHitsCut(), psetGsfEleSCEtaMultiRangeCut(), psetGsfEleTrkPtNoJetCoreIsoCut(), and psetMinPtCut().

499 def configureHEEPElectronID_V60_80XAOD(idName, wpEB, wpEE):
500  """
501  This function configures the full cms.PSet for a VID ID and returns it.
502  The inputs: two objects of the type HEEP_WorkingPoint_V1, one
503  containing the cuts for the Barrel (EB) and the other one for the Endcap (EE).
504  """
505  parameterSet = cms.PSet(
506  idName = cms.string(idName),
507  cutFlow = cms.VPSet(
508  psetMinPtCut(cutValue=35.), #0
510  psetGsfEleDEtaInSeedCut(wpEB,wpEE), #2
511  psetGsfEleDPhiInCut(wpEB,wpEE), #3
514  psetGsfEleHadronicOverEMLinearCut(wpEB,wpEE), #6
515  psetGsfEleTrkPtNoJetCoreIsoCut(wpEB,wpEE), #7
516  psetGsfEleEmHadD1IsoRhoCut(wpEB,wpEE), #8
517  psetGsfEleDxyCut(wpEB,wpEE), #9
518  psetGsfEleMissingHitsCut(wpEB,wpEE), #10,
519  psetGsfEleEcalDrivenCut(wpEB,wpEE) #11
520  )
521  )
522  return parameterSet
def heepElectronID_tools.configureHEEPElectronID_V61 (   wpEB,
  wpEE 
)
This function configures the full cms.PSet for a VID ID and returns it.
The inputs: two objects of the type HEEP_WorkingPoint_V2, one
containing the cuts for the Barrel (EB) and the other one for the Endcap (EE).

Definition at line 523 of file heepElectronID_tools.py.

References psetGsfEleDEtaInSeedCut(), psetGsfEleDPhiInCut(), psetGsfEleDxyCut(), psetGsfEleEcalDrivenCut(), psetGsfEleEmHadD1IsoRhoCut(), psetGsfEleFull5x5E2x5OverE5x5Cut(), psetGsfEleFull5x5SigmaIEtaIEtaCut(), psetGsfEleHadronicOverEMLinearCut(), psetGsfEleMissingHitsCut(), psetGsfEleSCEtaMultiRangeCut(), psetGsfEleTrkPtIsoRhoCut(), and psetMinPtCut().

524 def configureHEEPElectronID_V61(wpEB, wpEE):
525  """
526  This function configures the full cms.PSet for a VID ID and returns it.
527  The inputs: two objects of the type HEEP_WorkingPoint_V2, one
528  containing the cuts for the Barrel (EB) and the other one for the Endcap (EE).
529  """
530  parameterSet = cms.PSet(
531  idName = cms.string("heepElectronID-HEEPV61"),
532  cutFlow = cms.VPSet(
533  psetMinPtCut(cutValue=35.), #0
535  psetGsfEleDEtaInSeedCut(wpEB,wpEE), #2
536  psetGsfEleDPhiInCut(wpEB,wpEE), #3
537  psetGsfEleFull5x5SigmaIEtaIEtaCut(wpEB,wpEE), #4
538  psetGsfEleFull5x5E2x5OverE5x5Cut(wpEB,wpEE), #5
539  psetGsfEleHadronicOverEMLinearCut(wpEB,wpEE), #6
540  psetGsfEleTrkPtIsoRhoCut(wpEB,wpEE), #7
541  psetGsfEleEmHadD1IsoRhoCut(wpEB,wpEE), #8
542  psetGsfEleDxyCut(wpEB,wpEE), #9
543  psetGsfEleMissingHitsCut(wpEB,wpEE), #10,
544  psetGsfEleEcalDrivenCut(wpEB,wpEE) #11
545  )
546  )
547  return parameterSet
def heepElectronID_tools.configureHEEPElectronID_V70 (   idName,
  wpEB,
  wpEE 
)
This function configures the full cms.PSet for a VID ID and returns it.
The inputs: two objects of the type HEEP_WorkingPoint_V1, one
containing the cuts for the Barrel (EB) and the other one for the Endcap (EE).

Definition at line 447 of file heepElectronID_tools.py.

References psetGsfEleDEtaInSeedCut(), psetGsfEleDPhiInCut(), psetGsfEleDxyCut(), psetGsfEleEcalDrivenCut(), psetGsfEleEmHadD1IsoRhoCut(), psetGsfEleFull5x5E2x5OverE5x5WithSatCut(), psetGsfEleFull5x5SigmaIEtaIEtaWithSatCut(), psetGsfEleHadronicOverEMLinearCut(), psetGsfEleMissingHitsCut(), psetGsfEleSCEtaMultiRangeCut(), psetGsfEleTrkPtIsoCut(), and psetMinPtCut().

448 def configureHEEPElectronID_V70(idName, wpEB, wpEE):
449  """
450  This function configures the full cms.PSet for a VID ID and returns it.
451  The inputs: two objects of the type HEEP_WorkingPoint_V1, one
452  containing the cuts for the Barrel (EB) and the other one for the Endcap (EE).
453  """
454  parameterSet = cms.PSet(
455  idName = cms.string(idName),
456  cutFlow = cms.VPSet(
457  psetMinPtCut(cutValue=35.), #0
459  psetGsfEleDEtaInSeedCut(wpEB,wpEE), #2
460  psetGsfEleDPhiInCut(wpEB,wpEE), #3
463  psetGsfEleHadronicOverEMLinearCut(wpEB,wpEE), #6
464  psetGsfEleTrkPtIsoCut(wpEB,wpEE,useHEEPIso=True),#7
465  psetGsfEleEmHadD1IsoRhoCut(wpEB,wpEE), #8
466  psetGsfEleDxyCut(wpEB,wpEE), #9
467  psetGsfEleMissingHitsCut(wpEB,wpEE), #10,
468  psetGsfEleEcalDrivenCut(wpEB,wpEE) #11
469  )
470  )
471  return parameterSet
def heepElectronID_tools.configureHEEPElectronID_V71 (   idName,
  wpEB,
  wpEE,
  minEt 
)
This function configures the full cms.PSet for a VID ID and returns it.
The inputs: two objects of the type HEEP_WorkingPoint_V1, one
containing the cuts for the Barrel (EB) and the other one for the Endcap (EE).

Definition at line 472 of file heepElectronID_tools.py.

References psetGsfEleDEtaInSeedCut(), psetGsfEleDPhiInCut(), psetGsfEleDxyCut(), psetGsfEleEcalDrivenCut(), psetGsfEleEmHadD1IsoRhoCut(), psetGsfEleFull5x5E2x5OverE5x5WithSatCut(), psetGsfEleFull5x5SigmaIEtaIEtaWithSatCut(), psetGsfEleHadronicOverEMLinearCut(), psetGsfEleMissingHitsCut(), psetGsfEleSCEtaMultiRangeCut(), psetGsfEleTrkPtIsoCut(), and psetMinEcalEtCut().

473 def configureHEEPElectronID_V71(idName, wpEB, wpEE,minEt):
474  """
475  This function configures the full cms.PSet for a VID ID and returns it.
476  The inputs: two objects of the type HEEP_WorkingPoint_V1, one
477  containing the cuts for the Barrel (EB) and the other one for the Endcap (EE).
478  """
479  parameterSet = cms.PSet(
480  idName = cms.string(idName),
481  cutFlow = cms.VPSet(
482  psetMinEcalEtCut(cutValue=minEt), #0
484  psetGsfEleDEtaInSeedCut(wpEB,wpEE), #2
485  psetGsfEleDPhiInCut(wpEB,wpEE), #3
488  psetGsfEleHadronicOverEMLinearCut(wpEB,wpEE), #6
489  psetGsfEleTrkPtIsoCut(wpEB,wpEE,useHEEPIso=True),#7
490  psetGsfEleEmHadD1IsoRhoCut(wpEB,wpEE,energyType="Ecal"), #8
491  psetGsfEleDxyCut(wpEB,wpEE), #9
492  psetGsfEleMissingHitsCut(wpEB,wpEE), #10,
493  psetGsfEleEcalDrivenCut(wpEB,wpEE) #11
494  )
495  )
496  return parameterSet
497 
def heepElectronID_tools.psetGsfEleDEtaInSeedCut (   wpEB,
  wpEE 
)

Definition at line 172 of file heepElectronID_tools.py.

Referenced by configureHEEPElectronID_V51(), configureHEEPElectronID_V60(), configureHEEPElectronID_V60_80XAOD(), configureHEEPElectronID_V61(), configureHEEPElectronID_V70(), and configureHEEPElectronID_V71().

173 def psetGsfEleDEtaInSeedCut(wpEB, wpEE):
174  return cms.PSet(
175  cutName = cms.string('GsfEleEBEECut'),
176  cutString = cms.string("abs(deltaEtaSeedClusterTrackAtVtx)"),
177  cutValueEB = cms.double( wpEB.dEtaInSeedCut ),
178  cutValueEE = cms.double( wpEE.dEtaInSeedCut ),
179  needsAdditionalProducts = cms.bool(False),
180  isIgnored = cms.bool(False)
181  )
182 
# Configure the cut on the dPhiIn
def heepElectronID_tools.psetGsfEleDPhiInCut (   wpEB,
  wpEE 
)

Definition at line 183 of file heepElectronID_tools.py.

Referenced by configureHEEPElectronID_V51(), configureHEEPElectronID_V60(), configureHEEPElectronID_V60_80XAOD(), configureHEEPElectronID_V61(), configureHEEPElectronID_V70(), and configureHEEPElectronID_V71().

184 def psetGsfEleDPhiInCut(wpEB, wpEE):
185  return cms.PSet(
186  cutName = cms.string('GsfEleEBEECut'),
187  cutString = cms.string("abs(deltaPhiSuperClusterTrackAtVtx)"),
188  cutValueEB = cms.double( wpEB.dPhiInCut ),
189  cutValueEE = cms.double( wpEE.dPhiInCut ),
190  needsAdditionalProducts = cms.bool(False),
191  isIgnored = cms.bool(False)
192  )
193 
# Confugure the full 5x5 sigmaIEtaIEta cut
def heepElectronID_tools.psetGsfEleDxyCut (   wpEB,
  wpEE 
)

Definition at line 360 of file heepElectronID_tools.py.

Referenced by configureHEEPElectronID_V51(), configureHEEPElectronID_V60(), configureHEEPElectronID_V60_80XAOD(), configureHEEPElectronID_V61(), configureHEEPElectronID_V70(), and configureHEEPElectronID_V71().

361 def psetGsfEleDxyCut(wpEB, wpEE):
362  return cms.PSet(
363  cutName = cms.string('GsfEleDxyCut'),
364  dxyCutValueEB = cms.double( wpEB.dxyCut ),
365  dxyCutValueEE = cms.double( wpEE.dxyCut ),
366  vertexSrc = cms.InputTag("offlinePrimaryVertices"),
367  vertexSrcMiniAOD = cms.InputTag("offlineSlimmedPrimaryVertices"),
368  barrelCutOff = cms.double(ebCutOff),
369  needsAdditionalProducts = cms.bool(True),
370  isIgnored = cms.bool(False)
371  )
372 
# Configure the cut on missing hits
def heepElectronID_tools.psetGsfEleEcalDrivenCut (   wpEB,
  wpEE 
)

Definition at line 382 of file heepElectronID_tools.py.

Referenced by configureHEEPElectronID_V51(), configureHEEPElectronID_V60(), configureHEEPElectronID_V60_80XAOD(), configureHEEPElectronID_V61(), configureHEEPElectronID_V70(), and configureHEEPElectronID_V71().

383 def psetGsfEleEcalDrivenCut(wpEB, wpEE):
384  return cms.PSet(
385  cutName = cms.string('GsfEleEcalDrivenCut'),
386  ecalDrivenEB = cms.int32( wpEB.ecalDrivenCut ),
387  ecalDrivenEE = cms.int32( wpEE.ecalDrivenCut ),
388  barrelCutOff = cms.double(ebCutOff),
389  needsAdditionalProducts = cms.bool(False),
390  isIgnored = cms.bool(False)
391  )
392 
393 # ==============================================================
394 # Define the complete cut sets
395 # ==============================================================
def heepElectronID_tools.psetGsfEleEmHadD1IsoRhoCut (   wpEB,
  wpEE,
  energyType = "EcalTrk" 
)

Definition at line 343 of file heepElectronID_tools.py.

Referenced by configureHEEPElectronID_V51(), configureHEEPElectronID_V60(), configureHEEPElectronID_V60_80XAOD(), configureHEEPElectronID_V61(), configureHEEPElectronID_V70(), and configureHEEPElectronID_V71().

344 def psetGsfEleEmHadD1IsoRhoCut(wpEB, wpEE,energyType="EcalTrk"):
345  return cms.PSet(
346  cutName = cms.string('GsfEleEmHadD1IsoRhoCut'),
347  slopeTermEB = cms.double( wpEB.ehIsoSlopeTerm ),
348  slopeTermEE = cms.double( wpEE.ehIsoSlopeTerm ),
349  slopeStartEB = cms.double( wpEB.ehIsoSlopeStart ),
350  slopeStartEE = cms.double( wpEE.ehIsoSlopeStart ),
351  constTermEB = cms.double( wpEB.ehIsoConstTerm ),
352  constTermEE = cms.double( wpEE.ehIsoConstTerm ),
353  energyType = cms.string( energyType ),
354  rhoConstant = cms.double( wpEB.effAreaForEHIso), # expected to be the same for EB and EE
355  rho = cms.InputTag("fixedGridRhoFastjetAll"),
356  needsAdditionalProducts = cms.bool(True),
357  isIgnored = cms.bool(False)
358  )
359 
# Configure the dxy cut
def heepElectronID_tools.psetGsfEleFull5x5E2x5OverE5x5Cut (   wpEB,
  wpEE 
)

Definition at line 216 of file heepElectronID_tools.py.

Referenced by configureHEEPElectronID_V51(), configureHEEPElectronID_V60(), and configureHEEPElectronID_V61().

217 def psetGsfEleFull5x5E2x5OverE5x5Cut(wpEB, wpEE):
218  return cms.PSet(
219  cutName = cms.string('GsfEleFull5x5E2x5OverE5x5Cut'),
220  # E1x5 / E5x5
221  minE1x5OverE5x5EB = cms.double( wpEB.minE1x5OverE5x5Cut ),
222  minE1x5OverE5x5EE = cms.double( wpEE.minE1x5OverE5x5Cut ),
223  # E2x5 / E5x5
224  minE2x5OverE5x5EB = cms.double( wpEB.minE2x5OverE5x5Cut ),
225  minE2x5OverE5x5EE = cms.double( wpEE.minE2x5OverE5x5Cut ),
226  needsAdditionalProducts = cms.bool(False),
227  isIgnored = cms.bool(False)
228  )
# Configure XxX shower shape cuts
def heepElectronID_tools.psetGsfEleFull5x5E2x5OverE5x5WithSatCut (   wpEB,
  wpEE 
)

Definition at line 229 of file heepElectronID_tools.py.

Referenced by configureHEEPElectronID_V60_80XAOD(), configureHEEPElectronID_V70(), and configureHEEPElectronID_V71().

231  return cms.PSet(
232  cutName = cms.string('GsfEleFull5x5E2x5OverE5x5WithSatCut'),
233  # E1x5 / E5x5
234  minE1x5OverE5x5EB = cms.double( wpEB.minE1x5OverE5x5Cut ),
235  minE1x5OverE5x5EE = cms.double( wpEE.minE1x5OverE5x5Cut ),
236  # E2x5 / E5x5
237  minE2x5OverE5x5EB = cms.double( wpEB.minE2x5OverE5x5Cut ),
238  minE2x5OverE5x5EE = cms.double( wpEE.minE2x5OverE5x5Cut ),
239  maxNrSatCrysIn5x5EB =cms.int32( 0 ),
240  maxNrSatCrysIn5x5EE =cms.int32( 0 ),
241  needsAdditionalProducts = cms.bool(False),
242  isIgnored = cms.bool(False)
243  )
# Configure the cut of E/H
def heepElectronID_tools.psetGsfEleFull5x5SigmaIEtaIEtaCut (   wpEB,
  wpEE 
)

Definition at line 194 of file heepElectronID_tools.py.

Referenced by configureHEEPElectronID_V51(), configureHEEPElectronID_V60(), and configureHEEPElectronID_V61().

195 def psetGsfEleFull5x5SigmaIEtaIEtaCut(wpEB, wpEE):
196  return cms.PSet(
197  cutName = cms.string('GsfEleEBEECut'),
198  cutString = cms.string("full5x5_sigmaIetaIeta"),
199  cutValueEB = cms.double( wpEB.full5x5SigmaIEtaIEtaCut ),
200  cutValueEE = cms.double( wpEE.full5x5SigmaIEtaIEtaCut ),
201  needsAdditionalProducts = cms.bool(False),
202  isIgnored = cms.bool(False)
203  )
def heepElectronID_tools.psetGsfEleFull5x5SigmaIEtaIEtaWithSatCut (   wpEB,
  wpEE 
)

Definition at line 204 of file heepElectronID_tools.py.

Referenced by configureHEEPElectronID_V60_80XAOD(), configureHEEPElectronID_V70(), and configureHEEPElectronID_V71().

206  return cms.PSet(
207  cutName = cms.string('GsfEleFull5x5SigmaIEtaIEtaWithSatCut'),
208  maxSigmaIEtaIEtaEB = cms.double( wpEB.full5x5SigmaIEtaIEtaCut ),
209  maxSigmaIEtaIEtaEE = cms.double( wpEE.full5x5SigmaIEtaIEtaCut ),
210  maxNrSatCrysIn5x5EB =cms.int32( 0 ),
211  maxNrSatCrysIn5x5EE =cms.int32( 0 ),
212  needsAdditionalProducts = cms.bool(False),
213 
214  isIgnored = cms.bool(False)
215  )
# Configure XxX shower shape cuts
def heepElectronID_tools.psetGsfEleHadronicOverEMLinearCut (   wpEB,
  wpEE 
)

Definition at line 244 of file heepElectronID_tools.py.

Referenced by configureHEEPElectronID_V51(), configureHEEPElectronID_V60(), configureHEEPElectronID_V60_80XAOD(), configureHEEPElectronID_V61(), configureHEEPElectronID_V70(), and configureHEEPElectronID_V71().

245 def psetGsfEleHadronicOverEMLinearCut(wpEB, wpEE) :
246  return cms.PSet(
247  cutName = cms.string('GsfEleHadronicOverEMLinearCut'),
248  # Three constants for the GsfEleHadronicOverEMLinearCut
249  # cut = constTerm if value < slopeStart
250  # cut = slopeTerm * (value - slopeStart) + constTerm if value >= slopeStart
251  slopeTermEB = cms.double( wpEB.hOverESlopeTerm ),
252  slopeTermEE = cms.double( wpEE.hOverESlopeTerm ),
253  slopeStartEB = cms.double( wpEB.hOverESlopeStart ),
254  slopeStartEE = cms.double( wpEE.hOverESlopeStart ),
255  constTermEB = cms.double( wpEB.hOverEConstTerm ),
256  constTermEE = cms.double( wpEE.hOverEConstTerm ),
257  needsAdditionalProducts = cms.bool(False),
258  isIgnored = cms.bool(False)
259  )
260 
# Configure the cut on the tracker isolation
def heepElectronID_tools.psetGsfEleMissingHitsCut (   wpEB,
  wpEE 
)

Definition at line 373 of file heepElectronID_tools.py.

Referenced by configureHEEPElectronID_V51(), configureHEEPElectronID_V60(), configureHEEPElectronID_V60_80XAOD(), configureHEEPElectronID_V61(), configureHEEPElectronID_V70(), and configureHEEPElectronID_V71().

374 def psetGsfEleMissingHitsCut(wpEB, wpEE):
375  return cms.PSet(
376  cutName = cms.string('GsfEleMissingHitsCut'),
377  maxMissingHitsEB = cms.uint32( wpEB.maxMissingHitsCut ),
378  maxMissingHitsEE = cms.uint32( wpEE.maxMissingHitsCut ),
379  barrelCutOff = cms.double(ebCutOff),
380  needsAdditionalProducts = cms.bool(False),
381  isIgnored = cms.bool(False)
)
def heepElectronID_tools.psetGsfEleSCEtaMultiRangeCut ( )

Definition at line 157 of file heepElectronID_tools.py.

Referenced by configureHEEPElectronID_V51(), configureHEEPElectronID_V60(), configureHEEPElectronID_V60_80XAOD(), configureHEEPElectronID_V61(), configureHEEPElectronID_V70(), and configureHEEPElectronID_V71().

159  return cms.PSet(
160  cutName = cms.string("GsfEleSCEtaMultiRangeCut"),
161  useAbsEta = cms.bool(True),
162  allowedEtaRanges = cms.VPSet(
163  cms.PSet( minEta = cms.double(0.0),
164  maxEta = cms.double(ebMax) ),
165  cms.PSet( minEta = cms.double(eeMin),
166  maxEta = cms.double(2.5) )
167  ),
168  needsAdditionalProducts = cms.bool(False),
169  isIgnored = cms.bool(False)
170  )
171 
# Configure the cut on the dEtaIn for the seed
def heepElectronID_tools.psetGsfEleTrkPtFall16IsoCut (   wpEB,
  wpEE 
)

Definition at line 320 of file heepElectronID_tools.py.

321 def psetGsfEleTrkPtFall16IsoCut(wpEB, wpEE):
322  return cms.PSet(
323  cutName = cms.string('GsfEleValueMapIsoRhoCut'),
324  # Three constants for the GsfEleTrkPtIsoCut
325  # cut = constTerm if value < slopeStart
326  # cut = slopeTerm * (value - slopeStart) + constTerm if value >= slopeStart
327  slopeTermEB = cms.double( wpEB.trkIsoSlopeTerm ),
328  slopeTermEE = cms.double( wpEE.trkIsoSlopeTerm ),
329  slopeStartEB = cms.double( wpEB.trkIsoSlopeStart ),
330  slopeStartEE = cms.double( wpEE.trkIsoSlopeStart ),
331  constTermEB = cms.double( wpEB.trkIsoConstTerm ),
332  constTermEE = cms.double( wpEE.trkIsoConstTerm ),
333  #no rho so we zero it out, if the input tag is empty, its ignored anyways
334  rhoEtStartEB = cms.double( 999999.),
335  rhoEtStartEE = cms.double( 999999.),
336  rhoEAEB = cms.double( 0. ),
337  rhoEAEE = cms.double( 0. ),
338  rho = cms.InputTag(""),
339  value = cms.InputTag("heepIDVarValueMaps","eleTrkPtIso"),
340  needsAdditionalProducts = cms.bool(True),
341  isIgnored = cms.bool(False)
342  )
# Configure the cut on the EM + Had_depth_1 isolation with rho correction
def heepElectronID_tools.psetGsfEleTrkPtIsoCut (   wpEB,
  wpEE,
  useHEEPIso = False 
)

Definition at line 261 of file heepElectronID_tools.py.

Referenced by configureHEEPElectronID_V51(), configureHEEPElectronID_V60(), configureHEEPElectronID_V70(), and configureHEEPElectronID_V71().

262 def psetGsfEleTrkPtIsoCut(wpEB, wpEE, useHEEPIso = False):
263  return cms.PSet(
264  cutName = cms.string('GsfEleTrkPtIsoCut'),
265  # Three constants for the GsfEleTrkPtIsoCut
266  # cut = constTerm if value < slopeStart
267  # cut = slopeTerm * (value - slopeStart) + constTerm if value >= slopeStart
268  slopeTermEB = cms.double( wpEB.trkIsoSlopeTerm ),
269  slopeTermEE = cms.double( wpEE.trkIsoSlopeTerm ),
270  slopeStartEB = cms.double( wpEB.trkIsoSlopeStart ),
271  slopeStartEE = cms.double( wpEE.trkIsoSlopeStart ),
272  constTermEB = cms.double( wpEB.trkIsoConstTerm ),
273  constTermEE = cms.double( wpEE.trkIsoConstTerm ),
274  useHEEPIso = cms.bool(useHEEPIso),
275  needsAdditionalProducts = cms.bool(False),
276  isIgnored = cms.bool(False)
277  )
# Configure the cut on the tracker isolation with a rho correction (hack for 76X)
def heepElectronID_tools.psetGsfEleTrkPtIsoRhoCut (   wpEB,
  wpEE 
)

Definition at line 278 of file heepElectronID_tools.py.

Referenced by configureHEEPElectronID_V61().

279 def psetGsfEleTrkPtIsoRhoCut(wpEB, wpEE):
280  return cms.PSet(
281  cutName = cms.string('GsfEleTrkPtIsoRhoCut'),
282  # Three constants for the GsfEleTrkPtIsoCut
283  # cut = constTerm if value < slopeStart
284  # cut = slopeTerm * (value - slopeStart) + constTerm if value >= slopeStart
285  slopeTermEB = cms.double( wpEB.trkIsoSlopeTerm ),
286  slopeTermEE = cms.double( wpEE.trkIsoSlopeTerm ),
287  slopeStartEB = cms.double( wpEB.trkIsoSlopeStart ),
288  slopeStartEE = cms.double( wpEE.trkIsoSlopeStart ),
289  constTermEB = cms.double( wpEB.trkIsoConstTerm ),
290  constTermEE = cms.double( wpEE.trkIsoConstTerm ),
291  rhoEtStartEB = cms.double( wpEB.trkIsoRhoCorrStart),
292  rhoEtStartEE = cms.double( wpEE.trkIsoRhoCorrStart),
293  rhoEAEB = cms.double( wpEB.trkIsoEffArea),
294  rhoEAEE = cms.double( wpEE.trkIsoEffArea),
295  rho = cms.InputTag("fixedGridRhoFastjetAll"),
296  needsAdditionalProducts = cms.bool(True),
297  isIgnored = cms.bool(False)
)
def heepElectronID_tools.psetGsfEleTrkPtNoJetCoreIsoCut (   wpEB,
  wpEE 
)

Definition at line 298 of file heepElectronID_tools.py.

Referenced by configureHEEPElectronID_V60_80XAOD().

299 def psetGsfEleTrkPtNoJetCoreIsoCut(wpEB, wpEE):
300  return cms.PSet(
301  cutName = cms.string('GsfEleValueMapIsoRhoCut'),
302  # Three constants for the GsfEleTrkPtIsoCut
303  # cut = constTerm if value < slopeStart
304  # cut = slopeTerm * (value - slopeStart) + constTerm if value >= slopeStart
305  slopeTermEB = cms.double( wpEB.trkIsoSlopeTerm ),
306  slopeTermEE = cms.double( wpEE.trkIsoSlopeTerm ),
307  slopeStartEB = cms.double( wpEB.trkIsoSlopeStart ),
308  slopeStartEE = cms.double( wpEE.trkIsoSlopeStart ),
309  constTermEB = cms.double( wpEB.trkIsoConstTerm ),
310  constTermEE = cms.double( wpEE.trkIsoConstTerm ),
311  #no rho so we zero it out, if the input tag is empty, its ignored anyways
312  rhoEtStartEB = cms.double( 999999.),
313  rhoEtStartEE = cms.double( 999999.),
314  rhoEAEB = cms.double( 0. ),
315  rhoEAEE = cms.double( 0. ),
316  rho = cms.InputTag(""),
317  value = cms.InputTag("heepIDVarValueMaps","eleTrkPtIsoNoJetCore"),
318  needsAdditionalProducts = cms.bool(True),
319  isIgnored = cms.bool(False)
)
def heepElectronID_tools.psetMinEcalEtCut (   cutValue = 5.)

Definition at line 148 of file heepElectronID_tools.py.

Referenced by configureHEEPElectronID_V71().

149 def psetMinEcalEtCut(cutValue=5.):
150  return cms.PSet(
151  cutName = cms.string("GsfEleMinEcalEtCut"),
152  minEt = cms.double(cutValue),
153  needsAdditionalProducts = cms.bool(False),
154  isIgnored = cms.bool(False)
155  )
156 
# Take all particles in the eta ranges 0-ebMax and eeMin-2.5
def heepElectronID_tools.psetMinPtCut (   cutValue = 35.)

Definition at line 139 of file heepElectronID_tools.py.

Referenced by configureHEEPElectronID_V51(), configureHEEPElectronID_V60(), configureHEEPElectronID_V60_80XAOD(), configureHEEPElectronID_V61(), and configureHEEPElectronID_V70().

140 def psetMinPtCut(cutValue=35.):
141  return cms.PSet(
142  cutName = cms.string("MinPtCut"),
143  minPt = cms.double(cutValue),
144  needsAdditionalProducts = cms.bool(False),
145  isIgnored = cms.bool(False)
146  )
147 
# The mininum pt cut is set to 5 GeV

Variable Documentation

float heepElectronID_tools.ebCutOff = 1.479

Definition at line 6 of file heepElectronID_tools.py.

float heepElectronID_tools.ebMax = 1.4442

Definition at line 4 of file heepElectronID_tools.py.

float heepElectronID_tools.eeMin = 1.566

Definition at line 5 of file heepElectronID_tools.py.