Fork me on GitHub

source: svn/trunk/readers/DelphesLHEF.cpp@ 1157

Last change on this file since 1157 was 1130, checked in by Pavel Demin, 11 years ago

add MaxEvents and SkipEvents

File size: 5.5 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, maxEvents, skipEvents;
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 maxEvents = confReader->GetInt("::MaxEvents", 0);
92 skipEvents = confReader->GetInt("::SkipEvents", 0);
93
94 if(maxEvents < 0)
95 {
96 throw runtime_error("MaxEvents must be zero or positive");
97 }
98
99 if(skipEvents < 0)
100 {
101 throw runtime_error("SkipEvents must be zero or positive");
102 }
103
104 modularDelphes = new Delphes("Delphes");
105 modularDelphes->SetConfReader(confReader);
106 modularDelphes->SetTreeWriter(treeWriter);
107
108 factory = modularDelphes->GetFactory();
109 allParticleOutputArray = modularDelphes->ExportArray("allParticles");
110 stableParticleOutputArray = modularDelphes->ExportArray("stableParticles");
111 partonOutputArray = modularDelphes->ExportArray("partons");
112
113 reader = new DelphesLHEFReader;
114
115 modularDelphes->InitTask();
116
117 i = 3;
118 do
119 {
120 if(interrupted) break;
121
122 if(i == argc || strncmp(argv[i], "-", 2) == 0)
123 {
124 cout << "** Reading standard input" << endl;
125 inputFile = stdin;
126 length = -1;
127 }
128 else
129 {
130 cout << "** Reading " << argv[i] << endl;
131 inputFile = fopen(argv[i], "r");
132
133 if(inputFile == NULL)
134 {
135 message << "can't open " << argv[i];
136 throw runtime_error(message.str());
137 }
138
139 fseek(inputFile, 0L, SEEK_END);
140 length = ftello(inputFile);
141 fseek(inputFile, 0L, SEEK_SET);
142
143 if(length <= 0)
144 {
145 fclose(inputFile);
146 ++i;
147 continue;
148 }
149 }
150
151 reader->SetInputFile(inputFile);
152
153 ExRootProgressBar progressBar(length);
154
155 // Loop over all objects
156 eventCounter = 0;
157 treeWriter->Clear();
158 modularDelphes->Clear();
159 reader->Clear();
160 readStopWatch.Start();
161 while((maxEvents <= 0 || eventCounter - skipEvents < maxEvents) &&
162 reader->ReadBlock(factory, allParticleOutputArray,
163 stableParticleOutputArray, partonOutputArray) && !interrupted)
164 {
165 if(reader->EventReady())
166 {
167 ++eventCounter;
168
169 readStopWatch.Stop();
170
171 if(eventCounter > skipEvents)
172 {
173 readStopWatch.Stop();
174 procStopWatch.Start();
175 modularDelphes->ProcessTask();
176 procStopWatch.Stop();
177
178 reader->AnalyzeEvent(branchEvent, eventCounter, &readStopWatch, &procStopWatch);
179
180 treeWriter->Fill();
181
182 treeWriter->Clear();
183 }
184
185 modularDelphes->Clear();
186 reader->Clear();
187
188 readStopWatch.Start();
189 }
190 progressBar.Update(ftello(inputFile), eventCounter);
191 }
192
193 fseek(inputFile, 0L, SEEK_END);
194 progressBar.Update(ftello(inputFile), eventCounter, kTRUE);
195 progressBar.Finish();
196
197 if(inputFile != stdin) fclose(inputFile);
198
199 ++i;
200 }
201 while(i < argc);
202
203 modularDelphes->FinishTask();
204 treeWriter->Write();
205
206 cout << "** Exiting..." << endl;
207
208 delete reader;
209 delete modularDelphes;
210 delete confReader;
211 delete treeWriter;
212 delete outputFile;
213
214 return 0;
215 }
216 catch(runtime_error &e)
217 {
218 if(treeWriter) delete treeWriter;
219 if(outputFile) delete outputFile;
220 cerr << "** ERROR: " << e.what() << endl;
221 return 1;
222 }
223}
Note: See TracBrowser for help on using the repository browser.