Fork me on GitHub

source: git/examples/Example7.py@ d612dec

Last change on this file since d612dec was b03c0d0, checked in by Michele Selvaggi <michele.selvaggi@…>, 4 years ago

updated examples to include gen/LHE weights

  • Property mode set to 100644
File size: 2.5 KB
Line 
1#!/usr/bin/env python
2
3import sys
4
5import ROOT
6
7try:
8 input = raw_input
9except:
10 pass
11
12if len(sys.argv) < 2:
13 print(" Usage: Example1.py input_file")
14 sys.exit(1)
15
16ROOT.gSystem.Load("libDelphes")
17
18try:
19 ROOT.gInterpreter.Declare('#include "classes/DelphesClasses.h"')
20 ROOT.gInterpreter.Declare('#include "external/ExRootAnalysis/ExRootTreeReader.h"')
21except:
22 pass
23
24inputFile = sys.argv[1]
25
26# Create chain of root trees
27chain = ROOT.TChain("Delphes")
28chain.Add(inputFile)
29
30# Create object of class ExRootTreeReader
31treeReader = ROOT.ExRootTreeReader(chain)
32numberOfEntries = treeReader.GetEntries()
33
34# Get pointers to branches used in this analysis
35branchJet = treeReader.UseBranch("JetPUPPITight")
36branchElectron = treeReader.UseBranch("ElectronMedium")
37branchWeight = treeReader.UseBranch("Weight")
38branchEvent = treeReader.UseBranch("Event")
39
40# Book histograms
41histJetPT = ROOT.TH1F("jet_pt", "jet P_{T}", 100, 0.0, 1000.0)
42histElectronPT = ROOT.TH1F("Electron_pt", "electron P_{T}", 100, 0.0, 1000.0)
43
44# Loop over all events
45for entry in range(0, numberOfEntries):
46 # Load selected branches with data from specified event
47 treeReader.ReadEntry(entry)
48
49 ## main MC event weight
50 w = branchEvent[0].Weight
51
52 ## read lhe event weight
53 for weight in branchWeight:
54 lhe_weight = weight.Weight
55 ## do stuff ...
56 ## print lhe_weight
57
58 # If event contains at least 1 jet
59 if branchJet.GetEntries() > 0:
60 # Take first jet
61 jet = branchJet.At(0)
62
63 ## 0 - Loose , 1 - Medium, 2 - Tight
64 wp = 1
65
66 BtagOk = ( jet.BTag & (1 << wp) )
67 pt = jet.PT
68 eta = abs(jet.Eta)
69
70 # Plot jet transverse momentum
71 if (BtagOk and pt > 30. and eta < 5.):
72 histJetPT.Fill(jet.PT, w)
73
74
75 # If event contains at least 1 electron
76 if branchElectron.GetEntries() > 0:
77 # Take first electron
78 electron = branchElectron.At(0)
79
80 pt = electron.PT
81 eta = abs(electron.Eta)
82
83 ## looseCut = 0.3, mediumCut = 0.2, tightCut = 0.1
84 IsoCut = 0.2
85 IsoOk = electron.IsolationVar < IsoCut
86
87 # Plot electron transverse momentum
88 if (IsoOk and pt > 10. and eta < 5.):
89 histElectronPT.Fill(electron.PT, w)
90
91
92# Show resulting histograms
93cnv = ROOT.TCanvas("cnv", "cnv", 50, 50, 800, 500)
94cnv.Divide(2, 1)
95cnv.cd(1)
96ROOT.gStyle.SetOptStat(0)
97
98histJetPT.Draw()
99
100cnv.cd(2)
101ROOT.gStyle.SetOptStat(0)
102histElectronPT.Draw()
103
104input("Press Enter to continue...")
Note: See TracBrowser for help on using the repository browser.