Fork me on GitHub

source: git/readers/DelphesROOT.cpp@ 861ad5a

ImprovedOutputFile Timing
Last change on this file since 861ad5a 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: 7.2 KB
Line 
1/*
2 * Delphes: a framework for fast simulation of a generic collider experiment
3 * Copyright (C) 2012-2014 Universite catholique de Louvain (UCL), Belgium
4 *
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.
9 *
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.
14 *
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
19#include <algorithm>
20#include <iostream>
21#include <memory>
22#include <sstream>
23#include <stdexcept>
24
25#include <map>
26#include <vector>
27
28#include <signal.h>
29#include <stdio.h>
30#include <stdlib.h>
31
32#include "TApplication.h"
33#include "TROOT.h"
34
35#include "TClonesArray.h"
36#include "TDatabasePDG.h"
37#include "TFile.h"
38#include "TLorentzVector.h"
39#include "TObjArray.h"
40#include "TParticlePDG.h"
41#include "TStopwatch.h"
42
43#include "classes/DelphesClasses.h"
44#include "classes/DelphesFactory.h"
45#include "classes/DelphesStream.h"
46#include "modules/Delphes.h"
47
48#include "ExRootAnalysis/ExRootProgressBar.h"
49#include "ExRootAnalysis/ExRootTreeBranch.h"
50#include "ExRootAnalysis/ExRootTreeReader.h"
51#include "ExRootAnalysis/ExRootTreeWriter.h"
52
53using namespace std;
54
55//---------------------------------------------------------------------------
56
57//---------------------------------------------------------------------------
58
59static bool interrupted = false;
60
61void SignalHandler(int sig)
62{
63 interrupted = true;
64}
65
66//---------------------------------------------------------------------------
67
68int main(int argc, char *argv[])
69{
70 char appName[] = "DelphesROOT";
71 stringstream message;
72 TFile *inputFile = 0;
73 TFile *outputFile = 0;
74 TStopwatch eventStopWatch;
75 ExRootTreeWriter *treeWriter = 0;
76 ExRootTreeBranch *branchEvent = 0;
77 ExRootConfReader *confReader = 0;
78 Delphes *modularDelphes = 0;
79 DelphesFactory *factory = 0;
80 GenParticle *gen;
81 HepMCEvent *element, *eve;
82 Candidate *candidate;
83 Int_t pdgCode;
84
85 const Double_t c_light = 2.99792458E8;
86
87 TObjArray *allParticleOutputArray = 0, *stableParticleOutputArray = 0, *partonOutputArray = 0;
88 Int_t i;
89 Long64_t eventCounter, numberOfEvents;
90
91 if(argc < 4)
92 {
93 cout << " Usage: " << appName << " config_file"
94 << " output_file"
95 << " input_file(s)" << endl;
96 cout << " config_file - configuration file in Tcl format," << endl;
97 cout << " output_file - output file in ROOT format," << endl;
98 cout << " input_file(s) - input file(s) in ROOT format." << endl;
99 return 1;
100 }
101
102 signal(SIGINT, SignalHandler);
103
104 gROOT->SetBatch();
105
106 int appargc = 1;
107 char *appargv[] = {appName};
108 TApplication app(appName, &appargc, appargv);
109
110 try
111 {
112 outputFile = TFile::Open(argv[2], "CREATE");
113
114 if(outputFile == NULL)
115 {
116 message << "can't open " << argv[2] << endl;
117 throw runtime_error(message.str());
118 }
119
120 treeWriter = new ExRootTreeWriter(outputFile, "Delphes");
121
122 branchEvent = treeWriter->NewBranch("Event", HepMCEvent::Class());
123
124 confReader = new ExRootConfReader;
125 confReader->ReadFile(argv[1]);
126
127 modularDelphes = new Delphes("Delphes");
128 modularDelphes->SetConfReader(confReader);
129 modularDelphes->SetTreeWriter(treeWriter);
130
131 TChain *chain = new TChain("Delphes");
132
133 factory = modularDelphes->GetFactory();
134 allParticleOutputArray = modularDelphes->ExportArray("allParticles");
135 stableParticleOutputArray = modularDelphes->ExportArray("stableParticles");
136 partonOutputArray = modularDelphes->ExportArray("partons");
137
138 modularDelphes->InitTask();
139
140 for(i = 3; i < argc && !interrupted; ++i)
141 {
142 cout << "** Reading " << argv[i] << endl;
143
144 chain->Add(argv[i]);
145 ExRootTreeReader *treeReader = new ExRootTreeReader(chain);
146
147 inputFile = TFile::Open(argv[i]);
148
149 if(inputFile == NULL)
150 {
151 message << "can't open " << argv[i] << endl;
152 throw runtime_error(message.str());
153 }
154
155 numberOfEvents = treeReader->GetEntries();
156 TClonesArray *branchParticle = treeReader->UseBranch("Particle");
157 TClonesArray *branchHepMCEvent = treeReader->UseBranch("Event");
158
159 if(numberOfEvents <= 0) continue;
160
161 // ExRootProgressBar progressBar(numberOfEvents - 1);
162 ExRootProgressBar progressBar(-1);
163
164 // Loop over all objects
165 eventCounter = 0;
166 modularDelphes->Clear();
167 treeWriter->Clear();
168 for(Int_t entry = 0; entry < numberOfEvents && !interrupted; ++entry)
169 {
170
171 treeReader->ReadEntry(entry);
172
173 // -- TBC need also to include event weights --
174
175 eve = (HepMCEvent *)branchHepMCEvent->At(0);
176 element = static_cast<HepMCEvent *>(branchEvent->NewEntry());
177
178 element->Number = eventCounter;
179
180 element->ProcessID = eve->ProcessID;
181 element->MPI = eve->MPI;
182 element->Weight = eve->Weight;
183 element->Scale = eve->Scale;
184 element->AlphaQED = eve->AlphaQED;
185 element->AlphaQCD = eve->AlphaQCD;
186
187 element->ID1 = eve->ID1;
188 element->ID2 = eve->ID2;
189 element->X1 = eve->X1;
190 element->X2 = eve->X2;
191 element->ScalePDF = eve->ScalePDF;
192 element->PDF1 = eve->PDF1;
193 element->PDF2 = eve->PDF2;
194
195 element->ReadTime = eve->ReadTime;
196 element->ProcTime = eve->ProcTime;
197
198 for(Int_t j = 0; j < branchParticle->GetEntriesFast(); j++)
199 {
200
201 gen = (GenParticle *)branchParticle->At(j);
202 candidate = factory->NewCandidate();
203
204 candidate->Momentum = gen->P4();
205 candidate->Position.SetXYZT(gen->X, gen->Y, gen->Z, gen->T * 1.0E3 * c_light);
206
207 candidate->PID = gen->PID;
208 candidate->Status = gen->Status;
209
210 candidate->M1 = gen->M1;
211 candidate->M2 = gen->M2;
212
213 candidate->D1 = gen->D1;
214 candidate->D2 = gen->D2;
215
216 candidate->Charge = gen->Charge;
217 candidate->Mass = gen->Mass;
218
219 allParticleOutputArray->Add(candidate);
220
221 pdgCode = TMath::Abs(gen->PID);
222
223 if(gen->Status == 1)
224 {
225 stableParticleOutputArray->Add(candidate);
226 }
227 else if(pdgCode <= 5 || pdgCode == 21 || pdgCode == 15)
228 {
229 partonOutputArray->Add(candidate);
230 }
231 }
232
233 modularDelphes->ProcessTask();
234
235 treeWriter->Fill();
236
237 modularDelphes->Clear();
238 treeWriter->Clear();
239
240 progressBar.Update(eventCounter, eventCounter);
241 ++eventCounter;
242 }
243
244 progressBar.Update(eventCounter, eventCounter, kTRUE);
245 progressBar.Finish();
246
247 inputFile->Close();
248
249 delete treeReader;
250 }
251
252 modularDelphes->FinishTask();
253 treeWriter->Write();
254
255 cout << "** Exiting..." << endl;
256
257 delete modularDelphes;
258 delete confReader;
259 delete treeWriter;
260 delete outputFile;
261 delete chain;
262
263 return 0;
264 }
265 catch(runtime_error &e)
266 {
267 if(treeWriter) delete treeWriter;
268 if(outputFile) delete outputFile;
269 cerr << "** ERROR: " << e.what() << endl;
270 return 1;
271 }
272}
Note: See TracBrowser for help on using the repository browser.