Fork me on GitHub

source: git/readers/DelphesHepMC3.cpp@ 4564cba

Last change on this file since 4564cba was fc6bce2, checked in by Pavel Demin <pavel.demin@…>, 3 years ago

remove HepMC3 library

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