Fork me on GitHub

source: git/examples/Example2.C@ d244bc9

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

removed muon check in jet constituent loop. Added comments to example macros

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