CMS 3D CMS Logo

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