Fork me on GitHub

source: git/examples/Example3.C@ f2d7d0c

ImprovedOutputFile Timing dual_readout llp
Last change on this file since f2d7d0c was f2d7d0c, checked in by pavel <pavel@…>, 11 years ago

fix cards and Example3

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