Fork me on GitHub

source: git/converters/hepmc2pileup.cpp@ ededa33

ImprovedOutputFile Timing
Last change on this file since ededa33 was 341014c, checked in by Pavel Demin <pavel-demin@…>, 6 years ago

apply .clang-format to all .h, .cc and .cpp files

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