forked from delphes/delphes
-
Notifications
You must be signed in to change notification settings - Fork 0
/
Example2.C
188 lines (135 loc) · 4.73 KB
/
Example2.C
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
/*
Simple macro showing how to access branches from the delphes output root file,
loop over events, store histograms in a root file and print them as image files.
root -l examples/Example2.C'("delphes_output.root")'
*/
#include "TH1.h"
#include "TSystem.h"
#ifdef __CLING__
R__LOAD_LIBRARY(libDelphes)
#include "classes/DelphesClasses.h"
#include "external/ExRootAnalysis/ExRootTreeReader.h"
#include "external/ExRootAnalysis/ExRootResult.h"
#endif
//------------------------------------------------------------------------------
struct MyPlots
{
TH1 *fJetPT[2];
TH1 *fMissingET;
TH1 *fElectronPT;
};
//------------------------------------------------------------------------------
class ExRootResult;
class ExRootTreeReader;
//------------------------------------------------------------------------------
void BookHistograms(ExRootResult *result, MyPlots *plots)
{
THStack *stack;
TLegend *legend;
TPaveText *comment;
// book 2 histograms for PT of 1st and 2nd leading jets
plots->fJetPT[0] = result->AddHist1D(
"jet_pt_0", "leading jet P_{T}",
"jet P_{T}, GeV/c", "number of jets",
50, 0.0, 100.0);
plots->fJetPT[1] = result->AddHist1D(
"jet_pt_1", "2nd leading jet P_{T}",
"jet P_{T}, GeV/c", "number of jets",
50, 0.0, 100.0);
plots->fJetPT[0]->SetLineColor(kRed);
plots->fJetPT[1]->SetLineColor(kBlue);
// book 1 stack of 2 histograms
stack = result->AddHistStack("jet_pt_all", "1st and 2nd jets P_{T}");
stack->Add(plots->fJetPT[0]);
stack->Add(plots->fJetPT[1]);
// book legend for stack of 2 histograms
legend = result->AddLegend(0.72, 0.86, 0.98, 0.98);
legend->AddEntry(plots->fJetPT[0], "leading jet", "l");
legend->AddEntry(plots->fJetPT[1], "second jet", "l");
// attach legend to stack (legend will be printed over stack in .eps file)
result->Attach(stack, legend);
// book more histograms
plots->fElectronPT = result->AddHist1D(
"electron_pt", "electron P_{T}",
"electron P_{T}, GeV/c", "number of electrons",
50, 0.0, 100.0);
plots->fMissingET = result->AddHist1D(
"missing_et", "Missing E_{T}",
"Missing E_{T}, GeV", "number of events",
60, 0.0, 30.0);
// book general comment
comment = result->AddComment(0.64, 0.86, 0.98, 0.98);
comment->AddText("demonstration plot");
comment->AddText("produced by Example2.C");
// attach comment to single histograms
result->Attach(plots->fJetPT[0], comment);
result->Attach(plots->fJetPT[1], comment);
result->Attach(plots->fElectronPT, comment);
// show histogram statisics for MissingET
plots->fMissingET->SetStats();
}
//------------------------------------------------------------------------------
void AnalyseEvents(ExRootTreeReader *treeReader, MyPlots *plots)
{
TClonesArray *branchJet = treeReader->UseBranch("Jet");
TClonesArray *branchElectron = treeReader->UseBranch("Electron");
TClonesArray *branchMissingET = treeReader->UseBranch("MissingET");
Long64_t allEntries = treeReader->GetEntries();
cout << "** Chain contains " << allEntries << " events" << endl;
Jet *jet[2];
MissingET *met;
Electron *electron;
Long64_t entry;
Int_t i;
// Loop over all events
for(entry = 0; entry < allEntries; ++entry)
{
// Load selected branches with data from specified event
treeReader->ReadEntry(entry);
// Analyse two leading jets
if(branchJet->GetEntriesFast() >= 2)
{
jet[0] = (Jet*) branchJet->At(0);
jet[1] = (Jet*) branchJet->At(1);
plots->fJetPT[0]->Fill(jet[0]->PT);
plots->fJetPT[1]->Fill(jet[1]->PT);
}
// Analyse missing ET
if(branchMissingET->GetEntriesFast() > 0)
{
met = (MissingET*) branchMissingET->At(0);
plots->fMissingET->Fill(met->MET);
}
// Loop over all electrons in event
for(i = 0; i < branchElectron->GetEntriesFast(); ++i)
{
electron = (Electron*) branchElectron->At(i);
plots->fElectronPT->Fill(electron->PT);
}
}
}
//------------------------------------------------------------------------------
void PrintHistograms(ExRootResult *result, MyPlots *plots)
{
result->Print("png");
}
//------------------------------------------------------------------------------
void Example2(const char *inputFile)
{
gSystem->Load("libDelphes");
TChain *chain = new TChain("Delphes");
chain->Add(inputFile);
ExRootTreeReader *treeReader = new ExRootTreeReader(chain);
ExRootResult *result = new ExRootResult();
MyPlots *plots = new MyPlots;
BookHistograms(result, plots);
AnalyseEvents(treeReader, plots);
PrintHistograms(result, plots);
result->Write("results.root");
cout << "** Exiting..." << endl;
delete plots;
delete result;
delete treeReader;
delete chain;
}
//------------------------------------------------------------------------------