Fork me on GitHub

source: git/examples/Example3.C@ 239e1d0

ImprovedOutputFile Timing dual_readout llp 3.3.3pre11
Last change on this file since 239e1d0 was 4bf6061, checked in by Michele Selvaggi <michele.selvaggi@…>, 9 years ago

revert last commit

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