source: trunk/test/ExRootLHCOlympicsWriter.cpp@ 8

Last change on this file since 8 was 8, checked in by Pavel Demin, 16 years ago

first working version

  • Property svn:executable set to *
File size: 11.4 KB
RevLine 
[8]1
2#include <iostream>
3#include <fstream>
4#include <sstream>
5#include <map>
6
7#include "TROOT.h"
8#include "TApplication.h"
9
10#include "TFile.h"
11#include "TChain.h"
12#include "TString.h"
13
14#include "TH2.h"
15#include "THStack.h"
16#include "TLegend.h"
17#include "TPaveText.h"
18#include "TLorentzVector.h"
19
20#include "ExRootAnalysis/ExRootClasses.h"
21
22#include "ExRootAnalysis/ExRootTreeReader.h"
23
24#include "ExRootAnalysis/ExRootUtilities.h"
25#include "ExRootAnalysis/ExRootProgressBar.h"
26
27using namespace std;
28
29/*
30LHC Olympics format discription from http://www.jthaler.net/olympicswiki/doku.php?id=lhc_olympics:data_file_format
31
32 * The first column of each row is just a counter that labels the object.
33 * The event begins with a row labelled "0"; this row contains the event number and the triggering information. The last row of the event is always the missing transverse momentum (MET).
34 * The second column of each row gives the type of object being listed [0, 1, 2, 3, 4, 6 = photon, electron, muon, hadronically-decaying tau, jet, missing transverse energy].
35 * The next three columns give the pseudorapidity, the azimuthal angle, and the transverse momentum of the object.
36 * The sixth column gives the invariant mass of the object.
37 * The seventh column gives the number of tracks associated with the object; in the case of a lepton, this number is multiplied by the charge of the lepton.
38 * The eighth column is 1 or 2 for a jet that has been "tagged" as containing a b-quark (actually a heavy flavor tag that sometimes indicates c-quarks), otherwise it is 0. For muons, the integer part of this number is the identity of the jet (see column 1) that is closest ot this muon in Delta R.
39 * The ninth column is the ratio of the hadronic versus electromagnetic energy deposited in the calorimeter cells associated with the object. For muons to the left of the decimal point is the summed pT in a R=0.4 cone (excluding the muon). To the right of the decimal point is etrat, which is a percentage between .00 and .99. It is the ratio of the transverse energy in a 3x3 grid surrounding the muon to the pT of the muon.
40*/
41
42struct LHCOlympicsObject
43{
44 enum {maxIntParam = 2, maxDblParam = 9};
45
46 Int_t intParam[maxIntParam];
47 Double_t dblParam[maxDblParam];
48};
49
50//------------------------------------------------------------------------------
51
52class LHCOlympicsWriter
53{
54public:
55 LHCOlympicsWriter(ExRootTreeReader *treeReader,
56 ofstream *outputStream);
57 ~LHCOlympicsWriter();
58
59 void ProcessEvent();
60
61private:
62
63 void ResetLHCOlympicsObject();
64 void WriteLHCOlympicsObject();
65
66 void AnalyseEvent();
67
68 void AnalysePhotons();
69 void AnalyseElectrons();
70 void AnalyseMuons();
71 void AnalyseTaus();
72 void AnalyseJets();
73
74 void AnalyseMissingET();
75
76 LHCOlympicsObject fCurrentObject;
77
78 Long64_t fTriggerWord, fEventNumber;
79
80 ExRootTreeReader *fTreeReader;
81 ofstream *fOutputStream;
82
83 TClonesArray *fBranchEvent;
84 TClonesArray *fBranchPhoton;
85 TClonesArray *fBranchElectron;
86 TClonesArray *fBranchMuon;
87 TClonesArray *fBranchTau;
88 TClonesArray *fBranchJet;
89 TClonesArray *fBranchMissingET;
90
91 TIterator *fItPhoton;
92 TIterator *fItElectron;
93 TIterator *fItMuon;
94 TIterator *fItTau;
95 TIterator *fItJet;
96
97};
98
99//------------------------------------------------------------------------------
100
101LHCOlympicsWriter::LHCOlympicsWriter(ExRootTreeReader *treeReader,
102 ofstream *outputStream) :
103 fTriggerWord(0), fEventNumber(1), fTreeReader(0), fOutputStream(0)
104{
105 fTreeReader = treeReader;
106 fOutputStream = outputStream;
107
108 // information about reconstructed event
109 fBranchEvent = fTreeReader->UseBranch("Event");
110 // reconstructed photons
111 fBranchPhoton = fTreeReader->UseBranch("Photon");
112 fItPhoton = fBranchPhoton->MakeIterator();
113 // reconstructed electrons
114 fBranchElectron = fTreeReader->UseBranch("Electron");
115 fItElectron = fBranchElectron->MakeIterator();
116 // reconstructed muons
117 fBranchMuon = fTreeReader->UseBranch("Muon");
118 fItMuon = fBranchMuon->MakeIterator();
119 // reconstructed hadronically-decaying tau leptons
120 fBranchTau = fTreeReader->UseBranch("Tau");
121 fItTau = fBranchTau->MakeIterator();
122 // reconstructed jets
123 fBranchJet = fTreeReader->UseBranch("Jet");
124 fItJet = fBranchJet->MakeIterator();
125 // missing transverse energy
126 fBranchMissingET = fTreeReader->UseBranch("MissingET");
127}
128
129//------------------------------------------------------------------------------
130
131LHCOlympicsWriter::~LHCOlympicsWriter()
132{
133}
134
135//---------------------------------------------------------------------------
136
137void LHCOlympicsWriter::ProcessEvent()
138{
139 fCurrentObject.intParam[0] = 0;
140
141 AnalyseEvent();
142
143 AnalysePhotons();
144 AnalyseElectrons();
145 AnalyseMuons();
146 AnalyseTaus();
147 AnalyseJets();
148
149 AnalyseMissingET();
150}
151
152//---------------------------------------------------------------------------
153
154void LHCOlympicsWriter::ResetLHCOlympicsObject()
155{
156 int i;
157 for(i = 1; i < LHCOlympicsObject::maxIntParam; ++i)
158 {
159 fCurrentObject.intParam[i] = 0;
160 }
161
162 for(i = 0; i < LHCOlympicsObject::maxDblParam; ++i)
163 {
164 fCurrentObject.dblParam[i] = 0.0;
165 }
166}
167
168//---------------------------------------------------------------------------
169
170void LHCOlympicsWriter::WriteLHCOlympicsObject()
171{
172 int i;
173 for(i = 0; i < LHCOlympicsObject::maxIntParam; ++i)
174 {
175 (*fOutputStream) << fCurrentObject.intParam[i] << "\t";
176 }
177
178 for(i = 0; i < LHCOlympicsObject::maxDblParam; ++i)
179 {
180 (*fOutputStream) << fCurrentObject.dblParam[i] << "\t";
181 }
182
183 (*fOutputStream) << endl;
184
185 ++fCurrentObject.intParam[0];
186}
187
188//---------------------------------------------------------------------------
189
190void LHCOlympicsWriter::AnalyseEvent()
191{
192 ExRootEvent *element;
193
194 element = static_cast<ExRootEvent*>(fBranchEvent->At(0));
195
196 (*fOutputStream) << fCurrentObject.intParam[0] << "\t";
197 (*fOutputStream) << element->Number << "\t";
198 (*fOutputStream) << element->Trigger << endl;
199
200 ++fCurrentObject.intParam[0];
201}
202
203//---------------------------------------------------------------------------
204
205void LHCOlympicsWriter::AnalysePhotons()
206{
207 ExRootPhoton *element;
208
209 fItPhoton->Reset();
210 while((element = static_cast<ExRootPhoton*>(fItPhoton->Next())))
211 {
212 ResetLHCOlympicsObject();
213
214 fCurrentObject.intParam[1] = 0;
215
216 fCurrentObject.dblParam[0] = element->Eta;
217 fCurrentObject.dblParam[1] = element->Phi;
218 fCurrentObject.dblParam[2] = element->PT;
219
220 fCurrentObject.dblParam[6] = element->EhadOverEem;
221
222 WriteLHCOlympicsObject();
223 }
224}
225
226//---------------------------------------------------------------------------
227
228void LHCOlympicsWriter::AnalyseElectrons()
229{
230 ExRootElectron *element;
231
232 fItElectron->Reset();
233 while((element = static_cast<ExRootElectron*>(fItElectron->Next())))
234 {
235 ResetLHCOlympicsObject();
236
237 fCurrentObject.intParam[1] = 1;
238
239 fCurrentObject.dblParam[0] = element->Eta;
240 fCurrentObject.dblParam[1] = element->Phi;
241 fCurrentObject.dblParam[2] = element->PT;
242
243 fCurrentObject.dblParam[4] = element->Ntrk * element->Charge;
244
245 fCurrentObject.dblParam[6] = element->EhadOverEem;
246
247 WriteLHCOlympicsObject();
248 }
249}
250
251//---------------------------------------------------------------------------
252
253void LHCOlympicsWriter::AnalyseMuons()
254{
255 ExRootMuon *element;
256
257 fItMuon->Reset();
258 while((element = static_cast<ExRootMuon*>(fItMuon->Next())))
259 {
260 ResetLHCOlympicsObject();
261
262 fCurrentObject.intParam[1] = 2;
263
264 fCurrentObject.dblParam[0] = element->Eta;
265 fCurrentObject.dblParam[1] = element->Phi;
266 fCurrentObject.dblParam[2] = element->PT;
267
268 fCurrentObject.dblParam[4] = element->Ntrk * element->Charge;
269
270 fCurrentObject.dblParam[5] = element->JetIndex;
271
272 fCurrentObject.dblParam[6] = element->PTiso + element->ETiso;
273
274 WriteLHCOlympicsObject();
275 }
276}
277
278//---------------------------------------------------------------------------
279
280void LHCOlympicsWriter::AnalyseTaus()
281{
282 ExRootTau *element;
283
284 fItTau->Reset();
285 while((element = static_cast<ExRootTau*>(fItTau->Next())))
286 {
287 ResetLHCOlympicsObject();
288
289 fCurrentObject.intParam[1] = 3;
290
291 fCurrentObject.dblParam[0] = element->Eta;
292 fCurrentObject.dblParam[1] = element->Phi;
293 fCurrentObject.dblParam[2] = element->PT;
294
295 fCurrentObject.dblParam[4] = element->Ntrk * element->Charge;
296
297 fCurrentObject.dblParam[6] = element->EhadOverEem;
298
299 WriteLHCOlympicsObject();
300 }
301}
302
303//---------------------------------------------------------------------------
304
305void LHCOlympicsWriter::AnalyseJets()
306{
307 ExRootJet *element;
308
309 fItJet->Reset();
310 while((element = static_cast<ExRootJet*>(fItJet->Next())))
311 {
312 ResetLHCOlympicsObject();
313
314 fCurrentObject.intParam[1] = 4;
315
316 fCurrentObject.dblParam[0] = element->Eta;
317 fCurrentObject.dblParam[1] = element->Phi;
318 fCurrentObject.dblParam[2] = element->PT;
319
320 fCurrentObject.dblParam[3] = element->Mass;
321
322 fCurrentObject.dblParam[4] = element->Ntrk;
323
324 fCurrentObject.dblParam[5] = element->BTag;
325
326 fCurrentObject.dblParam[6] = element->EhadOverEem;
327
328 WriteLHCOlympicsObject();
329 }
330}
331
332//---------------------------------------------------------------------------
333
334void LHCOlympicsWriter::AnalyseMissingET()
335{
336 ExRootMissingET *element;
337
338 element = static_cast<ExRootMissingET*>(fBranchMissingET->At(0));
339
340 ResetLHCOlympicsObject();
341
342 fCurrentObject.intParam[1] = 6;
343
344 fCurrentObject.dblParam[1] = element->Phi;
345 fCurrentObject.dblParam[2] = element->MET;
346
347 WriteLHCOlympicsObject();
348}
349
350//---------------------------------------------------------------------------
351
352int main(int argc, char *argv[])
353{
354 char *appName = "ExRootLHCOlympicsWriter";
355
356 if(argc != 3)
357 {
358 cout << " Usage: " << appName << " input_file" << " output_file" << endl;
359 cout << " input_file - input file in ROOT format." << endl;
360 cout << " output_file - output file in LHEF format," << endl;
361 return 1;
362 }
363
364 int appargc = 1;
365 char *appargv[] = {appName};
366 TApplication app(appName, &appargc, appargv);
367
368 // Open a stream connected to an output file:
369 ofstream outputFileStream(argv[2]);
370
371 if(!outputFileStream.is_open())
372 {
373 cerr << "** ERROR: Can't open '" << argv[2] << "' for output" << endl;
374 return 1;
375 }
376
377 TChain *chain = new TChain("LHCO");
378 chain->Add(argv[1]);
379
380 ExRootTreeReader *treeReader = new ExRootTreeReader(chain);
381
382 cout << "** Calculating number of events to process. Please wait..." << endl;
383 Long64_t allEntries = treeReader->GetEntries();
384 cout << "** Input file contains " << allEntries << " events" << endl;
385
386 Long64_t entry;
387
388 if(allEntries > 0)
389 {
390 // Create LHC Olympics converter:
391 LHCOlympicsWriter *writer = new LHCOlympicsWriter(treeReader, &outputFileStream);
392
393 ExRootProgressBar progressBar(allEntries);
394 // Loop over all events
395 for(entry = 0; entry < allEntries; ++entry)
396 {
397 if(!treeReader->ReadEntry(entry))
398 {
399 cout << "** ERROR: cannot read event " << entry << endl;
400 break;
401 }
402
403 writer->ProcessEvent();
404
405 progressBar.Update(entry);
406 }
407 progressBar.Finish();
408
409 delete writer;
410 }
411
412 cout << "** Exiting..." << endl;
413
414 delete treeReader;
415 delete chain;
416}
417
418
Note: See TracBrowser for help on using the repository browser.