Fork me on GitHub

source: git/converters/lhco2root.cpp@ a6db2d6

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

move branches/ModularDelphes to trunk

  • Property mode set to 100644
File size: 11.9 KB
RevLine 
[d7d2da3]1#include <stdexcept>
2#include <iostream>
3#include <sstream>
4
5#include <stdlib.h>
6#include <signal.h>
7#include <stdio.h>
8
9#include "TROOT.h"
10#include "TApplication.h"
11
12#include "TFile.h"
13#include "TObjArray.h"
14#include "TStopwatch.h"
15#include "TDatabasePDG.h"
16#include "TParticlePDG.h"
17#include "TLorentzVector.h"
18
19#include "modules/Delphes.h"
20#include "classes/DelphesStream.h"
21#include "classes/DelphesClasses.h"
22#include "classes/DelphesFactory.h"
23
24#include "ExRootAnalysis/ExRootTreeWriter.h"
25#include "ExRootAnalysis/ExRootTreeBranch.h"
26#include "ExRootAnalysis/ExRootProgressBar.h"
27
28using namespace std;
29
30static const int kBufferSize = 1024;
31
32/*
33LHC Olympics format discription from http://www.jthaler.net/olympicswiki/doku.php?id=lhc_olympics:data_file_format
34
35 * The first column of each row is just a counter that labels the object.
36 * 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).
37 * 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].
38 * The next three columns give the pseudorapidity, the azimuthal angle, and the transverse momentum of the object.
39 * The sixth column gives the invariant mass of the object.
40 * 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.
41 * 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.
42 * 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.
43*/
44
45//------------------------------------------------------------------------------
46
47class LHCOConverter
48{
49public:
50 LHCOConverter(TFile *outputFile);
51 ~LHCOConverter();
52
53 void Write();
54
55 Bool_t ReadLine(FILE *inputFile);
56
57private:
58
59 void AddMissingEvents();
60
61 void AnalyseEvent(ExRootTreeBranch *branch);
62
63 void AnalysePhoton(ExRootTreeBranch *branch);
64 void AnalyseElectron(ExRootTreeBranch *branch);
65 void AnalyseMuon(ExRootTreeBranch *branch);
66 void AnalyseTau(ExRootTreeBranch *branch);
67 void AnalyseJet(ExRootTreeBranch *branch);
68 void AnalyseMissingET(ExRootTreeBranch *branch);
69
70 enum {kIntParamSize = 2, kDblParamSize = 7};
71 Int_t fIntParam[kIntParamSize];
72 Double_t fDblParam[kDblParamSize];
73
74 Bool_t fIsReadyToFill;
75
76 Int_t fTriggerWord, fEventNumber;
77
78 char *fBuffer;
79
80 ExRootTreeWriter *fTreeWriter;
81
82 ExRootTreeBranch *fBranchEvent;
83 ExRootTreeBranch *fBranchTrack;
84 ExRootTreeBranch *fBranchTower;
85 ExRootTreeBranch *fBranchPhoton;
86 ExRootTreeBranch *fBranchElectron;
87 ExRootTreeBranch *fBranchMuon;
88 ExRootTreeBranch *fBranchJet;
89 ExRootTreeBranch *fBranchMissingET;
90
91};
92
93//------------------------------------------------------------------------------
94
95LHCOConverter::LHCOConverter(TFile *outputFile) :
96 fIsReadyToFill(kFALSE),
97 fTriggerWord(0), fEventNumber(1),
98 fBuffer(0), fTreeWriter(0)
99{
100 fBuffer = new char[kBufferSize];
101 fTreeWriter = new ExRootTreeWriter(outputFile, "Delphes");
102
103 // information about reconstructed event
104 fBranchEvent = fTreeWriter->NewBranch("Event", LHCOEvent::Class());
105 // reconstructed tracks
106 fBranchTrack = fTreeWriter->NewBranch("Track", Track::Class());
107 /// calorimeter towers
108 fBranchTower = fTreeWriter->NewBranch("Tower", Tower::Class());
109 // reconstructed photons
110 fBranchPhoton = fTreeWriter->NewBranch("Photon", Photon::Class());
111 // reconstructed electrons
112 fBranchElectron = fTreeWriter->NewBranch("Electron", Electron::Class());
113 // reconstructed muons
114 fBranchMuon = fTreeWriter->NewBranch("Muon", Muon::Class());
115 // reconstructed jets
116 fBranchJet = fTreeWriter->NewBranch("Jet", Jet::Class());
117 // missing transverse energy
118 fBranchMissingET = fTreeWriter->NewBranch("MissingET", MissingET::Class());
119}
120
121//------------------------------------------------------------------------------
122
123LHCOConverter::~LHCOConverter()
124{
125 if(fTreeWriter) delete fTreeWriter;
126 if(fBuffer) delete[] fBuffer;
127}
128
129//------------------------------------------------------------------------------
130
131Bool_t LHCOConverter::ReadLine(FILE *inputFile)
132{
133 int rc;
134
135 if(!fgets(fBuffer, kBufferSize, inputFile)) return kFALSE;
136
137 DelphesStream bufferStream(fBuffer);
138
139 rc = bufferStream.ReadInt(fIntParam[0]);
140
141 if(!rc)
142 {
143 return kTRUE;
144 }
145
146 if(fIntParam[0] == 0)
147 {
148 rc = bufferStream.ReadInt(fEventNumber)
149 && bufferStream.ReadInt(fTriggerWord);
150
151 if(!rc)
152 {
153 cerr << "** ERROR: " << "invalid event format" << endl;
154 return kFALSE;
155 }
156
157 if(fIsReadyToFill && fTreeWriter)
158 {
159 fTreeWriter->Fill();
160 fTreeWriter->Clear();
161 }
162
163 AnalyseEvent(fBranchEvent);
164 fIsReadyToFill = kTRUE;
165 }
166 else
167 {
168 rc = bufferStream.ReadInt(fIntParam[1])
169 && bufferStream.ReadDbl(fDblParam[0])
170 && bufferStream.ReadDbl(fDblParam[1])
171 && bufferStream.ReadDbl(fDblParam[2])
172 && bufferStream.ReadDbl(fDblParam[3])
173 && bufferStream.ReadDbl(fDblParam[4])
174 && bufferStream.ReadDbl(fDblParam[5])
175 && bufferStream.ReadDbl(fDblParam[6]);
176
177 if(!rc)
178 {
179 cerr << "** ERROR: " << "invalid object format" << endl;
180 return kFALSE;
181 }
182
183 switch(fIntParam[1])
184 {
185 case 0: AnalysePhoton(fBranchPhoton); break;
186 case 1: AnalyseElectron(fBranchElectron); break;
187 case 2: AnalyseMuon(fBranchMuon); break;
188 case 3: AnalyseTau(fBranchJet); break;
189 case 4: AnalyseJet(fBranchJet); break;
190 case 6: AnalyseMissingET(fBranchMissingET); break;
191 }
192 }
193
194 return kTRUE;
195}
196
197//---------------------------------------------------------------------------
198
199void LHCOConverter::Write()
200{
201 if(fIsReadyToFill && fTreeWriter) fTreeWriter->Fill();
202 if(fTreeWriter) fTreeWriter->Write();
203 fIsReadyToFill = kFALSE;
204}
205
206//---------------------------------------------------------------------------
207
208void LHCOConverter::AnalyseEvent(ExRootTreeBranch *branch)
209{
210 LHCOEvent *element;
211
212 element = static_cast<LHCOEvent*>(branch->NewEntry());
213
214 element->Number = fEventNumber;
215 element->Trigger = fTriggerWord;
216}
217
218//---------------------------------------------------------------------------
219
220void LHCOConverter::AnalysePhoton(ExRootTreeBranch *branch)
221{
222 Photon *element;
223
224 element = static_cast<Photon*>(branch->NewEntry());
225
226 element->Eta = fDblParam[0];
227 element->Phi = fDblParam[1];
228 element->PT = fDblParam[2];
229 element->EhadOverEem = fDblParam[6];
230}
231
232//---------------------------------------------------------------------------
233
234void LHCOConverter::AnalyseElectron(ExRootTreeBranch *branch)
235{
236 Electron *element;
237
238 element = static_cast<Electron*>(branch->NewEntry());
239
240 element->Eta = fDblParam[0];
241 element->Phi = fDblParam[1];
242 element->PT = fDblParam[2];
243
244 element->Charge = fDblParam[4] < 0.0 ? -1 : 1;
245/*
246 element->Ntrk = TMath::Abs(fDblParam[4]);
247*/
248 element->EhadOverEem = fDblParam[6];
249}
250
251//---------------------------------------------------------------------------
252
253void LHCOConverter::AnalyseMuon(ExRootTreeBranch *branch)
254{
255 Muon *element;
256
257 element = static_cast<Muon*>(branch->NewEntry());
258
259 element->Eta = fDblParam[0];
260 element->Phi = fDblParam[1];
261 element->PT = fDblParam[2];
262
263 element->Charge = fDblParam[4] < 0.0 ? -1 : 1;
264/*
265 element->Ntrk = TMath::Abs(fDblParam[4]);
266
267 element->JetIndex = Int_t(fDblParam[5]);
268
269 element->PTiso = Int_t(fDblParam[6]);
270 element->ETiso = fDblParam[6] - element->PTiso;
271*/
272}
273
274//---------------------------------------------------------------------------
275
276void LHCOConverter::AnalyseTau(ExRootTreeBranch *branch)
277{
278 Jet *element;
279
280 element = static_cast<Jet*>(branch->NewEntry());
281
282 element->Eta = fDblParam[0];
283 element->Phi = fDblParam[1];
284 element->PT = fDblParam[2];
285
286 element->Mass = fDblParam[3];
287
288 element->BTag = 0;
289 element->TauTag = 1;
290
291 element->Charge = fDblParam[4] < 0 ? -1 : 1;
292/*
293 element->Ntrk = TMath::Abs(fDblParam[4]);
294*/
295 element->EhadOverEem = fDblParam[6];
296}
297
298//---------------------------------------------------------------------------
299
300void LHCOConverter::AnalyseJet(ExRootTreeBranch *branch)
301{
302 Jet *element;
303
304 element = static_cast<Jet*>(branch->NewEntry());
305
306 element->Eta = fDblParam[0];
307 element->Phi = fDblParam[1];
308 element->PT = fDblParam[2];
309
310 element->Mass = fDblParam[3];
311/*
312 element->Ntrk = TMath::Abs(Int_t(fDblParam[4]));
313*/
314 element->BTag = Int_t(fDblParam[5]);
315 element->TauTag = 0;
316
317 element->Charge = 0;
318
319 element->EhadOverEem = fDblParam[6];
320/*
321 element->Index = fIntParam[0];
322*/
323}
324
325//---------------------------------------------------------------------------
326
327void LHCOConverter::AnalyseMissingET(ExRootTreeBranch *branch)
328{
329 MissingET *element;
330
331 element = static_cast<MissingET*>(branch->NewEntry());
332
333 element->Phi = fDblParam[1];
334 element->MET = fDblParam[2];
335}
336
337//---------------------------------------------------------------------------
338
339static bool interrupted = false;
340
341void SignalHandler(int sig)
342{
343 interrupted = true;
344}
345
346//---------------------------------------------------------------------------
347
348int main(int argc, char *argv[])
349{
350 char appName[] = "lhco2root";
351 stringstream message;
352 FILE *inputFile = 0;
353 TFile *outputFile = 0;
354 LHCOConverter *converter = 0;
355 Int_t i;
356 Long64_t length, eventCounter;
357
358 if(argc < 2)
359 {
360 cout << " Usage: " << appName << " output_file" << " [input_file(s)]" << endl;
361 cout << " output_file - output file in ROOT format," << endl;
362 cout << " input_file(s) - input file(s) in LHCO format," << endl;
363 cout << " with no input_file, or when input_file is -, read standard input." << endl;
364 return 1;
365 }
366
367 signal(SIGINT, SignalHandler);
368
369 gROOT->SetBatch();
370
371 int appargc = 1;
372 char *appargv[] = {appName};
373 TApplication app(appName, &appargc, appargv);
374
375 try
376 {
377 outputFile = TFile::Open(argv[1], "CREATE");
378
379 if(outputFile == NULL)
380 {
381 message << "can't open " << argv[1];
382 throw runtime_error(message.str());
383 }
384
385 converter = new LHCOConverter(outputFile);
386
387 i = 2;
388 do
389 {
390 if(interrupted) break;
391
392 if(i == argc || strcmp(argv[i], "-") == 0)
393 {
394 cout << "** Reading standard input" << endl;
395 inputFile = stdin;
396 length = -1;
397 }
398 else
399 {
400 cout << "** Reading " << argv[i] << endl;
401 inputFile = fopen(argv[i], "r");
402
403 if(inputFile == NULL)
404 {
405 message << "can't open " << argv[i];
406 throw runtime_error(message.str());
407 }
408
409 fseek(inputFile, 0L, SEEK_END);
410 length = ftello(inputFile);
411 fseek(inputFile, 0L, SEEK_SET);
412
413 if(length <= 0)
414 {
415 fclose(inputFile);
416 ++i;
417 continue;
418 }
419 }
420
421 eventCounter = 0;
422 ExRootProgressBar progressBar(length);
423
424 // Loop over all objects
425 while(converter->ReadLine(inputFile) && !interrupted)
426 {
427 ++eventCounter;
428
429 progressBar.Update(ftello(inputFile), eventCounter);
430 }
431 converter->Write();
432
433 fseek(inputFile, 0L, SEEK_END);
434 progressBar.Update(ftello(inputFile), eventCounter, kTRUE);
435 progressBar.Finish();
436
437 if(inputFile != stdin) fclose(inputFile);
438
439 ++i;
440 }
441 while(i < argc);
442
443 cout << "** Exiting..." << endl;
444
445 delete converter;
446 delete outputFile;
447
448 return 0;
449 }
450 catch(runtime_error &e)
451 {
452 if(converter) delete converter;
453 if(outputFile) delete outputFile;
454 cerr << "** ERROR: " << e.what() << endl;
455 return 1;
456 }
457}
458
459
Note: See TracBrowser for help on using the repository browser.