cepgen
is hosted by
Hepforge
,
IPPP Durham
CepGen
1.2.5
Central exclusive processes event generator
Loading...
Searching...
No Matches
LHEFHandler.cpp
Go to the documentation of this file.
1
/*
2
* CepGen: a central exclusive processes event generator
3
* Copyright (C) 2016-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 <sstream>
20
21
#include "
CepGen/Core/Exception.h
"
22
#include "
CepGen/Event/Event.h
"
23
#include "
CepGen/EventFilter/EventExporter.h
"
24
#include "
CepGen/Modules/EventExporterFactory.h
"
25
#include "
CepGen/Utils/Caller.h
"
26
#include "
CepGen/Utils/Filesystem.h
"
27
#include "
CepGen/Utils/String.h
"
28
#include "
CepGen/Utils/Value.h
"
29
#include "
CepGenAddOns/Pythia8Wrapper/PythiaEventInterface.h
"
30
31
namespace
cepgen
{
32
namespace
pythia8
{
36
class
LHEFHandler
:
public
EventExporter
{
37
public
:
38
explicit
LHEFHandler
(
const
ParametersList
& params)
39
:
EventExporter
(params),
40
pythia_(new
Pythia8
::Pythia),
41
lhaevt_(new
Pythia8
::CepGenEvent),
42
compress_event_(
steer
<bool>(
"compress"
)),
43
filename_(
steer
<std::string>(
"filename"
)) {
44
if
(
utils::fileExtension
(filename_) ==
".gz"
) {
45
#ifdef GZIP_BIN
46
utils::replaceAll
(filename_,
".gz"
,
""
);
47
#else
48
CG_WARNING
(
"pythia8:LHEFHandler"
)
49
<<
"gzip compression requested, but the executable was not linked at Pythia8 wrapper compile time."
;
50
#endif
51
gzip_ =
true
;
52
}
53
{
54
auto
file_tmp = std::ofstream(filename_);
55
if
(!file_tmp.is_open())
56
throw
CG_FATAL
(
"pythia8:LHEFHandler"
)
57
<<
"Failed to open output filename \""
<< filename_ <<
"\" for writing!"
;
58
}
59
lhaevt_->openLHEF(filename_);
60
}
61
inline
~LHEFHandler
() {
62
if
(lhaevt_)
63
lhaevt_->closeLHEF(
false
);
// we do not want to rewrite the init block
64
if
(gzip_)
65
#ifdef GZIP_BIN
66
utils::Caller::call
({GZIP_BIN,
"-f"
, filename_});
67
#endif
68
}
69
70
inline
static
ParametersDescription
description
() {
71
auto
desc =
EventExporter::description
();
72
desc.setDescription(
"Pythia 8-based LHEF output module"
);
73
desc.add<
bool
>(
"compress"
,
true
);
74
desc.add<std::string>(
"filename"
,
"output.lhe"
).setDescription(
"Output filename"
);
75
return
desc;
76
}
77
78
inline
void
initialise
()
override
{
79
std::ostringstream oss_init;
80
oss_init <<
"<!--\n"
<<
banner
() <<
"\n-->"
;
81
oss_init << std::endl;
// LHEF is usually not as beautifully parsed as a standard XML...
82
// we're physicists, what do you expect?
83
lhaevt_->addComments(oss_init.str());
84
lhaevt_->initialise(
runParameters
());
85
#if PYTHIA_VERSION_INTEGER < 8300
86
pythia_->setLHAupPtr(lhaevt_.get());
87
#else
88
pythia_->setLHAupPtr(lhaevt_);
89
#endif
90
pythia_->settings.flag(
"ProcessLevel:all"
,
false
);
// we do not want Pythia to interfere...
91
pythia_->settings.flag(
"PartonLevel:all"
,
false
);
// we do not want Pythia to interfere...
92
pythia_->settings.flag(
"HadronLevel:all"
,
false
);
// we do not want Pythia to interfere...
93
pythia_->settings.mode(
"Beams:frameType"
, 5);
// LHEF event readout
94
pythia_->settings.mode(
"Next:numberCount"
, 0);
// remove some of the Pythia output
95
pythia_->init();
96
lhaevt_->initLHEF();
97
}
98
99
inline
bool
operator<<
(
const
Event
& ev)
override
{
100
lhaevt_->feedEvent(compress_event_ ? ev : ev.
compress
(),
101
Pythia8::CepGenEvent::Type::centralAndFullBeamRemnants
);
102
pythia_->next();
103
lhaevt_->eventLHEF();
104
return
true
;
105
}
106
inline
void
setCrossSection
(
const
Value
& cross_section)
override
{
107
lhaevt_->setCrossSection(0, cross_section, cross_section.
uncertainty
());
108
}
109
110
private
:
111
const
std::unique_ptr<Pythia8::Pythia> pythia_;
112
const
std::shared_ptr<Pythia8::CepGenEvent> lhaevt_;
113
const
bool
compress_event_;
114
std::string filename_;
115
bool
gzip_{
false
};
116
};
117
}
// namespace pythia8
118
}
// namespace cepgen
119
using
cepgen::pythia8::LHEFHandler
;
120
REGISTER_EXPORTER
(
"lhef"
,
LHEFHandler
);
Caller.h
EventExporterFactory.h
REGISTER_EXPORTER
#define REGISTER_EXPORTER(name, obj)
Add a generic export module definition to the factory.
Definition
EventExporterFactory.h:25
EventExporter.h
Event.h
Exception.h
CG_FATAL
#define CG_FATAL(mod)
Definition
Exception.h:61
Filesystem.h
CG_WARNING
#define CG_WARNING(mod)
Definition
Message.h:228
PythiaEventInterface.h
String.h
Value.h
Pythia8::CepGenEvent::Type::centralAndFullBeamRemnants
@ centralAndFullBeamRemnants
include dissociated beam remnants and central system
cepgen::EventExporter
Output format handler for events export.
Definition
EventExporter.h:30
cepgen::EventExporter::banner
std::string banner(const std::string &prep="") const
Print a banner containing runtime parameters.
Definition
EventExporter.cpp:32
cepgen::EventHandler::runParameters
const RunParameters & runParameters() const
List of run parameters.
Definition
EventHandler.cpp:44
cepgen::Event
Container for the information on the in- and outgoing particles' kinematics.
Definition
Event.h:28
cepgen::Event::compress
Event compress() const
Compress the event record.
Definition
Event.cpp:115
cepgen::ParametersDescription
A description object for parameters collection.
Definition
ParametersDescription.h:26
cepgen::ParametersList
Definition
ParametersList.h:52
cepgen::Steerable::description
static ParametersDescription description()
Description of all object parameters.
Definition
Steerable.cpp:42
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::Value::uncertainty
double uncertainty() const
Absolute uncertainty around the central value.
Definition
Value.h:35
cepgen::pythia8::LHEFHandler
Pythia8 handler for the LHE file output.
Definition
LHEFHandler.cpp:36
cepgen::pythia8::LHEFHandler::LHEFHandler
LHEFHandler(const ParametersList ¶ms)
Definition
LHEFHandler.cpp:38
cepgen::pythia8::LHEFHandler::~LHEFHandler
~LHEFHandler()
Definition
LHEFHandler.cpp:61
cepgen::pythia8::LHEFHandler::description
static ParametersDescription description()
Definition
LHEFHandler.cpp:70
cepgen::pythia8::LHEFHandler::setCrossSection
void setCrossSection(const Value &cross_section) override
Specify the cross section value, in pb.
Definition
LHEFHandler.cpp:106
cepgen::pythia8::LHEFHandler::operator<<
bool operator<<(const Event &ev) override
Writer operator.
Definition
LHEFHandler.cpp:99
cepgen::pythia8::LHEFHandler::initialise
void initialise() override
Definition
LHEFHandler.cpp:78
cepgen::utils::Caller::call
static std::string call(const std::vector< std::string > &commands)
Start a logged call command.
Definition
Caller.cpp:45
Pythia8
Definition
PythiaEventInterface.cpp:26
cepgen::utils::replaceAll
size_t replaceAll(std::string &str, const std::string &from, const std::string &to)
Replace all occurrences of a text by another.
Definition
String.cpp:118
cepgen::utils::fileExtension
std::string fileExtension(const std::string &file)
Small utility to retrieve the extension of a filename.
Definition
Filesystem.cpp:29
cepgen
Common namespace for this Monte Carlo generator.
Definition
CommandLineHandler.cpp:36
python.Config.Hadronisation.pythia8_cfi.pythia8
pythia8
Definition
pythia8_cfi.py:13
CepGenAddOns
Pythia8Wrapper
LHEFHandler.cpp
Generated on Mon Jul 29 2024 for CepGen by
1.9.7