Fork me on GitHub

source: git/readers/DelphesCMSFWLite.cpp@ b1e03e5

ImprovedOutputFile Timing dual_readout llp
Last change on this file since b1e03e5 was 2886328, checked in by Pavel Demin <pavel-demin@…>, 7 years ago

replace branchRwgt with branchWeight in DelphesCMSFWLite.cpp

  • Property mode set to 100644
File size: 13.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
[975405a]19#include <algorithm>
[d7d2da3]20#include <stdexcept>
21#include <iostream>
22#include <sstream>
23#include <memory>
24
25#include <map>
[975405a]26#include <vector>
[d7d2da3]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 "TObjArray.h"
37#include "TStopwatch.h"
38#include "TDatabasePDG.h"
39#include "TParticlePDG.h"
40#include "TLorentzVector.h"
41
42#include "modules/Delphes.h"
43#include "classes/DelphesStream.h"
44#include "classes/DelphesClasses.h"
45#include "classes/DelphesFactory.h"
46
47#include "ExRootAnalysis/ExRootTreeWriter.h"
48#include "ExRootAnalysis/ExRootTreeBranch.h"
49#include "ExRootAnalysis/ExRootProgressBar.h"
50
[cd699d0]51#include "FWCore/FWLite/interface/FWLiteEnabler.h"
[d7d2da3]52#include "DataFormats/FWLite/interface/Event.h"
53#include "DataFormats/FWLite/interface/Handle.h"
54#include "DataFormats/HepMCCandidate/interface/GenParticle.h"
[386f526]55#include "DataFormats/PatCandidates/interface/PackedGenParticle.h"
[f29758e]56#include "SimDataFormats/GeneratorProducts/interface/HepMCProduct.h"
57#include "SimDataFormats/GeneratorProducts/interface/GenEventInfoProduct.h"
[3241a0e]58#include "SimDataFormats/GeneratorProducts/interface/LHEEventProduct.h"
59#include "SimDataFormats/GeneratorProducts/interface/WeightsInfo.h"
[d7d2da3]60
61using namespace std;
62
63//---------------------------------------------------------------------------
64
[f29758e]65void ConvertInput(fwlite::Event &event, Long64_t eventCounter,
[2886328]66 ExRootTreeBranch *branchEvent, ExRootTreeBranch *branchWeight,
[5df97ac]67 DelphesFactory *factory, TObjArray *allParticleOutputArray,
68 TObjArray *stableParticleOutputArray, TObjArray *partonOutputArray, Bool_t firstEvent)
[d7d2da3]69{
[f29758e]70
[480f9ed]71 fwlite::Handle< GenEventInfoProduct > handleGenEventInfo;
[3241a0e]72 fwlite::Handle< LHEEventProduct > handleLHEEvent;
[d7d2da3]73 fwlite::Handle< vector< reco::GenParticle > > handleParticle;
[386f526]74 fwlite::Handle< vector< pat::PackedGenParticle > > handlePackedParticle;
[480f9ed]75
[975405a]76 vector< reco::GenParticle >::const_iterator itParticle;
[386f526]77 vector< pat::PackedGenParticle >::const_iterator itPackedParticle;
[975405a]78
79 vector< const reco::Candidate * > vectorCandidate;
80 vector< const reco::Candidate * >::iterator itCandidate;
81
[f29758e]82 handleGenEventInfo.getByLabel(event, "generator");
[480f9ed]83
[5df97ac]84 if(!((handleLHEEvent.getBranchNameFor(event, "source")).empty()))
85 {
[480f9ed]86 handleLHEEvent.getByLabel(event, "source");
87 }
[5df97ac]88 else if(!((handleLHEEvent.getBranchNameFor(event, "externalLHEProducer")).empty()))
[480f9ed]89 {
90 handleLHEEvent.getByLabel(event, "externalLHEProducer");
91 }
[5df97ac]92 else if(firstEvent)
[480f9ed]93 {
[5df97ac]94 cout << "Wrong LHEEvent Label! Please, check the input file." << endl;
[480f9ed]95 }
96
[5df97ac]97 if(!((handleParticle.getBranchNameFor(event, "genParticles")).empty()))
[480f9ed]98 {
99 handleParticle.getByLabel(event, "genParticles");
100 }
[5df97ac]101 else if(!((handlePackedParticle.getBranchNameFor(event, "packedGenParticles")).empty()) && !((handleParticle.getBranchNameFor(event,"prunedGenParticles")).empty()))
[480f9ed]102 {
[837fa70]103 handleParticle.getByLabel(event, "prunedGenParticles");
[386f526]104 handlePackedParticle.getByLabel(event, "packedGenParticles");
[480f9ed]105 }
106 else
107 {
108 std::cout<<"Wrong GenParticle Label! Please, check the input file."<<std::endl;
109 exit(-1);
110 }
[d7d2da3]111
[7611cb9]112 Bool_t foundLHE = !((handleLHEEvent.getBranchNameFor(event, "source")).empty()) || !((handleLHEEvent.getBranchNameFor(event, "externalLHEProducer")).empty());
[386f526]113 Bool_t isMiniAOD = !((handlePackedParticle.getBranchNameFor(event, "packedGenParticles")).empty()) && ((handleParticle.getBranchNameFor(event, "genParticles")).empty()) ;
[d7d2da3]114
[f29758e]115 HepMCEvent *element;
[3241a0e]116 Weight *weight;
[d7d2da3]117 Candidate *candidate;
118 TDatabasePDG *pdg;
119 TParticlePDG *pdgParticle;
120 Int_t pdgCode;
121
122 Int_t pid, status;
[80d4a34]123 Double_t px, py, pz, e, mass;
[d7d2da3]124 Double_t x, y, z;
125
[f29758e]126 element = static_cast<HepMCEvent *>(branchEvent->NewEntry());
127
128 element->Number = eventCounter;
129
130 element->ProcessID = handleGenEventInfo->signalProcessID();
131 element->MPI = 1;
[529fe78]132 element->Weight = handleGenEventInfo->weight();
[f29758e]133 element->Scale = handleGenEventInfo->qScale();
134 element->AlphaQED = handleGenEventInfo->alphaQED();
135 element->AlphaQCD = handleGenEventInfo->alphaQCD();
136
137 element->ID1 = 0;
138 element->ID2 = 0;
139 element->X1 = 0.0;
140 element->X2 = 0.0;
141 element->ScalePDF = 0.0;
142 element->PDF1 = 0.0;
143 element->PDF2 = 0.0;
144
145 element->ReadTime = 0.0;
146 element->ProcTime = 0.0;
147
[7611cb9]148 if(foundLHE)
[3241a0e]149 {
[7611cb9]150 const vector< gen::WeightsInfo > &vectorWeightsInfo = handleLHEEvent->weights();
151 vector< gen::WeightsInfo >::const_iterator itWeightsInfo;
[386f526]152
[7611cb9]153 for(itWeightsInfo = vectorWeightsInfo.begin(); itWeightsInfo != vectorWeightsInfo.end(); ++itWeightsInfo)
154 {
[2886328]155 weight = static_cast<Weight *>(branchWeight->NewEntry());
[7611cb9]156 weight->Weight = itWeightsInfo->wgt;
[5df97ac]157 }
[3241a0e]158 }
159
[d7d2da3]160 pdg = TDatabasePDG::Instance();
161
[837fa70]162 for(itParticle = handleParticle->begin(); itParticle != handleParticle->end(); ++itParticle)
163 {
164 const reco::GenParticle &particle = *itParticle;
[5df97ac]165 if(!isMiniAOD || particle.status() != 1) vectorCandidate.push_back(&*itParticle);
[837fa70]166 }
[d7d2da3]167
[837fa70]168 for(itParticle = handleParticle->begin(); itParticle != handleParticle->end(); ++itParticle)
169 {
170 const reco::GenParticle &particle = *itParticle;
[d7d2da3]171
[837fa70]172 pid = particle.pdgId();
173 status = particle.status();
[5df97ac]174 if(isMiniAOD && particle.status() == 1) continue;
[837fa70]175 px = particle.px(); py = particle.py(); pz = particle.pz(); e = particle.energy(); mass = particle.mass();
176 x = particle.vx(); y = particle.vy(); z = particle.vz();
[d7d2da3]177
[837fa70]178 candidate = factory->NewCandidate();
[d7d2da3]179
[837fa70]180 candidate->PID = pid;
181 pdgCode = TMath::Abs(candidate->PID);
[d7d2da3]182
[837fa70]183 candidate->Status = status;
[f29758e]184
[837fa70]185 if(particle.mother())
186 {
187 itCandidate = find(vectorCandidate.begin(), vectorCandidate.end(), particle.mother());
188 if(itCandidate != vectorCandidate.end()) candidate->M1 = distance(vectorCandidate.begin(), itCandidate);
189 }
[975405a]190
[837fa70]191 itCandidate = find(vectorCandidate.begin(), vectorCandidate.end(), particle.daughter(0));
192 if(itCandidate != vectorCandidate.end()) candidate->D1 = distance(vectorCandidate.begin(), itCandidate);
[975405a]193
[837fa70]194 itCandidate = find(vectorCandidate.begin(), vectorCandidate.end(), particle.daughter(particle.numberOfDaughters() - 1));
195 if(itCandidate != vectorCandidate.end()) candidate->D2 = distance(vectorCandidate.begin(), itCandidate);
[d7d2da3]196
[837fa70]197 pdgParticle = pdg->GetParticle(pid);
198 candidate->Charge = pdgParticle ? Int_t(pdgParticle->Charge()/3.0) : -999;
199 candidate->Mass = mass;
[d7d2da3]200
[837fa70]201 candidate->Momentum.SetPxPyPzE(px, py, pz, e);
[d7d2da3]202
[837fa70]203 candidate->Position.SetXYZT(x*10.0, y*10.0, z*10.0, 0.0);
[d7d2da3]204
[837fa70]205 allParticleOutputArray->Add(candidate);
[d7d2da3]206
[837fa70]207 if(!pdgParticle) continue;
208
[5df97ac]209 if(status == 1)
[837fa70]210 {
211 // Prevent duplicated particle.
[5df97ac]212 if(!isMiniAOD) stableParticleOutputArray->Add(candidate);
[386f526]213 }
[837fa70]214 else if(pdgCode <= 5 || pdgCode == 21 || pdgCode == 15)
[1e1f73f]215 {
[837fa70]216 partonOutputArray->Add(candidate);
[1e1f73f]217 }
[837fa70]218 }
[5df97ac]219
220 if(!isMiniAOD) return;
[837fa70]221 // For MiniAOD sample,
222 // Only status==1 particles are saved to packedGenParticles.
223 for(itPackedParticle = handlePackedParticle->begin(); itPackedParticle != handlePackedParticle->end(); ++itPackedParticle)
224 {
225 vectorCandidate.push_back(&*itPackedParticle);
226 }
[386f526]227
[837fa70]228 for(itPackedParticle = handlePackedParticle->begin(); itPackedParticle != handlePackedParticle->end(); ++itPackedParticle)
229 {
230 const pat::PackedGenParticle &particle = *itPackedParticle;
[386f526]231
[837fa70]232 pid = particle.pdgId();
233 status = particle.status();
234 px = particle.px(); py = particle.py(); pz = particle.pz(); e = particle.energy(); mass = particle.mass();
235 x = particle.vx(); y = particle.vy(); z = particle.vz();
[386f526]236
[837fa70]237 candidate = factory->NewCandidate();
[386f526]238
[837fa70]239 candidate->PID = pid;
240 pdgCode = TMath::Abs(candidate->PID);
[386f526]241
[837fa70]242 candidate->Status = status;
[386f526]243
[837fa70]244 if(particle.mother(0))
245 {
246 itCandidate = find(vectorCandidate.begin(), vectorCandidate.end(), particle.mother(0));
247 if(itCandidate != vectorCandidate.end()) candidate->M1 = distance(vectorCandidate.begin(), itCandidate);
248 }
[386f526]249
[837fa70]250 itCandidate = find(vectorCandidate.begin(), vectorCandidate.end(), particle.daughter(0));
251 if(itCandidate != vectorCandidate.end()) candidate->D1 = distance(vectorCandidate.begin(), itCandidate);
[386f526]252
[837fa70]253 itCandidate = find(vectorCandidate.begin(), vectorCandidate.end(), particle.daughter(particle.numberOfDaughters() - 1));
254 if(itCandidate != vectorCandidate.end()) candidate->D2 = distance(vectorCandidate.begin(), itCandidate);
[386f526]255
[837fa70]256 pdgParticle = pdg->GetParticle(pid);
257 candidate->Charge = pdgParticle ? Int_t(pdgParticle->Charge()/3.0) : -999;
258 candidate->Mass = mass;
[386f526]259
[837fa70]260 candidate->Momentum.SetPxPyPzE(px, py, pz, e);
[386f526]261
[837fa70]262 candidate->Position.SetXYZT(x*10.0, y*10.0, z*10.0, 0.0);
[386f526]263
[837fa70]264 allParticleOutputArray->Add(candidate);
[386f526]265
[837fa70]266 if(!pdgParticle) continue;
267
268 if(status == 1)
269 {
270 stableParticleOutputArray->Add(candidate);
[d7d2da3]271 }
272 }
273}
274
275//---------------------------------------------------------------------------
276
277static bool interrupted = false;
278
279void SignalHandler(int sig)
280{
281 interrupted = true;
282}
283
284//---------------------------------------------------------------------------
285
286int main(int argc, char *argv[])
287{
288 char appName[] = "DelphesCMSFWLite";
289 stringstream message;
290 TFile *inputFile = 0;
291 TFile *outputFile = 0;
292 TStopwatch eventStopWatch;
293 ExRootTreeWriter *treeWriter = 0;
[2886328]294 ExRootTreeBranch *branchEvent = 0, *branchWeight = 0;
[d7d2da3]295 ExRootConfReader *confReader = 0;
296 Delphes *modularDelphes = 0;
297 DelphesFactory *factory = 0;
298 TObjArray *allParticleOutputArray = 0, *stableParticleOutputArray = 0, *partonOutputArray = 0;
299 Int_t i;
[5ca3d52]300 Long64_t eventCounter, numberOfEvents;
[7611cb9]301 Bool_t firstEvent = kTRUE;
[d7d2da3]302
303 if(argc < 4)
304 {
305 cout << " Usage: " << appName << " config_file" << " output_file" << " input_file(s)" << endl;
306 cout << " config_file - configuration file in Tcl format," << endl;
307 cout << " output_file - output file in ROOT format," << endl;
308 cout << " input_file(s) - input file(s) in ROOT format." << endl;
309 return 1;
310 }
311
312 signal(SIGINT, SignalHandler);
313
314 gROOT->SetBatch();
315
316 int appargc = 1;
317 char *appargv[] = {appName};
318 TApplication app(appName, &appargc, appargv);
319
[cd699d0]320 FWLiteEnabler::enable();
[386f526]321
[d7d2da3]322 try
323 {
324 outputFile = TFile::Open(argv[2], "CREATE");
325
326 if(outputFile == NULL)
327 {
328 message << "can't open " << argv[2] << endl;
329 throw runtime_error(message.str());
330 }
331
332 treeWriter = new ExRootTreeWriter(outputFile, "Delphes");
333
[f29758e]334 branchEvent = treeWriter->NewBranch("Event", HepMCEvent::Class());
[2886328]335 branchWeight = treeWriter->NewBranch("Weight", Weight::Class());
[f29758e]336
[d7d2da3]337 confReader = new ExRootConfReader;
338 confReader->ReadFile(argv[1]);
339
340 modularDelphes = new Delphes("Delphes");
341 modularDelphes->SetConfReader(confReader);
342 modularDelphes->SetTreeWriter(treeWriter);
343
344 factory = modularDelphes->GetFactory();
345 allParticleOutputArray = modularDelphes->ExportArray("allParticles");
346 stableParticleOutputArray = modularDelphes->ExportArray("stableParticles");
347 partonOutputArray = modularDelphes->ExportArray("partons");
348
349 modularDelphes->InitTask();
350
351 for(i = 3; i < argc && !interrupted; ++i)
352 {
353 cout << "** Reading " << argv[i] << endl;
354
355 inputFile = TFile::Open(argv[i]);
356
357 if(inputFile == NULL)
358 {
359 message << "can't open " << argv[i] << endl;
360 throw runtime_error(message.str());
361 }
362
363 fwlite::Event event(inputFile);
364
[5ca3d52]365 numberOfEvents = event.size();
[d7d2da3]366
[5ca3d52]367 if(numberOfEvents <= 0) continue;
[d7d2da3]368
[a0538b9]369 // ExRootProgressBar progressBar(numberOfEvents - 1);
370 ExRootProgressBar progressBar(-1);
[d7d2da3]371
372 // Loop over all objects
[5ca3d52]373 eventCounter = 0;
[d7d2da3]374 modularDelphes->Clear();
375 treeWriter->Clear();
[7611cb9]376
[d7d2da3]377 for(event.toBegin(); !event.atEnd() && !interrupted; ++event)
378 {
[2886328]379 ConvertInput(event, eventCounter, branchEvent, branchWeight, factory,
[5df97ac]380 allParticleOutputArray, stableParticleOutputArray, partonOutputArray, firstEvent);
[d7d2da3]381 modularDelphes->ProcessTask();
[386f526]382
[7611cb9]383 firstEvent = kFALSE;
[d7d2da3]384
385 treeWriter->Fill();
386
387 modularDelphes->Clear();
388 treeWriter->Clear();
389
[a0538b9]390 progressBar.Update(eventCounter, eventCounter);
[5ca3d52]391 ++eventCounter;
[d7d2da3]392 }
[a0538b9]393
394 progressBar.Update(eventCounter, eventCounter, kTRUE);
[d7d2da3]395 progressBar.Finish();
396
397 inputFile->Close();
398 }
399
400 modularDelphes->FinishTask();
401 treeWriter->Write();
402
403 cout << "** Exiting..." << endl;
404
405 delete modularDelphes;
406 delete confReader;
407 delete treeWriter;
408 delete outputFile;
409
410 return 0;
411 }
412 catch(runtime_error &e)
413 {
414 if(treeWriter) delete treeWriter;
415 if(outputFile) delete outputFile;
416 cerr << "** ERROR: " << e.what() << endl;
417 return 1;
418 }
419}
Note: See TracBrowser for help on using the repository browser.