cepgen is hosted by Hepforge, IPPP Durham
CepGen 1.2.5
Central exclusive processes event generator
Loading...
Searching...
No Matches
EventExporter.cpp
Go to the documentation of this file.
1/*
2 * CepGen: a central exclusive processes event generator
3 * Copyright (C) 2014-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 <TFile.h>
20
21#include <sstream>
22
25#include "CepGen/Event/Event.h"
30#include "CepGen/Utils/String.h"
31#include "CepGen/Utils/Value.h"
32#include "CepGen/Version.h"
34
35namespace cepgen {
36 namespace root {
40 class EventExporter final : public cepgen::EventExporter {
41 public:
42 explicit EventExporter(const ParametersList& params)
43 : cepgen::EventExporter(params),
44 filename_(steer<std::string>("filename")),
45 compress_(steer<bool>("compress")),
46 file_(TFile::Open(filename_.data(), "recreate")) {
47 if (!file_->IsOpen())
48 throw CG_FATAL("root:EventExporter") << "Failed to create the output file!";
49 }
51 run_tree_.fill();
52 file_->Write();
53 }
54
57 desc.setDescription("ROOT TTree storage module");
58 desc.add<std::string>("filename", "output.root").setDescription("Output filename");
59 desc.add<bool>("compress", false).setDescription("Compress the event content? (merge down two-parton system)");
60 desc.add<bool>("autoFilename", false).setDescription("automatically generate the output filename");
61 return desc;
62 }
63
64 bool operator<<(const Event& ev) override {
65 evt_tree_.fill(ev, compress_);
66 run_tree_.num_events += 1;
67 return true;
68 }
69 void setCrossSection(const Value& cross_section) override {
70 run_tree_.xsect = cross_section;
71 run_tree_.errxsect = cross_section.uncertainty();
72 }
73
74 private:
75 void initialise() override {
76 if (steer<bool>("autoFilename")) {
77 auto filename = generateFilename();
78 CG_INFO("root:EventExporter") << "Output ROOT filename automatically set to '" << filename << "'.";
79 if (file_.reset(TFile::Open(filename.data(), "recreate")); !file_->IsOpen())
80 throw CG_FATAL("root:EventExporter") << "Failed to create the output file!";
81 }
82 run_tree_.create();
83 evt_tree_.create();
84 run_tree_.litigious_events = 0;
85 if (runParameters().hasProcess()) {
89 }
90 }
91 std::string generateFilename() const {
92 std::string evt_mods, proc_mode;
93 for (const auto& mod : runParameters().eventModifiersSequence())
94 evt_mods += (evt_mods.empty() ? "" : "-") + mod->name();
95 const auto symm = runParameters().process().parameters().get<bool>("symmetrise");
96 const auto sf_info =
97 utils::sanitise(runParameters().process().kinematics().incomingBeams().structureFunctions().serialise());
98 switch (runParameters().process().kinematics().incomingBeams().mode()) {
100 proc_mode = "el";
101 break;
103 proc_mode = symm ? "sd" : "sdie_" + sf_info;
104 break;
106 proc_mode = symm ? "sd" : "sdei_" + sf_info;
107 break;
109 proc_mode = "dd_" + sf_info;
110 break;
112 break;
113 }
114 return utils::format("cepgen%s_%s_%s_%gTeV%s.root",
116 runParameters().processName().data(),
117 proc_mode.data(),
118 runParameters().kinematics().incomingBeams().sqrtS() / 1000.,
119 evt_mods.data());
120 }
121
122 const std::string filename_;
123 const bool compress_;
124 std::unique_ptr<TFile> file_;
125 ROOT::CepGenRun run_tree_;
126 ROOT::CepGenEvent evt_tree_;
127 };
128 } // namespace root
129} // namespace cepgen
#define REGISTER_EXPORTER(name, obj)
Add a generic export module definition to the factory.
#define CG_FATAL(mod)
Definition Exception.h:61
#define CG_INFO(mod)
Definition Message.h:216
All useful information about a generated event.
void fill(const cepgen::Event &, bool compress=false)
Fill the tree with a new event.
void create()
Populate the tree and all associated branches.
All useful information about a generation run.
std::string process_name
Unique name of the process generated in this run.
unsigned int litigious_events
Number of litigious events in run.
double xsect
Process cross section, in pb.
double errxsect
Uncertainty on process cross section, in pb.
std::string process_parameters
Serialised process parameters.
unsigned int num_events
Number of events generated in run.
void create()
Populate the run tree.
void fill()
Fill the run tree.
double sqrt_s
Centre of mass energy for beam particles.
Output format handler for events export.
const RunParameters & runParameters() const
List of run parameters.
Container for the information on the in- and outgoing particles' kinematics.
Definition Event.h:28
double sqrtS() const
Incoming beams centre of mass energy (in GeV)
IncomingBeams & incomingBeams()
Beam/primary particle kinematics.
Definition Kinematics.h:36
const std::string & name() const
Module unique indexing name.
Definition NamedModule.h:42
A description object for parameters collection.
std::string serialise() const
Serialise a parameters collection into a parseable string.
T get(const std::string &key, const T &def=default_arg< T >::get()) const
Get a parameter value.
const Kinematics & kinematics() const
Events kinematics for phase space definition.
proc::Process & process()
Process object for cross-section computation/events generation.
std::string processName() const
Name of the process considered.
bool hasProcess() const
Is this parameters collection holding any physics process?
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
const ParametersList & parameters() const override
Module user-defined parameters.
A scalar value with its uncertainty.
Definition Value.h:26
double uncertainty() const
Absolute uncertainty around the central value.
Definition Value.h:35
Handler for the storage of events in a ROOT format.
EventExporter(const ParametersList &params)
static ParametersDescription description()
void setCrossSection(const Value &cross_section) override
Specify the cross section value, in pb.
bool operator<<(const Event &ev) override
Writer operator.
@ ElasticElastic
proton-proton elastic case
@ InelasticElastic
proton-proton single-dissociative (or elastic-inelastic) case
@ InelasticInelastic
proton-proton double-dissociative case
@ ElasticInelastic
proton-proton single-dissociative (or inelastic-elastic) case
std::string format(const std::string &fmt, Args... args)
Format a string using a printf style format descriptor.
Definition String.h:61
std::string sanitise(const std::string &str)
Replace all unsafe characters to build a computer-readable (and filename-safe) string.
Definition String.cpp:105
Common namespace for this Monte Carlo generator.
static const std::string tag
CepGen version.
Definition Version.h:28