Fork me on GitHub

source: git/examples/Example7.py@ 7e6c201

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

added snowmass examples

  • Property mode set to 100644
File size: 2.3 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")
38
39# Book histograms
40histJetPT = ROOT.TH1F("jet_pt", "jet P_{T}", 100, 0.0, 1000.0)
41histElectronPT = ROOT.TH1F("Electron_pt", "electron P_{T}", 100, 0.0, 1000.0)
42
43# Loop over all events
44for entry in range(0, numberOfEntries):
45 # Load selected branches with data from specified event
46 treeReader.ReadEntry(entry)
47
48 w = 1.0
49 ## read MC event weight
50 if branchWeight.GetEntries() > 0:
51 weight = branchWeight.At(0).Weight
52
53 # If event contains at least 1 jet
54 if branchJet.GetEntries() > 0:
55 # Take first jet
56 jet = branchJet.At(0)
57
58 ## 0 - Loose , 1 - Medium, 2 - Tight
59 wp = 1
60
61 BtagOk = ( jet.BTag & (1 << wp) )
62 pt = jet.PT
63 eta = abs(jet.Eta)
64
65 # Plot jet transverse momentum
66 if (BtagOk and pt > 30. and eta < 5.):
67 histJetPT.Fill(jet.PT, w)
68
69
70 # If event contains at least 1 electron
71 if branchElectron.GetEntries() > 0:
72 # Take first electron
73 electron = branchElectron.At(0)
74
75 pt = electron.PT
76 eta = abs(electron.Eta)
77
78 IsoCut = 0.2
79 IsoOk = electron.IsolationVar < IsoCut
80
81 # Plot electron transverse momentum
82 if (IsoOk and pt > 10. and eta < 5.):
83 histElectronPT.Fill(electron.PT, w)
84
85
86# Show resulting histograms
87cnv = ROOT.TCanvas("cnv", "cnv", 50, 50, 800, 500)
88cnv.Divide(2, 1)
89cnv.cd(1)
90ROOT.gStyle.SetOptStat(0)
91
92histJetPT.Draw()
93
94cnv.cd(2)
95ROOT.gStyle.SetOptStat(0)
96histElectronPT.Draw()
97
98input("Press Enter to continue...")
Note: See TracBrowser for help on using the repository browser.