CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
lowPtElectrons_cff.py
Go to the documentation of this file.
1 import FWCore.ParameterSet.Config as cms
4 
5 ################################################################################
6 # Modules
7 ################################################################################
8 
9 from RecoEgamma.EgammaTools.lowPtElectronModifier_cfi import lowPtElectronModifier
10 from RecoEgamma.EgammaElectronProducers.lowPtGsfElectrons_cfi import lowPtRegressionModifier
11 modifiedLowPtElectrons = cms.EDProducer(
12  "ModifiedElectronProducer",
13  src = cms.InputTag("slimmedLowPtElectrons"),
14  modifierConfig = cms.PSet(
15  modifications = cms.VPSet(lowPtElectronModifier,lowPtRegressionModifier)
16  )
17 )
18 
20 updatedLowPtElectrons = cms.EDProducer(
21  "PATElectronUpdater",
22  src = cms.InputTag("modifiedLowPtElectrons"),
23  vertices = cms.InputTag("offlineSlimmedPrimaryVertices"),
24  computeMiniIso = cms.bool(True),
25  fixDxySign = cms.bool(False),
26  pfCandsForMiniIso = cms.InputTag("packedPFCandidates"),
27  miniIsoParamsB = PhysicsTools.PatAlgos.producersLayer1.electronProducer_cfi.patElectrons.miniIsoParamsB,
28  miniIsoParamsE = PhysicsTools.PatAlgos.producersLayer1.electronProducer_cfi.patElectrons.miniIsoParamsE,
29 )
30 
31 from RecoEgamma.EgammaElectronProducers.lowPtGsfElectronID_cfi import lowPtGsfElectronID
32 lowPtPATElectronID = lowPtGsfElectronID.clone(
33  usePAT = True,
34  electrons = "updatedLowPtElectrons",
35  unbiased = "",
36  ModelWeights = [
37  'RecoEgamma/ElectronIdentification/data/LowPtElectrons/LowPtElectrons_ID_2020Nov28.root',
38  ],
39 )
40 
41 isoForLowPtEle = cms.EDProducer(
42  "EleIsoValueMapProducer",
43  src = cms.InputTag("updatedLowPtElectrons"),
44  relative = cms.bool(True),
45  rho_MiniIso = cms.InputTag("fixedGridRhoFastjetAll"),
46  rho_PFIso = cms.InputTag("fixedGridRhoFastjetAll"),
47  EAFile_MiniIso = cms.FileInPath("RecoEgamma/ElectronIdentification/data/Fall17/effAreaElectrons_cone03_pfNeuHadronsAndPhotons_94X.txt"),
48  EAFile_PFIso = cms.FileInPath("RecoEgamma/ElectronIdentification/data/Fall17/effAreaElectrons_cone03_pfNeuHadronsAndPhotons_94X.txt"),
49 )
50 
51 updatedLowPtElectronsWithUserData = cms.EDProducer(
52  "PATElectronUserDataEmbedder",
53  src = cms.InputTag("updatedLowPtElectrons"),
54  userFloats = cms.PSet(
55  ID = cms.InputTag("lowPtPATElectronID"),
56  miniIsoChg = cms.InputTag("isoForLowPtEle:miniIsoChg"),
57  miniIsoAll = cms.InputTag("isoForLowPtEle:miniIsoAll"),
58  ),
59  userIntFromBools = cms.PSet(),
60  userInts = cms.PSet(),
61  userCands = cms.PSet(),
62 )
63 
64 finalLowPtElectrons = cms.EDFilter(
65  "PATElectronRefSelector",
66  src = cms.InputTag("updatedLowPtElectronsWithUserData"),
67  cut = cms.string("pt > 1. && userFloat('ID') > -0.25"),
68 )
69 
70 ################################################################################
71 # electronTable
72 ################################################################################
73 
74 lowPtElectronTable = cms.EDProducer(
75  "SimpleCandidateFlatTableProducer",
76  src = cms.InputTag("finalLowPtElectrons"),
77  cut = cms.string(""),
78  name= cms.string("LowPtElectron"),
79  doc = cms.string("slimmedLowPtElectrons after basic selection (" + finalLowPtElectrons.cut.value()+")"),
80  singleton = cms.bool(False), # the number of entries is variable
81  extension = cms.bool(False), # this is the main table for the electrons
82  variables = cms.PSet(
83  # Basic variables
84  CandVars,
85  # BDT scores and WPs
86  embeddedID = Var("electronID('ID')",float,doc="ID, BDT (raw) score"),
87  ID = Var("userFloat('ID')",float,doc="New ID, BDT (raw) score"),
88  unbiased = Var("electronID('unbiased')",float,doc="ElectronSeed, pT- and dxy- agnostic BDT (raw) score"),
89  ptbiased = Var("electronID('ptbiased')",float,doc="ElectronSeed, pT- and dxy- dependent BDT (raw) score"),
90  # Isolation
91  miniPFRelIso_chg = Var("userFloat('miniIsoChg')",float,
92  doc="mini PF relative isolation, charged component"),
93  miniPFRelIso_all = Var("userFloat('miniIsoAll')",float,
94  doc="mini PF relative isolation, total (with scaled rho*EA PU corrections)"),
95  # Conversions
96  convVeto = Var("passConversionVeto()",bool,doc="pass conversion veto"),
97  convWP = Var("userInt('convOpen')*1 + userInt('convLoose')*2 + userInt('convTight')*4",
98  int,doc="conversion flag bit map: 1=Veto, 2=Loose, 3=Tight"),
99  convVtxRadius = Var("userFloat('convVtxRadius')",float,doc="conversion vertex radius (cm)",precision=7),
100  # Tracking
101  lostHits = Var("gsfTrack.hitPattern.numberOfLostHits('MISSING_INNER_HITS')","uint8",doc="number of missing inner hits"),
102  # Cluster-related
103  energyErr = Var("p4Error('P4_COMBINATION')",float,doc="energy error of the cluster-track combination",precision=6),
104  deltaEtaSC = Var("superCluster().eta()-eta()",float,doc="delta eta (SC,ele) with sign",precision=10),
105  r9 = Var("full5x5_r9()",float,doc="R9 of the SC, calculated with full 5x5 region",precision=10),
106  sieie = Var("full5x5_sigmaIetaIeta()",float,doc="sigma_IetaIeta of the SC, calculated with full 5x5 region",precision=10),
107  eInvMinusPInv = Var("(1-eSuperClusterOverP())/ecalEnergy()",float,doc="1/E_SC - 1/p_trk",precision=10),
108  scEtOverPt = Var("(superCluster().energy()/(pt*cosh(superCluster().eta())))-1",float,doc="(SC energy)/pt-1",precision=8),
109  hoe = Var("hadronicOverEm()",float,doc="H over E",precision=8),
110  # Displacement
111  dxy = Var("dB('PV2D')",float,doc="dxy (with sign) wrt first PV, in cm",precision=10),
112  dxyErr = Var("edB('PV2D')",float,doc="dxy uncertainty, in cm",precision=6),
113  dz = Var("dB('PVDZ')",float,doc="dz (with sign) wrt first PV, in cm",precision=10),
114  dzErr = Var("abs(edB('PVDZ'))",float,doc="dz uncertainty, in cm",precision=6),
115  # Cross-referencing
116  #jetIdx
117  #photonIdx
118  ),
119 )
120 
121 ################################################################################
122 # electronTable (MC)
123 ################################################################################
124 
125 # Depends on tautaggerForMatching being run in electrons_cff
126 matchingLowPtElecPhoton = cms.EDProducer(
127  "GenJetGenPartMerger",
128  srcJet =cms.InputTag("particleLevel:leptons"),
129  srcPart=cms.InputTag("particleLevel:photons"),
130  cut = cms.string(""),
131  hasTauAnc=cms.InputTag("tautaggerForMatching"),
132 )
133 
134 lowPtElectronsMCMatchForTableAlt = cms.EDProducer(
135  "GenJetMatcherDRPtByDR", # cut on deltaR, deltaPt/Pt; pick best by deltaR
136  src = lowPtElectronTable.src, # final reco collection
137  matched = cms.InputTag("matchingLowPtElecPhoton:merged"), # final mc-truth particle collection
138  mcPdgId = cms.vint32(11,22), # one or more PDG ID (11 = el, 22 = pho); absolute values (see below)
139  checkCharge = cms.bool(False), # True = require RECO and MC objects to have the same charge
140  mcStatus = cms.vint32(),
141  maxDeltaR = cms.double(0.3), # Minimum deltaR for the match
142  maxDPtRel = cms.double(0.5), # Minimum deltaPt/Pt for the match
143  resolveAmbiguities = cms.bool(True), # Forbid two RECO objects to match to the same GEN object
144  resolveByMatchQuality = cms.bool(True), # False = just match input in order; True = pick lowest deltaR pair first
145 )
146 
147 lowPtElectronsMCMatchForTable = cms.EDProducer(
148  "MCMatcher", # cut on deltaR, deltaPt/Pt; pick best by deltaR
149  src = lowPtElectronTable.src, # final reco collection
150  matched = cms.InputTag("finalGenParticles"), # final mc-truth particle collection
151  mcPdgId = cms.vint32(11), # one or more PDG ID (11 = ele); absolute values (see below)
152  checkCharge = cms.bool(False), # True = require RECO and MC objects to have the same charge
153  mcStatus = cms.vint32(1), # PYTHIA status code (1 = stable, 2 = shower, 3 = hard scattering)
154  maxDeltaR = cms.double(0.3), # Minimum deltaR for the match
155  maxDPtRel = cms.double(0.5), # Minimum deltaPt/Pt for the match
156  resolveAmbiguities = cms.bool(True), # Forbid two RECO objects to match to the same GEN object
157  resolveByMatchQuality = cms.bool(True), # False = just match input in order; True = pick lowest deltaR pair first
158 )
159 
160 from PhysicsTools.NanoAOD.electrons_cff import electronMCTable
161 lowPtElectronMCTable = cms.EDProducer(
162  "CandMCMatchTableProducer",
163  src = lowPtElectronTable.src,
164  mcMapDressedLep = cms.InputTag("lowPtElectronsMCMatchForTableAlt"),
165  mcMap = cms.InputTag("lowPtElectronsMCMatchForTable"),
166  mapTauAnc = cms.InputTag("matchingLowPtElecPhoton:hasTauAnc"),
167  objName = lowPtElectronTable.name,
168  objType = electronMCTable.objType,
169  branchName = cms.string("genPart"),
170  docString = cms.string("MC matching to status==1 electrons or photons"),
171  genparticles = cms.InputTag("finalGenParticles"),
172 )
173 
174 ################################################################################
175 # Tasks
176 ################################################################################
177 
178 lowPtElectronTask = cms.Task(modifiedLowPtElectrons,
179  updatedLowPtElectrons,
180  lowPtPATElectronID,
181  isoForLowPtEle,
182  updatedLowPtElectronsWithUserData,
183  finalLowPtElectrons)
184 
185 lowPtElectronTablesTask = cms.Task(lowPtElectronTable)
186 
187 lowPtElectronMCTask = cms.Task(
188  matchingLowPtElecPhoton,
189  lowPtElectronsMCMatchForTable,
190  lowPtElectronsMCMatchForTableAlt,
191  lowPtElectronMCTable)
192 
193 ################################################################################
194 # Modifiers
195 ################################################################################
196 
197 _modifiers = ( run2_miniAOD_80XLegacy |
198  run2_nanoAOD_94XMiniAODv1 |
199  run2_nanoAOD_94XMiniAODv2 |
200  run2_nanoAOD_94X2016 |
201  run2_nanoAOD_102Xv1 |
202  run2_nanoAOD_106Xv1 )
203 (_modifiers).toReplaceWith(lowPtElectronTask,cms.Task())
204 (_modifiers).toReplaceWith(lowPtElectronTablesTask,cms.Task())
205 (_modifiers).toReplaceWith(lowPtElectronMCTask,cms.Task())