CMS 3D CMS Logo

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