Fork me on GitHub

source: git/readers/DelphesROOT.cpp@ 28077b2

ImprovedOutputFile Timing dual_readout llp
Last change on this file since 28077b2 was 94083e7, checked in by Pavel Demin <pavel-demin@…>, 6 years ago

remove maxEvents from DelphesROOT.cpp

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