source: trunk/test/ExRootLHCOlympicsConverter.cpp@ 22

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

first commit

File size: 13.6 KB
RevLine 
[2]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/ExRootTreeWriter.h"
23#include "ExRootAnalysis/ExRootTreeBranch.h"
24
25#include "ExRootAnalysis/ExRootUtilities.h"
26#include "ExRootAnalysis/ExRootProgressBar.h"
27
28using namespace std;
29
30/*
31LHC Olympics format discription from http://www.jthaler.net/olympicswiki/doku.php?id=lhc_olympics:data_file_format
32
33 * The first column of each row is just a counter that labels the object.
34 * 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).
35 * 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].
36 * The next three columns give the pseudorapidity, the azimuthal angle, and the transverse momentum of the object.
37 * The sixth column gives the invariant mass of the object.
38 * 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.
39 * 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.
40 * 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.
41*/
42
43struct LHCOlympicsObject
44{
45 enum {maxIntParam = 2, maxDblParam = 7};
46
47 Int_t intParam[maxIntParam];
48 Double_t dblParam[maxDblParam];
49};
50
51//------------------------------------------------------------------------------
52
53class LHCOlympicsConverter
54{
55public:
56 LHCOlympicsConverter(const char *outputFileName);
57 ~LHCOlympicsConverter();
58
59 void ProcessObject();
60 void Write();
61
62 Long64_t GetNumberOfObjects(ifstream &inputFileStream);
63 Bool_t ReadObject(ifstream &inputFileStream);
64
65private:
66
67 void AddMissingEvents();
68
69 void AnalyseEvent(ExRootTreeBranch *branch,
70 Long64_t eventNumber, Int_t triggerWord);
71
72 void AnalysePhoton(ExRootTreeBranch *branch);
73 void AnalyseElectron(ExRootTreeBranch *branch);
74 void AnalyseMuon(ExRootTreeBranch *branch);
75 void AnalyseTau(ExRootTreeBranch *branch);
76 void AnalyseJet(ExRootTreeBranch *branch);
77 void AnalyseMissingET(ExRootTreeBranch *branch);
78
79 istringstream fBufferStream;
80 string fBuffer;
81
82 LHCOlympicsObject fCurrentObject;
83
84 Bool_t fIsFirstEvent, fIsNewEvent, fIsReadyToFill;
85
86 Long64_t fPreviousObjectNumber, fTriggerWord, fEventNumber, fRecordNumber;
87
88 TFile *fOutputFile;
89 ExRootTreeWriter *fTreeWriter;
90
91 ExRootTreeBranch *fBranchEvent;
92 ExRootTreeBranch *fBranchPhoton;
93 ExRootTreeBranch *fBranchElectron;
94 ExRootTreeBranch *fBranchMuon;
95 ExRootTreeBranch *fBranchTau;
96 ExRootTreeBranch *fBranchJet;
97 ExRootTreeBranch *fBranchMissingET;
98
99};
100
101//------------------------------------------------------------------------------
102
103LHCOlympicsConverter::LHCOlympicsConverter(const char *outputFileName) :
104 fIsFirstEvent(kTRUE), fIsNewEvent(kFALSE), fIsReadyToFill(kFALSE),
105 fPreviousObjectNumber(0), fTriggerWord(0), fEventNumber(1), fRecordNumber(1),
106 fOutputFile(0), fTreeWriter(0)
107{
108 fOutputFile = TFile::Open(outputFileName, "RECREATE");
109 fTreeWriter = new ExRootTreeWriter(fOutputFile, "LHCO");
110
111 // information about reconstructed event
112 fBranchEvent = fTreeWriter->NewBranch("Event", ExRootEvent::Class());
113 // reconstructed photons
114 fBranchPhoton = fTreeWriter->NewBranch("Photon", ExRootPhoton::Class());
115 // reconstructed electrons
116 fBranchElectron = fTreeWriter->NewBranch("Electron", ExRootElectron::Class());
117 // reconstructed muons
118 fBranchMuon = fTreeWriter->NewBranch("Muon", ExRootMuon::Class());
119 // reconstructed hadronically-decaying tau leptons
120 fBranchTau = fTreeWriter->NewBranch("Tau", ExRootTau::Class());
121 // reconstructed jets
122 fBranchJet = fTreeWriter->NewBranch("Jet", ExRootJet::Class());
123 // missing transverse energy
124 fBranchMissingET = fTreeWriter->NewBranch("MissingET", ExRootMissingET::Class());
125}
126
127//------------------------------------------------------------------------------
128
129LHCOlympicsConverter::~LHCOlympicsConverter()
130{
131 if(fTreeWriter) delete fTreeWriter;
132 if(fOutputFile) delete fOutputFile;
133}
134
135//------------------------------------------------------------------------------
136
137Long64_t LHCOlympicsConverter::GetNumberOfObjects(ifstream &inputFileStream)
138{
139 Long64_t counter = 0;
140 Bool_t canReadNumber, canReadFile = kTRUE;
141 Int_t number;
142 int position = inputFileStream.tellg();
143 inputFileStream.seekg(0, std::ios::beg);
144
145 inputFileStream.clear();
146
147 while(canReadFile)
148 {
149 do
150 {
151 getline(inputFileStream, fBuffer);
152
153 if(!inputFileStream.good())
154 {
155 canReadFile = kFALSE;
156 break;
157 }
158
159 fBufferStream.clear();
160 fBufferStream.str(fBuffer);
161
162 canReadNumber = (fBufferStream >> number);
163 }
164 while(!canReadNumber);
165
166 ++counter;
167 }
168
169 inputFileStream.clear();
170
171 inputFileStream.seekg(position, std::ios::beg);
172
173 return (counter - 1);
174}
175
176//------------------------------------------------------------------------------
177
178Bool_t LHCOlympicsConverter::ReadObject(ifstream &inputFileStream)
179{
180 Int_t i;
181 Bool_t canReadNumber;
182
183 do
184 {
185 getline(inputFileStream, fBuffer);
186
187 if(!inputFileStream.good()) return kFALSE;
188
189 fBufferStream.clear();
190 fBufferStream.str(fBuffer);
191
192 canReadNumber = kTRUE;
193
194 for(i = 0; canReadNumber && i < LHCOlympicsObject::maxIntParam; ++i)
195 {
196 canReadNumber = (fBufferStream >> fCurrentObject.intParam[i]);
197 }
198
199 if(canReadNumber && fCurrentObject.intParam[0] == 0)
200 {
201 fEventNumber = fCurrentObject.intParam[1];
202 canReadNumber = (fBufferStream >> fTriggerWord);
203 }
204 else
205 {
206 for(i = 0; canReadNumber && i < LHCOlympicsObject::maxDblParam; ++i)
207 {
208 canReadNumber = (fBufferStream >> fCurrentObject.dblParam[i]);
209 }
210 }
211 }
212 while(!canReadNumber);
213
214 return kTRUE;
215}
216
217//---------------------------------------------------------------------------
218
219void LHCOlympicsConverter::Write()
220{
221 if(fIsReadyToFill && fTreeWriter) fTreeWriter->Fill();
222 if(fTreeWriter) fTreeWriter->Write();
223 fIsReadyToFill = kFALSE;
224}
225
226//---------------------------------------------------------------------------
227// add empty events for missing event numbers
228
229void LHCOlympicsConverter::AddMissingEvents()
230{
231 while(fRecordNumber < fEventNumber)
232 {
233 fTreeWriter->Clear();
234 AnalyseEvent(fBranchEvent, fRecordNumber, 0);
235 fTreeWriter->Fill();
236
237 ++fRecordNumber;
238 }
239
240 fTreeWriter->Clear();
241
242 fIsReadyToFill = kFALSE;
243}
244
245//---------------------------------------------------------------------------
246
247void LHCOlympicsConverter::ProcessObject()
248{
249 fIsNewEvent = (fCurrentObject.intParam[0] <= fPreviousObjectNumber);
250
251 fPreviousObjectNumber = fCurrentObject.intParam[0];
252
253 if(fIsNewEvent && fIsFirstEvent && fTreeWriter)
254 {
255 fIsFirstEvent = kFALSE;
256
257 AddMissingEvents();
258 }
259
260 if(fIsNewEvent && fIsReadyToFill && fTreeWriter)
261 {
262 fIsReadyToFill = kFALSE;
263
264 fTreeWriter->Fill();
265 fTreeWriter->Clear();
266
267 ++fRecordNumber;
268
269 AddMissingEvents();
270 }
271
272 if(fCurrentObject.intParam[0] == 0)
273 {
274 AnalyseEvent(fBranchEvent, fEventNumber, fTriggerWord);
275 }
276 else
277 {
278 switch(fCurrentObject.intParam[1])
279 {
280 case 0: AnalysePhoton(fBranchPhoton); break;
281 case 1: AnalyseElectron(fBranchElectron); break;
282 case 2: AnalyseMuon(fBranchMuon); break;
283 case 3: AnalyseTau(fBranchTau); break;
284 case 4: AnalyseJet(fBranchJet); break;
285 case 6: AnalyseMissingET(fBranchMissingET); break;
286 }
287 }
288
289 fIsReadyToFill = kTRUE;
290}
291
292//---------------------------------------------------------------------------
293
294void LHCOlympicsConverter::AnalyseEvent(ExRootTreeBranch *branch,
295 Long64_t eventNumber, Int_t triggerWord)
296{
297 ExRootEvent *element;
298
299 element = static_cast<ExRootEvent*>(branch->NewEntry());
300
301 element->Number = eventNumber;
302 element->Trigger = triggerWord;
303}
304
305//---------------------------------------------------------------------------
306
307void LHCOlympicsConverter::AnalysePhoton(ExRootTreeBranch *branch)
308{
309 ExRootPhoton *element;
310
311 element = static_cast<ExRootPhoton*>(branch->NewEntry());
312
313 element->Eta = fCurrentObject.dblParam[0];
314 element->Phi = fCurrentObject.dblParam[1];
315 element->PT = fCurrentObject.dblParam[2];
316
317 element->EhadOverEem = fCurrentObject.dblParam[6];
318}
319
320//---------------------------------------------------------------------------
321
322void LHCOlympicsConverter::AnalyseElectron(ExRootTreeBranch *branch)
323{
324 ExRootElectron *element;
325
326 element = static_cast<ExRootElectron*>(branch->NewEntry());
327
328 element->Eta = fCurrentObject.dblParam[0];
329 element->Phi = fCurrentObject.dblParam[1];
330 element->PT = fCurrentObject.dblParam[2];
331
332 element->Ntrk = TMath::Abs(fCurrentObject.dblParam[4]);
333
334 element->Charge = fCurrentObject.dblParam[4] < 0.0 ? -1.0 : 1.0;
335
336 element->EhadOverEem = fCurrentObject.dblParam[6];
337}
338
339//---------------------------------------------------------------------------
340
341void LHCOlympicsConverter::AnalyseMuon(ExRootTreeBranch *branch)
342{
343 ExRootMuon *element;
344
345 element = static_cast<ExRootMuon*>(branch->NewEntry());
346
347 element->Eta = fCurrentObject.dblParam[0];
348 element->Phi = fCurrentObject.dblParam[1];
349 element->PT = fCurrentObject.dblParam[2];
350
351 element->Ntrk = TMath::Abs(fCurrentObject.dblParam[4]);
352
353 element->Charge = fCurrentObject.dblParam[4] < 0.0 ? -1.0 : 1.0;
354
355 element->JetIndex = Int_t(fCurrentObject.dblParam[5]);
356
357 element->PTiso = Int_t(fCurrentObject.dblParam[6]);
358 element->ETiso = fCurrentObject.dblParam[6] - element->PTiso;
359}
360
361//---------------------------------------------------------------------------
362
363void LHCOlympicsConverter::AnalyseTau(ExRootTreeBranch *branch)
364{
365 ExRootTau *element;
366
367 element = static_cast<ExRootTau*>(branch->NewEntry());
368
369 element->Eta = fCurrentObject.dblParam[0];
370 element->Phi = fCurrentObject.dblParam[1];
371 element->PT = fCurrentObject.dblParam[2];
372
373 element->Ntrk = TMath::Abs(fCurrentObject.dblParam[4]);
374
375 element->Charge = fCurrentObject.dblParam[4] < 0 ? -1.0 : 1.0;
376
377 element->EhadOverEem = fCurrentObject.dblParam[6];
378}
379
380//---------------------------------------------------------------------------
381
382void LHCOlympicsConverter::AnalyseJet(ExRootTreeBranch *branch)
383{
384 ExRootJet *element;
385
386 element = static_cast<ExRootJet*>(branch->NewEntry());
387
388 element->Eta = fCurrentObject.dblParam[0];
389 element->Phi = fCurrentObject.dblParam[1];
390 element->PT = fCurrentObject.dblParam[2];
391
392 element->Mass = fCurrentObject.dblParam[3];
393
394 element->Ntrk = TMath::Abs(fCurrentObject.dblParam[4]);
395
396 element->BTag = fCurrentObject.dblParam[5];
397
398 element->EhadOverEem = fCurrentObject.dblParam[6];
399
400 element->Index = fCurrentObject.intParam[0];
401
402}
403
404//---------------------------------------------------------------------------
405
406void LHCOlympicsConverter::AnalyseMissingET(ExRootTreeBranch *branch)
407{
408 ExRootMissingET *element;
409
410 element = static_cast<ExRootMissingET*>(branch->NewEntry());
411
412 element->Phi = fCurrentObject.dblParam[1];
413 element->MET = fCurrentObject.dblParam[2];
414}
415
416//---------------------------------------------------------------------------
417
418int main(int argc, char *argv[])
419{
420 char *appName = "ExRootLHEFConverter";
421
422 if(argc != 3)
423 {
424 cout << " Usage: " << appName << " input_file" << " output_file" << endl;
425 cout << " input_file - input file in LHEF format," << endl;
426 cout << " output_file - output file in ROOT format." << endl;
427 return 1;
428 }
429
430 gROOT->SetBatch();
431
432 int appargc = 1;
433 char *appargv[] = {appName};
434 TApplication app(appName, &appargc, appargv);
435
436 // Open a stream connected to an event file:
437 ifstream inputFileStream(argv[1]);
438
439 if(!inputFileStream.is_open())
440 {
441 cerr << "** ERROR: Can't open '" << argv[1] << "' for input" << endl;
442 return 1;
443 }
444
445 // Create LHC Olympics converter:
446 LHCOlympicsConverter *converter = new LHCOlympicsConverter(argv[2]);
447
448 cout << "** Calculating number of objects to process. Please wait..." << endl;
449 Long64_t allEntries = converter->GetNumberOfObjects(inputFileStream);
450 cout << "** Input file contains " << allEntries << " objects" << endl;
451
452 if(allEntries > 0)
453 {
454 ExRootProgressBar progressBar(allEntries);
455
456 // Loop over all objects
457 Long64_t entry = 0;
458 while(converter->ReadObject(inputFileStream))
459 {
460 converter->ProcessObject();
461
462 progressBar.Update(entry);
463
464 ++entry;
465 }
466 progressBar.Finish();
467
468 converter->Write();
469 }
470
471 cout << "** Exiting..." << endl;
472
473 delete converter;
474}
475
476
Note: See TracBrowser for help on using the repository browser.