-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathjetmet_analyzer.py
431 lines (347 loc) · 14.3 KB
/
jetmet_analyzer.py
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
# import ROOT in batch mode
import sys
import time
oldargv = sys.argv[:]
sys.argv = [ '-b-' ]
import ROOT
ROOT.gROOT.SetBatch(True)
sys.argv = oldargv
PLOTRAW=True #Set True if you want only raw pt from met objects
from FWCore.ParameterSet.VarParsing import VarParsing
options = VarParsing('python')
#default options
options.inputFiles="../step3_inMINIAODSIM.root"
options.outputFile="jetmetNtuples.root"
options.maxEvents=-1
#overwrite if given any command line arguments
options.parseArguments()
#in case of txt input file, read the information from txt
li_f=[]
for iF,F in enumerate(options.inputFiles):
print F
if F.split('.')[-1] == "txt":
options.inputFiles.pop(iF)
with open(F) as f:
li_f+=f.read().splitlines()
options.inputFiles+=li_f
print "analyzing files:"
for F in options.inputFiles: print F
# define deltaR
from math import hypot, pi, sqrt, fabs, cos, sin
import numpy as n
from jetmet_tree import *
from functions import *
# create an oput tree.
#fout = ROOT.TFile ("jetmet.root","recreate")
fout = ROOT.TFile (options.outputFile,"recreate")
t = ROOT.TTree ("events","events")
declare_branches(t)
# load FWLite C++ libraries
ROOT.gSystem.Load("libFWCoreFWLite.so");
ROOT.gSystem.Load("libDataFormatsFWLite.so");
ROOT.FWLiteEnabler.enable()
# load FWlite python libraries
from DataFormats.FWLite import Handle, Events
pfs, pfLabel = Handle("std::vector<pat::PackedCandidate>"), "packedPFCandidates"
muons, muonLabel = Handle("std::vector<pat::Muon>"), "slimmedMuons"
jets, jetLabel = Handle("std::vector<pat::Jet>"), "slimmedJets"
pjets, pjetLabel = Handle("std::vector<pat::Jet>"), "slimmedJetsPuppi"
mets, metLabel = Handle("std::vector<pat::MET>"), "slimmedMETs"
pmets, pmetLabel = Handle("std::vector<pat::MET>"), "slimmedMETsPuppi"
vertex, vertexLabel = Handle("std::vector<reco::Vertex>"),"offlineSlimmedPrimaryVertices"
rhoall_, rhoallLabel = Handle("double"), "fixedGridRhoFastjetAll"
rhocentral_, rhocentralLabel = Handle("double"), "fixedGridRhoFastjetCentral"
rhoneutral_, rhoneutralLabel = Handle("double"), "fixedGridRhoFastjetCentralNeutral"
rhochargedpileup_, rhochargedpileupLabel = Handle("double"), "fixedGridRhoFastjetCentralChargedPileUp"
#in order to print out the progress
def print_same_line(s):
sys.stdout.write(s) # just print
sys.stdout.flush() # needed for flush when using \x08
sys.stdout.write((b'\x08' * len(s)).decode())# back n chars
#time.sleep(0.2)
#prepare some variable distribution plots:
h_qt = ROOT.TH1F("qt","qt",40,0,200)
h_Zmass = ROOT.TH1F("Zmass","Zmass",25,50,150)
h_qeta = ROOT.TH1F("Zeta","Zeta",40,-5,5)
h_qrapidity = ROOT.TH1F("Zrapidity","Zrapidity",60,-3,3)
h_upll = ROOT.TH1F("upll","upll",40,-200,200)
h_upllqt = ROOT.TH1F("upllqt","upll+qt",40,-200,200)
h_uprp = ROOT.TH1F("uprp","uprp",40,-200,200)
# open file (you can use 'edmFileUtil -d /store/whatever.root' to get the physical file name)
#events = Events("file:/eos/cms/store/relval/CMSSW_9_4_0_pre3/RelValTTbar_13/MINIAODSIM/PU25ns_94X_mc2017_realistic_PixFailScenario_Run305081_FIXED_HS_AVE50-v1/10000/02B605A1-86C2-E711-A445-4C79BA28012B.root")
events = Events(options)
nevents = int(events.size())
print "total events: ", events.size()
for ievent,event in enumerate(events):
if options.maxEvents is not -1:
if ievent > options.maxEvents: continue
event.getByLabel(pfLabel, pfs)
event.getByLabel(muonLabel, muons)
event.getByLabel(jetLabel, jets)
event.getByLabel(pjetLabel, pjets)
event.getByLabel(metLabel, mets)
event.getByLabel(pmetLabel, pmets)
event.getByLabel(vertexLabel, vertex)
event.getByLabel(rhoallLabel,rhoall_)
event.getByLabel(rhocentralLabel,rhocentral_)
event.getByLabel(rhoneutralLabel,rhoneutral_)
event.getByLabel(rhochargedpileupLabel,rhochargedpileup_)
rhoall[0] = rhoall_.product()[0]
rhocentral[0] = rhocentral_.product()[0]
rhoneutral[0] = rhoneutral_.product()[0]
rhochargedpileup[0] = rhochargedpileup_.product()[0]
#print "\nEvent %d: run %6d, lumi %4d, event %12d" % (ievent,event.eventAuxiliary().run(), event.eventAuxiliary().luminosityBlock(),event.eventAuxiliary().event())
print_same_line(str(round(100.*ievent/nevents,2))+'%')
nrun[0] = event.eventAuxiliary().run()
nlumi[0] = event.eventAuxiliary().luminosityBlock()
nevent[0] = event.eventAuxiliary().event()
npv[0] = vertex.product().size()
njet[0] = jets.product().size()
npjet[0] = pjets.product().size()
###CHS JETS
for i,j in enumerate(jets.product()):
if i>=maxjet: break
jet_pt[i] = j.pt()
jet_eta[i] = j.eta()
jet_phi[i] = j.phi()
jet_energy[i] = j.energy()
rawjet_pt[i] = j.pt()*j.jecFactor("Uncorrected")
rawjet_eta[i] = j.eta()
rawjet_phi[i] = j.phi()
rawjet_energy[i] = j.energy()*j.jecFactor("Uncorrected")
NHF[i] = j.neutralHadronEnergyFraction()
NEMF[i] = j.neutralEmEnergyFraction();
CHF[i] = j.chargedHadronEnergyFraction();
MUF[i] = j.muonEnergyFraction();
CEMF[i] = j.chargedEmEnergyFraction();
NumConst[i] = j.chargedMultiplicity()+j.neutralMultiplicity();
NumNeutralParticle[i] = j.neutralMultiplicity();
CHM[i] = j.chargedMultiplicity();
eta = j.eta();
jet_loose[i] = (NHF[i]<0.99 and NEMF[i]<0.99 and NumConst[i]>1) and ((abs(eta)<=2.4 and CHF[i]>0 and CHM[i]>0 and CEMF[i]<0.99) or abs(eta)>2.4) and abs(eta)<=2.7
jet_loose[i] = jet_loose[i] or (NHF[i]<0.98 and NEMF[i]>0.01 and NumNeutralParticle[i]>2 and abs(eta)>2.7 and abs(eta)<=3.0 )
jet_loose[i] = jet_loose[i] or (NEMF[i]<0.90 and NumNeutralParticle[i]>10 and abs(eta)>3.0)
if not (j.genJet() == None):
genjet_pt[i] = j.genJet().pt()
genjet_eta[i] = j.genJet().eta()
genjet_phi[i] = j.genJet().phi()
genjet_energy[i] = j.genJet().energy()
else:
genjet_pt[i] = -999.
genjet_eta[i] = -999.
genjet_phi[i] = -999.
genjet_energy[i] = -999.
sourceCandidate = set()
# now get a list of the PF candidates used to build this jet
for isource in xrange(j.numberOfSourceCandidatePtrs()):
sourceCandidate.add(j.sourceCandidatePtr(isource).key()) # the key is the index in the pf collection
neutral[i], charged[i], photon[i], muon[i], electron[i], hhf[i], ehf[i], other[i] = 0, 0, 0, 0, 0, 0, 0, 0
neutral_e[i], charged_e[i], photon_e[i],muon_e[i],electron_e[i],hhf_e[i], ehf_e[i], other_e[i] = 0, 0, 0, 0, 0, 0, 0, 0
neutral_n[i], charged_n[i], photon_n[i],muon_n[i],electron_n[i],hhf_n[i], ehf_n[i], other_n[i] = 0, 0, 0, 0, 0, 0, 0, 0
for ipf,pf in enumerate(pfs.product()):
# if pf candidate is part of the jet, lets categorize
if ipf in sourceCandidate:
if ( abs(pf.pdgId()) == 211 ):
charged[i] += pf.pt();
charged_e[i] += pf.energy();
charged_n[i] += 1;
elif abs(pf.pdgId()) == 130:
neutral[i] += pf.pt();
neutral_e[i] += pf.energy();
neutral_n[i] += 1;
elif abs(pf.pdgId()) == 22:
photon[i] += pf.pt();
photon_e[i] += pf.energy();
photon_n[i] += 1;
elif abs(pf.pdgId()) == 13:
muon[i] += pf.pt();
muon_e[i] += pf.energy();
muon_n[i] += 1;
elif abs(pf.pdgId()) == 11:
electron[i] += pf.pt();
electron_e[i] += pf.energy();
electron_n[i] += 1;
elif abs(pf.pdgId()) == 1:
hhf[i] += pf.pt();
hhf_e[i] += pf.energy();
hhf_n[i] += 1;
elif abs(pf.pdgId()) == 2:
ehf[i] += pf.pt();
ehf_e[i] += pf.energy();
ehf_n[i] += 1;
else:
other[i] += pf.pt();
other_e[i] += pf.energy();
other_n[i] += 1;
# store the hcal depth information
for depth in range(1,8):
if abs(pf.pdgId())==130 and abs(pf.eta())>1.3 and abs(pf.eta())<5.0:
jet_depth[i][depth-1] += pf.hcalDepthEnergyFraction(depth)*pf.energy()*pf.hcalFraction()
jet_depth_uncorrected[i][depth-1] += pf.hcalDepthEnergyFraction(depth)*pf.energy()
###Puppi JETS
for i,j in enumerate(pjets.product()):
if i>=maxjet: break
pjet_pt[i] = j.pt()
pjet_eta[i] = j.eta()
pjet_phi[i] = j.phi()
pjet_energy[i] = j.energy()
rawpjet_pt[i] = j.pt()*j.jecFactor("Uncorrected")
rawpjet_eta[i] = j.eta()
rawpjet_phi[i] = j.phi()
rawpjet_energy[i] = j.energy()*j.jecFactor("Uncorrected")
NHF[i] = j.neutralHadronEnergyFraction()
NEMF[i] = j.neutralEmEnergyFraction();
CHF[i] = j.chargedHadronEnergyFraction();
MUF[i] = j.muonEnergyFraction();
CEMF[i] = j.chargedEmEnergyFraction();
NumConst[i] = j.chargedMultiplicity()+j.neutralMultiplicity();
NumNeutralParticle[i] = j.neutralMultiplicity();
CHM[i] = j.chargedMultiplicity();
eta = j.eta();
pjet_loose[i] = (NHF[i]<0.99 and NEMF[i]<0.99 and NumConst[i]>1) and ((abs(eta)<=2.4 and CHF[i]>0 and CHM[i]>0 and CEMF[i]<0.99) or abs(eta)>2.4) and abs(eta)<=2.7
pjet_loose[i] = jet_loose[i] or (NHF[i]<0.98 and NEMF[i]>0.01 and NumNeutralParticle[i]>2 and abs(eta)>2.7 and abs(eta)<=3.0 )
pjet_loose[i] = jet_loose[i] or (NEMF[i]<0.90 and NumNeutralParticle[i]>10 and abs(eta)>3.0)
if not (j.genJet() == None):
genpjet_pt[i] = j.genJet().pt()
genpjet_eta[i] = j.genJet().eta()
genpjet_phi[i] = j.genJet().phi()
genpjet_energy[i] = j.genJet().energy()
else:
genpjet_pt[i] = -999.
genpjet_eta[i] = -999.
genpjet_phi[i] = -999.
genpjet_energy[i] = -999.
sourceCandidate = set()
# now get a list of the PF candidates used to build this jet
for isource in xrange(j.numberOfSourceCandidatePtrs()):
sourceCandidate.add(j.sourceCandidatePtr(isource).key()) # the key is the index in the pf collection
for ipf,pf in enumerate(pfs.product()):
if ipf in sourceCandidate:
if ( abs(pf.pdgId()) == 211 ):
charged_pjet[i] += pf.pt()*pf.puppiWeight();
charged_e_pjet[i] += pf.energy()*pf.puppiWeight();
charged_n_pjet[i] += 1*pf.puppiWeight();
elif abs(pf.pdgId()) == 130:
neutral_pjet[i] += pf.pt()*pf.puppiWeight();
neutral_e_pjet[i] += pf.energy()*pf.puppiWeight();
neutral_n_pjet[i] += 1*pf.puppiWeight();
elif abs(pf.pdgId()) == 22:
photon_pjet[i] += pf.pt()*pf.puppiWeight();
photon_e_pjet[i] += pf.energy()*pf.puppiWeight();
photon_n_pjet[i] += 1*pf.puppiWeight();
elif abs(pf.pdgId()) == 13:
muon_pjet[i] += pf.pt()*pf.puppiWeight();
muon_e_pjet[i] += pf.energy()*pf.puppiWeight();
muon_n_pjet[i] += 1*pf.puppiWeight();
elif abs(pf.pdgId()) == 11:
electron_pjet[i] += pf.pt()*pf.puppiWeight();
electron_e_pjet[i] += pf.energy()*pf.puppiWeight();
electron_n_pjet[i] += 1*pf.puppiWeight();
elif abs(pf.pdgId()) == 1:
hhf_pjet[i] += pf.pt()*pf.puppiWeight();
hhf_e_pjet[i] += pf.energy()*pf.puppiWeight();
hhf_n_pjet[i] += 1*pf.puppiWeight();
elif abs(pf.pdgId()) == 2:
ehf_pjet[i] += pf.pt()*pf.puppiWeight();
ehf_e_pjet[i] += pf.energy()*pf.puppiWeight();
ehf_n_pjet[i] += 1*pf.puppiWeight();
else:
other_pjet[i] += pf.pt()*pf.puppiWeight();
other_e_pjet[i] += pf.energy()*pf.puppiWeight();
other_n_pjet[i] += 1*pf.puppiWeight();
# store the hcal depth information
for depth in range(1,8):
if abs(pf.pdgId())==130 and abs(pf.eta())>1.3 and abs(pf.eta())<5.0:
pjet_depth[i][depth-1] += pf.hcalDepthEnergyFraction(depth)*pf.energy()*pf.hcalFraction()*pf.puppiWeight()
pjet_depth_uncorrected[i][depth-1] += pf.hcalDepthEnergyFraction(depth)*pf.energy()*pf.puppiWeight()
#met information
metprod = mets.product().front()
met[0] = metprod.pt()
genmet[0] = metprod.genMET().pt()
rawmet[0] = metprod.uncorPt()
for ipf,pf in enumerate(pfs.product()):
if ( abs(pf.pdgId()) == 211 ):
charged_met[0] += pf.pt();
#charged_n_met[0] += 1;
elif abs(pf.pdgId()) == 130:
neutral_met[0] += pf.pt();
#neutral_n_met[0] += 1;
elif abs(pf.pdgId()) == 22:
photon_met[0] += pf.pt();
#photon_n_met[0] += 1;
elif abs(pf.pdgId()) == 13:
muon_met[0] += pf.pt();
#muon_n_met[0] += 1;
elif abs(pf.pdgId()) == 11:
electron_met[0] += pf.pt();
#electron_n_met[0] += 1;
elif abs(pf.pdgId()) == 1:
hhf_met[0] += pf.pt();
#hhf_n_met[0] += 1;
elif abs(pf.pdgId()) == 2:
ehf_met[0] += pf.pt();
#ehf_n_met[0] += 1;
else:
other_met[0] += pf.pt();
#other_n_met[0] += 1;
#puppi met information
pmetprod = pmets.product().front()
pmet[0] = pmetprod.pt()
genpmet[0] = pmetprod.genMET().pt()
rawpmet[0] = pmetprod.uncorPt()
###Reconstruct Z boson:
nmuons = len(muons.product())
Zvector = None
if nmuons <2: continue
muonCollection=[]
for imu,mu in enumerate(muons.product()):
if len(muonCollection)>1: continue
if len(muonCollection)==0 and mu.pt() < 20: continue
if len(muonCollection)==1 and mu.pt() < 10: continue
if not mu.isLooseMuon(): continue
if abs(mu.eta())>2.4: continue
if mu.trackIso()/mu.pt()>0.10: continue #Loose working point ;
#if mu.particleIso()>0.25: always -1, why?
muonCollection.append(mu.p4())
if len(muonCollection)<2: continue
mu1_pt[0]=muonCollection[0].pt()
mu1_eta[0]=muonCollection[0].eta()
mu1_phi[0]=muonCollection[0].phi()
mu2_pt[0]=muonCollection[1].pt()
mu2_eta[0]=muonCollection[1].eta()
mu2_phi[0]=muonCollection[1].phi()
Zvector=muonCollection[0]+muonCollection[1] #Sum of two muon Lerentz vectors
h_Zmass.Fill(Zvector.M())
if Zvector.M()<75 or Zvector.M()>105: continue
Zmass[0]=Zvector.M()
if PLOTRAW:
pfmetvector = mets.product().front().uncorP4()
puppimetvector = pmets.product().front().uncorP4()
else:
pfmetvector = mets.product().front().corP4()
puppimetvector = pmets.product().front().corP4()
pfrecoil = pfmetvector + Zvector
puppirecoil = puppimetvector + Zvector
qt[0] = Zvector.pt()
q_eta[0] = Zvector.eta()
u_pll_pf[0] = - pfrecoil.pt() * cos(Zvector.phi()-pfrecoil.phi())
u_pll_puppi[0] = - puppirecoil.pt() * cos(Zvector.phi()-puppirecoil.phi())
u_prp_pf[0] = pfmetvector.pt() * sin(Zvector.phi()-pfmetvector.phi())
u_prp_puppi[0] = puppimetvector.pt() * sin(Zvector.phi()-puppimetvector.phi())
h_qeta.Fill(Zvector.eta())
h_qrapidity.Fill(Zvector.Rapidity())
h_qt.Fill(qt[0])
h_upll.Fill(u_pll_pf[0])
h_upllqt.Fill(u_pll_pf[0]+qt[0])
h_uprp.Fill(u_prp_pf[0])
t.Fill()
canv = ROOT.TCanvas("canv","canv",1000,1000);
h_Zmass.Draw();canv.Print("result/zmass.png");
h_qeta.Draw();canv.Print("result/qeta.png");
h_qrapidity.Draw();canv.Print("result/qrapidity.png");
h_qt.Draw();canv.Print("result/qt.png");
h_upll.Draw();canv.Print("result/upll.png");
h_upllqt.Draw();canv.Print("result/upllqt.png");
h_uprp.Draw();canv.Print("result/uprp.png");
fout.Write()
fout.Close()