Fork me on GitHub

source: git/readers/DelphesSTDHEP.cpp@ bdf9e2e

Last change on this file since bdf9e2e was 1fa50c2, checked in by Pavel Demin <pavel.demin@…>, 10 years ago

fix GPLv3 header

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