CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
METAnalyzer.py
Go to the documentation of this file.
1 from PhysicsTools.Heppy.analyzers.core.Analyzer import Analyzer
2 from PhysicsTools.Heppy.analyzers.core.AutoHandle import AutoHandle
3 from PhysicsTools.Heppy.physicsobjects.PhysicsObjects import Jet
5 from PhysicsTools.HeppyCore.statistics.counter import Counter, Counters
6 from PhysicsTools.Heppy.physicsutils.JetReCalibrator import Type1METCorrector, setFakeRawMETOnOldMiniAODs
7 import PhysicsTools.HeppyCore.framework.config as cfg
8 
9 import copy
10 import ROOT
11 from math import hypot
12 
13 from copy import deepcopy
14 
15 class METAnalyzer( Analyzer ):
16  def __init__(self, cfg_ana, cfg_comp, looperName ):
17  super(METAnalyzer,self).__init__(cfg_ana,cfg_comp,looperName)
18  self.recalibrateMET = cfg_ana.recalibrate
19  self.applyJetSmearing = cfg_comp.isMC and cfg_ana.applyJetSmearing
20  self.old74XMiniAODs = cfg_ana.old74XMiniAODs
21  self.jetAnalyzerPostFix = getattr(cfg_ana, 'jetAnalyzerPostFix', '')
22  if self.recalibrateMET in [ "type1", True ]:
23  if self.recalibrateMET == "type1":
24  self.type1METCorrector = Type1METCorrector(self.old74XMiniAODs)
25  elif self.recalibrateMET != False:
26  raise RuntimeError("Unsupported value %r for option 'recalibrate': allowed are True, False, 'type1'" % self.recalibrateMET)
27 
28  def declareHandles(self):
29  super(METAnalyzer, self).declareHandles()
30  self.handles['met'] = AutoHandle( self.cfg_ana.metCollection, 'std::vector<pat::MET>' )
31  if self.cfg_ana.doMetNoPU:
32  self.handles['nopumet'] = AutoHandle( self.cfg_ana.noPUMetCollection, 'std::vector<pat::MET>' )
33  if self.cfg_ana.doTkMet:
34  self.handles['cmgCand'] = AutoHandle( self.cfg_ana.candidates, self.cfg_ana.candidatesTypes )
35  #self.handles['vertices'] = AutoHandle( "offlineSlimmedPrimaryVertices", 'std::vector<reco::Vertex>', fallbackLabel="offlinePrimaryVertices" )
36  self.mchandles['packedGen'] = AutoHandle( 'packedGenParticles', 'std::vector<pat::PackedGenParticle>' )
37 
38  def beginLoop(self, setup):
39  super(METAnalyzer,self).beginLoop(setup)
40  self.counters.addCounter('events')
41  count = self.counters.counter('events')
42  count.register('all events')
43 
44  def applyDeltaMet(self, met, deltaMet):
45  px,py = self.met.px()+deltaMet[0], self.met.py()+deltaMet[1]
46  met.setP4(ROOT.reco.Particle.LorentzVector(px,py, 0, hypot(px,py)))
47 
48  def adduParaPerp(self, met, boson, postfix):
49 
50  upara = 0
51  uperp = 0
52  uX = - met.px() - boson.px()
53  uY = - met.py() - boson.py()
54  u1 = (uX*boson.px() + uY*boson.py())/boson.pt()
55  u2 = (uX*boson.py() - uY*boson.px())/boson.pt()
56 
57  setattr(met, "upara"+postfix, u1)
58  setattr(met, "uperp"+postfix, u2)
59 
60  def makeTkMETs(self, event):
61  charged = []
62  chargedchs = []
63  chargedPVLoose = []
64  chargedPVTight = []
65  dochs=getattr(self.cfg_ana,"includeTkMetCHS",True)
66  dotight=getattr(self.cfg_ana,"includeTkMetPVTight",True)
67  doloose=getattr(self.cfg_ana,"includeTkMetPVLoose",True)
68  pfcands = self.handles['cmgCand'].product()
69 
70  for pfcand in pfcands:
71 
72 ## ===> require the Track Candidate charge and with a minimum dz
73  if (pfcand.charge()!=0):
74 
75  pvflag = pfcand.fromPV()
76  pxy = pfcand.px(), pfcand.py()
77 
78  if abs(pfcand.dz())<=self.cfg_ana.dzMax:
79  charged.append(pxy)
80 
81  if dochs and pvflag>0:
82  chargedchs.append(pxy)
83 
84  if doloose and pvflag>1:
85  chargedPVLoose.append(pxy)
86 
87  if dotight and pvflag>2:
88  chargedPVTight.append(pxy)
89 
90  def sumXY(pxys):
91  px, py = sum(x[0] for x in pxys), sum(x[1] for x in pxys)
92  return ROOT.reco.Particle.LorentzVector(-px, -py, 0, hypot(px,py))
93  setattr(event, "tkMet"+self.cfg_ana.collectionPostFix, sumXY(charged))
94  setattr(event, "tkMetPVchs"+self.cfg_ana.collectionPostFix, sumXY(chargedchs))
95  setattr(event, "tkMetPVLoose"+self.cfg_ana.collectionPostFix, sumXY(chargedPVLoose))
96  setattr(event, "tkMetPVTight"+self.cfg_ana.collectionPostFix, sumXY(chargedPVTight))
97  getattr(event,"tkMet"+self.cfg_ana.collectionPostFix).sumEt = sum([hypot(x[0],x[1]) for x in charged])
98  getattr(event,"tkMetPVchs"+self.cfg_ana.collectionPostFix).sumEt = sum([hypot(x[0],x[1]) for x in chargedchs])
99  getattr(event,"tkMetPVLoose"+self.cfg_ana.collectionPostFix).sumEt = sum([hypot(x[0],x[1]) for x in chargedPVLoose])
100  getattr(event,"tkMetPVTight"+self.cfg_ana.collectionPostFix).sumEt = sum([hypot(x[0],x[1]) for x in chargedPVTight])
101 
102  if hasattr(event,'zll_p4'):
103  self.adduParaPerp(getattr(event,"tkMet"+self.cfg_ana.collectionPostFix), event.zll_p4,"_zll")
104  self.adduParaPerp(getattr(event,"tkMetPVchs"+self.cfg_ana.collectionPostFix), event.zll_p4,"_zll")
105  self.adduParaPerp(getattr(event,"tkMetPVLoose"+self.cfg_ana.collectionPostFix), event.zll_p4,"_zll")
106  self.adduParaPerp(getattr(event,"tkMetPVTight"+self.cfg_ana.collectionPostFix), event.zll_p4,"_zll")
107 
108  def makeGenTkMet(self, event):
109  genCharged = [ (x.px(),x.py()) for x in self.mchandles['packedGen'].product() if x.charge() != 0 and abs(x.eta()) < 2.4 ]
110  px, py = sum(x[0] for x in genCharged), sum(x[1] for x in genCharged)
111  setattr(event,"tkGenMet"+self.cfg_ana.collectionPostFix, ROOT.reco.Particle.LorentzVector(-px , -py, 0, hypot(px,py)))
112 
113  def makeMETNoMu(self, event):
114  self.metNoMu = copy.deepcopy(self.met)
115  if self.cfg_ana.doMetNoPU: self.metNoMuNoPU = copy.deepcopy(self.metNoPU)
116 
117  mupx = 0
118  mupy = 0
119  #sum muon momentum
120  for mu in event.selectedMuons:
121  mupx += mu.px()
122  mupy += mu.py()
123 
124  #subtract muon momentum and construct met
125  px,py = self.metNoMu.px()+mupx, self.metNoMu.py()+mupy
126  self.metNoMu.setP4(ROOT.reco.Particle.LorentzVector(px,py, 0, hypot(px,py)))
127  px,py = self.metNoMuNoPU.px()+mupx, self.metNoMuNoPU.py()+mupy
128  self.metNoMuNoPU.setP4(ROOT.reco.Particle.LorentzVector(px,py, 0, hypot(px,py)))
129  setattr(event, "metNoMu"+self.cfg_ana.collectionPostFix, self.metNoMu)
130  if self.cfg_ana.doMetNoPU: setattr(event, "metNoMuNoPU"+self.cfg_ana.collectionPostFix, self.metNoMuNoPU)
131 
132 
133  def makeMETNoEle(self, event):
134  self.metNoEle = copy.deepcopy(self.met)
135  if self.cfg_ana.doMetNoPU: self.metNoEleNoPU = copy.deepcopy(self.metNoPU)
136 
137  elepx = 0
138  elepy = 0
139  #sum electron momentum
140  for ele in event.selectedElectrons:
141  elepx += ele.px()
142  elepy += ele.py()
143 
144  #subtract electron momentum and construct met
145  px,py = self.metNoEle.px()+elepx, self.metNoEle.py()+elepy
146  self.metNoEle.setP4(ROOT.reco.Particle.LorentzVector(px,py, 0, hypot(px,py)))
147 
148  px,py = self.metNoEleNoPU.px()+elepx, self.metNoEleNoPU.py()+elepy
149  self.metNoEleNoPU.setP4(ROOT.reco.Particle.LorentzVector(px,py, 0, hypot(px,py)))
150  setattr(event, "metNoEle"+self.cfg_ana.collectionPostFix, self.metNoEle)
151  if self.cfg_ana.doMetNoPU: setattr(event, "metNoEleNoPU"+self.cfg_ana.collectionPostFix, self.metNoEleNoPU)
152 
153  def makeMETNoPhoton(self, event):
154  self.metNoPhoton = copy.deepcopy(self.met)
155 
156  phopx = 0
157  phopy = 0
158  #sum photon momentum
159  for pho in event.selectedPhotons:
160  phopx += pho.px()
161  phopy += pho.py()
162 
163  #subtract photon momentum and construct met
164  px,py = self.metNoPhoton.px()+phopx, self.metNoPhoton.py()+phopy
165  self.metNoPhoton.setP4(ROOT.reco.Particle.LorentzVector(px,py, 0, hypot(px,py)))
166  setattr(event, "metNoPhoton"+self.cfg_ana.collectionPostFix, self.metNoPhoton)
167  if self.cfg_ana.doMetNoPU:
168  self.metNoPhotonNoPU = copy.deepcopy(self.metNoPU)
169  px,py = self.metNoPhotonNoPU.px()+phopx, self.metNoPhotonNoPU.py()+phopy
170  self.metNoPhotonNoPU.setP4(ROOT.reco.Particle.LorentzVector(px,py, 0, hypot(px,py)))
171  setattr(event, "metNoPhotonNoPU"+self.cfg_ana.collectionPostFix, self.metNoPhotonNoPU)
172 
173 
174  def makeMETs(self, event):
175  import ROOT
176  if self.cfg_ana.copyMETsByValue:
177  self.met = ROOT.pat.MET(self.handles['met'].product()[0])
178  if self.cfg_ana.doMetNoPU: self.metNoPU = ROOT.pat.MET(self.handles['nopumet'].product()[0])
179  else:
180  self.met = self.handles['met'].product()[0]
181  if self.cfg_ana.doMetNoPU: self.metNoPU = self.handles['nopumet'].product()[0]
182 
183  if self.recalibrateMET == "type1":
184  type1METCorr = getattr(event, 'type1METCorr'+self.jetAnalyzerPostFix)
185  self.type1METCorrector.correct(self.met, type1METCorr)
186  elif self.recalibrateMET == True:
187  deltaMetJEC = getattr(event, 'deltaMetFromJEC'+self.jetAnalyzerPostFix)
188  self.applyDeltaMet(self.met, deltaMetJEC)
189  if self.applyJetSmearing:
190  deltaMetSmear = getattr(event, 'deltaMetFromJetSmearing'+self.jetAnalyzerPostFix)
191  self.applyDeltaMet(self.met, deltaMetSmear)
192 
193 
194  if (not self.cfg_ana.copyMETsByValue) and getattr(self.cfg_ana, 'makeShiftedMETs', True):
195  shifts = []
196  for obj in 'JetEn', 'JetRes', 'MuonEn', 'ElectronEn', 'PhotonEn', 'TauEn', 'UnclusteredEn':
197  for sh in 'Up', 'Down':
198  shifts.append( (obj+sh, getattr(self.met,obj+sh)) )
199  shifts.append( ('NoShift', self.met.NoShift) )
200  for name,i in shifts:
201  key = i
202  m = ROOT.pat.MET(self.met)
203  if self.old74XMiniAODs:
204  if key > 12: key = 12
205  elif key <= 3: key = { 'JetEnUp':0, 'JetEnDown':1, 'JetResUp':2, 'JetResDown':3 }[name]
206  m.setP4(self.met.shiftedP4(key))
207  setattr(event, "met{0}_shifted_{1}".format(self.cfg_ana.collectionPostFix, i),m)
208  setattr(event, "met{0}_shifted_{1}".format(self.cfg_ana.collectionPostFix, name),m)
209 
210  self.met_sig = self.met.significance()
211  self.met_sumet = self.met.sumEt()
212 
213  if self.old74XMiniAODs and self.recalibrateMET != "type1":
214  oldraw = self.met.shiftedP2_74x(12,0);
215  setFakeRawMETOnOldMiniAODs( self.met, oldraw.px, oldraw.py, self.met.shiftedSumEt_74x(12,0) )
216  px, py = oldraw.px, oldraw.py
217  else:
218  px, py = self.met.uncorPx(), self.met.uncorPy()
219  self.met_raw = ROOT.reco.Particle.LorentzVector(px,py,0,hypot(px,py))
220 
221  if hasattr(event,'zll_p4'):
222  self.adduParaPerp(self.met,event.zll_p4,"_zll")
223  self.adduParaPerp(self.met_raw, event.zll_p4,"_zll")
224  setattr(event,"met_raw"+self.cfg_ana.collectionPostFix, self.met_raw)
225  setattr(event,"met_raw.upara_zll"+self.cfg_ana.collectionPostFix, self.met_raw.upara_zll)
226  setattr(event,"met_raw.uperp_zll"+self.cfg_ana.collectionPostFix, self.met_raw.uperp_zll)
227 
228  if hasattr(event,"met"+self.cfg_ana.collectionPostFix): raise RuntimeError("Event already contains met with the following postfix: "+self.cfg_ana.collectionPostFix)
229  setattr(event, "met"+self.cfg_ana.collectionPostFix, self.met)
230  if self.cfg_ana.doMetNoPU: setattr(event, "metNoPU"+self.cfg_ana.collectionPostFix, self.metNoPU)
231  setattr(event, "met_sig"+self.cfg_ana.collectionPostFix, self.met_sig)
232  setattr(event, "met_sumet"+self.cfg_ana.collectionPostFix, self.met_sumet)
233 
234  genMET = self.met.genMET()
235  if genMET:
236  setattr(event, "met_genPt"+self.cfg_ana.collectionPostFix, genMET.pt())
237  setattr(event, "met_genPhi"+self.cfg_ana.collectionPostFix, genMET.phi())
238  else:
239  setattr(event, "met_genPt"+self.cfg_ana.collectionPostFix, float('nan'))
240  setattr(event, "met_genPhi"+self.cfg_ana.collectionPostFix, float('nan'))
241 
242  if self.cfg_ana.doMetNoMu and hasattr(event, 'selectedMuons'):
243  self.makeMETNoMu(event)
244 
245  if self.cfg_ana.doMetNoEle and hasattr(event, 'selectedElectrons'):
246  self.makeMETNoEle(event)
247 
248  if self.cfg_ana.doMetNoPhoton and hasattr(event, 'selectedPhotons'):
249  self.makeMETNoPhoton(event)
250 
251  def process(self, event):
252  self.readCollections( event.input)
253  self.counters.counter('events').inc('all events')
254 
255  self.makeMETs(event)
256 
257  if self.cfg_ana.doTkMet:
258  self.makeTkMETs(event);
259 
260  if getattr(self.cfg_ana,"doTkGenMet",self.cfg_ana.doTkMet) and self.cfg_comp.isMC and hasattr(event, 'genParticles'):
261  self.makeGenTkMet(event)
262 
263  return True
264 
265 
266 setattr(METAnalyzer,"defaultConfig", cfg.Analyzer(
267  class_object = METAnalyzer,
268  metCollection = "slimmedMETs",
269  noPUMetCollection = "slimmedMETs",
270  copyMETsByValue = False,
271  recalibrate = True,
272  applyJetSmearing = True,
273  jetAnalyzerPostFix = "",
274  old74XMiniAODs = False, # need to set to True to get proper Raw MET on plain 74X MC produced with CMSSW <= 7_4_12
275  makeShiftedMETs = True,
276  doTkMet = False,
277  includeTkMetCHS = True,
278  includeTkMetPVLoose = True,
279  includeTkMetPVTight = True,
280  doMetNoPU = False, # Not existing in MiniAOD at the moment
281  doMetNoMu = False,
282  doMetNoEle = False,
283  doMetNoPhoton = False,
284  candidates='packedPFCandidates',
285  candidatesTypes='std::vector<pat::PackedCandidate>',
286  dzMax = 0.1,
287  collectionPostFix = "",
288  )
289 )
Abs< T >::type abs(const T &t)
Definition: Abs.h:22