Fork me on GitHub

source: git/converters/hepmc2pileup.cpp@ 0a0a5ef

ImprovedOutputFile Timing dual_readout llp
Last change on this file since 0a0a5ef was d7d2da3, checked in by pavel <pavel@…>, 12 years ago

move branches/ModularDelphes to trunk

  • Property mode set to 100644
File size: 4.3 KB
RevLine 
[d7d2da3]1#include <stdexcept>
2#include <iostream>
3#include <sstream>
4
5#include <signal.h>
6
7#include "TROOT.h"
8#include "TApplication.h"
9
10#include "TFile.h"
11#include "TObjArray.h"
12#include "TStopwatch.h"
13#include "TDatabasePDG.h"
14#include "TParticlePDG.h"
15#include "TLorentzVector.h"
16
17#include "classes/DelphesClasses.h"
18#include "classes/DelphesFactory.h"
19#include "classes/DelphesHepMCReader.h"
20#include "classes/DelphesPileUpWriter.h"
21
22#include "ExRootAnalysis/ExRootTreeWriter.h"
23#include "ExRootAnalysis/ExRootTreeBranch.h"
24#include "ExRootAnalysis/ExRootProgressBar.h"
25
26using namespace std;
27
28//---------------------------------------------------------------------------
29
30static bool interrupted = false;
31
32void SignalHandler(int sig)
33{
34 interrupted = true;
35}
36
37//---------------------------------------------------------------------------
38
39int main(int argc, char *argv[])
40{
41 char appName[] = "hepmc2pileup";
42 stringstream message;
43 FILE *inputFile = 0;
44 DelphesFactory *factory = 0;
45 TObjArray *stableParticleOutputArray = 0, *allParticleOutputArray = 0, *partonOutputArray = 0;
46 TIterator *itParticle = 0;
47 Candidate *candidate = 0;
48 DelphesPileUpWriter *writer = 0;
49 DelphesHepMCReader *reader = 0;
50 Int_t i;
51 Long64_t length, eventCounter;
52
53 if(argc < 2)
54 {
55 cout << " Usage: " << appName << " output_file" << " [input_file(s)]" << endl;
56 cout << " output_file - output binary pile-up file," << endl;
57 cout << " input_file(s) - input file(s) in HepMC format," << endl;
58 cout << " with no input_file, or when input_file is -, read standard input." << endl;
59 return 1;
60 }
61
62 signal(SIGINT, SignalHandler);
63
64 gROOT->SetBatch();
65
66 int appargc = 1;
67 char *appargv[] = {appName};
68 TApplication app(appName, &appargc, appargv);
69
70 try
71 {
72 writer = new DelphesPileUpWriter(argv[1]);
73
74 factory = new DelphesFactory("ObjectFactory");
75 allParticleOutputArray = factory->NewPermanentArray();
76 stableParticleOutputArray = factory->NewPermanentArray();
77 partonOutputArray = factory->NewPermanentArray();
78
79 itParticle = stableParticleOutputArray->MakeIterator();
80
81 reader = new DelphesHepMCReader;
82
83 i = 2;
84 do
85 {
86 if(interrupted) break;
87
88 if(i == argc || strncmp(argv[i], "-", 2) == 0)
89 {
90 cout << "** Reading standard input" << endl;
91 inputFile = stdin;
92 length = -1;
93 }
94 else
95 {
96 cout << "** Reading " << argv[i] << endl;
97 inputFile = fopen(argv[i], "r");
98
99 if(inputFile == NULL)
100 {
101 message << "can't open " << argv[i];
102 throw runtime_error(message.str());
103 }
104
105 fseek(inputFile, 0L, SEEK_END);
106 length = ftello(inputFile);
107 fseek(inputFile, 0L, SEEK_SET);
108
109 if(length <= 0)
110 {
111 fclose(inputFile);
112 ++i;
113 continue;
114 }
115 }
116
117 reader->SetInputFile(inputFile);
118
119 ExRootProgressBar progressBar(length);
120
121 // Loop over all objects
122 eventCounter = 0;
123 factory->Clear();
124 reader->Clear();
125 while(reader->ReadBlock(factory, allParticleOutputArray,
126 stableParticleOutputArray, partonOutputArray) && !interrupted)
127 {
128 if(reader->EventReady())
129 {
130 ++eventCounter;
131
132 itParticle->Reset();
133 while((candidate = static_cast<Candidate*>(itParticle->Next())))
134 {
135 const TLorentzVector &position = candidate->Position;
136 const TLorentzVector &momentum = candidate->Momentum;
137 writer->WriteParticle(candidate->PID,
138 position.X(), position.Y(), position.Z(), position.T(),
139 momentum.Px(), momentum.Py(), momentum.Pz(), momentum.E());
140 }
141
142 writer->WriteEntry();
143
144 factory->Clear();
145 reader->Clear();
146 }
147 progressBar.Update(ftello(inputFile), eventCounter);
148 }
149
150 fseek(inputFile, 0L, SEEK_END);
151 progressBar.Update(ftello(inputFile), eventCounter, kTRUE);
152 progressBar.Finish();
153
154 if(inputFile != stdin) fclose(inputFile);
155
156 ++i;
157 }
158 while(i < argc);
159
160 writer->WriteIndex();
161
162 cout << "** Exiting..." << endl;
163
164 delete reader;
165 delete factory;
166 delete writer;
167
168 return 0;
169 }
170 catch(runtime_error &e)
171 {
172 if(writer) delete writer;
173 cerr << "** ERROR: " << e.what() << endl;
174 return 1;
175 }
176}
Note: See TracBrowser for help on using the repository browser.