Fork me on GitHub

source: git/readers/DelphesProIO.cpp@ 1292d20

ImprovedOutputFile Timing llp
Last change on this file since 1292d20 was 538edc4, checked in by Chekanov <chekanov@…>, 6 years ago

Corrected ProIO reader to make it consistent ProIO.

  • Property mode set to 100644
File size: 11.9 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 * @author S.Chekanov (ANL)
18 */
19
20#include <stdexcept>
21#include <iostream>
22#include <sstream>
23#include <memory>
24
25#include <map>
26
27#include <stdlib.h>
28#include <signal.h>
29#include <stdio.h>
30
31#include "TROOT.h"
32#include "TApplication.h"
33
34#include "TFile.h"
35#include "TObjArray.h"
36#include "TStopwatch.h"
37#include "TDatabasePDG.h"
38#include "TParticlePDG.h"
39#include "TLorentzVector.h"
40
41#include "modules/Delphes.h"
42#include "classes/DelphesStream.h"
43#include "classes/DelphesClasses.h"
44#include "classes/DelphesFactory.h"
45
46#include "ExRootAnalysis/ExRootTreeWriter.h"
47#include "ExRootAnalysis/ExRootTreeBranch.h"
48#include "ExRootAnalysis/ExRootProgressBar.h"
49
50#include <proio/reader.h>
51#include <proio/event.h>
52#include <proio/model/mc.pb.h>
53namespace model=proio::model::mc;
54
55
56using namespace std;
57
58//---------------------------------------------------------------------------
59// This method dynamically checks the message type (varint or not) depending on
60// non-zero value of units momentumUnit and positionUnit.
61void ConvertInput(proio::Event *event, double momentumUnit, double positionUnit,
62 ExRootTreeBranch *branch, DelphesFactory *factory,
63 TObjArray *allParticleOutputArray, TObjArray *stableParticleOutputArray,
64 TObjArray *partonOutputArray, TStopwatch *readStopWatch, TStopwatch *procStopWatch)
65{
66
67 HepMCEvent *element;
68 Candidate *candidate;
69 TDatabasePDG *pdg;
70 TParticlePDG *pdgParticle;
71 Int_t pdgCode;
72
73 Int_t pid, status;
74 Double_t px, py, pz, mass;
75 Double_t x, y, z, t;
76
77 pdg = TDatabasePDG::Instance();
78
79 // event information
80 element = static_cast<HepMCEvent *>(branch->NewEntry());
81
82 int nID=0;
83 double weight=0;
84 int process_id=0;
85 auto mcentries = event->TaggedEntries("MCParameters");
86 for (uint64_t mcentryID : mcentries) {
87 auto mcpar = dynamic_cast<proio::model::mc::MCParameters *>(event->GetEntry(mcentryID));
88 nID=mcpar->number();
89 weight=mcpar->weight();
90 process_id=mcpar->processid();
91 break; // consider only most generic from 1st entry
92 };
93
94 element->Number = nID;
95 element->ProcessID =process_id;
96 element->Weight = weight;
97
98/*
99 // Pythia8 specific
100 element->MPI = mutableEvent->mpi();
101 element->Scale = mutableEvent->scale();
102 element->AlphaQED = mutableEvent->alpha_qed();
103 element->AlphaQCD = mutableEvent->alpha_qcd();
104 element->ID1 = mutableEvent->id1();
105 element->ID2 = mutableEvent->id2();
106 element->X1 = mutableEvent->x1();
107 element->X2 = mutableEvent->x2();
108 element->ScalePDF = mutableEvent->scale_pdf();
109 element->PDF1 = mutableEvent->pdf1();
110 element->PDF2 = mutableEvent->pdf2();
111*/
112
113 element->ReadTime = readStopWatch->RealTime();
114 element->ProcTime = procStopWatch->RealTime();
115
116
117
118 if ( momentumUnit >0 && positionUnit>0) {
119
120
121 auto entries = event->TaggedEntries("VarintPackedParticles");
122
123 for (uint64_t entryID : entries) {
124
125 auto mutableParticles = dynamic_cast<model::VarintPackedParticles *>(event->GetEntry(entryID));
126
127 for(int i = 0; i < mutableParticles->pdg_size(); ++i)
128 {
129 pid = mutableParticles->pdg(i);
130 status = mutableParticles->status(i);
131 px = mutableParticles->px(i)/momentumUnit;
132 py = mutableParticles->py(i)/momentumUnit;
133 pz = mutableParticles->pz(i)/momentumUnit;
134 mass = mutableParticles->mass(i)/momentumUnit;
135 x = mutableParticles->x(i)/positionUnit;
136 y = mutableParticles->y(i)/positionUnit;
137 z = mutableParticles->z(i)/positionUnit;
138 t = mutableParticles->t(i)/positionUnit;
139
140 candidate = factory->NewCandidate();
141 candidate->PID = pid;
142 pdgCode = TMath::Abs(candidate->PID);
143 candidate->Status = status;
144 candidate->M1 = mutableParticles->parent1(i);
145 candidate->M2 = mutableParticles->parent2(i);
146 candidate->D1 = mutableParticles->child1(i);
147 candidate->D2 = mutableParticles->child2(i);
148 pdgParticle = pdg->GetParticle(pid);
149 candidate->Charge = mutableParticles->charge(i)/3.0;
150 //candidate->Charge = pdgParticle ? Int_t(pdgParticle->Charge()/3.0) : -999;
151 candidate->Mass = mass;
152 candidate->Momentum.SetXYZM(px, py, pz, mass);
153 candidate->Position.SetXYZT(x, y, z, t);
154 allParticleOutputArray->Add(candidate);
155 if(!pdgParticle) continue;
156
157 if(status == 1)
158 {
159 stableParticleOutputArray->Add(candidate);
160 }
161 else if(pdgCode <= 5 || pdgCode == 21 || pdgCode == 15)
162 {
163 partonOutputArray->Add(candidate);
164 }
165 }
166
167 }
168
169 } else {
170
171
172 auto entries = event->TaggedEntries("Particle");
173
174 for (uint64_t entryID : entries) {
175 auto part = dynamic_cast<model::Particle *>(event->GetEntry(entryID));
176 pid = part->pdg();
177 status = part->status();
178 model::XYZF pvector=part->p();
179 px=pvector.x();
180 py=pvector.y();
181 pz=pvector.z();
182 mass = part->mass();
183 auto v =part->vertex();
184 x=v.x();
185 y=v.y();
186 z=v.z();
187 t=v.t();
188
189 candidate = factory->NewCandidate();
190 candidate->PID = pid;
191 pdgCode = TMath::Abs(candidate->PID);
192 candidate->Status = status;
193
194 int M1=0;
195 int M2=0;
196 for (int k1=0; k1<part->parent_size(); k1++){
197 if (k1==0) {
198 auto mother = dynamic_cast<model::Particle *>(event->GetEntry(part->parent(0)));
199 M1=mother->barcode();
200 }
201 if (k1==1) {
202 auto mother = dynamic_cast<model::Particle *>(event->GetEntry(part->parent(1)));
203 M2=mother->barcode();
204 }
205 }
206
207
208 int D1=0;
209 int D2=0;
210 for (int k1=0; k1<part->child_size(); k1++){
211 if (k1==0) {
212 auto child = dynamic_cast<model::Particle *>(event->GetEntry(part->child(0)));
213 D1=child->barcode();
214 }
215
216 if (k1==1) {
217 auto child = dynamic_cast<model::Particle *>(event->GetEntry(part->child(1)));
218 D2=child->barcode();
219 };
220 }
221
222
223 candidate->M1 = M1;
224 candidate->M2 = M2;
225 candidate->D1 = D1;
226 candidate->D2 = D2;
227
228 pdgParticle = pdg->GetParticle(pid);
229 candidate->Charge = pdgParticle ? Int_t(pdgParticle->Charge()/3.0) : -999;
230 candidate->Mass = mass;
231
232 candidate->Momentum.SetXYZM(px, py, pz, mass);
233
234 candidate->Position.SetXYZT(x, y, z, t);
235
236 allParticleOutputArray->Add(candidate);
237
238 if(!pdgParticle) continue;
239
240 if(status == 1)
241 {
242 stableParticleOutputArray->Add(candidate);
243 }
244 else if(pdgCode <= 5 || pdgCode == 21 || pdgCode == 15)
245 {
246 partonOutputArray->Add(candidate);
247 }
248
249 }
250
251
252 } // end particle type
253
254}
255
256//---------------------------------------------------------------------------
257
258static bool interrupted = false;
259
260void SignalHandler(int sig)
261{
262 interrupted = true;
263}
264
265//---------------------------------------------------------------------------
266
267int main(int argc, char *argv[])
268{
269 char appName[] = "DelphesProIO";
270 stringstream message;
271 proio::Reader *inputFile = 0;
272 TFile *outputFile = 0;
273 TStopwatch readStopWatch, procStopWatch;
274 ExRootTreeWriter *treeWriter = 0;
275 ExRootTreeBranch *branchEvent = 0;
276 ExRootConfReader *confReader = 0;
277 Delphes *modularDelphes = 0;
278 DelphesFactory *factory = 0;
279 TObjArray *allParticleOutputArray = 0, *stableParticleOutputArray = 0, *partonOutputArray = 0;
280 Int_t i;
281 Long64_t eventCounter, numberOfEvents;
282
283 if(argc < 4)
284 {
285 cout << " Usage: " << appName << " config_file" << " output_file" << " input_file(s)" << endl;
286 cout << " config_file - configuration file in Tcl format," << endl;
287 cout << " output_file - output file in ROOT format," << endl;
288 cout << " input_file(s) - input file(s) in ProIO format." << endl;
289 return 1;
290 }
291
292 signal(SIGINT, SignalHandler);
293
294 gROOT->SetBatch();
295
296 int appargc = 1;
297 char *appargv[] = {appName};
298 TApplication app(appName, &appargc, appargv);
299
300 try
301 {
302 outputFile = TFile::Open(argv[2], "CREATE");
303
304 if(outputFile == NULL)
305 {
306 message << "can't open " << argv[2] << endl;
307 throw runtime_error(message.str());
308 }
309
310 treeWriter = new ExRootTreeWriter(outputFile, "Delphes");
311
312 branchEvent = treeWriter->NewBranch("Event", HepMCEvent::Class());
313
314 confReader = new ExRootConfReader;
315 confReader->ReadFile(argv[1]);
316
317 modularDelphes = new Delphes("Delphes");
318 modularDelphes->SetConfReader(confReader);
319 modularDelphes->SetTreeWriter(treeWriter);
320
321 factory = modularDelphes->GetFactory();
322 allParticleOutputArray = modularDelphes->ExportArray("allParticles");
323 stableParticleOutputArray = modularDelphes->ExportArray("stableParticles");
324 partonOutputArray = modularDelphes->ExportArray("partons");
325
326 modularDelphes->InitTask();
327
328 for(i = 3; i < argc && !interrupted; ++i)
329 {
330 cout << "** INFO: Reading " << argv[i] << endl;
331
332 inputFile=new proio::Reader( argv[i] );
333
334 if(inputFile == NULL)
335 {
336 message << "can't open " << argv[i] << endl;
337 throw runtime_error(message.str());
338 }
339
340
341/*
342// this is slow method, but general
343 inputFile->SeekToStart();
344 int nn=0;
345 auto event = new proio::Event();
346 while(inputFile->Next(event)){
347 auto entries = event->TaggedEntries("Particle");
348 if (entries.size()>0) nn++;
349 }
350*/
351
352
353 auto event = new proio::Event();
354
355
356 double varint_energy=0;
357 double varint_length=0;
358
359
360 auto max_n_events = std::numeric_limits<uint64_t>::max();
361 auto nn = inputFile->Skip(max_n_events);
362 cout << "** INFO: " << nn-1 << " events found in ProIO file" << endl;
363 inputFile->SeekToStart();
364 numberOfEvents = nn-1; // last event has only metadata (no particle record)
365
366 if(numberOfEvents <= 0) continue;
367
368 ExRootProgressBar progressBar(numberOfEvents - 1);
369
370
371 // Loop over all objects
372 modularDelphes->Clear();
373 treeWriter->Clear();
374 readStopWatch.Start();
375
376 for(eventCounter = 0; eventCounter < numberOfEvents && !interrupted; ++eventCounter)
377 {
378 inputFile->Next(event);
379 if(event == 0) continue;
380
381 // get metadata
382 if (eventCounter == 0) {
383 auto metadata = event->Metadata();
384 std::cout << "** INFO: ProIO file metadata:" << std::endl;
385 for (auto element : metadata) {
386 string key=(string)element.first;
387 string value=(string)(*element.second);
388 std::cout << "** INFO: " << key << " = " << value << std::endl;
389 if (key=="info:varint_energy") varint_energy=std::stod(value);
390 if (key=="info:varint_length") varint_length=std::stod(value);
391 }
392 }
393
394
395 readStopWatch.Stop();
396
397 procStopWatch.Start();
398
399 ConvertInput(event, varint_energy, varint_length,
400 branchEvent, factory,
401 allParticleOutputArray, stableParticleOutputArray,
402 partonOutputArray, &readStopWatch, &procStopWatch);
403
404 modularDelphes->ProcessTask();
405 procStopWatch.Stop();
406
407 treeWriter->Fill();
408
409 modularDelphes->Clear();
410 treeWriter->Clear();
411
412 readStopWatch.Start();
413 progressBar.Update(eventCounter);
414 }
415
416
417 progressBar.Update(eventCounter, eventCounter, kTRUE);
418 progressBar.Finish();
419
420 delete inputFile;
421 }
422
423
424 modularDelphes->FinishTask();
425 treeWriter->Write();
426
427 cout << "** Exiting..." << endl;
428
429 delete modularDelphes;
430 delete confReader;
431 delete treeWriter;
432 return 0;
433 }
434 catch(runtime_error &e)
435 {
436 if(treeWriter) delete treeWriter;
437 if(outputFile) delete outputFile;
438 cerr << "** ERROR: " << e.what() << endl;
439 return 1;
440 }
441}
Note: See TracBrowser for help on using the repository browser.