Skip to content

Commit db6a2ac

Browse files
authored
Merge pull request #48989 from sihyunjeon/cmssw_merge/hepmc_writer
add hepmc writer
2 parents 5ad997d + 0d12c38 commit db6a2ac

File tree

2 files changed

+28
-3
lines changed

2 files changed

+28
-3
lines changed

GeneratorInterface/RivetInterface/plugins/GenParticles2HepMCConverter.cc

Lines changed: 25 additions & 3 deletions
Original file line numberDiff line numberDiff line change
@@ -22,6 +22,7 @@
2222
#include "HepMC3/GenParticle.h"
2323
#include "HepMC3/GenVertex.h"
2424
#include "HepMC3/Print.h"
25+
#include "HepMC3/WriterAscii.h"
2526

2627
#include <iostream>
2728
#include <map>
@@ -31,7 +32,7 @@ using namespace std;
3132
class GenParticles2HepMCConverter : public edm::stream::EDProducer<> {
3233
public:
3334
explicit GenParticles2HepMCConverter(const edm::ParameterSet& pset);
34-
~GenParticles2HepMCConverter() override {}
35+
~GenParticles2HepMCConverter() override;
3536

3637
void beginRun(edm::Run const& iRun, edm::EventSetup const&) override;
3738
void produce(edm::Event& event, const edm::EventSetup& eventSetup) override;
@@ -43,7 +44,10 @@ class GenParticles2HepMCConverter : public edm::stream::EDProducer<> {
4344
edm::ESGetToken<HepPDT::ParticleDataTable, PDTRecord> pTable_;
4445

4546
const double cmEnergy_;
47+
const bool writeHepMC_;
48+
const std::string outputFile_;
4649
HepMC3::GenCrossSectionPtr xsec_;
50+
std::shared_ptr<HepMC3::Writer> writer_;
4751

4852
private:
4953
inline HepMC3::FourVector FourVector(const reco::Candidate::Point& point) {
@@ -57,14 +61,27 @@ class GenParticles2HepMCConverter : public edm::stream::EDProducer<> {
5761
};
5862

5963
GenParticles2HepMCConverter::GenParticles2HepMCConverter(const edm::ParameterSet& pset)
60-
// dummy value to set incident proton pz for particle gun samples
61-
: cmEnergy_(pset.getUntrackedParameter<double>("cmEnergy", 13000)) {
64+
: cmEnergy_(pset.getUntrackedParameter<double>("cmEnergy", 13000)),
65+
writeHepMC_(pset.getUntrackedParameter<bool>("writeHepMC", false)),
66+
outputFile_(pset.getUntrackedParameter<std::string>("outputFile", "events.hepmc")) {
6267
genParticlesToken_ = consumes<reco::CandidateView>(pset.getParameter<edm::InputTag>("genParticles"));
6368
genEventInfoToken_ = consumes<GenEventInfoProduct>(pset.getParameter<edm::InputTag>("genEventInfo"));
6469
genRunInfoToken_ = consumes<GenRunInfoProduct, edm::InRun>(pset.getParameter<edm::InputTag>("genEventInfo"));
6570
pTable_ = esConsumes<HepPDT::ParticleDataTable, PDTRecord>();
6671

6772
produces<edm::HepMC3Product>("unsmeared");
73+
74+
// Initialize HepMC3 writer
75+
if (writeHepMC_) {
76+
writer_ = std::make_shared<HepMC3::WriterAscii>(outputFile_);
77+
}
78+
}
79+
80+
GenParticles2HepMCConverter::~GenParticles2HepMCConverter() {
81+
if (writeHepMC_) {
82+
writer_->close();
83+
writer_.reset();
84+
}
6885
}
6986

7087
void GenParticles2HepMCConverter::beginRun(edm::Run const& iRun, edm::EventSetup const&) {
@@ -208,6 +225,11 @@ void GenParticles2HepMCConverter::produce(edm::Event& event, const edm::EventSet
208225
}
209226
}
210227

228+
// Write event to HepMC file
229+
if (writeHepMC_) {
230+
writer_->write_event(hepmc_event);
231+
}
232+
211233
// Finalize HepMC event record
212234
auto hepmc_product = std::make_unique<edm::HepMC3Product>(&hepmc_event);
213235
event.put(std::move(hepmc_product), "unsmeared");

GeneratorInterface/RivetInterface/test/genParticles2HepMC_cfg.py

Lines changed: 3 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -26,5 +26,8 @@
2626
process.load("GeneratorInterface.RivetInterface.genParticles2HepMC_cfi")
2727
process.genParticles2HepMC.genParticles = cms.InputTag("mergedGenParticles")
2828

29+
# Turn on to write HepMC file
30+
#process.genParticles2HepMC.writeHepMC = cms.untracked.bool(True)
31+
2932
process.path = cms.Path(process.mergedGenParticles*process.genParticles2HepMC)
3033

0 commit comments

Comments
 (0)