Fork me on GitHub

source: svn/trunk/readers/DelphesSTDHEP.cpp@ 1253

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

add MaxEvents and SkipEvents

File size: 5.5 KB
RevLine 
[618]1#include <stdexcept>
2#include <iostream>
3#include <sstream>
4
[779]5#include <signal.h>
[879]6
[618]7#include "TROOT.h"
8#include "TApplication.h"
9
10#include "TFile.h"
[619]11#include "TObjArray.h"
[705]12#include "TStopwatch.h"
[679]13#include "TDatabasePDG.h"
[678]14#include "TParticlePDG.h"
[618]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/DelphesSTDHEPReader.h"
[618]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
[618]39int main(int argc, char *argv[])
40{
[779]41 char appName[] = "DelphesSTDHEP";
[621]42 stringstream message;
[655]43 FILE *inputFile = 0;
[621]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;
[1052]51 TObjArray *stableParticleOutputArray = 0, *allParticleOutputArray = 0, *partonOutputArray = 0;
52 DelphesSTDHEPReader *reader = 0;
[1130]53 Int_t i, maxEvents, skipEvents;
[795]54 Long64_t length, eventCounter;
[618]55
[795]56 if(argc < 3)
[618]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 STDHEP format," << endl;
62 cout << " with no input_file, or when input_file is -, read standard input." << endl;
[618]63 return 1;
64 }
65
[779]66 signal(SIGINT, SignalHandler);
67
[618]68 gROOT->SetBatch();
69
70 int appargc = 1;
71 char *appargv[] = {appName};
72 TApplication app(appName, &appargc, appargv);
73
[621]74 try
[618]75 {
[917]76 outputFile = TFile::Open(argv[2], "CREATE");
[618]77
[621]78 if(outputFile == NULL)
79 {
[917]80 message << "can't create output file " << argv[2];
[621]81 throw runtime_error(message.str());
82 }
[619]83
[649]84 treeWriter = new ExRootTreeWriter(outputFile, "Delphes");
[618]85
[866]86 branchEvent = treeWriter->NewBranch("Event", LHEFEvent::Class());
[619]87
[621]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");
[621]105 modularDelphes->SetConfReader(confReader);
106 modularDelphes->SetTreeWriter(treeWriter);
[960]107
[621]108 factory = modularDelphes->GetFactory();
[905]109 allParticleOutputArray = modularDelphes->ExportArray("allParticles");
110 stableParticleOutputArray = modularDelphes->ExportArray("stableParticles");
[689]111 partonOutputArray = modularDelphes->ExportArray("partons");
[618]112
[1052]113 reader = new DelphesSTDHEPReader;
[619]114
[621]115 modularDelphes->InitTask();
[618]116
[795]117 i = 3;
118 do
[618]119 {
[795]120 if(interrupted) break;
121
[952]122 if(i == argc || strncmp(argv[i], "-", 2) == 0)
[621]123 {
[780]124 cout << "** Reading standard input" << endl;
125 inputFile = stdin;
[797]126 length = -1;
[621]127 }
[780]128 else
129 {
130 cout << "** Reading " << argv[i] << endl;
131 inputFile = fopen(argv[i], "r");
[619]132
[780]133 if(inputFile == NULL)
134 {
[863]135 message << "can't open " << argv[i];
[780]136 throw runtime_error(message.str());
137 }
[649]138
[780]139 fseek(inputFile, 0L, SEEK_END);
[1038]140 length = ftello(inputFile);
[780]141 fseek(inputFile, 0L, SEEK_SET);
[1038]142
143 if(length <= 0)
144 {
145 fclose(inputFile);
146 ++i;
147 continue;
148 }
[780]149 }
150
[1052]151 reader->SetInputFile(inputFile);
[621]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)
[618]164 {
[1052]165 if(reader->EventReady())
[618]166 {
[795]167 ++eventCounter;
[618]168
[843]169 readStopWatch.Stop();
[619]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);
[843]178
[1130]179 treeWriter->Fill();
180
181 treeWriter->Clear();
182 }
183
[843]184 modularDelphes->Clear();
[1052]185 reader->Clear();
[843]186
187 readStopWatch.Start();
[618]188 }
[1042]189 progressBar.Update(ftello(inputFile), eventCounter);
[618]190 }
[779]191
[795]192 fseek(inputFile, 0L, SEEK_END);
[1042]193 progressBar.Update(ftello(inputFile), eventCounter, kTRUE);
[795]194 progressBar.Finish();
[780]195
196 if(inputFile != stdin) fclose(inputFile);
[795]197
198 ++i;
[621]199 }
[795]200 while(i < argc);
201
[621]202 modularDelphes->FinishTask();
203 treeWriter->Write();
[619]204
[621]205 cout << "** Exiting..." << endl;
[619]206
[1052]207 delete reader;
[621]208 delete modularDelphes;
[625]209 delete confReader;
[621]210 delete treeWriter;
211 delete outputFile;
[960]212
[621]213 return 0;
214 }
215 catch(runtime_error &e)
216 {
[649]217 if(treeWriter) delete treeWriter;
[621]218 if(outputFile) delete outputFile;
219 cerr << "** ERROR: " << e.what() << endl;
220 return 1;
221 }
[618]222}
Note: See TracBrowser for help on using the repository browser.