cepgen is hosted by Hepforge, IPPP Durham
CepGen 1.2.3
A generic central exclusive processes event generator
Loading...
Searching...
No Matches
RunParameters.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 <iomanip>
20
24#include "CepGen/Event/Event.h"
33#include "CepGen/Physics/PDG.h"
36#include "CepGen/Utils/String.h"
38
39namespace cepgen {
42 integrator_(steer<ParametersList>("integrator")),
43 generation_(steer<ParametersList>("generation")) {}
44
46 : SteeredObject(param),
47 process_(std::move(param.process_)),
48 evt_modifiers_(std::move(param.evt_modifiers_)),
49 evt_exporters_(std::move(param.evt_exporters_)),
50 taming_functions_(std::move(param.taming_functions_)),
51 total_gen_time_(param.total_gen_time_),
52 num_gen_events_(param.num_gen_events_),
53 integrator_(param.integrator_),
54 generation_(param.generation_),
55 tmr_(std::move(param.tmr_)) {}
56
58 : SteeredObject(param),
59 total_gen_time_(param.total_gen_time_),
60 num_gen_events_(param.num_gen_events_),
61 integrator_(param.integrator_),
62 generation_(param.generation_) {}
63
64 RunParameters::~RunParameters() {} // required for unique_ptr initialisation!
65
67 process_ = std::move(param.process_);
68 evt_modifiers_ = std::move(param.evt_modifiers_);
69 evt_exporters_ = std::move(param.evt_exporters_);
70 taming_functions_ = std::move(param.taming_functions_);
71 total_gen_time_ = param.total_gen_time_;
72 num_gen_events_ = param.num_gen_events_;
73 integrator_ = param.integrator_;
74 generation_ = param.generation_;
75 tmr_ = std::move(param.tmr_);
76 return *this;
77 }
78
80 // prepare the event modifications algorithms for event generation
81 for (auto& mod : evt_modifiers_)
82 mod->initialise(*this);
83 // prepare the output modules for event generation
84 for (auto& mod : evt_exporters_)
85 mod->initialise(*this);
86 }
87
89 if (tmr_)
90 tmr_->clear();
91 CG_TICKER(tmr_.get());
92
93 //--- clear the run statistics
94 total_gen_time_ = 0.;
95 num_gen_events_ = 0ul;
96 }
97
98 void RunParameters::setTimeKeeper(utils::TimeKeeper* kpr) { tmr_.reset(kpr); }
99
100 void RunParameters::addGenerationTime(double gen_time) {
101 total_gen_time_ += gen_time;
102 num_gen_events_++;
103 }
104
105 proc::Process& RunParameters::process() { return *process_.get(); }
106
107 const proc::Process& RunParameters::process() const { return *process_.get(); }
108
109 std::string RunParameters::processName() const {
110 if (!process_)
111 return "no process";
112 return process_->name();
113 }
114
115 void RunParameters::clearProcess() { process_.release(); }
116
117 void RunParameters::setProcess(std::unique_ptr<proc::Process> proc) { process_ = std::move(proc); }
118
120 if (!proc)
121 throw CG_FATAL("RunParameters") << "Trying to clone an invalid process!";
122 process_.reset(proc);
123 }
124
126 if (!process_)
127 throw CG_FATAL("RunParameters") << "Process must be defined before its kinematics is retrieved!";
128 return process_->kinematics();
129 }
130
131 EventModifier& RunParameters::eventModifier(size_t i) { return *evt_modifiers_.at(i); }
132
133 void RunParameters::clearEventModifiersSequence() { evt_modifiers_.clear(); }
134
135 void RunParameters::addModifier(std::unique_ptr<EventModifier> mod) { evt_modifiers_.emplace_back(std::move(mod)); }
136
138 evt_modifiers_.emplace_back(std::move(std::unique_ptr<EventModifier>(mod)));
139 }
140
141 EventExporter& RunParameters::eventExporter(size_t i) { return *evt_exporters_.at(i); }
142
143 void RunParameters::clearEventExportersSequence() { evt_exporters_.clear(); }
144
145 void RunParameters::addEventExporter(std::unique_ptr<EventExporter> mod) {
146 evt_exporters_.emplace_back(std::move(mod));
147 }
148
150 evt_exporters_.emplace_back(std::unique_ptr<EventExporter>(mod));
151 }
152
153 void RunParameters::addTamingFunction(std::unique_ptr<utils::Functional> fct) {
154 taming_functions_.emplace_back(std::move(fct));
155 }
156
157 std::ostream& operator<<(std::ostream& os, const RunParameters& param) {
158 const int wb = 90, wt = 33;
159
160 os << std::left << "\n"
161 << std::setfill('_') << std::setw(wb + 3) << "_/¯ RUN INFORMATION ¯\\_" << std::setfill(' ') << "\n\n";
162 if (param.process_) {
163 const auto& proc_params = param.process().parameters();
164 os << std::setw(wt) << "Process to generate:"
165 << utils::boldify(ProcessFactory::get().describeParameters(proc_params).description()) << "\n";
166 for (const auto& key : proc_params.keys(false)) {
167 if (key == "kinematics" || key == "partonFluxes" || key == "ktFluxes") // these are shown below
168 continue;
169 os << std::setw(wt) << "" << key << ": "
170 << (proc_params.has<ParametersList>(key) ? proc_params.get<ParametersList>(key).print(true)
171 : proc_params.getString(key))
172 << "\n";
173 }
174 }
175 if (!param.evt_modifiers_.empty() || !param.evt_exporters_.empty() || !param.taming_functions_.empty())
176 os << "\n"
177 << std::setfill('-') << std::setw(wb + 6) << utils::boldify(" Event treatment ") << std::setfill(' ')
178 << "\n\n";
179 if (!param.evt_modifiers_.empty()) {
180 std::string mod_name = utils::s("Event modifier", param.evt_modifiers_.size(), false), sep;
181 for (const auto& mod : param.evt_modifiers_)
182 os << std::setw(wt) << mod_name << sep << utils::boldify(mod->name()) << "\n", sep = "+ ", mod_name.clear();
183 os << "\n";
184 }
185 if (!param.evt_exporters_.empty()) {
186 os << utils::s("Output module", param.evt_exporters_.size(), false);
187 for (const auto& mod : param.evt_exporters_)
188 os << "\n\t*) " << EventExporterFactory::get().describeParameters(mod->name(), mod->parameters()).describe(1);
189 }
190 if (!param.taming_functions_.empty()) {
191 os << std::setw(wt) << utils::s("Taming function", param.taming_functions_.size(), false) << "\n";
192 for (const auto& tf : param.taming_functions_)
193 os << std::setw(wt) << "" << tf->variables().at(0) << ": " << tf->expression() << "\n";
194 }
195 os << "\n\n"
196 << std::setfill('-') << std::setw(wb + 6) << utils::boldify(" Integration/generation parameters ")
197 << std::setfill(' ') << "\n\n"
198 << std::setw(wt) << "Integration" << utils::boldify(param.integrator_.name("N/A")) << "\n";
199 for (const auto& key : param.integrator_.keys(false))
200 os << std::setw(wt) << "" << key << ": " << param.integrator_.getString(key) << "\n";
201 os << std::setw(wt) << "Event generation? " << utils::yesno(param.generation_.enabled()) << "\n"
202 << std::setw(wt) << "Number of events to generate" << utils::boldify(param.generation_.maxGen()) << "\n"
203 << std::setw(wt) << "Generator worker"
204 << param.generation_.parameters().get<ParametersList>("worker").print(true) << "\n";
205 if (param.generation_.numThreads() > 1)
206 os << std::setw(wt) << "Number of threads" << param.generation_.numThreads() << "\n";
207 os << std::setw(wt) << "Number of points to try per bin" << param.generation_.numPoints() << "\n"
208 << std::setw(wt) << "Verbosity level " << utils::Logger::get().level() << "\n";
209 const auto& kin = param.process().kinematics();
210 const auto& beams = kin.incomingBeams();
211 os << "\n"
212 << std::setfill('_') << std::setw(wb + 3) << "_/¯ EVENTS KINEMATICS ¯\\_" << std::setfill(' ') << "\n\n"
213 << std::setw(wt) << "Incoming particles" << beams.positive() << ",\n"
214 << std::setw(wt) << "" << beams.negative() << "\n"
215 << std::setw(wt) << "C.m. energy (GeV)" << utils::format("%g", beams.sqrtS()) << "\n";
216 if (beams.mode() != mode::Kinematics::ElasticElastic)
217 os << std::setw(wt) << "Structure functions"
219 StructureFunctionsFactory::get().describeParameters(beams.structureFunctions()).description())
220 << "\n"
221 << std::setw(wt) << "" << beams.structureFunctions().print(true) << "\n ";
222
223 // helper to print cuts list for a given category
224 const auto dump_cuts = [&os](const auto& obj) {
225 for (const auto& lim : obj.parameters().template keysOf<Limits>())
226 if (const auto& limit = obj.parameters().template get<Limits>(lim); limit.valid() && obj.description().has(lim))
227 os << std::setw(wt) << obj.description().get(lim).description() << limit << "\n";
228 for (const auto& vlim : obj.parameters().template keysOf<std::vector<Limits> >())
229 if (const auto& limit = obj.parameters().template get<std::vector<Limits> >(vlim); obj.description().has(vlim))
230 os << std::setw(wt) << obj.description().get(vlim).description() << utils::repr(limit, " and ") << "\n";
231 };
232
233 os << "\n"
234 << std::setfill('-') << std::setw(wb + 6) << utils::boldify(" Incoming partons ") << std::setfill(' ') << "\n\n";
235 const auto& cuts = kin.cuts();
236 dump_cuts(cuts.initial);
237 os << "\n"
238 << std::setfill('-') << std::setw(wb + 6) << utils::boldify(" Outgoing central system ") << std::setfill(' ')
239 << "\n\n";
240 if (!kin.minimumFinalState().empty()) {
241 os << std::setw(wt) << "Minimum final state";
242 std::string sep;
243 for (const auto& part : kin.minimumFinalState())
244 os << sep << (PDG::Id)part, sep = ", ";
245 os << "\n";
246 }
247 dump_cuts(cuts.central);
248 if (cuts.central_particles.size() > 0) {
249 os << std::setw(wt) << utils::boldify(">>> per-particle cuts:") << "\n";
250 for (const auto& part_per_lim : cuts.central_particles) {
251 os << " * all single " << std::setw(wt - 3) << (PDG::Id)part_per_lim.first << "\n";
252 for (const auto& lim : part_per_lim.second.parameters().keysOf<Limits>()) {
253 const auto& limit = part_per_lim.second.parameters().get<Limits>(lim);
254 if (limit.valid())
255 os << " - " << std::setw(wt - 5) << cuts::Central::description().get(lim).description() << limit << "\n";
256 }
257 }
258 }
259 os << "\n";
260 os << std::setfill('-') << std::setw(wb + 6) << utils::boldify(" Proton / remnants ") << std::setfill(' ')
261 << "\n\n";
262 dump_cuts(cuts.remnants);
263 return os << "\n"
264 << std::setfill('_') << std::setw(wb) << ""
265 << "\n";
266 }
267
269 auto desc = ParametersDescription();
270 desc.add<ParametersDescription>("integrator", IntegratorFactory::get().describeParameters("Vegas"));
272 return desc;
273 }
274
275 //-----------------------------------------------------------------------------------------------
276
278 (*this)
279 .add("maxgen", max_gen_)
280 .add("printEvery", gen_print_every_)
281 .add("targetLumi", target_lumi_)
282 .add("symmetrise", symmetrise_)
283 .add("numThreads", num_threads_)
284 .add("numPoints", num_points_);
285 }
286
288 auto desc = ParametersDescription();
289 desc.add<ParametersDescription>("worker", GeneratorWorkerFactory::get().describeParameters("grid_optimised"))
290 .setDescription("type of generator worker to use for event generation");
291 desc.add<int>("maxgen", 0).setDescription("Number of events to generate");
292 desc.add<int>("printEvery", 10000).setDescription("Printing frequency for the events content");
293 desc.add<double>("targetLumi", -1.).setDescription("Target luminosity (in pb-1) to reach for this run");
294 desc.add<bool>("symmetrise", false).setDescription("Are events to be symmetrised wrt beam collinear axis");
295 desc.add<int>("numThreads", 1).setDescription("Number of threads to use for event generation");
296 desc.add<int>("numPoints", 100);
297 return desc;
298 }
299} // namespace cepgen
#define CG_FATAL(mod)
Definition Exception.h:61
#define CG_TICKER(tmr)
Definition TimeKeeper.h:30
Output format handler for events export.
Class template to interface (external/internal) events modification algorithms.
const Beam & positive() const
Reference to the positive-z beam information.
List of kinematic constraints to apply on the process phase space.
Definition Kinematics.h:27
IncomingBeams & incomingBeams()
Beam/primary particle kinematics.
Definition Kinematics.h:36
Validity interval for a variable.
Definition Limits.h:28
A class-in-the-middle PDG identifier for printout operations.
Definition PDG.h:55
A description object for parameters collection.
ParametersDescription & add(const std::string &name, const T &def)
Add the description to a new parameter.
const std::string & description() const
Description of this parameter (or parameters collection)
const ParametersDescription & get(const std::string &) const
Get the description of a sub-object.
std::string name(const std::string &def="") const
Retrieve the module name if any.
const ParametersList & print(std::ostream &) const
Debugging-like printout of a parameters container.
std::vector< std::string > keys(bool name_key=true) const
T get(const std::string &key, const T &def=default_arg< T >::get()) const
Get a parameter value.
std::string getString(const std::string &key, bool wrap=false) const
Get a string-converted version of a value.
size_t numPoints() const
Number of points to "shoot" in each integration bin.
size_t numThreads() const
Number of threads to perform the events generation.
size_t maxGen() const
Maximal number of events to generate.
static ParametersDescription description()
Generation(const ParametersList &=ParametersList())
Build a generation parameters collection from a user input.
bool enabled() const
Are we generating events?
List of parameters used to start and run the simulation job.
const Kinematics & kinematics() const
Events kinematics for phase space definition.
void initialiseModules()
Initialise the event handling modules for an event generation.
proc::Process & process()
Process object for cross-section computation/events generation.
EventExporter & eventExporter(size_t)
Output module.
void clearEventModifiersSequence()
Remove all event modifiers from sequence.
void addModifier(std::unique_ptr< EventModifier >)
Add a new event modification algorithm to the sequence.
EventModifier & eventModifier(size_t)
Event modification algorithm.
RunParameters & operator=(RunParameters)
Assignment operator.
void prepareRun()
Reset total generation time and number of events generated for this run, prepare kinematics.
void clearEventExportersSequence()
Remove all output modules from sequence.
void setProcess(std::unique_ptr< proc::Process >)
Set a process configuration.
std::string processName() const
Name of the process considered.
void addTamingFunction(std::unique_ptr< utils::Functional >)
Set a new taming function definition.
void addGenerationTime(double gen_time)
Add a new timing into the total generation time.
void clearProcess()
Remove the process pointer.
static ParametersDescription description()
void addEventExporter(std::unique_ptr< EventExporter >)
Set a new output module definition.
void setTimeKeeper(utils::TimeKeeper *)
Initialise the timekeeper instance.
Base user-steerable object.
const ParametersList & parameters() const override
Module user-defined parameters.
Class template to define any process to compute using this MC integrator/events generator.
Definition Process.h:34
const Kinematics & kinematics() const
Constant reference to the process kinematics.
Definition Process.h:50
static Logger & get(std::ostream *=nullptr)
Retrieve the running instance of the logger.
Definition Logger.cpp:31
Level level() const
Logging threshold.
Definition Logger.h:50
Collection of clocks to benchmark execution blocks.
Definition TimeKeeper.h:35
@ ElasticElastic
proton-proton elastic case
std::string format(const std::string &fmt, Args... args)
Format a string using a printf style format descriptor.
Definition String.h:59
std::string yesno(bool test)
Human-readable boolean printout Boldify a string for TTY-type output streams.
Definition String.cpp:45
std::string s(const std::string &word, float num, bool show_number)
Add a trailing "s" when needed.
Definition String.cpp:219
std::string boldify(std::string str)
String implementation of the boldification procedure.
Definition String.cpp:49
std::string repr(const std::vector< T > &vec, const std::function< std::string(const T &)> &printer, const std::string &sep=",")
Helper to print a vector.
Definition String.h:154
Common namespace for this Monte Carlo generator.
std::ostream & operator<<(std::ostream &os, const Exception::Type &type)
Definition Exception.cpp:59
static ParametersDescription description()
Definition Cuts.cpp:77