Fork me on GitHub

source: git/examples/Example2.C@ 3a105e5

Last change on this file since 3a105e5 was c9803f7, checked in by Pavel Demin <pavel.demin@…>, 9 years ago

adapt examples to ROOT 6.04

  • Property mode set to 100644
File size: 4.7 KB
Line 
1/*
2Simple macro showing how to access branches from the delphes output root file,
3loop over events, store histograms in a root file and print them as image files.
4
5root -l examples/Example2.C'("delphes_output.root")'
6*/
7
8#include "TH1.h"
9#include "TSystem.h"
10
11#ifdef __CLING__
12R__LOAD_LIBRARY(libDelphes)
13#include "classes/DelphesClasses.h"
14#include "external/ExRootAnalysis/ExRootTreeReader.h"
15#include "external/ExRootAnalysis/ExRootResult.h"
16#endif
17
18//------------------------------------------------------------------------------
19
20struct MyPlots
21{
22 TH1 *fJetPT[2];
23 TH1 *fMissingET;
24 TH1 *fElectronPT;
25};
26
27//------------------------------------------------------------------------------
28
29class ExRootResult;
30class ExRootTreeReader;
31
32//------------------------------------------------------------------------------
33
34void BookHistograms(ExRootResult *result, MyPlots *plots)
35{
36 THStack *stack;
37 TLegend *legend;
38 TPaveText *comment;
39
40 // book 2 histograms for PT of 1st and 2nd leading jets
41
42 plots->fJetPT[0] = result->AddHist1D(
43 "jet_pt_0", "leading jet P_{T}",
44 "jet P_{T}, GeV/c", "number of jets",
45 50, 0.0, 100.0);
46
47 plots->fJetPT[1] = result->AddHist1D(
48 "jet_pt_1", "2nd leading jet P_{T}",
49 "jet P_{T}, GeV/c", "number of jets",
50 50, 0.0, 100.0);
51
52 plots->fJetPT[0]->SetLineColor(kRed);
53 plots->fJetPT[1]->SetLineColor(kBlue);
54
55 // book 1 stack of 2 histograms
56
57 stack = result->AddHistStack("jet_pt_all", "1st and 2nd jets P_{T}");
58 stack->Add(plots->fJetPT[0]);
59 stack->Add(plots->fJetPT[1]);
60
61 // book legend for stack of 2 histograms
62
63 legend = result->AddLegend(0.72, 0.86, 0.98, 0.98);
64 legend->AddEntry(plots->fJetPT[0], "leading jet", "l");
65 legend->AddEntry(plots->fJetPT[1], "second jet", "l");
66
67 // attach legend to stack (legend will be printed over stack in .eps file)
68
69 result->Attach(stack, legend);
70
71 // book more histograms
72
73 plots->fElectronPT = result->AddHist1D(
74 "electron_pt", "electron P_{T}",
75 "electron P_{T}, GeV/c", "number of electrons",
76 50, 0.0, 100.0);
77
78 plots->fMissingET = result->AddHist1D(
79 "missing_et", "Missing E_{T}",
80 "Missing E_{T}, GeV", "number of events",
81 60, 0.0, 30.0);
82
83 // book general comment
84
85 comment = result->AddComment(0.64, 0.86, 0.98, 0.98);
86 comment->AddText("demonstration plot");
87 comment->AddText("produced by Example2.C");
88
89 // attach comment to single histograms
90
91 result->Attach(plots->fJetPT[0], comment);
92 result->Attach(plots->fJetPT[1], comment);
93 result->Attach(plots->fElectronPT, comment);
94
95 // show histogram statisics for MissingET
96 plots->fMissingET->SetStats();
97}
98
99//------------------------------------------------------------------------------
100
101void AnalyseEvents(ExRootTreeReader *treeReader, MyPlots *plots)
102{
103 TClonesArray *branchJet = treeReader->UseBranch("Jet");
104 TClonesArray *branchElectron = treeReader->UseBranch("Electron");
105 TClonesArray *branchMissingET = treeReader->UseBranch("MissingET");
106
107 Long64_t allEntries = treeReader->GetEntries();
108
109 cout << "** Chain contains " << allEntries << " events" << endl;
110
111 Jet *jet[2];
112 MissingET *met;
113 Electron *electron;
114
115 Long64_t entry;
116
117 Int_t i;
118
119 // Loop over all events
120 for(entry = 0; entry < allEntries; ++entry)
121 {
122 // Load selected branches with data from specified event
123 treeReader->ReadEntry(entry);
124
125 // Analyse two leading jets
126 if(branchJet->GetEntriesFast() >= 2)
127 {
128 jet[0] = (Jet*) branchJet->At(0);
129 jet[1] = (Jet*) branchJet->At(1);
130
131 plots->fJetPT[0]->Fill(jet[0]->PT);
132 plots->fJetPT[1]->Fill(jet[1]->PT);
133 }
134
135 // Analyse missing ET
136 if(branchMissingET->GetEntriesFast() > 0)
137 {
138 met = (MissingET*) branchMissingET->At(0);
139 plots->fMissingET->Fill(met->MET);
140 }
141
142 // Loop over all electrons in event
143 for(i = 0; i < branchElectron->GetEntriesFast(); ++i)
144 {
145 electron = (Electron*) branchElectron->At(i);
146 plots->fElectronPT->Fill(electron->PT);
147 }
148 }
149}
150
151//------------------------------------------------------------------------------
152
153void PrintHistograms(ExRootResult *result, MyPlots *plots)
154{
155 result->Print("png");
156}
157
158//------------------------------------------------------------------------------
159
160void Example2(const char *inputFile)
161{
162 gSystem->Load("libDelphes");
163
164 TChain *chain = new TChain("Delphes");
165 chain->Add(inputFile);
166
167 ExRootTreeReader *treeReader = new ExRootTreeReader(chain);
168 ExRootResult *result = new ExRootResult();
169
170 MyPlots *plots = new MyPlots;
171
172 BookHistograms(result, plots);
173
174 AnalyseEvents(treeReader, plots);
175
176 PrintHistograms(result, plots);
177
178 result->Write("results.root");
179
180 cout << "** Exiting..." << endl;
181
182 delete plots;
183 delete result;
184 delete treeReader;
185 delete chain;
186}
187
188//------------------------------------------------------------------------------
Note: See TracBrowser for help on using the repository browser.