Fork me on GitHub

source: git/examples/Example3.C@ 6fb1a5d

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