Fork me on GitHub

source: svn/trunk/readers/DelphesLHEF.cpp@ 1324

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

read reweighting information from LHEF

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