cepgen is hosted by Hepforge, IPPP Durham
CepGen 1.2.5
Central exclusive processes event generator
Loading...
Searching...
No Matches
HepMC3Handler.cpp
Go to the documentation of this file.
1/*
2 * CepGen: a central exclusive processes event generator
3 * Copyright (C) 2013-2024 Laurent Forthomme
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 * 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 <HepMC3/GenCrossSection.h>
20#include <HepMC3/Version.h>
21
22#include <memory>
23
24#include "CepGen/Event/Event.h"
28#include "CepGen/Utils/Value.h"
30
31using namespace HepMC3;
32
33namespace cepgen {
38 template <typename T>
40 public:
41 explicit HepMC3Handler(const ParametersList& params)
42 : EventExporter(params),
43 output_(new T(steer<std::string>("filename").c_str())),
44 xs_(new GenCrossSection),
45 run_info_(new GenRunInfo) {
46 output_->set_run_info(run_info_);
47 run_info_->set_weight_names({"Default"});
48 CG_INFO("HepMC") << "Interfacing module initialised "
49 << "for HepMC version " << HEPMC3_VERSION << ".";
50 }
51 virtual ~HepMC3Handler() { output_->close(); }
52
54 auto desc = EventExporter::description();
55 desc.setDescription("HepMC3 ASCII file output module");
56 desc.add<std::string>("filename", "output.hepmc").setDescription("Output filename");
57 return desc;
58 }
59
60 bool operator<<(const Event& cg_event) override {
61 CepGenEvent event(cg_event);
62 event.set_cross_section(xs_);
63 event.set_run_info(run_info_);
64 event.set_event_number(event_num_++);
65 output_->write_event(event);
66 return !output_->failed();
67 }
68 void setCrossSection(const Value& cross_section) override {
69 xs_->set_cross_section(cross_section, cross_section.uncertainty());
70 }
71
72 private:
73 void initialise() override {}
74
75 std::unique_ptr<T> output_;
76 std::shared_ptr<GenCrossSection> xs_;
77 std::shared_ptr<GenRunInfo> run_info_;
78 };
79} // namespace cepgen
80
81//----------------------------------------------------------------------
82// Defining the various templated plugins made available by this
83// specific version of HepMC
84//----------------------------------------------------------------------
85#include <HepMC3/WriterAscii.h>
86#include <HepMC3/WriterHEPEVT.h>
91
92#if HEPMC3_VERSION_CODE >= 3001000
93#include <HepMC3/WriterAsciiHepMC2.h>
94typedef cepgen::HepMC3Handler<WriterAsciiHepMC2> HepMC3HepMC2Handler;
95REGISTER_EXPORTER("hepmc3_hepmc2", HepMC3HepMC2Handler);
96#if HEPMC3_VERSION_CODE >= 3002005 && HEPMC3_USE_COMPRESSION
97#include <HepMC3/WriterGZ.h>
104REGISTER_EXPORTER("hepmc_z", HepMC3AsciiZHandler);
105REGISTER_EXPORTER("hepevt_z", HepMC3HEPEVTZHandler);
106REGISTER_EXPORTER("hepmc_lzma", HepMC3AsciiLZMAHandler);
107REGISTER_EXPORTER("hepevt_lzma", HepMC3HEPEVTLZMAHandler);
108REGISTER_EXPORTER("hepmc_bz2", HepMC3AsciiBZ2Handler);
109REGISTER_EXPORTER("hepevt_bz2", HepMC3HEPEVTBZ2Handler);
110#endif
111#endif
112
113#ifdef HEPMC3_ROOTIO
114#include <HepMC3/WriterRoot.h>
115#include <HepMC3/WriterRootTree.h>
116typedef cepgen::HepMC3Handler<WriterRoot> HepMC3RootHandler;
117typedef cepgen::HepMC3Handler<WriterRootTree> HepMC3RootTreeHandler;
118REGISTER_EXPORTER("hepmc_root", HepMC3RootHandler);
119REGISTER_EXPORTER("hepmc_root_tree", HepMC3RootTreeHandler);
120#endif
121
122#ifdef HEPMC3_EXTRA_PLUGINS
123#include <ConvertExample/include/WriterDOT.h>
124#include <ConvertExample/include/WriterRootTreeOPAL.h>
125typedef cepgen::HepMC3Handler<WriterDOT> HepMC3DOTHandler;
126typedef cepgen::HepMC3Handler<WriterRootTreeOPAL> HepMC3RootTreeOPALHandler;
127REGISTER_EXPORTER("hepmc_dot", HepMC3DOTHandler);
128REGISTER_EXPORTER("hepmc_root_tree_opal", HepMC3RootTreeOPALHandler);
129#endif
#define REGISTER_EXPORTER(name, obj)
Add a generic export module definition to the factory.
cepgen::HepMC3Handler< WriterAscii > HepMC3AsciiHandler
cepgen::HepMC3Handler< WriterHEPEVT > HepMC3HEPEVTHandler
#define CG_INFO(mod)
Definition Message.h:216
Interfacing between CepGen and HepMC event definitions.
Output format handler for events export.
unsigned long long event_num_
Event index.
Container for the information on the in- and outgoing particles' kinematics.
Definition Event.h:28
Handler for the HepMC3 file output.
HepMC3Handler(const ParametersList &params)
bool operator<<(const Event &cg_event) override
Writer operator.
static ParametersDescription description()
void setCrossSection(const Value &cross_section) override
Specify the cross section value, in pb.
A description object for parameters collection.
static ParametersDescription description()
Description of all object parameters.
Definition Steerable.cpp:42
T steer(const std::string &key) const
Retrieve a parameters as previously steered.
Definition Steerable.h:39
A scalar value with its uncertainty.
Definition Value.h:26
double uncertainty() const
Absolute uncertainty around the central value.
Definition Value.h:35
Common namespace for this Monte Carlo generator.