Fork me on GitHub

source: git/readers/DelphesLHEF.cpp@ 84a1f7d

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

move branches/ModularDelphes to trunk

  • Property mode set to 100644
File size: 5.0 KB
Line 
1#include <stdexcept>
2#include <iostream>
3#include <sstream>
4
5#include <signal.h>
6
7#include "TROOT.h"
8#include "TApplication.h"
9
10#include "TFile.h"
11#include "TObjArray.h"
12#include "TStopwatch.h"
13#include "TDatabasePDG.h"
14#include "TParticlePDG.h"
15#include "TLorentzVector.h"
16
17#include "modules/Delphes.h"
18#include "classes/DelphesClasses.h"
19#include "classes/DelphesFactory.h"
20#include "classes/DelphesLHEFReader.h"
21
22#include "ExRootAnalysis/ExRootTreeWriter.h"
23#include "ExRootAnalysis/ExRootTreeBranch.h"
24#include "ExRootAnalysis/ExRootProgressBar.h"
25
26using namespace std;
27
28//---------------------------------------------------------------------------
29
30static bool interrupted = false;
31
32void SignalHandler(int sig)
33{
34 interrupted = true;
35}
36
37//---------------------------------------------------------------------------
38
39int main(int argc, char *argv[])
40{
41 char appName[] = "DelphesLHEF";
42 stringstream message;
43 FILE *inputFile = 0;
44 TFile *outputFile = 0;
45 TStopwatch readStopWatch, procStopWatch;
46 ExRootTreeWriter *treeWriter = 0;
47 ExRootTreeBranch *branchEvent = 0;
48 ExRootConfReader *confReader = 0;
49 Delphes *modularDelphes = 0;
50 DelphesFactory *factory = 0;
51 TObjArray *stableParticleOutputArray = 0, *allParticleOutputArray = 0, *partonOutputArray = 0;
52 DelphesLHEFReader *reader = 0;
53 Int_t i;
54 Long64_t length, eventCounter;
55
56 if(argc < 3)
57 {
58 cout << " Usage: " << appName << " config_file" << " output_file" << " [input_file(s)]" << endl;
59 cout << " config_file - configuration file in Tcl format," << endl;
60 cout << " output_file - output file in ROOT format," << endl;
61 cout << " input_file(s) - input file(s) in LHEF format," << endl;
62 cout << " with no input_file, or when input_file is -, read standard input." << endl;
63 return 1;
64 }
65
66 signal(SIGINT, SignalHandler);
67
68 gROOT->SetBatch();
69
70 int appargc = 1;
71 char *appargv[] = {appName};
72 TApplication app(appName, &appargc, appargv);
73
74 try
75 {
76 outputFile = TFile::Open(argv[2], "CREATE");
77
78 if(outputFile == NULL)
79 {
80 message << "can't create output file " << argv[2];
81 throw runtime_error(message.str());
82 }
83
84 treeWriter = new ExRootTreeWriter(outputFile, "Delphes");
85
86 branchEvent = treeWriter->NewBranch("Event", LHEFEvent::Class());
87
88 confReader = new ExRootConfReader;
89 confReader->ReadFile(argv[1]);
90
91 modularDelphes = new Delphes("Delphes");
92 modularDelphes->SetConfReader(confReader);
93 modularDelphes->SetTreeWriter(treeWriter);
94
95 factory = modularDelphes->GetFactory();
96 allParticleOutputArray = modularDelphes->ExportArray("allParticles");
97 stableParticleOutputArray = modularDelphes->ExportArray("stableParticles");
98 partonOutputArray = modularDelphes->ExportArray("partons");
99
100 reader = new DelphesLHEFReader;
101
102 modularDelphes->InitTask();
103
104 i = 3;
105 do
106 {
107 if(interrupted) break;
108
109 if(i == argc || strncmp(argv[i], "-", 2) == 0)
110 {
111 cout << "** Reading standard input" << endl;
112 inputFile = stdin;
113 length = -1;
114 }
115 else
116 {
117 cout << "** Reading " << argv[i] << endl;
118 inputFile = fopen(argv[i], "r");
119
120 if(inputFile == NULL)
121 {
122 message << "can't open " << argv[i];
123 throw runtime_error(message.str());
124 }
125
126 fseek(inputFile, 0L, SEEK_END);
127 length = ftello(inputFile);
128 fseek(inputFile, 0L, SEEK_SET);
129
130 if(length <= 0)
131 {
132 fclose(inputFile);
133 ++i;
134 continue;
135 }
136 }
137
138 reader->SetInputFile(inputFile);
139
140 ExRootProgressBar progressBar(length);
141
142 // Loop over all objects
143 eventCounter = 0;
144 treeWriter->Clear();
145 modularDelphes->Clear();
146 reader->Clear();
147 readStopWatch.Start();
148 while(reader->ReadBlock(factory, allParticleOutputArray,
149 stableParticleOutputArray, partonOutputArray) && !interrupted)
150 {
151 if(reader->EventReady())
152 {
153 ++eventCounter;
154
155 readStopWatch.Stop();
156 procStopWatch.Start();
157 modularDelphes->ProcessTask();
158 procStopWatch.Stop();
159
160 reader->AnalyzeEvent(branchEvent, eventCounter, &readStopWatch, &procStopWatch);
161
162 treeWriter->Fill();
163
164 treeWriter->Clear();
165 modularDelphes->Clear();
166 reader->Clear();
167
168 readStopWatch.Start();
169 }
170 progressBar.Update(ftello(inputFile), eventCounter);
171 }
172
173 fseek(inputFile, 0L, SEEK_END);
174 progressBar.Update(ftello(inputFile), eventCounter, kTRUE);
175 progressBar.Finish();
176
177 if(inputFile != stdin) fclose(inputFile);
178
179 ++i;
180 }
181 while(i < argc);
182
183 modularDelphes->FinishTask();
184 treeWriter->Write();
185
186 cout << "** Exiting..." << endl;
187
188 delete reader;
189 delete modularDelphes;
190 delete confReader;
191 delete treeWriter;
192 delete outputFile;
193
194 return 0;
195 }
196 catch(runtime_error &e)
197 {
198 if(treeWriter) delete treeWriter;
199 if(outputFile) delete outputFile;
200 cerr << "** ERROR: " << e.what() << endl;
201 return 1;
202 }
203}
Note: See TracBrowser for help on using the repository browser.