Fork me on GitHub

source: git/examples/Example2.C@ 14ae668

ImprovedOutputFile Timing dual_readout llp
Last change on this file since 14ae668 was e7e90df, checked in by Michele <michele.selvaggi@…>, 10 years ago

Merge branch 'master' into TestFastJet310b1

Conflicts:

Makefile
examples/Example1.C
examples/Example2.C
examples/Example3.C

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