test
CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
JetAnalyzer.py
Go to the documentation of this file.
1 import math, os
2 from PhysicsTools.Heppy.analyzers.core.Analyzer import Analyzer
3 from PhysicsTools.Heppy.analyzers.core.AutoHandle import AutoHandle
4 from PhysicsTools.Heppy.physicsobjects.PhysicsObjects import Jet
5 from PhysicsTools.HeppyCore.utils.deltar import deltaR2, deltaPhi, matchObjectCollection, matchObjectCollection2, bestMatch,matchObjectCollection3
6 from PhysicsTools.Heppy.physicsutils.JetReCalibrator import JetReCalibrator
7 import PhysicsTools.HeppyCore.framework.config as cfg
8 
9 from PhysicsTools.Heppy.physicsutils.QGLikelihoodCalculator import QGLikelihoodCalculator
10 
11 import copy
12 def cleanNearestJetOnly(jets,leptons,deltaR):
13  dr2 = deltaR**2
14  good = [ True for j in jets ]
15  for l in leptons:
16  ibest, d2m = -1, dr2
17  for i,j in enumerate(jets):
18  d2i = deltaR2(l.eta(),l.phi(), j.eta(),j.phi())
19  if d2i < d2m:
20  ibest, d2m = i, d2i
21  if ibest != -1: good[ibest] = False
22  return [ j for (i,j) in enumerate(jets) if good[i] == True ]
23 
24 def cleanJetsAndLeptons(jets,leptons,deltaR,arbitration):
25  dr2 = deltaR**2
26  goodjet = [ True for j in jets ]
27  goodlep = [ True for l in leptons ]
28  for il, l in enumerate(leptons):
29  ibest, d2m = -1, dr2
30  for i,j in enumerate(jets):
31  d2i = deltaR2(l.eta(),l.phi(), j.eta(),j.phi())
32  if d2i < dr2:
33  choice = arbitration(j,l)
34  if choice == j:
35  # if the two match, and we prefer the jet, then drop the lepton and be done
36  goodlep[il] = False
37  break
38  elif choice == (j,l) or choice == (l,j):
39  # asked to keep both, so we don't consider this match
40  continue
41  if d2i < d2m:
42  ibest, d2m = i, d2i
43  # this lepton has been killed by a jet, then we clean the jet that best matches it
44  if not goodlep[il]: continue
45  if ibest != -1: goodjet[ibest] = False
46  return ( [ j for (i ,j) in enumerate(jets) if goodjet[i ] == True ],
47  [ l for (il,l) in enumerate(leptons) if goodlep[il] == True ] )
48 
49 
50 
51 def shiftJERfactor(JERShift, aeta):
52  factor = 1.079 + JERShift*0.026
53  if aeta > 3.2: factor = 1.056 + JERShift * 0.191
54  elif aeta > 2.8: factor = 1.395 + JERShift * 0.063
55  elif aeta > 2.3: factor = 1.254 + JERShift * 0.062
56  elif aeta > 1.7: factor = 1.208 + JERShift * 0.046
57  elif aeta > 1.1: factor = 1.121 + JERShift * 0.029
58  elif aeta > 0.5: factor = 1.099 + JERShift * 0.028
59  return factor
60 
61 
62 
63 
64 
65 class JetAnalyzer( Analyzer ):
66  """Taken from RootTools.JetAnalyzer, simplified, modified, added corrections """
67  def __init__(self, cfg_ana, cfg_comp, looperName):
68  super(JetAnalyzer,self).__init__(cfg_ana, cfg_comp, looperName)
69  mcGT = cfg_ana.mcGT if hasattr(cfg_ana,'mcGT') else "PHYS14_25_V2"
70  dataGT = cfg_ana.dataGT if hasattr(cfg_ana,'dataGT') else "GR_70_V2_AN1"
71  self.shiftJEC = self.cfg_ana.shiftJEC if hasattr(self.cfg_ana, 'shiftJEC') else 0
72  self.recalibrateJets = self.cfg_ana.recalibrateJets
73  self.addJECShifts = self.cfg_ana.addJECShifts if hasattr(self.cfg_ana, 'addJECShifts') else 0
74  if self.recalibrateJets == "MC" : self.recalibrateJets = self.cfg_comp.isMC
75  elif self.recalibrateJets == "Data": self.recalibrateJets = not self.cfg_comp.isMC
76  elif self.recalibrateJets not in [True,False]: raise RuntimeError, "recalibrateJets must be any of { True, False, 'MC', 'Data' }, while it is %r " % self.recalibrateJets
77 
78  calculateSeparateCorrections = getattr(cfg_ana,"calculateSeparateCorrections", False);
79  calculateType1METCorrection = getattr(cfg_ana,"calculateType1METCorrection", False);
80  self.doJEC = self.recalibrateJets or (self.shiftJEC != 0) or self.addJECShifts or calculateSeparateCorrections or calculateType1METCorrection
81  if self.doJEC:
82  doResidual = getattr(cfg_ana, 'applyL2L3Residual', 'Data')
83  if doResidual == "MC": doResidual = self.cfg_comp.isMC
84  elif doResidual == "Data": doResidual = not self.cfg_comp.isMC
85  elif doResidual not in [True,False]: raise RuntimeError, "If specified, applyL2L3Residual must be any of { True, False, 'MC', 'Data'(default)}"
86  GT = mcGT if self.cfg_comp.isMC else dataGT
87  # Now take care of the optional arguments
88  kwargs = { 'calculateSeparateCorrections':calculateSeparateCorrections,
89  'calculateType1METCorrection' :calculateType1METCorrection, }
90  if kwargs['calculateType1METCorrection']: kwargs['type1METParams'] = cfg_ana.type1METParams
91  # instantiate the jet re-calibrator
92  self.jetReCalibrator = JetReCalibrator(GT, cfg_ana.recalibrationType, doResidual, cfg_ana.jecPath, **kwargs)
93  self.doPuId = getattr(self.cfg_ana, 'doPuId', True)
94  self.jetLepDR = getattr(self.cfg_ana, 'jetLepDR', 0.4)
95  self.jetLepArbitration = getattr(self.cfg_ana, 'jetLepArbitration', lambda jet,lepton: lepton)
96  self.lepPtMin = getattr(self.cfg_ana, 'minLepPt', -1)
97  self.lepSelCut = getattr(self.cfg_ana, 'lepSelCut', lambda lep : True)
98  self.jetGammaDR = getattr(self.cfg_ana, 'jetGammaDR', 0.4)
99  if(self.cfg_ana.doQG):
100  qgdefname="{CMSSW_BASE}/src/PhysicsTools/Heppy/data/pdfQG_AK4chs_13TeV_v2b.root"
101  self.qglcalc = QGLikelihoodCalculator(getattr(self.cfg_ana,"QGpath",qgdefname).format(CMSSW_BASE= os.environ['CMSSW_BASE']))
102  if not hasattr(self.cfg_ana ,"collectionPostFix"):self.cfg_ana.collectionPostFix=""
103 
104  def declareHandles(self):
105  super(JetAnalyzer, self).declareHandles()
106  self.handles['jets'] = AutoHandle( self.cfg_ana.jetCol, 'std::vector<pat::Jet>' )
107  self.handles['genJet'] = AutoHandle( self.cfg_ana.genJetCol, 'vector<reco::GenJet>' )
108  self.shiftJER = self.cfg_ana.shiftJER if hasattr(self.cfg_ana, 'shiftJER') else 0
109  self.addJERShifts = self.cfg_ana.addJERShifts if hasattr(self.cfg_ana, 'addJERShifts') else 0
110  self.handles['rho'] = AutoHandle( self.cfg_ana.rho, 'double' )
111 
112  def beginLoop(self, setup):
113  super(JetAnalyzer,self).beginLoop(setup)
114 
115  def process(self, event):
116  self.readCollections( event.input )
117  rho = float(self.handles['rho'].product()[0])
118  self.rho = rho
119 
120  ## Read jets, if necessary recalibrate and shift MET
121  if self.cfg_ana.copyJetsByValue:
122  import ROOT
123  allJets = map(lambda j:Jet(ROOT.pat.Jet(ROOT.edm.Ptr(ROOT.pat.Jet)(ROOT.edm.ProductID(),j,0))), self.handles['jets'].product()) #copy-by-value is safe if JetAnalyzer is ran more than once
124  else:
125  allJets = map(Jet, self.handles['jets'].product())
126 
127  #set dummy MC flavour for all jets in case we want to ntuplize discarded jets later
128  for jet in allJets:
129  jet.mcFlavour = 0
130 
131  self.deltaMetFromJEC = [0.,0.]
132  self.type1METCorr = [0.,0.,0.]
133 # print "before. rho",self.rho,self.cfg_ana.collectionPostFix,'allJets len ',len(allJets),'pt', [j.pt() for j in allJets]
134  if self.doJEC:
135  if not self.recalibrateJets: # check point that things won't change
136  jetsBefore = [ (j.pt(),j.eta(),j.phi(),j.rawFactor()) for j in allJets ]
137  self.jetReCalibrator.correctAll(allJets, rho, delta=self.shiftJEC,
138  addCorr=True, addShifts=self.addJECShifts,
139  metShift=self.deltaMetFromJEC, type1METCorr=self.type1METCorr )
140  if not self.recalibrateJets:
141  jetsAfter = [ (j.pt(),j.eta(),j.phi(),j.rawFactor()) for j in allJets ]
142  if len(jetsBefore) != len(jetsAfter):
143  print "ERROR: I had to recompute jet corrections, and they rejected some of the jets:\nold = %s\n new = %s\n" % (jetsBefore,jetsAfter)
144  else:
145  for (told, tnew) in zip(jetsBefore,jetsAfter):
146  if (deltaR2(told[1],told[2],tnew[1],tnew[2])) > 0.0001:
147  print "ERROR: I had to recompute jet corrections, and one jet direction moved: old = %s, new = %s\n" % (told, tnew)
148  elif abs(told[0]-tnew[0])/(told[0]+tnew[0]) > 0.5e-3 or abs(told[3]-tnew[3]) > 0.5e-3:
149  print "ERROR: I had to recompute jet corrections, and one jet pt or corr changed: old = %s, new = %s\n" % (told, tnew)
150  self.allJetsUsedForMET = allJets
151 # print "after. rho",self.rho,self.cfg_ana.collectionPostFix,'allJets len ',len(allJets),'pt', [j.pt() for j in allJets]
152 
153  if self.cfg_comp.isMC:
154  self.genJets = [ x for x in self.handles['genJet'].product() ]
155  if self.cfg_ana.do_mc_match:
156  for igj, gj in enumerate(self.genJets):
157  gj.index = igj
158  self.matchJets(event, allJets)
159  if getattr(self.cfg_ana, 'smearJets', False):
160  self.smearJets(event, allJets)
161 
162 
163 
164 
165  ##Sort Jets by pT
166  allJets.sort(key = lambda j : j.pt(), reverse = True)
167 
168  leptons = []
169  if hasattr(event, 'selectedLeptons'):
170  leptons = [ l for l in event.selectedLeptons if l.pt() > self.lepPtMin and self.lepSelCut(l) ]
171  if self.cfg_ana.cleanJetsFromTaus and hasattr(event, 'selectedTaus'):
172  leptons = leptons[:] + event.selectedTaus
173  if self.cfg_ana.cleanJetsFromIsoTracks and hasattr(event, 'selectedIsoCleanTrack'):
174  leptons = leptons[:] + event.selectedIsoCleanTrack
175 
176  ## Apply jet selection
177  self.jets = []
178  self.jetsFailId = []
179  self.jetsAllNoID = []
180  self.jetsIdOnly = []
181  for jet in allJets:
182  #Check if lepton and jet have overlapping PF candidates
183  leps_with_overlaps = []
184  if getattr(self.cfg_ana, 'checkLeptonPFOverlap', True):
185  for i in range(jet.numberOfSourceCandidatePtrs()):
186  p1 = jet.sourceCandidatePtr(i) #Ptr<Candidate> p1
187  for lep in leptons:
188  for j in range(lep.numberOfSourceCandidatePtrs()):
189  p2 = lep.sourceCandidatePtr(j)
190  has_overlaps = p1.key() == p2.key() and p1.refCore().id().productIndex() == p2.refCore().id().productIndex() and p1.refCore().id().processIndex() == p2.refCore().id().processIndex()
191  if has_overlaps:
192  leps_with_overlaps += [lep]
193  if len(leps_with_overlaps)>0:
194  for lep in leps_with_overlaps:
195  lep.jetOverlap = jet
196  if self.testJetNoID( jet ):
197  self.jetsAllNoID.append(jet)
198  if(self.cfg_ana.doQG):
199  jet.qgl_calc = self.qglcalc.computeQGLikelihood
200  jet.qgl_rho = rho
201  if self.testJetID( jet ):
202  self.jets.append(jet)
203  self.jetsIdOnly.append(jet)
204  else:
205  self.jetsFailId.append(jet)
206  elif self.testJetID (jet ):
207  self.jetsIdOnly.append(jet)
208 
209  jetsEtaCut = [j for j in self.jets if abs(j.eta()) < self.cfg_ana.jetEta ]
210  self.cleanJetsAll, cleanLeptons = cleanJetsAndLeptons(jetsEtaCut, leptons, self.jetLepDR, self.jetLepArbitration)
211 
212  self.cleanJets = [j for j in self.cleanJetsAll if abs(j.eta()) < self.cfg_ana.jetEtaCentral ]
213  self.cleanJetsFwd = [j for j in self.cleanJetsAll if abs(j.eta()) >= self.cfg_ana.jetEtaCentral ]
214  self.discardedJets = [j for j in self.jets if j not in self.cleanJetsAll]
215  if hasattr(event, 'selectedLeptons') and self.cfg_ana.cleanSelectedLeptons:
216  event.discardedLeptons = [ l for l in leptons if l not in cleanLeptons ]
217  event.selectedLeptons = [ l for l in event.selectedLeptons if l not in event.discardedLeptons ]
218  for lep in leptons:
219  if hasattr(lep, "jetOverlap"):
220  if lep.jetOverlap in self.cleanJetsAll:
221  #print "overlap reco", lep.p4().pt(), lep.p4().eta(), lep.p4().phi(), lep.jetOverlap.p4().pt(), lep.jetOverlap.p4().eta(), lep.jetOverlap.p4().phi()
222  lep.jetOverlapIdx = self.cleanJetsAll.index(lep.jetOverlap)
223  elif lep.jetOverlap in self.discardedJets:
224  #print "overlap discarded", lep.p4().pt(), lep.p4().eta(), lep.p4().phi(), lep.jetOverlap.p4().pt(), lep.jetOverlap.p4().eta(), lep.jetOverlap.p4().phi()
225  lep.jetOverlapIdx = 1000 + self.discardedJets.index(lep.jetOverlap)
226 
227  ## First cleaning, then Jet Id
228  self.noIdCleanJetsAll, cleanLeptons = cleanJetsAndLeptons(self.jetsAllNoID, leptons, self.jetLepDR, self.jetLepArbitration)
229  self.noIdCleanJets = [j for j in self.noIdCleanJetsAll if abs(j.eta()) < self.cfg_ana.jetEtaCentral ]
230  self.noIdCleanJetsFwd = [j for j in self.noIdCleanJetsAll if abs(j.eta()) >= self.cfg_ana.jetEtaCentral ]
231  self.noIdDiscardedJets = [j for j in self.jetsAllNoID if j not in self.noIdCleanJetsAll]
232 
233  ## Clean Jets from photons (first cleaning, then Jet Id)
234  photons = []
235  if hasattr(event, 'selectedPhotons'):
236  if self.cfg_ana.cleanJetsFromFirstPhoton:
237  photons = event.selectedPhotons[:1]
238  else:
239  photons = [ g for g in event.selectedPhotons ]
240 
241  self.gamma_cleanJetsAll = cleanNearestJetOnly(self.cleanJetsAll, photons, self.jetGammaDR)
242  self.gamma_cleanJets = [j for j in self.gamma_cleanJetsAll if abs(j.eta()) < self.cfg_ana.jetEtaCentral ]
243  self.gamma_cleanJetsFwd = [j for j in self.gamma_cleanJetsAll if abs(j.eta()) >= self.cfg_ana.jetEtaCentral ]
244  self.gamma_noIdCleanJetsAll = cleanNearestJetOnly(self.noIdCleanJetsAll, photons, self.jetGammaDR)
245  self.gamma_noIdCleanJets = [j for j in self.gamma_noIdCleanJetsAll if abs(j.eta()) < self.cfg_ana.jetEtaCentral ]
246  self.gamma_noIdCleanJetsFwd = [j for j in self.gamma_noIdCleanJetsAll if abs(j.eta()) >= self.cfg_ana.jetEtaCentral ]
247  ###
248 
249  if self.cfg_ana.alwaysCleanPhotons:
250  self.cleanJets = self.gamma_cleanJets
251  self.cleanJetsAll = self.gamma_cleanJetsAll
252  self.cleanJetsFwd = self.gamma_cleanJetsFwd
253  #
254  self.noIdCleanJets = self.gamma_noIdCleanJets
255  self.noIdCleanJetsAll = self.gamma_noIdCleanJetsAll
256  self.noIdCleanJetsFwd = self.gamma_noIdCleanJetsFwd
257 
258  ## Jet Id, after jet/lepton cleaning
259  self.cleanJetsFailIdAll = []
260  for jet in self.noIdCleanJetsAll:
261  if not self.testJetID( jet ):
262  self.cleanJetsFailIdAll.append(jet)
263 
264  self.cleanJetsFailId = [j for j in self.cleanJetsFailIdAll if abs(j.eta()) < self.cfg_ana.jetEtaCentral ]
265 
266  ## Jet Id, after jet/photon cleaning
267  self.gamma_cleanJetsFailIdAll = []
268  for jet in self.gamma_noIdCleanJetsAll:
269  if not self.testJetID( jet ):
270  self.gamma_cleanJetsFailIdAll.append(jet)
271 
272  self.gamma_cleanJetsFailId = [j for j in self.gamma_cleanJetsFailIdAll if abs(j.eta()) < self.cfg_ana.jetEtaCentral ]
273 
274  ## Associate jets to leptons
275  incleptons = event.inclusiveLeptons if hasattr(event, 'inclusiveLeptons') else event.selectedLeptons
276  jlpairs = matchObjectCollection(incleptons, allJets, self.jetLepDR**2)
277 
278  for jet in allJets:
279  jet.leptons = [l for l in jlpairs if jlpairs[l] == jet ]
280  for lep in incleptons:
281  jet = jlpairs[lep]
282  if jet is None:
283  setattr(lep,"jet"+self.cfg_ana.collectionPostFix,lep)
284  else:
285  setattr(lep,"jet"+self.cfg_ana.collectionPostFix,jet)
286  ## Associate jets to taus
287  taus = getattr(event,'selectedTaus',[])
288  jtaupairs = matchObjectCollection( taus, allJets, self.jetLepDR**2)
289 
290  for jet in allJets:
291  jet.taus = [l for l in jtaupairs if jtaupairs[l] == jet ]
292  for tau in taus:
293  setattr(tau,"jet"+self.cfg_ana.collectionPostFix,jtaupairs[tau])
294 
295  #MC stuff
296  if self.cfg_comp.isMC:
297  self.deltaMetFromJetSmearing = [0, 0]
298  for j in self.cleanJetsAll:
299  if hasattr(j, 'deltaMetFromJetSmearing'):
300  self.deltaMetFromJetSmearing[0] += j.deltaMetFromJetSmearing[0]
301  self.deltaMetFromJetSmearing[1] += j.deltaMetFromJetSmearing[1]
302 
303  self.cleanGenJets = cleanNearestJetOnly(self.genJets, leptons, self.jetLepDR)
304 
305  if self.cfg_ana.cleanGenJetsFromPhoton:
306  self.cleanGenJets = cleanNearestJetOnly(self.cleanGenJets, photons, self.jetLepDR)
307 
308  if getattr(self.cfg_ana, 'attachNeutrinos', True) and hasattr(self.cfg_ana,"genNuSelection") :
309  jetNus=[x for x in event.genParticles if abs(x.pdgId()) in [12,14,16] and self.cfg_ana.genNuSelection(x) ]
310  pairs= matchObjectCollection (jetNus, self.genJets, 0.4**2)
311 
312  for (nu,genJet) in pairs.iteritems() :
313  if genJet is not None :
314  if not hasattr(genJet,"nu") :
315  genJet.nu=nu.p4()
316  else :
317  genJet.nu+=nu.p4()
318 
319 
320  if self.cfg_ana.do_mc_match:
321  self.jetFlavour(event)
322 
323  if hasattr(event,"jets"+self.cfg_ana.collectionPostFix): raise RuntimeError, "Event already contains a jet collection with the following postfix: "+self.cfg_ana.collectionPostFix
324  setattr(event,"rho" +self.cfg_ana.collectionPostFix, self.rho )
325  setattr(event,"deltaMetFromJEC" +self.cfg_ana.collectionPostFix, self.deltaMetFromJEC )
326  setattr(event,"type1METCorr" +self.cfg_ana.collectionPostFix, self.type1METCorr )
327  setattr(event,"allJetsUsedForMET" +self.cfg_ana.collectionPostFix, self.allJetsUsedForMET )
328  setattr(event,"jets" +self.cfg_ana.collectionPostFix, self.jets )
329  setattr(event,"jetsFailId" +self.cfg_ana.collectionPostFix, self.jetsFailId )
330  setattr(event,"jetsAllNoID" +self.cfg_ana.collectionPostFix, self.jetsAllNoID )
331  setattr(event,"jetsIdOnly" +self.cfg_ana.collectionPostFix, self.jetsIdOnly )
332  setattr(event,"cleanJetsAll" +self.cfg_ana.collectionPostFix, self.cleanJetsAll )
333  setattr(event,"cleanJets" +self.cfg_ana.collectionPostFix, self.cleanJets )
334  setattr(event,"cleanJetsFwd" +self.cfg_ana.collectionPostFix, self.cleanJetsFwd )
335  setattr(event,"cleanJetsFailIdAll" +self.cfg_ana.collectionPostFix, self.cleanJetsFailIdAll )
336  setattr(event,"cleanJetsFailId" +self.cfg_ana.collectionPostFix, self.cleanJetsFailId )
337  setattr(event,"discardedJets" +self.cfg_ana.collectionPostFix, self.discardedJets )
338  setattr(event,"gamma_cleanJetsAll" +self.cfg_ana.collectionPostFix, self.gamma_cleanJetsAll )
339  setattr(event,"gamma_cleanJets" +self.cfg_ana.collectionPostFix, self.gamma_cleanJets )
340  setattr(event,"gamma_cleanJetsFwd" +self.cfg_ana.collectionPostFix, self.gamma_cleanJetsFwd )
341  setattr(event,"gamma_cleanJetsFailIdAll" +self.cfg_ana.collectionPostFix, self.gamma_cleanJetsFailIdAll )
342  setattr(event,"gamma_cleanJetsFailId" +self.cfg_ana.collectionPostFix, self.gamma_cleanJetsFailId )
343 
344 
345  if self.cfg_comp.isMC:
346  setattr(event,"deltaMetFromJetSmearing"+self.cfg_ana.collectionPostFix, self.deltaMetFromJetSmearing)
347  setattr(event,"cleanGenJets" +self.cfg_ana.collectionPostFix, self.cleanGenJets )
348  setattr(event,"genJets" +self.cfg_ana.collectionPostFix, self.genJets )
349  if self.cfg_ana.do_mc_match:
350  setattr(event,"bqObjects" +self.cfg_ana.collectionPostFix, self.bqObjects )
351  setattr(event,"cqObjects" +self.cfg_ana.collectionPostFix, self.cqObjects )
352  setattr(event,"partons" +self.cfg_ana.collectionPostFix, self.partons )
353  setattr(event,"heaviestQCDFlavour" +self.cfg_ana.collectionPostFix, self.heaviestQCDFlavour )
354 
355 
356  return True
357 
358 
359 
360  def testJetID(self, jet):
361  jet.puJetIdPassed = jet.puJetId()
362  jet.pfJetIdPassed = jet.jetID('POG_PFID_Loose')
363  if self.cfg_ana.relaxJetId:
364  return True
365  else:
366  return jet.pfJetIdPassed and (jet.puJetIdPassed or not(self.doPuId))
367 
368  def testJetNoID( self, jet ):
369  # 2 is loose pile-up jet id
370  return jet.pt() > self.cfg_ana.jetPt and \
371  abs( jet.eta() ) < self.cfg_ana.jetEta;
372 
373  def jetFlavour(self,event):
374  def isFlavour(x,f):
375  id = abs(x.pdgId())
376  if id > 999: return (id/1000)%10 == f
377  if id > 99: return (id/100)%10 == f
378  return id % 100 == f
379 
380 
381 
382  self.bqObjects = [ p for p in event.genParticles if (p.status() == 2 and isFlavour(p,5)) ]
383  self.cqObjects = [ p for p in event.genParticles if (p.status() == 2 and isFlavour(p,4)) ]
384 
385  self.partons = [ p for p in event.genParticles if ((p.status() == 23 or p.status() == 3) and abs(p.pdgId())>0 and (abs(p.pdgId()) in [1,2,3,4,5,21]) ) ]
386  match = matchObjectCollection2(self.cleanJetsAll,
387  self.partons,
388  deltaRMax = 0.3)
389 
390  for jet in self.cleanJetsAll:
391  parton = match[jet]
392  jet.partonId = (parton.pdgId() if parton != None else 0)
393  jet.partonMotherId = (parton.mother(0).pdgId() if parton != None and parton.numberOfMothers()>0 else 0)
394 
395  for jet in self.jets:
396  (bmatch, dr) = bestMatch(jet, self.bqObjects)
397  if dr < 0.4:
398  jet.mcFlavour = 5
399  else:
400  (cmatch, dr) = bestMatch(jet, self.cqObjects)
401  if dr < 0.4:
402  jet.mcFlavour = 4
403  else:
404  jet.mcFlavour = 0
405 
406  self.heaviestQCDFlavour = 5 if len(self.bqObjects) else (4 if len(self.cqObjects) else 1);
407 
408  def matchJets(self, event, jets):
409  match = matchObjectCollection2(jets,
410  event.genbquarks + event.genwzquarks,
411  deltaRMax = 0.3)
412  for jet in jets:
413  gen = match[jet]
414  jet.mcParton = gen
415  jet.mcMatchId = (gen.sourceId if gen != None else 0)
416  jet.mcMatchFlav = (abs(gen.pdgId()) if gen != None else 0)
417 
418  match = matchObjectCollection2(jets,
419  self.genJets,
420  deltaRMax = 0.3)
421  for jet in jets:
422  jet.mcJet = match[jet]
423 
424 
425 
426  def smearJets(self, event, jets):
427  # https://twiki.cern.ch/twiki/bin/viewauth/CMS/TWikiTopRefSyst#Jet_energy_resolution
428  for jet in jets:
429  gen = jet.mcJet
430  if gen != None:
431  genpt, jetpt, aeta = gen.pt(), jet.pt(), abs(jet.eta())
432  # from https://twiki.cern.ch/twiki/bin/view/CMS/JetResolution
433  #8 TeV tables
434  factor = shiftJERfactor(self.shiftJER, aeta)
435  ptscale = max(0.0, (jetpt + (factor-1)*(jetpt-genpt))/jetpt)
436  #print "get with pt %.1f (gen pt %.1f, ptscale = %.3f)" % (jetpt,genpt,ptscale)
437  jet.deltaMetFromJetSmearing = [ -(ptscale-1)*jet.rawFactor()*jet.px(), -(ptscale-1)*jet.rawFactor()*jet.py() ]
438  if ptscale != 0:
439  jet.setP4(jet.p4()*ptscale)
440  # leave the uncorrected unchanged for sync
441  jet.setRawFactor(jet.rawFactor()/ptscale)
442  #else: print "jet with pt %.1d, eta %.2f is unmatched" % (jet.pt(), jet.eta())
443  if (self.shiftJER==0) and (self.addJERShifts):
444  setattr(jet, "corrJER", ptscale )
445  factorJERUp= shiftJERfactor(1, aeta)
446  ptscaleJERUp = max(0.0, (jetpt + (factorJERUp-1)*(jetpt-genpt))/jetpt)
447  setattr(jet, "corrJERUp", ptscaleJERUp)
448  factorJERDown= shiftJERfactor(-1, aeta)
449  ptscaleJERDown = max(0.0, (jetpt + (factorJERDown-1)*(jetpt-genpt))/jetpt)
450  setattr(jet, "corrJERDown", ptscaleJERDown)
451 
452 
453 
454 
455 
456 setattr(JetAnalyzer,"defaultConfig", cfg.Analyzer(
457  class_object = JetAnalyzer,
458  jetCol = 'slimmedJets',
459  copyJetsByValue = False, #Whether or not to copy the input jets or to work with references (should be 'True' if JetAnalyzer is run more than once)
460  genJetCol = 'slimmedGenJets',
461  rho = ('fixedGridRhoFastjetAll','',''),
462  jetPt = 25.,
463  jetEta = 4.7,
464  jetEtaCentral = 2.4,
465  jetLepDR = 0.4,
466  jetLepArbitration = (lambda jet,lepton : lepton), # you can decide which to keep in case of overlaps; e.g. if the jet is b-tagged you might want to keep the jet
467  cleanSelectedLeptons = True, #Whether to clean 'selectedLeptons' after disambiguation. Treat with care (= 'False') if running Jetanalyzer more than once
468  minLepPt = 10,
469  lepSelCut = lambda lep : True,
470  relaxJetId = False,
471  doPuId = False, # Not commissioned in 7.0.X
472  doQG = False,
473  checkLeptonPFOverlap = True,
474  recalibrateJets = False,
475  applyL2L3Residual = 'Data', # if recalibrateJets, apply L2L3Residual to Data only
476  recalibrationType = "AK4PFchs",
477  shiftJEC = 0, # set to +1 or -1 to apply +/-1 sigma shift to the nominal jet energies
478  addJECShifts = False, # if true, add "corr", "corrJECUp", and "corrJECDown" for each jet (requires uncertainties to be available!)
479  smearJets = True,
480  shiftJER = 0, # set to +1 or -1 to get +/-1 sigma shifts
481  jecPath = "",
482  calculateSeparateCorrections = False,
483  calculateType1METCorrection = False,
484  type1METParams = { 'jetPtThreshold':15., 'skipEMfractionThreshold':0.9, 'skipMuons':True },
485  addJERShifts = 0, # add +/-1 sigma shifts to jets, intended to be used with shiftJER=0
486  cleanJetsFromFirstPhoton = False,
487  cleanJetsFromTaus = False,
488  cleanJetsFromIsoTracks = False,
489  alwaysCleanPhotons = False,
490  do_mc_match=True,
491  cleanGenJetsFromPhoton = False,
492  attachNeutrinos = True,
493  genNuSelection = lambda nu : True, #FIXME: add here check for ispromptfinalstate
494  collectionPostFix = ""
495  )
496 )
def bestMatch
Definition: deltar.py:139
deltaMetFromJEC
Read jets, if necessary recalibrate and shift MET.
Definition: JetAnalyzer.py:131
def matchObjectCollection
Definition: deltar.py:152
Definition: Jet.py:1
def matchObjectCollection2
Definition: deltar.py:167
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
if(c.getParameter< edm::InputTag >("puppiValueMap").label().size()!=0)
double deltaR2(const T1 &t1, const T2 &t2)
Definition: deltaR.h:36