CMS 3D CMS Logo

JetReCalibrator.py
Go to the documentation of this file.
1 import ROOT
2 import os, types
3 from math import *
5 
7  def __init__(self,globalTag,jetFlavour,doResidualJECs,jecPath,upToLevel=3,
8  calculateSeparateCorrections=False,
9  calculateType1METCorrection=False, type1METParams={'jetPtThreshold':15., 'skipEMfractionThreshold':0.9, 'skipMuons':True} ):
10  """Create a corrector object that reads the payloads from the text dumps of a global tag under
11  CMGTools/RootTools/data/jec (see the getJec.py there to make the dumps).
12  It will apply the L1,L2,L3 and possibly the residual corrections to the jets.
13  If configured to do so, it will also compute the type1 MET corrections."""
14  self.globalTag = globalTag
15  self.jetFlavour = jetFlavour
16  self.doResidualJECs = doResidualJECs
17  self.jecPath = jecPath
18  self.upToLevel = upToLevel
19  self.calculateType1METCorr = calculateType1METCorrection
20  self.type1METParams = type1METParams
21  # Make base corrections
22  path = os.path.expandvars(jecPath) #"%s/src/CMGTools/RootTools/data/jec" % os.environ['CMSSW_BASE'];
23  self.L1JetPar = ROOT.JetCorrectorParameters("%s/%s_L1FastJet_%s.txt" % (path,globalTag,jetFlavour),"");
24  self.L2JetPar = ROOT.JetCorrectorParameters("%s/%s_L2Relative_%s.txt" % (path,globalTag,jetFlavour),"");
25  self.L3JetPar = ROOT.JetCorrectorParameters("%s/%s_L3Absolute_%s.txt" % (path,globalTag,jetFlavour),"");
26  self.vPar = ROOT.vector(ROOT.JetCorrectorParameters)()
27  self.vPar.push_back(self.L1JetPar);
28  if upToLevel >= 2: self.vPar.push_back(self.L2JetPar);
29  if upToLevel >= 3: self.vPar.push_back(self.L3JetPar);
30  # Add residuals if needed
31  if doResidualJECs :
32  self.ResJetPar = ROOT.JetCorrectorParameters("%s/%s_L2L3Residual_%s.txt" % (path,globalTag,jetFlavour))
33  self.vPar.push_back(self.ResJetPar);
34  #Step3 (Construct a FactorizedJetCorrector object)
35  self.JetCorrector = ROOT.FactorizedJetCorrector(self.vPar)
36  if os.path.exists("%s/%s_Uncertainty_%s.txt" % (path,globalTag,jetFlavour)):
37  self.JetUncertainty = ROOT.JetCorrectionUncertainty("%s/%s_Uncertainty_%s.txt" % (path,globalTag,jetFlavour));
38  elif os.path.exists("%s/Uncertainty_FAKE.txt" % path):
39  self.JetUncertainty = ROOT.JetCorrectionUncertainty("%s/Uncertainty_FAKE.txt" % path);
40  else:
41  print 'Missing JEC uncertainty file "%s/%s_Uncertainty_%s.txt", so jet energy uncertainties will not be available' % (path,globalTag,jetFlavour)
42  self.JetUncertainty = None
44  if calculateSeparateCorrections or calculateType1METCorrection:
45  self.vParL1 = ROOT.vector(ROOT.JetCorrectorParameters)()
46  self.vParL1.push_back(self.L1JetPar)
47  self.separateJetCorrectors["L1"] = ROOT.FactorizedJetCorrector(self.vParL1)
48  if upToLevel >= 2 and calculateSeparateCorrections:
49  self.vParL2 = ROOT.vector(ROOT.JetCorrectorParameters)()
50  for i in [self.L1JetPar,self.L2JetPar]: self.vParL2.push_back(i)
51  self.separateJetCorrectors["L1L2"] = ROOT.FactorizedJetCorrector(self.vParL2)
52  if upToLevel >= 3 and calculateSeparateCorrections:
53  self.vParL3 = ROOT.vector(ROOT.JetCorrectorParameters)()
54  for i in [self.L1JetPar,self.L2JetPar,self.L3JetPar]: self.vParL3.push_back(i)
55  self.separateJetCorrectors["L1L2L3"] = ROOT.FactorizedJetCorrector(self.vParL3)
56  if doResidualJECs and calculateSeparateCorrections:
57  self.vParL3Res = ROOT.vector(ROOT.JetCorrectorParameters)()
58  for i in [self.L1JetPar,self.L2JetPar,self.L3JetPar,self.ResJetPar]: self.vParL3Res.push_back(i)
59  self.separateJetCorrectors["L1L2L3Res"] = ROOT.FactorizedJetCorrector(self.vParL3Res)
60 
61  def getCorrection(self,jet,rho,delta=0,corrector=None):
62  if not corrector: corrector = self.JetCorrector
63  if corrector != self.JetCorrector and delta!=0: raise RuntimeError('Configuration not supported')
64  corrector.setJetEta(jet.eta())
65  corrector.setJetPt(jet.pt()*jet.rawFactor())
66  corrector.setJetA(jet.jetArea())
67  corrector.setRho(rho)
68  corr = corrector.getCorrection()
69  if delta != 0:
70  if not self.JetUncertainty: raise RuntimeError("Jet energy scale uncertainty shifts requested, but not available")
71  self.JetUncertainty.setJetEta(jet.eta())
72  self.JetUncertainty.setJetPt(corr * jet.pt() * jet.rawFactor())
73  try:
74  jet.jetEnergyCorrUncertainty = self.JetUncertainty.getUncertainty(True)
75  except RuntimeError as r:
76  print "Caught %s when getting uncertainty for jet of pt %.1f, eta %.2f\n" % (r,corr * jet.pt() * jet.rawFactor(),jet.eta())
77  jet.jetEnergyCorrUncertainty = 0.5
78  #print " jet with corr pt %6.2f has an uncertainty %.2f " % (jet.pt()*jet.rawFactor()*corr, jet.jetEnergyCorrUncertainty)
79  corr *= max(0, 1+delta*jet.jetEnergyCorrUncertainty)
80  return corr
81 
82  def rawP4forType1MET_(self, jet):
83  """Return the raw 4-vector, after subtracting the muons (if requested),
84  or None if the jet fails the EMF cut."""
85  p4 = jet.p4() * jet.rawFactor()
86  emf = ( jet.physObj.neutralEmEnergy() + jet.physObj.chargedEmEnergy() )/p4.E()
87  if emf > self.type1METParams['skipEMfractionThreshold']:
88  return None
89  if self.type1METParams['skipMuons']:
90  for idau in xrange(jet.numberOfDaughters()):
91  pfcand = jet.daughter(idau)
92  if pfcand.isGlobalMuon() or pfcand.isStandAloneMuon():
93  p4 -= pfcand.p4()
94  return p4
95 
96  def correct(self,jet,rho,delta=0,addCorr=False,addShifts=False, metShift=[0,0],type1METCorr=[0,0,0]):
97  """Corrects a jet energy (optionally shifting it also by delta times the JEC uncertainty)
98 
99  If addCorr, set jet.corr to the correction.
100  If addShifts, set also the +1 and -1 jet shifts
101 
102  The metShift vector will accumulate the x and y changes to the MET from the JEC, i.e. the
103  negative difference between the new and old jet momentum, for jets eligible for type1 MET
104  corrections, and after subtracting muons. The pt cut is applied to the new corrected pt.
105  This shift can be applied on top of the *OLD TYPE1 MET*, but only if there was no change
106  in the L1 corrections nor in the definition of the type1 MET (e.g. jet pt cuts).
107 
108  The type1METCorr vector, will accumulate the x, y, sumEt type1 MET corrections, to be
109  applied to the *RAW MET* (if the feature was turned on in the constructor of the class).
110  """
111  raw = jet.rawFactor()
112  corr = self.getCorrection(jet,rho,delta)
113  if addCorr:
114  jet.corr = corr
115  for sepcorr in self.separateJetCorrectors.keys():
116  setattr(jet,"CorrFactor_"+sepcorr,self.getCorrection(jet,rho,delta=0,corrector=self.separateJetCorrectors[sepcorr]))
117  if addShifts:
118  for cdelta,shift in [(1.0, "JECUp"), (-1.0, "JECDown")]:
119  cshift = self.getCorrection(jet,rho,delta+cdelta)
120  setattr(jet, "corr"+shift, cshift)
121  if corr <= 0:
122  return False
123  newpt = jet.pt()*raw*corr
124  if newpt > self.type1METParams['jetPtThreshold']:
125  rawP4forT1 = self.rawP4forType1MET_(jet)
126  if rawP4forT1 and rawP4forT1.Pt()*corr > self.type1METParams['jetPtThreshold']:
127  metShift[0] -= rawP4forT1.Px() * (corr - 1.0/raw)
128  metShift[1] -= rawP4forT1.Py() * (corr - 1.0/raw)
129  if self.calculateType1METCorr:
130  l1corr = self.getCorrection(jet,rho,delta=0,corrector=self.separateJetCorrectors["L1"])
131  #print "\tfor jet with raw pt %.5g, eta %.5g, dpx = %.5g, dpy = %.5g" % (
132  # jet.pt()*raw, jet.eta(),
133  # rawP4forT1.Px()*(corr - l1corr),
134  # rawP4forT1.Py()*(corr - l1corr))
135  type1METCorr[0] -= rawP4forT1.Px() * (corr - l1corr)
136  type1METCorr[1] -= rawP4forT1.Py() * (corr - l1corr)
137  type1METCorr[2] += rawP4forT1.Et() * (corr - l1corr)
138  jet.setCorrP4(jet.p4() * (corr * raw))
139  return True
140 
141  def correctAll(self,jets,rho,delta=0, addCorr=False, addShifts=False, metShift=[0.,0.], type1METCorr=[0.,0.,0.]):
142  """Applies 'correct' to all the jets, discard the ones that have bad corrections (corrected pt <= 0)"""
143  badJets = []
144  if metShift != [0.,0. ]: raise RuntimeError("input metShift tuple is not initialized to zeros")
145  if type1METCorr != [0.,0.,0.]: raise RuntimeError("input type1METCorr tuple is not initialized to zeros")
146  for j in jets:
147  ok = self.correct(j,rho,delta,addCorr=addCorr,addShifts=addShifts,metShift=metShift,type1METCorr=type1METCorr)
148  if not ok: badJets.append(j)
149  if len(badJets) > 0:
150  print "Warning: %d out of %d jets flagged bad by JEC." % (len(badJets), len(jets))
151  for bj in badJets:
152  jets.remove(bj)
153 
155  def __init__(self, old74XMiniAODs):
156  """Object to apply type1 corrections computed by the JetReCalibrator to the MET.
157  old74XMiniAODs should be True if using inputs produced with CMSSW_7_4_11 or earlier."""
158  self.oldMC = old74XMiniAODs
159  def correct(self,met,type1METCorrections):
160  oldpx, oldpy = met.px(), met.py()
161  #print "old met: px %+10.5f, py %+10.5f" % (oldpx, oldpy)
162  if self.oldMC:
163  raw = met.shiftedP2_74x(12,0);
164  rawsumet = met.shiftedSumEt_74x(12,0);
165  else:
166  raw = met.uncorP2()
167  rawsumet = met.uncorSumEt();
168  rawpx, rawpy = raw.px, raw.py
169  #print "raw met: px %+10.5f, py %+10.5f" % (rawpx, rawpy)
170  corrpx = rawpx + type1METCorrections[0]
171  corrpy = rawpy + type1METCorrections[1]
172  corrsumet = rawsumet + type1METCorrections[2]
173  #print "done met: px %+10.5f, py %+10.5f\n" % (corrpx,corrpy)
174  met.setP4(ROOT.reco.Particle.LorentzVector(corrpx,corrpy,0,hypot(corrpx,corrpy)))
175  ## workaround for missing setSumEt in reco::MET and pat::MET
176  met._sumEt = corrsumet
177  met.sumEt = types.MethodType(lambda myself : myself._sumEt, met, met.__class__)
178  if not self.oldMC:
179  met.setCorShift(rawpx, rawpy, rawsumet, met.Raw)
180  else:
181  # to avoid segfauls in pat::MET, I need a more ugly solution
182  setFakeRawMETOnOldMiniAODs(met, rawpx,rawpy, rawsumet)
183 
184 def setFakeRawMETOnOldMiniAODs(met, rawpx, rawpy, rawsumet):
185  met._rawSumEt = rawsumet
186  met._rawP4 = ROOT.reco.Particle.LorentzVector(rawpx,rawpy,0,hypot(rawpx,rawpy))
187  met.uncorPt = types.MethodType(lambda myself : myself._rawP4.Pt(), met, met.__class__)
188  met.uncorPx = types.MethodType(lambda myself : myself._rawP4.Px(), met, met.__class__)
189  met.uncorPy = types.MethodType(lambda myself : myself._rawP4.Py(), met, met.__class__)
190  met.uncorPhi = types.MethodType(lambda myself : myself._rawP4.Phi(), met, met.__class__)
191  met.uncorP4 = types.MethodType(lambda myself : myself._rawP4, met, met.__class__)
192  met.uncorSumEt = types.MethodType(lambda myself : myself._rawSumEt, met, met.__class__)
193  # the two below are a bit more tricky, but probably less needed, but something dummy
194  met.uncorP2 = types.MethodType(lambda myself : None, met, met.__class__)
195  met.uncorP3 = types.MethodType(lambda myself : None, met, met.__class__)
196 
def correctAll(self, jets, rho, delta=0, addCorr=False, addShifts=False, metShift=[0., type1METCorr=[0.)
def correct(self, jet, rho, delta=0, addCorr=False, addShifts=False, metShift=[0, type1METCorr=[0)
def getCorrection(self, jet, rho, delta=0, corrector=None)
def correct(self, met, type1METCorrections)
def __init__(self, old74XMiniAODs)
def setFakeRawMETOnOldMiniAODs(met, rawpx, rawpy, rawsumet)