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
21
#include "
CepGen/Core/Exception.h
"
22
#include "
CepGen/Core/RunParameters.h
"
23
#include "
CepGen/Event/Event.h
"
24
#include "
CepGen/EventFilter/EventImporter.h
"
25
#include "
CepGen/Modules/EventImporterFactory.h
"
26
#include "
CepGen/Utils/Value.h
"
27
#include "
CepGenAddOns/ROOTWrapper/ROOTTreeInfo.h
"
28
29
namespace
cepgen
{
30
namespace
root
{
34
class
EventImporter
:
public
cepgen::EventImporter
{
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
45
static
ParametersDescription
description
() {
46
auto
desc =
cepgen::EventImporter::description
();
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
63
using
ROOTEventImporter
=
cepgen::root::EventImporter
;
64
REGISTER_EVENT_IMPORTER
(
"root_tree"
,
ROOTEventImporter
);
EventImporterFactory.h
REGISTER_EVENT_IMPORTER
#define REGISTER_EVENT_IMPORTER(name, obj)
Add a generic import module definition to the factory.
Definition
EventImporterFactory.h:25
EventImporter.h
Event.h
Exception.h
CG_FATAL
#define CG_FATAL(mod)
Definition
Exception.h:61
ROOTTreeInfo.h
RunParameters.h
Value.h
ROOT::CepGenEvent
All useful information about a generated event.
Definition
ROOTTreeInfo.h:60
ROOT::CepGenEvent::next
bool next(cepgen::Event &)
Read the next event in the file.
Definition
ROOTTreeInfo.cpp:203
ROOT::CepGenEvent::attach
void attach()
Attach the event tree reader to a tree.
Definition
ROOTTreeInfo.cpp:135
ROOT::CepGenRun
All useful information about a generation run.
Definition
ROOTTreeInfo.h:29
ROOT::CepGenRun::xsect
double xsect
Process cross section, in pb.
Definition
ROOTTreeInfo.h:37
ROOT::CepGenRun::errxsect
double errxsect
Uncertainty on process cross section, in pb.
Definition
ROOTTreeInfo.h:38
ROOT::CepGenRun::attach
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.
Definition
ROOTTreeInfo.cpp:64
cepgen::EventImporter
Base event importer module.
Definition
EventImporter.h:30
cepgen::EventImporter::description
static ParametersDescription description()
Definition
EventImporter.h:34
cepgen::EventImporter::setCrossSection
void setCrossSection(const Value &xsec)
Specify the process cross section and uncertainty, in pb.
Definition
EventImporter.h:45
cepgen::Event
Container for the information on the in- and outgoing particles' kinematics.
Definition
Event.h:28
cepgen::ParametersDescription
A description object for parameters collection.
Definition
ParametersDescription.h:26
cepgen::ParametersList
Definition
ParametersList.h:52
cepgen::Steerable::steer
T steer(const std::string &key) const
Retrieve a parameters as previously steered.
Definition
Steerable.h:39
cepgen::Value
A scalar value with its uncertainty.
Definition
Value.h:26
cepgen::root::EventImporter
ROOT handler for an event tree import.
Definition
EventImporter.cpp:34
cepgen::root::EventImporter::EventImporter
EventImporter(const ParametersList ¶ms)
Definition
EventImporter.cpp:36
cepgen::root::EventImporter::description
static ParametersDescription description()
Definition
EventImporter.cpp:45
cepgen::root::EventImporter::operator>>
bool operator>>(Event &evt) override
Read the next event.
Definition
EventImporter.cpp:52
cepgen
Common namespace for this Monte Carlo generator.
Definition
CommandLineHandler.cpp:36
python_modules.Integrators.root_cfi.root
root
Definition
root_cfi.py:3
CepGenAddOns
ROOTWrapper
EventImporter.cpp
Generated on Mon Jul 29 2024 for CepGen by
1.9.7