Fork me on GitHub

source: svn/trunk/readers/DelphesHepMC.cpp@ 1179

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

add MaxEvents and SkipEvents

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