Fork me on GitHub

source: git/converters/lhco2root.cpp@ 8dd735c

Last change on this file since 8dd735c was 341014c, checked in by Pavel Demin <pavel-demin@…>, 5 years ago

apply .clang-format to all .h, .cc and .cpp files

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