Fork me on GitHub

source: git/examples/Example3.C@ ff9f67fb

ImprovedOutputFile Timing dual_readout llp
Last change on this file since ff9f67fb 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: 8.4 KB
RevLine 
[d7d2da3]1/*
[fc4e460]2This macro shows how to access the particle-level reference for reconstructed objects.
[fbb617b]3It is also shown how to loop over the jet constituents.
[fc4e460]4
[21f3c04]5root -l examples/Example3.C'("delphes_output.root")'
[d7d2da3]6*/
7
[c9803f7]8#ifdef __CLING__
9R__LOAD_LIBRARY(libDelphes)
10#include "classes/DelphesClasses.h"
11#include "external/ExRootAnalysis/ExRootTreeReader.h"
12#include "external/ExRootAnalysis/ExRootResult.h"
13#else
14class ExRootTreeReader;
15class ExRootResult;
16#endif
17
[d7d2da3]18//------------------------------------------------------------------------------
19
20struct TestPlots
21{
22 TH1 *fElectronDeltaPT;
23 TH1 *fElectronDeltaEta;
24
25 TH1 *fPhotonDeltaPT;
26 TH1 *fPhotonDeltaEta;
27 TH1 *fPhotonDeltaE;
28
29 TH1 *fMuonDeltaPT;
30 TH1 *fMuonDeltaEta;
31
32 TH1 *fTrackDeltaPT;
33 TH1 *fTrackDeltaEta;
34
35 TH1 *fJetDeltaPT;
36};
37
38//------------------------------------------------------------------------------
39
40void BookHistograms(ExRootResult *result, TestPlots *plots)
41{
42 TLegend *legend;
43 TPaveText *comment;
44
45 plots->fElectronDeltaPT = result->AddHist1D(
[fc4e460]46 "electron_delta_pt", "(p_{T}^{particle} - p_{T}^{electron})/p_{T}^{particle}",
[d7d2da3]47 "(p_{T}^{particle} - p_{T}^{electron})/p_{T}^{particle}", "number of electrons",
48 100, -0.1, 0.1);
49
50 plots->fElectronDeltaEta = result->AddHist1D(
[fc4e460]51 "electron_delta_eta", "(#eta^{particle} - #eta^{electron})/#eta^{particle}",
[d7d2da3]52 "(#eta^{particle} - #eta^{electron})/#eta^{particle}", "number of electrons",
53 100, -0.1, 0.1);
54
55 plots->fPhotonDeltaPT = result->AddHist1D(
[fc4e460]56 "photon_delta_pt", "(p_{T}^{particle} - p_{T}^{photon})/p_{T}^{particle}",
[d7d2da3]57 "(p_{T}^{particle} - p_{T}^{photon})/p_{T}^{particle}", "number of photons",
58 100, -0.1, 0.1);
59
60 plots->fPhotonDeltaEta = result->AddHist1D(
[fc4e460]61 "photon_delta_eta", "(#eta^{particle} - #eta^{photon})/#eta^{particle}",
[d7d2da3]62 "(#eta^{particle} - #eta^{photon})/#eta^{particle}", "number of photons",
63 100, -0.1, 0.1);
64
65 plots->fPhotonDeltaE = result->AddHist1D(
[fc4e460]66 "photon_delta_energy", "(E^{particle} - E^{photon})/E^{particle}",
[d7d2da3]67 "(E^{particle} - E^{photon})/E^{particle}", "number of photons",
68 100, -0.1, 0.1);
69
70 plots->fMuonDeltaPT = result->AddHist1D(
[fc4e460]71 "muon_delta_pt", "(p_{T}^{particle} - p_{T}^{muon})/p_{T}^{particle}",
[d7d2da3]72 "(p_{T}^{particle} - p_{T}^{muon})/p_{T}^{particle}", "number of muons",
73 100, -0.1, 0.1);
74
75 plots->fMuonDeltaEta = result->AddHist1D(
[fc4e460]76 "muon_delta_eta", "(#eta^{particle} - #eta^{muon})/#eta^{particle}",
[d7d2da3]77 "(#eta^{particle} - #eta^{muon})/#eta^{particle}", "number of muons",
78 100, -0.1, 0.1);
79
80 plots->fTrackDeltaPT = result->AddHist1D(
[fc4e460]81 "track_delta_pt", "(p_{T}^{particle} - p_{T}^{track})/p_{T}^{particle}",
[d7d2da3]82 "(p_{T}^{particle} - p_{T}^{track})/p_{T}^{particle}", "number of tracks",
83 100, -0.1, 0.1);
84
85 plots->fTrackDeltaEta = result->AddHist1D(
[fc4e460]86 "track_delta_eta", "(#eta^{particle} - #eta^{track})/#eta^{particle}",
[d7d2da3]87 "(#eta^{particle} - #eta^{track})/#eta^{particle}", "number of tracks",
88 100, -0.1, 0.1);
89
90 plots->fJetDeltaPT = result->AddHist1D(
[fc4e460]91 "jet_delta_pt", "(p_{T}^{jet} - p_{T}^{constituents})/p_{T}^{jet}",
[d7d2da3]92 "(p_{T}^{jet} - p_{T}^{constituents})/p_{T}^{jet}", "number of jets",
[8555f6d]93 100, -1.0e-1, 1.0e-1);
[d7d2da3]94
95}
96
97//------------------------------------------------------------------------------
98
99void AnalyseEvents(ExRootTreeReader *treeReader, TestPlots *plots)
100{
101 TClonesArray *branchParticle = treeReader->UseBranch("Particle");
102 TClonesArray *branchElectron = treeReader->UseBranch("Electron");
103 TClonesArray *branchPhoton = treeReader->UseBranch("Photon");
104 TClonesArray *branchMuon = treeReader->UseBranch("Muon");
105
[f2d7d0c]106 TClonesArray *branchTrack = treeReader->UseBranch("Track");
107 TClonesArray *branchTower = treeReader->UseBranch("Tower");
108
[d7d2da3]109 TClonesArray *branchEFlowTrack = treeReader->UseBranch("EFlowTrack");
[27bf162]110 TClonesArray *branchEFlowPhoton = treeReader->UseBranch("EFlowPhoton");
111 TClonesArray *branchEFlowNeutralHadron = treeReader->UseBranch("EFlowNeutralHadron");
[d7d2da3]112 TClonesArray *branchJet = treeReader->UseBranch("Jet");
113
114 Long64_t allEntries = treeReader->GetEntries();
115
116 cout << "** Chain contains " << allEntries << " events" << endl;
117
118 GenParticle *particle;
119 Electron *electron;
120 Photon *photon;
121 Muon *muon;
122
123 Track *track;
124 Tower *tower;
125
126 Jet *jet;
127 TObject *object;
128
129 TLorentzVector momentum;
130
131 Float_t Eem, Ehad;
132 Bool_t skip;
133
134 Long64_t entry;
135
136 Int_t i, j, pdgCode;
137
138 // Loop over all events
139 for(entry = 0; entry < allEntries; ++entry)
140 {
141 // Load selected branches with data from specified event
142 treeReader->ReadEntry(entry);
143
144 // Loop over all electrons in event
145 for(i = 0; i < branchElectron->GetEntriesFast(); ++i)
146 {
147 electron = (Electron*) branchElectron->At(i);
148 particle = (GenParticle*) electron->Particle.GetObject();
149
150 plots->fElectronDeltaPT->Fill((particle->PT - electron->PT)/particle->PT);
151 plots->fElectronDeltaEta->Fill((particle->Eta - electron->Eta)/particle->Eta);
152 }
153
154 // Loop over all photons in event
155 for(i = 0; i < branchPhoton->GetEntriesFast(); ++i)
156 {
157 photon = (Photon*) branchPhoton->At(i);
158
159 // skip photons with references to multiple particles
160 if(photon->Particles.GetEntriesFast() != 1) continue;
161
162 particle = (GenParticle*) photon->Particles.At(0);
163
164 plots->fPhotonDeltaPT->Fill((particle->PT - photon->PT)/particle->PT);
165 plots->fPhotonDeltaEta->Fill((particle->Eta - photon->Eta)/particle->Eta);
166 plots->fPhotonDeltaE->Fill((particle->E - photon->E)/particle->E);
167 }
168
169 // Loop over all muons in event
170 for(i = 0; i < branchMuon->GetEntriesFast(); ++i)
171 {
172 muon = (Muon*) branchMuon->At(i);
173 particle = (GenParticle*) muon->Particle.GetObject();
174
175 plots->fMuonDeltaPT->Fill((particle->PT - muon->PT)/particle->PT);
176 plots->fMuonDeltaEta->Fill((particle->Eta - muon->Eta)/particle->Eta);
177 }
178
179 // Loop over all tracks in event
[f2d7d0c]180 for(i = 0; i < branchTrack->GetEntriesFast(); ++i)
[d7d2da3]181 {
[f2d7d0c]182 track = (Track*) branchTrack->At(i);
[d7d2da3]183 particle = (GenParticle*) track->Particle.GetObject();
184
185 plots->fTrackDeltaPT->Fill((particle->PT - track->PT)/particle->PT);
186 plots->fTrackDeltaEta->Fill((particle->Eta - track->Eta)/particle->Eta);
187 }
188
[fbb617b]189 // cout << "-- New event -- " << endl;
[d7d2da3]190
191 // Loop over all jets in event
192 for(i = 0; i < branchJet->GetEntriesFast(); ++i)
193 {
194 jet = (Jet*) branchJet->At(i);
195
196 momentum.SetPxPyPzE(0.0, 0.0, 0.0, 0.0);
197
[fbb617b]198 // cout<<"Looping over jet constituents. Jet pt: "<<jet->PT<<", eta: "<<jet->Eta<<", phi: "<<jet->Phi<<endl;
[27bf162]199
[d7d2da3]200 // Loop over all jet's constituents
201 for(j = 0; j < jet->Constituents.GetEntriesFast(); ++j)
202 {
203 object = jet->Constituents.At(j);
204
205 // Check if the constituent is accessible
206 if(object == 0) continue;
207
208 if(object->IsA() == GenParticle::Class())
209 {
[f2d7d0c]210 particle = (GenParticle*) object;
[fbb617b]211 // cout << " GenPart pt: " << particle->PT << ", eta: " << particle->Eta << ", phi: " << particle->Phi << endl;
[f2d7d0c]212 momentum += particle->P4();
[d7d2da3]213 }
214 else if(object->IsA() == Track::Class())
215 {
[f2d7d0c]216 track = (Track*) object;
[fbb617b]217 // cout << " Track pt: " << track->PT << ", eta: " << track->Eta << ", phi: " << track->Phi << endl;
[f2d7d0c]218 momentum += track->P4();
[d7d2da3]219 }
220 else if(object->IsA() == Tower::Class())
221 {
[f2d7d0c]222 tower = (Tower*) object;
[fbb617b]223 // cout << " Tower pt: " << tower->ET << ", eta: " << tower->Eta << ", phi: " << tower->Phi << endl;
[f2d7d0c]224 momentum += tower->P4();
[d7d2da3]225 }
226 }
[8555f6d]227 plots->fJetDeltaPT->Fill((jet->PT - momentum.Pt())/jet->PT);
[d7d2da3]228 }
229 }
230}
231
232//------------------------------------------------------------------------------
233
234void PrintHistograms(ExRootResult *result, TestPlots *plots)
235{
236 result->Print("png");
237}
238
239//------------------------------------------------------------------------------
240
241void Example3(const char *inputFile)
242{
243 gSystem->Load("libDelphes");
244
245 TChain *chain = new TChain("Delphes");
246 chain->Add(inputFile);
247
248 ExRootTreeReader *treeReader = new ExRootTreeReader(chain);
249 ExRootResult *result = new ExRootResult();
250
251 TestPlots *plots = new TestPlots;
252
253 BookHistograms(result, plots);
254
255 AnalyseEvents(treeReader, plots);
256
257 PrintHistograms(result, plots);
258
259 result->Write("results.root");
260
261 cout << "** Exiting..." << endl;
262
263 delete plots;
264 delete result;
265 delete treeReader;
266 delete chain;
267}
268
269//------------------------------------------------------------------------------
Note: See TracBrowser for help on using the repository browser.