Fork me on GitHub

source: git/readers/DelphesLHEF.cpp@ 93af22b

ImprovedOutputFile Timing dual_readout llp
Last change on this file since 93af22b was 1e8afcc, checked in by Pavel Demin <pavel.demin@…>, 9 years ago

replace Rwgt with Weight in DelphesPythia8 and DelphesLHEF

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