cepgen is hosted by Hepforge, IPPP Durham
CepGen 1.2.5
Central exclusive processes event generator
Loading...
Searching...
No Matches
EventImporter.cpp
Go to the documentation of this file.
1/*
2 * CepGen: a central exclusive processes event generator
3 * Copyright (C) 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
23#include "CepGen/Event/Event.h"
26#include "CepGen/Utils/Value.h"
28
29namespace cepgen {
30 namespace root {
35 public:
36 explicit EventImporter(const ParametersList& params)
37 : cepgen::EventImporter(params), file_(TFile::Open(steer<std::string>("filename").data())) {
38 if (!file_)
39 throw CG_FATAL("root::EventImporter")
40 << "Failed to load the ROOT file '" << steer<std::string>("filename") << "'.";
41 run_tree_.attach(file_.get());
42 evt_tree_.attach(file_.get());
43 }
44
47 desc.setDescription("ROOT TTree importer module");
48 desc.add<std::string>("filename", "output.root").setDescription("Input filename");
49 return desc;
50 }
51
52 bool operator>>(Event& evt) override { return evt_tree_.next(evt); }
53
54 private:
55 void initialise() override { setCrossSection(Value{run_tree_.xsect, run_tree_.errxsect}); }
56
57 const std::unique_ptr<TFile> file_;
58 ROOT::CepGenRun run_tree_;
59 ROOT::CepGenEvent evt_tree_;
60 };
61 } // namespace root
62} // namespace cepgen
#define REGISTER_EVENT_IMPORTER(name, obj)
Add a generic import module definition to the factory.
#define CG_FATAL(mod)
Definition Exception.h:61
All useful information about a generated event.
bool next(cepgen::Event &)
Read the next event in the file.
void attach()
Attach the event tree reader to a tree.
All useful information about a generation run.
double xsect
Process cross section, in pb.
double errxsect
Uncertainty on process cross section, in pb.
void attach(TFile *file, const std::string &run_tree=TREE_NAME)
Attach the run tree reader to a given tree Attach the run tree reader to a given file.
Base event importer module.
static ParametersDescription description()
void setCrossSection(const Value &xsec)
Specify the process cross section and uncertainty, in pb.
Container for the information on the in- and outgoing particles' kinematics.
Definition Event.h:28
A description object for parameters collection.
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
ROOT handler for an event tree import.
EventImporter(const ParametersList &params)
static ParametersDescription description()
bool operator>>(Event &evt) override
Read the next event.
Common namespace for this Monte Carlo generator.