Fork me on GitHub

source: git/converters/root2lhco.cpp@ d7d2da3

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

move branches/ModularDelphes to trunk

  • Property mode set to 100644
File size: 12.0 KB
RevLine 
[d7d2da3]1#include <stdexcept>
2#include <iostream>
3#include <fstream>
4#include <sstream>
5#include <string>
6
7#include <stdlib.h>
8#include <signal.h>
9#include <stdio.h>
10
11#include "TROOT.h"
12#include "TApplication.h"
13
14#include "TFile.h"
15#include "TClonesArray.h"
16
17#include "classes/DelphesClasses.h"
18
19#include "ExRootAnalysis/ExRootTreeReader.h"
20#include "ExRootAnalysis/ExRootProgressBar.h"
21
22using namespace std;
23
24/*
25LHC Olympics format discription from http://www.jthaler.net/olympicswiki/doku.php?id=lhc_olympics:data_file_format
26
27 * The first column of each row is just a counter that labels the object.
28 * 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).
29 * 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].
30 * The next three columns give the pseudorapidity, the azimuthal angle, and the transverse momentum of the object.
31 * The sixth column gives the invariant mass of the object.
32 * 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.
33 * 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.
34 * 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.
35*/
36
37class LHCOWriter
38{
39public:
40 LHCOWriter(ExRootTreeReader *treeReader, FILE *outputFile);
41 ~LHCOWriter();
42
43 void ProcessEvent();
44
45private:
46
47 void Reset();
48 void Write();
49
50 void AnalyseEvent();
51
52 void AnalysePhotons();
53 void AnalyseElectrons();
54 void AnalyseMuons();
55 void AnalyseTauJets();
56 void AnalyseJets();
57
58 void AnalyseMissingET();
59
60 enum {kIntParamSize = 2, kDblParamSize = 9};
61 Int_t fIntParam[kIntParamSize];
62 Double_t fDblParam[kDblParamSize];
63
64 Long64_t fTriggerWord, fEventNumber;
65
66 ExRootTreeReader *fTreeReader;
67 FILE *fOutputFile;
68
69 TClonesArray *fBranchEvent;
70
71 TClonesArray *fBranchTrack;
72 TClonesArray *fBranchTower;
73
74 TClonesArray *fBranchPhoton;
75 TClonesArray *fBranchElectron;
76 TClonesArray *fBranchMuon;
77 TClonesArray *fBranchJet;
78 TClonesArray *fBranchMissingET;
79
80 TIterator *fItTrack;
81 TIterator *fItTower;
82
83 TIterator *fItPhoton;
84 TIterator *fItElectron;
85 TIterator *fItMuon;
86 TIterator *fItJet;
87};
88
89//------------------------------------------------------------------------------
90
91LHCOWriter::LHCOWriter(ExRootTreeReader *treeReader, FILE *outputFile) :
92 fTriggerWord(0), fEventNumber(1), fTreeReader(0), fOutputFile(0)
93{
94 fTreeReader = treeReader;
95 fOutputFile = outputFile;
96
97 // information about reconstructed event
98 fBranchEvent = fTreeReader->UseBranch("Event");
99 // reconstructed tracks
100 fBranchTrack = fTreeReader->UseBranch("Track");
101 fItTrack = fBranchTrack->MakeIterator();
102 // calorimeter towers
103 fBranchTower = fTreeReader->UseBranch("Tower");
104 fItTower = fBranchTower->MakeIterator();
105 // reconstructed photons
106 fBranchPhoton = fTreeReader->UseBranch("Photon");
107 fItPhoton = fBranchPhoton->MakeIterator();
108 // reconstructed electrons
109 fBranchElectron = fTreeReader->UseBranch("Electron");
110 fItElectron = fBranchElectron->MakeIterator();
111 // reconstructed muons
112 fBranchMuon = fTreeReader->UseBranch("Muon");
113 fItMuon = fBranchMuon->MakeIterator();
114 // reconstructed jets
115 fBranchJet = fTreeReader->UseBranch("Jet");
116 fItJet = fBranchJet->MakeIterator();
117 // missing transverse energy
118 fBranchMissingET = fTreeReader->UseBranch("MissingET");
119}
120
121//------------------------------------------------------------------------------
122
123LHCOWriter::~LHCOWriter()
124{
125}
126
127//---------------------------------------------------------------------------
128
129void LHCOWriter::ProcessEvent()
130{
131 fIntParam[0] = 0;
132
133 AnalyseEvent();
134
135 AnalysePhotons();
136 AnalyseElectrons();
137 AnalyseMuons();
138 AnalyseTauJets();
139 AnalyseJets();
140
141 AnalyseMissingET();
142}
143
144//---------------------------------------------------------------------------
145
146void LHCOWriter::Reset()
147{
148 int i;
149 for(i = 1; i < kIntParamSize; ++i)
150 {
151 fIntParam[i] = 0;
152 }
153
154 for(i = 0; i < kDblParamSize; ++i)
155 {
156 fDblParam[i] = 0.0;
157 }
158}
159
160//---------------------------------------------------------------------------
161
162void LHCOWriter::Write()
163{
164 fprintf(fOutputFile, "%4d %4d %8.3f %8.3f %7.2f %7.2f %6.1f %6.1f %7.2f %6.1f %6.1f\n",
165 fIntParam[0], fIntParam[1], fDblParam[0], fDblParam[1], fDblParam[2],
166 fDblParam[3], fDblParam[4], fDblParam[5], fDblParam[6], fDblParam[7], fDblParam[8]);
167
168 ++fIntParam[0];
169}
170
171//---------------------------------------------------------------------------
172
173void LHCOWriter::AnalyseEvent()
174{
175 Event *element;
176
177 element = static_cast<Event*>(fBranchEvent->At(0));
178
179 fprintf(fOutputFile, "%4d %13lld %8d\n", 0, element->Number, 0);
180
181 ++fIntParam[0];
182}
183
184//---------------------------------------------------------------------------
185
186void LHCOWriter::AnalysePhotons()
187{
188 Photon *element;
189
190 fItPhoton->Reset();
191 while((element = static_cast<Photon*>(fItPhoton->Next())))
192 {
193 Reset();
194
195 fIntParam[1] = 0;
196
197 fDblParam[0] = element->Eta;
198 fDblParam[1] = element->Phi;
199 fDblParam[2] = element->PT;
200
201 fDblParam[6] = element->EhadOverEem;
202
203 Write();
204 }
205}
206
207//---------------------------------------------------------------------------
208
209void LHCOWriter::AnalyseElectrons()
210{
211 Electron *element;
212
213 fItElectron->Reset();
214 while((element = static_cast<Electron*>(fItElectron->Next())))
215 {
216 Reset();
217
218 fIntParam[1] = 1;
219
220 fDblParam[0] = element->Eta;
221 fDblParam[1] = element->Phi;
222 fDblParam[2] = element->PT;
223
224 fDblParam[4] = element->Charge;
225
226 fDblParam[6] = element->EhadOverEem;
227
228 Write();
229 }
230}
231
232//---------------------------------------------------------------------------
233
234void LHCOWriter::AnalyseMuons()
235{
236 Muon *element;
237 Track *track;
238 Tower *tower;
239 Jet *jet;
240 Int_t muonCounter, tauCounter, jetCounter, minIndex;
241 Float_t sumPT, sumET, ratET, jetDR, minDR;
242
243 muonCounter = 0;
244 fItMuon->Reset();
245 while((element = static_cast<Muon*>(fItMuon->Next())))
246 {
247 Reset();
248
249 sumPT = 0.0;
250 fItTrack->Reset();
251 while((track = static_cast<Track*>(fItTrack->Next())))
252 {
253 if(element->P4().DeltaR(track->P4()) < 0.5) sumPT += track->PT;
254 }
255
256 sumET = 0.0;
257 fItTower->Reset();
258 while((tower = static_cast<Tower*>(fItTower->Next())))
259 {
260 if(element->P4().DeltaR(tower->P4()) < 0.5) sumET += tower->ET;
261 }
262
263 tauCounter = 0;
264 jetCounter = 0;
265 minIndex = -1;
266 minDR = 1.0E9;
267 fItJet->Reset();
268 while((jet = static_cast<Jet*>(fItJet->Next())))
269 {
270 if(jet->TauTag != 0)
271 {
272 ++tauCounter;
273 continue;
274 }
275
276 jetDR = element->P4().DeltaR(jet->P4());
277 if(jetDR < minDR)
278 {
279 minIndex = jetCounter;
280 minDR = jetDR;
281 }
282 ++jetCounter;
283 }
284
285 fIntParam[1] = 2;
286
287 fDblParam[0] = element->Eta;
288 fDblParam[1] = element->Phi;
289 fDblParam[2] = element->PT;
290
291 fDblParam[3] = 0.11;
292
293 fDblParam[4] = element->Charge;
294
295 if(minIndex >= 0)
296 {
297 fDblParam[5] = fIntParam[0] + fBranchMuon->GetEntriesFast() - muonCounter + tauCounter + minIndex;
298 }
299
300 ratET = sumET/element->PT;
301 fDblParam[6] = Float_t(TMath::Nint(sumPT)) + (ratET < 1.0 ? ratET : 0.99);
302
303 Write();
304 ++muonCounter;
305 }
306}
307
308//---------------------------------------------------------------------------
309
310void LHCOWriter::AnalyseTauJets()
311{
312 Jet *element;
313 Track *track;
314 Int_t counter;
315
316 fItJet->Reset();
317 while((element = static_cast<Jet*>(fItJet->Next())))
318 {
319 if(element->TauTag == 0) continue;
320
321 Reset();
322
323 counter = 0;
324 fItTrack->Reset();
325 while((track = static_cast<Track*>(fItTrack->Next())))
326 {
327 if(element->P4().DeltaR(track->P4()) < 0.5) ++counter;
328 }
329
330 fIntParam[1] = 3;
331
332 fDblParam[0] = element->Eta;
333 fDblParam[1] = element->Phi;
334 fDblParam[2] = element->PT;
335 fDblParam[3] = element->Mass;
336 fDblParam[4] = counter * element->Charge;
337
338 fDblParam[6] = element->EhadOverEem;
339
340 Write();
341 }
342}
343
344//---------------------------------------------------------------------------
345
346void LHCOWriter::AnalyseJets()
347{
348 Jet *element;
349 Track *track;
350 Int_t counter;
351
352 fItJet->Reset();
353 while((element = static_cast<Jet*>(fItJet->Next())))
354 {
355 if(element->TauTag != 0) continue;
356
357 Reset();
358
359 counter = 0;
360 fItTrack->Reset();
361 while((track = static_cast<Track*>(fItTrack->Next())))
362 {
363 if(element->P4().DeltaR(track->P4()) < 0.5) ++counter;
364 }
365
366 fIntParam[1] = 4;
367
368 fDblParam[0] = element->Eta;
369 fDblParam[1] = element->Phi;
370 fDblParam[2] = element->PT;
371 fDblParam[3] = element->Mass;
372 fDblParam[4] = counter;
373 fDblParam[5] = element->BTag;
374 fDblParam[6] = element->EhadOverEem;
375
376 Write();
377 }
378}
379
380//---------------------------------------------------------------------------
381
382void LHCOWriter::AnalyseMissingET()
383{
384 MissingET *element;
385
386 element = static_cast<MissingET*>(fBranchMissingET->At(0));
387
388 Reset();
389
390 fIntParam[1] = 6;
391
392 fDblParam[1] = element->Phi;
393 fDblParam[2] = element->MET;
394
395 Write();
396}
397
398//---------------------------------------------------------------------------
399
400static bool interrupted = false;
401
402void SignalHandler(int sig)
403{
404 interrupted = true;
405}
406
407//---------------------------------------------------------------------------
408
409int main(int argc, char *argv[])
410{
411 char appName[] = "root2lhco";
412 stringstream message;
413 FILE *outputFile = 0;
414 TChain *inputChain = 0;
415 LHCOWriter *writer = 0;
416 ExRootTreeReader *treeReader = 0;
417 Long64_t entry, allEntries;
418
419 if(argc < 2 || argc > 3)
420 {
421 cout << " Usage: " << appName << " input_file" << " [output_file]" << endl;
422 cout << " input_file - input file in ROOT format," << endl;
423 cout << " output_file - output file in LHCO format," << endl;
424 cout << " with no output_file, or when output_file is -, write to standard output." << endl;
425 return 1;
426 }
427
428 signal(SIGINT, SignalHandler);
429
430 gROOT->SetBatch();
431
432 int appargc = 1;
433 char *appargv[] = {appName};
434 TApplication app(appName, &appargc, appargv);
435
436 try
437 {
438 cout << "** Reading " << argv[1] << endl;
439 inputChain = new TChain("Delphes");
440 inputChain->Add(argv[1]);
441
442 ExRootTreeReader *treeReader = new ExRootTreeReader(inputChain);
443
444 if(argc == 2 || strcmp(argv[2], "-") == 0)
445 {
446 outputFile = stdout;
447 }
448 else
449 {
450 outputFile = fopen(argv[2], "w");
451
452 if(outputFile == NULL)
453 {
454 message << "can't open " << argv[2];
455 throw runtime_error(message.str());
456 }
457 }
458
459 allEntries = treeReader->GetEntries();
460 cout << "** Input file contains " << allEntries << " events" << endl;
461
462 if(allEntries > 0)
463 {
464 // Create LHC Olympics converter:
465 writer = new LHCOWriter(treeReader, outputFile);
466
467 ExRootProgressBar progressBar(allEntries - 1);
468 // Loop over all events
469 for(entry = 0; entry < allEntries && !interrupted; ++entry)
470 {
471 if(!treeReader->ReadEntry(entry))
472 {
473 cerr << "** ERROR: cannot read event " << entry << endl;
474 break;
475 }
476
477 writer->ProcessEvent();
478
479 progressBar.Update(entry);
480 }
481 progressBar.Finish();
482
483 delete writer;
484 }
485
486 cout << "** Exiting..." << endl;
487
488 if(outputFile != stdout) fclose(outputFile);
489 delete treeReader;
490 delete inputChain;
491 return 0;
492 }
493 catch(runtime_error &e)
494 {
495 if(writer) delete writer;
496 if(treeReader) delete treeReader;
497 if(inputChain) delete inputChain;
498 cerr << "** ERROR: " << e.what() << endl;
499 return 1;
500 }
501}
502
503
Note: See TracBrowser for help on using the repository browser.