1 from __future__
import print_function
2 from builtins
import range
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."""
24 path = os.path.expandvars(jecPath)
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)()
30 if upToLevel >= 2: self.vPar.push_back(self.
L2JetPar);
31 if upToLevel >= 3: self.vPar.push_back(self.
L3JetPar);
34 self.
ResJetPar = ROOT.JetCorrectorParameters(
"%s/%s_L2L3Residual_%s.txt" % (path,globalTag,jetFlavour))
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);
43 print(
'Missing JEC uncertainty file "%s/%s_Uncertainty_%s.txt", so jet energy uncertainties will not be available' % (path,globalTag,jetFlavour))
46 if calculateSeparateCorrections
or calculateType1METCorrection:
47 self.
vParL1 = ROOT.vector(ROOT.JetCorrectorParameters)()
50 if upToLevel >= 2
and calculateSeparateCorrections:
51 self.
vParL2 = ROOT.vector(ROOT.JetCorrectorParameters)()
54 if upToLevel >= 3
and calculateSeparateCorrections:
55 self.
vParL3 = ROOT.vector(ROOT.JetCorrectorParameters)()
58 if doResidualJECs
and calculateSeparateCorrections:
59 self.
vParL3Res = ROOT.vector(ROOT.JetCorrectorParameters)()
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())
70 corr = corrector.getCorrection()
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())
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
81 corr *=
max(0, 1+delta*jet.jetEnergyCorrUncertainty)
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()
92 for idau
in range(jet.numberOfDaughters()):
93 pfcand = jet.daughter(idau)
94 if pfcand.isGlobalMuon()
or pfcand.isStandAloneMuon():
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)
101 If addCorr, set jet.corr to the correction.
102 If addShifts, set also the +1 and -1 jet shifts
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).
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).
113 raw = jet.rawFactor()
117 for sepcorr
in self.separateJetCorrectors.keys():
120 for cdelta,shift
in [(1.0,
"JECUp"), (-1.0,
"JECDown")]:
122 setattr(jet,
"corr"+shift, cshift)
125 newpt = jet.pt()*raw*corr
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)
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))
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)"""
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")
149 ok = self.
correct(j,rho,delta,addCorr=addCorr,addShifts=addShifts,metShift=metShift,type1METCorr=type1METCorr)
150 if not ok: badJets.append(j)
152 print(
"Warning: %d out of %d jets flagged bad by JEC." % (len(badJets), len(jets)))
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."""
162 oldpx, oldpy = met.px(), met.py()
165 raw = met.shiftedP2_74x(12,0);
166 rawsumet = met.shiftedSumEt_74x(12,0);
169 rawsumet = met.uncorSumEt();
170 rawpx, rawpy = raw.px, raw.py
172 corrpx = rawpx + type1METCorrections[0]
173 corrpy = rawpy + type1METCorrections[1]
174 corrsumet = rawsumet + type1METCorrections[2]
176 met.setP4(ROOT.reco.Particle.LorentzVector(corrpx,corrpy,0,hypot(corrpx,corrpy)))
178 met._sumEt = corrsumet
179 met.sumEt = types.MethodType(
lambda myself : myself._sumEt, met, met.__class__)
181 met.setCorShift(rawpx, rawpy, rawsumet, met.Raw)
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__)
196 met.uncorP2 = types.MethodType(
lambda myself :
None, met, met.__class__)
197 met.uncorP3 = types.MethodType(
lambda myself :
None, met, met.__class__)
const uint16_t range(const Frame &aFrame)
void print(TMatrixD &m, const char *label=nullptr, bool mathematicaFormat=false)
def setFakeRawMETOnOldMiniAODs