CMS 3D CMS Logo

PhotonAnalyzer.py
Go to the documentation of this file.
1 from __future__ import print_function
2 import operator
3 import itertools
4 import copy
5 import types
6 import re
7 
8 from ROOT import TLorentzVector
9 from ROOT import heppy
10 
11 from PhysicsTools.Heppy.analyzers.core.Analyzer import Analyzer
12 from PhysicsTools.HeppyCore.framework.event import Event
13 from PhysicsTools.HeppyCore.statistics.counter import Counter, Counters
14 from PhysicsTools.Heppy.analyzers.core.AutoHandle import AutoHandle
15 from PhysicsTools.Heppy.physicsobjects.Photon import Photon
16 from PhysicsTools.Heppy.physicsutils.PhotonCalibrator import Run2PhotonCalibrator
17 
18 from PhysicsTools.HeppyCore.utils.deltar import deltaR, deltaPhi, bestMatch, matchObjectCollection3
19 
20 import PhysicsTools.HeppyCore.framework.config as cfg
21 
22 
24 
25 
26  def __init__(self, cfg_ana, cfg_comp, looperName ):
27  super(PhotonAnalyzer,self).__init__(cfg_ana,cfg_comp,looperName)
28  self.etaCentral = self.cfg_ana.etaCentral if hasattr(self.cfg_ana, 'etaCentral') else 9999
29  # footprint removed isolation
30  self.doFootprintRemovedIsolation = getattr(cfg_ana, 'doFootprintRemovedIsolation', False)
32  self.footprintRemovedIsolationPUCorr = self.cfg_ana.footprintRemovedIsolationPUCorr
34  #FIXME: only Embedded works
35  if self.cfg_ana.doPhotonScaleCorrections:
36  conf = cfg_ana.doPhotonScaleCorrections
37  self.photonEnergyCalibrator = Run2PhotonCalibrator(
38  conf['data'],
39  cfg_comp.isMC,
40  conf['isSync'] if 'isSync' in conf else False,
41  )
42 
43  def declareHandles(self):
44  super(PhotonAnalyzer, self).declareHandles()
45  self.handles['rhoPhoton'] = AutoHandle( self.cfg_ana.rhoPhoton, 'double')
46 
47  #----------------------------------------
48  # DECLARATION OF HANDLES OF PHOTONS STUFF
49  #----------------------------------------
50 
51  self.handles['photons'] = AutoHandle( self.cfg_ana.photons,'std::vector<pat::Photon>')
52  self.mchandles['packedGen'] = AutoHandle( 'packedGenParticles', 'std::vector<pat::PackedGenParticle>' )
53  self.mchandles['prunedGen'] = AutoHandle( 'prunedGenParticles', 'std::vector<reco::GenParticle>' )
54 
55  self.handles['packedCandidates'] = AutoHandle( 'packedPFCandidates', 'std::vector<pat::PackedCandidate>')
56  self.handles['jets'] = AutoHandle( "slimmedJets", 'std::vector<pat::Jet>' )
57 
58 
59  def beginLoop(self, setup):
60  super(PhotonAnalyzer,self).beginLoop(setup)
61  self.counters.addCounter('events')
62  count = self.counters.counter('events')
63  count.register('all events')
64  count.register('has >=1 gamma at preselection')
65  count.register('has >=1 selected gamma')
66 
67  def makePhotons(self, event):
68  event.allphotons = map( Photon, self.handles['photons'].product() )
69  event.allphotons.sort(key = lambda l : l.pt(), reverse = True)
70 
71  event.selectedPhotons = []
72  event.selectedPhotonsCentral = []
73 
75  # values are taken from EGamma implementation: https://github.com/cms-sw/cmssw/blob/CMSSW_7_6_X/RecoEgamma/PhotonIdentification/plugins/PhotonIDValueMapProducer.cc#L198-L199
76  self.IsolationComputer.setPackedCandidates(self.handles['packedCandidates'].product(), -1, 0.1, 0.2)
77 
78  # Photon scale calibrations
79  if self.cfg_ana.doPhotonScaleCorrections:
80  for gamma in event.allphotons:
81  self.photonEnergyCalibrator.correct(gamma, event.run)
82 
83  foundPhoton = False
84  for gamma in event.allphotons:
85  if gamma.pt() < self.cfg_ana.ptMin: continue
86  if abs(gamma.eta()) > self.cfg_ana.etaMax: continue
87  foundPhoton = True
88 
89  gamma.rho = float(self.handles['rhoPhoton'].product()[0])
90  # https://twiki.cern.ch/twiki/bin/view/CMS/CutBasedPhotonIdentificationRun2#Selection_implementation_details
91  if abs(gamma.eta()) < 1.0: gamma.EffectiveArea03 = [ 0.0234, 0.0053, 0.078 ]
92  elif abs(gamma.eta()) < 1.479: gamma.EffectiveArea03 = [ 0.0189, 0.0103, 0.0629 ]
93  elif abs(gamma.eta()) < 2.0: gamma.EffectiveArea03 = [ 0.0171, 0.0057, 0.0264 ]
94  elif abs(gamma.eta()) < 2.2: gamma.EffectiveArea03 = [ 0.0129, 0.0070, 0.0462 ]
95  elif abs(gamma.eta()) < 2.3: gamma.EffectiveArea03 = [ 0.0110, 0.0152, 0.0740 ]
96  elif abs(gamma.eta()) < 2.4: gamma.EffectiveArea03 = [ 0.0074, 0.0232, 0.0924 ]
97  else: gamma.EffectiveArea03 = [ 0.0035, 0.1709, 0.1484 ]
98 
101 
102  gamma.relIso = (max(gamma.chargedHadronIso()-gamma.rho*gamma.EffectiveArea03[0],0) + max(gamma.neutralHadronIso()-gamma.rho*gamma.EffectiveArea03[1],0) + max(gamma.photonIso() - gamma.rho*gamma.EffectiveArea03[2],0))/gamma.pt()
103 
104  def idWP(gamma,X):
105  """Create an integer equal to 1-2-3 for (loose,medium,tight)"""
106 
107  id=0
108  if gamma.photonID(X%"Loose"):
109  id=1
110  #if gamma.photonID(X%"Medium"):
111  # id=2
112  if gamma.photonID(X%"Tight"):
113  id=3
114  return id
115 
116  gamma.idCutBased = idWP(gamma, "PhotonCutBasedID%s")
117 
118 
119  keepThisPhoton = True
120 
121  if self.cfg_ana.gammaID=="PhotonCutBasedIDLoose_CSA14" or self.cfg_ana.gammaID=="PhotonCutBasedIDLoose_PHYS14" :
122  gamma.idCutBased = gamma.photonIDCSA14(self.cfg_ana.gammaID)
123  # we're keeing sigmaietaieta sidebands:
124  keepThisPhoton = gamma.photonIDCSA14(self.cfg_ana.gammaID, True)
125 
126  if gamma.hasPixelSeed():
127  keepThisPhoton = False
128  gamma.idCutBased = 0
129  elif "NoIso" in self.cfg_ana.gammaID:
130  idName = re.split('_NoIso',self.cfg_ana.gammaID)
131  keepThisPhoton = gamma.passPhotonID(idName[0],self.cfg_ana.conversionSafe_eleVeto)
132  basenameID = re.split('_looseSieie',idName[0])
133  gamma.idCutBased = gamma.passPhotonID(basenameID[0],self.cfg_ana.conversionSafe_eleVeto)
134  else:
135  # Reading from miniAOD directly
136  #keepThisPhoton = gamma.photonID(self.cfg_ana.gammaID)
137 
138  # implement cut based ID with CMGTools
139  keepThisPhoton = gamma.passPhotonID(self.cfg_ana.gammaID,self.cfg_ana.conversionSafe_eleVeto) and gamma.passPhotonIso(self.cfg_ana.gammaID,self.cfg_ana.gamma_isoCorr)
140 
141  if keepThisPhoton:
142  event.selectedPhotons.append(gamma)
143 
144  if keepThisPhoton and abs(gamma.eta()) < self.etaCentral:
145  event.selectedPhotonsCentral.append(gamma)
146 
147  event.selectedPhotons.sort(key = lambda l : l.pt(), reverse = True)
148  event.selectedPhotonsCentral.sort(key = lambda l : l.pt(), reverse = True)
149 
150  self.counters.counter('events').inc('all events')
151  if foundPhoton: self.counters.counter('events').inc('has >=1 gamma at preselection')
152  if len(event.selectedPhotons): self.counters.counter('events').inc('has >=1 selected gamma')
153 
154  def matchPhotons(self, event):
155  event.genPhotons = [ x for x in event.genParticles if x.status() == 1 and abs(x.pdgId()) == 22 ]
156  event.genPhotonsWithMom = [ x for x in event.genPhotons if x.numberOfMothers()>0 ]
157  event.genPhotonsWithoutMom = [ x for x in event.genPhotons if x.numberOfMothers()==0 ]
158  event.genPhotonsMatched = [ x for x in event.genPhotonsWithMom if abs(x.mother(0).pdgId())<23 or x.mother(0).pdgId()==2212 ]
159  match = matchObjectCollection3(event.allphotons, event.genPhotonsMatched, deltaRMax = 0.1)
160  matchNoMom = matchObjectCollection3(event.allphotons, event.genPhotonsWithoutMom, deltaRMax = 0.1)
161  packedGenParts = [ p for p in self.mchandles['packedGen'].product() if abs(p.eta()) < 3.1 ]
162  partons = [ p for p in self.mchandles['prunedGen'].product() if (p.status()==23 or p.status()==22) and abs(p.pdgId())<22 ]
163  for gamma in event.allphotons:
164  gen = match[gamma]
165  gamma.mcGamma = gen
166  if gen and gen.pt()>=0.5*gamma.pt() and gen.pt()<=2.*gamma.pt():
167  gamma.mcMatchId = 22
168  sumPt03 = 0.;
169  sumPt04 = 0.;
170  for part in packedGenParts:
171  if abs(part.pdgId())==12: continue # exclude neutrinos
172  if abs(part.pdgId())==14: continue
173  if abs(part.pdgId())==16: continue
174  if abs(part.pdgId())==18: continue
175  deltar = deltaR(gen.eta(), gen.phi(), part.eta(), part.phi())
176  if deltar <= 0.3:
177  sumPt03 += part.pt()
178  if deltar <= 0.4:
179  sumPt04 += part.pt()
180  sumPt03 -= gen.pt()
181  sumPt04 -= gen.pt()
182  if sumPt03<0. : sumPt03=0.
183  if sumPt04<0. : sumPt04=0.
184  gamma.genIso03 = sumPt03
185  gamma.genIso04 = sumPt04
186  # match to parton
187  deltaRmin = 999.
188  for p in partons:
189  deltar = deltaR(gen.eta(), gen.phi(), p.eta(), p.phi())
190  if deltar < deltaRmin:
191  deltaRmin = deltar
192  gamma.drMinParton = deltaRmin
193  else:
194  genNoMom = matchNoMom[gamma]
195  if genNoMom:
196  gamma.mcMatchId = 7
197  sumPt03 = 0.;
198  sumPt04 = 0.;
199  for part in packedGenParts:
200  if abs(part.pdgId())==12: continue # exclude neutrinos
201  if abs(part.pdgId())==14: continue
202  if abs(part.pdgId())==16: continue
203  if abs(part.pdgId())==18: continue
204  deltar = deltaR(genNoMom.eta(), genNoMom.phi(), part.eta(), part.phi())
205  if deltar <= 0.3:
206  sumPt03 += part.pt()
207  if deltar <= 0.4:
208  sumPt04 += part.pt()
209  sumPt03 -= genNoMom.pt()
210  sumPt04 -= genNoMom.pt()
211  if sumPt03<0. : sumPt03=0.
212  if sumPt04<0. : sumPt04=0.
213  gamma.genIso03 = sumPt03
214  gamma.genIso04 = sumPt04
215  # match to parton
216  deltaRmin = 999.
217  for p in partons:
218  deltar = deltaR(genNoMom.eta(), genNoMom.phi(), p.eta(), p.phi())
219  if deltar < deltaRmin:
220  deltaRmin = deltar
221  gamma.drMinParton = deltaRmin
222  else:
223  gamma.mcMatchId = 0
224  gamma.genIso03 = -1.
225  gamma.genIso04 = -1.
226  gamma.drMinParton = -1.
227 
228 
229 
230 
231 
232  def checkMatch( self, eta, phi, particles, deltar ):
233 
234  for part in particles:
235  if deltaR(eta, phi, part.eta(), part.phi()) < deltar:
236  return True
237 
238  return False
239 
240 
241 
242 
243 
244  def computeRandomCone( self, event, eta, phi, deltarmax, charged, jets, photons ):
245 
246  if self.checkMatch( eta, phi, jets, 2.*deltarmax ):
247  return -1.
248 
249  if self.checkMatch( eta, phi, photons, 2.*deltarmax ):
250  return -1.
251 
252  if self.checkMatch( eta, phi, event.selectedLeptons, deltarmax ):
253  return -1.
254 
255  iso = 0.
256 
257  for part in charged:
258  if deltaR(eta, phi, part.eta(), part.phi()) > deltarmax : continue
259  #if deltaR(eta, phi, part.eta(), part.phi()) < 0.02: continue
260  iso += part.pt()
261 
262  return iso
263 
264 
265 
266 
267 
268 
269  def randomCone( self, event ):
270 
271  patcands = self.handles['packedCandidates'].product()
272  jets = self.handles['jets'].product()
273 
274  charged = [ p for p in patcands if ( p.charge() != 0 and abs(p.pdgId())>20 and abs(p.dz())<=0.1 and p.fromPV()>1 and p.trackHighPurity() ) ]
275  photons10 = [ p for p in patcands if ( p.pdgId() == 22 and p.pt()>10. ) ]
276  jets20 = [ j for j in jets if j.pt() > 20 and abs(j.eta())<2.5 ]
277 
278  for gamma in event.allphotons:
279 
280  etaPhot = gamma.eta()
281  phiPhot = gamma.eta()
282  pi = 3.14159
283  phiRC = phiPhot + 0.5*pi
284  while phiRC>pi:
285  phiRC -= 2.*pi
286 
287 
288  gamma.chHadIsoRC03 = self.computeRandomCone( event, etaPhot, phiRC, 0.3, charged, jets20, photons10 )
289  gamma.chHadIsoRC04 = self.computeRandomCone( event, etaPhot, phiRC, 0.4, charged, jets20, photons10 )
290 
291 
292  #try other side
293  phiRC = phiPhot - 0.5*pi
294  while phiRC<-pi:
295  phiRC += 2.*pi
296 
297  if gamma.chHadIsoRC03<0. : gamma.chHadIsoRC03 = self.computeRandomCone( event, etaPhot, phiRC, 0.3, charged, jets20, photons10 )
298  if gamma.chHadIsoRC04<0. : gamma.chHadIsoRC04 = self.computeRandomCone( event, etaPhot, phiRC, 0.4, charged, jets20, photons10 )
299 
300 
302  # cone deltar=0.3
303  gamma.ftprAbsIsoCharged03 = self.IsolationComputer.chargedAbsIso(gamma.physObj, 0.3, 0, 0.0);
304  gamma.ftprAbsIsoPho03 = self.IsolationComputer.photonAbsIsoRaw( gamma.physObj, 0.3, 0, 0.0);
305  gamma.ftprAbsIsoNHad03 = self.IsolationComputer.neutralHadAbsIsoRaw( gamma.physObj, 0.3, 0, 0.0);
306  gamma.ftprAbsIso03 = gamma.ftprAbsIsoCharged03 + gamma.ftprAbsIsoPho03 + gamma.ftprAbsIsoNHad03
307  if self.cfg_ana.gamma_isoCorr == "rhoArea":
308  gamma.ftprAbsIso03 = (max(gamma.ftprAbsIsoCharged03-gamma.rho*gamma.EffectiveArea03[0],0) + max(gamma.ftprAbsIsoPho03-gamma.rho*gamma.EffectiveArea03[1],0) + max(gamma.ftprAbsIsoNHad03 - gamma.rho*gamma.EffectiveArea03[2],0))
309  elif self.cfg_ana.gamma_isoCorr != 'raw':
310  raise RuntimeError("Unsupported gamma_isoCorr name '" + str(self.cfg_ana.gamma_isoCorr) + "'! For now only 'rhoArea', 'raw' are supported.")
311  gamma.ftprRelIso03 = gamma.ftprAbsIso03/gamma.pt()
312 
313  def printInfo(self, event):
314  print('----------------')
315  if len(event.selectedPhotons)>0:
316  print('lenght: ',len(event.selectedPhotons))
317  print('gamma candidate pt: ',event.selectedPhotons[0].pt())
318  print('gamma candidate eta: ',event.selectedPhotons[0].eta())
319  print('gamma candidate phi: ',event.selectedPhotons[0].phi())
320  print('gamma candidate mass: ',event.selectedPhotons[0].mass())
321  print('gamma candidate HoE: ',event.selectedPhotons[0].hOVERe())
322  print('gamma candidate r9: ',event.selectedPhotons[0].full5x5_r9())
323  print('gamma candidate sigmaIetaIeta: ',event.selectedPhotons[0].full5x5_sigmaIetaIeta())
324  print('gamma candidate had iso: ',event.selectedPhotons[0].chargedHadronIso())
325  print('gamma candidate neu iso: ',event.selectedPhotons[0].neutralHadronIso())
326  print('gamma candidate gamma iso: ',event.selectedPhotons[0].photonIso())
327  print('gamma idCutBased',event.selectedPhotons[0].idCutBased)
328 
329 
330  def process(self, event):
331  self.readCollections( event.input )
332  self.makePhotons(event)
333 # self.printInfo(event)
334 
335  if self.cfg_ana.do_randomCone:
336  self.randomCone(event)
337 
338  if not self.cfg_comp.isMC:
339  return True
340 
341  if self.cfg_ana.do_mc_match and hasattr(event, 'genParticles'):
342  self.matchPhotons(event)
343 
344 
345  return True
346 
347 
348 setattr(PhotonAnalyzer,"defaultConfig",cfg.Analyzer(
349  class_object=PhotonAnalyzer,
350  photons='slimmedPhotons',
351  ptMin = 20,
352  etaMax = 2.5,
353  # energy scale corrections (off by default)
354  doPhotonScaleCorrections=False,
355  gammaID = "PhotonCutBasedIDLoose_CSA14",
356  rhoPhoton = 'fixedGridRhoFastjetAll',
357  gamma_isoCorr = 'rhoArea',
358  # Footprint-removed isolation, removing all the footprint of the photon
359  doFootprintRemovedIsolation = False, # off by default since it requires access to all PFCandidates
360  packedCandidates = 'packedPFCandidates',
361  footprintRemovedIsolationPUCorr = 'rhoArea', # Allowed options: 'rhoArea', 'raw' (uncorrected)
362  conversionSafe_eleVeto = False,
363  do_mc_match = True,
364  do_randomCone = False,
365  )
366 )
367 
def attachFootprintRemovedIsolation(self, gamma)
void print(TMatrixD &m, const char *label=nullptr, bool mathematicaFormat=false)
Definition: Utilities.cc:47
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
def __init__(self, cfg_ana, cfg_comp, looperName)
def computeRandomCone(self, event, eta, phi, deltarmax, charged, jets, photons)
def checkMatch(self, eta, phi, particles, deltar)
#define str(s)
def matchObjectCollection3
Definition: deltar.py:41