cepgen is hosted by Hepforge, IPPP Durham
CepGen 1.2.5
Central exclusive processes event generator
Loading...
Searching...
No Matches
CommandLineHandler.cpp
Go to the documentation of this file.
1/*
2 * CepGen: a central exclusive processes event generator
3 * Copyright (C) 2020-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
22#include "CepGen/Event/Event.h"
30#include "CepGen/Physics/PDG.h"
33#include "CepGen/Utils/String.h"
35
36namespace cepgen {
37 namespace card {
39 class CommandLineHandler final : public Handler {
40 public:
42 explicit CommandLineHandler(const ParametersList& params) : Handler(params) {
43 if (const auto argv = steer<std::vector<std::string> >("args"); !argv.empty())
44 parseCommands(argv);
45 }
46
48 auto desc = Handler::description();
49 desc.setDescription("Command line configuration parser");
50 desc.add<std::vector<std::string> >("args", {}).setDescription("Collection of arguments to be parsed");
51 return desc;
52 }
53
54 inline CommandLineHandler& parseFile(const std::string& filename) override {
55 if (filename.empty())
56 throw CG_FATAL("CommandLineHandler") << "Empty filename to be parsed! Aborting.";
57 return parseCommands(utils::split(utils::readFile(filename), '\n'));
58 }
59
60 inline CommandLineHandler& parseCommands(const std::vector<std::string>& commands) override {
61 if (commands.empty())
62 return *this;
63 ParametersList pars;
64 for (const auto& cmd : commands)
65 pars.feed(cmd);
66 CG_INFO("CommandLineHandler") << "Arguments list: " << commands << " unpacked to:\n\t" << pars << ".";
67
68 //----- timer definition
69 if (pars.get<bool>("timer", false))
71
72 //----- logging definition
73 if (pars.has<int>("logging"))
75 else if (pars.has<ParametersList>("logging")) {
76 const auto& log = pars.get<ParametersList>("logging");
77 if (log.has<int>("level"))
79 if (log.has<std::string>("modules"))
80 utils::Logger::get().addExceptionRule(log.get<std::string>("modules"));
81 else if (log.has<std::vector<std::string> >("modules"))
82 for (const auto& mod : log.get<std::vector<std::string> >("modules"))
84 utils::Logger::get().setExtended(log.get<bool>("extended", false));
85 }
86
87 //----- PDG definition
88 auto pars_pdg = pars.get<ParametersList>("pdg");
89 for (const auto& id : pars_pdg.keys())
90 PDG::get().define(pars_pdg.get<ParticleProperties>(id));
91
92 //----- phase space definition
93 auto pars_kin = pars.get<ParametersList>("kinematics");
94
95 //----- process definition
96 if (auto proc = pars.get<ParametersList>("process"); !proc.empty()) {
98 proc = ParametersList(runParameters()->process().parameters()) + proc;
99 runParameters()->setProcess(ProcessFactory::get().build(proc));
100 if (proc.has<int>("mode"))
101 pars_kin.set<int>("mode", proc.get<int>("mode"));
102 }
103
104 if (!pars_kin.empty()) {
105 //----- set auxiliary information for phase space definition
106 if (pars_kin.has<int>("strfun"))
107 pars_kin
108 .set<ParametersList>(
109 "structureFunctions",
110 StructureFunctionsFactory::get().describeParameters(pars_kin.get<int>("strfun")).parameters())
111 .erase("strfun");
112 else if (pars_kin.has<std::string>("strfun"))
113 pars_kin
114 .set<ParametersList>("structureFunctions",
115 StructureFunctionsFactory::get()
116 .describeParameters(pars_kin.get<std::string>("strfun"))
117 .parameters())
118 .erase("strfun");
119 else if (pars_kin.has<ParametersList>("strfun"))
120 pars_kin.rename("strfun", "structureFunctions");
121 pars_kin.rename("formfac", "formFactors");
122
123 //----- get the kinematics as already defined in the process object and modify it accordingly
124 pars_kin = runParameters()->process().kinematics().parameters() + pars_kin;
126 }
127
128 //----- integration
129 pars.fill<ParametersList>("integrator", runParameters()->integrator());
130
131 //----- events generation
132 const auto& gen = pars.get<ParametersList>("generation");
133 runParameters()->generation().setMaxGen(gen.get<int>("ngen", runParameters()->generation().maxGen()));
134 if (gen.has<int>("nthreads"))
135 runParameters()->generation().setNumThreads(gen.get<int>("nthreads"));
136 if (gen.has<int>("nprn"))
137 runParameters()->generation().setPrintEvery(gen.get<int>("nprn"));
138 if (gen.has<int>("seed"))
139 runParameters()->integrator().set<int>("seed", gen.get<int>("seed"));
140
141 //----- event modification modules
142 if (const auto& mod = pars.get<ParametersList>("eventmod"); !mod.keys(true).empty()) {
143 runParameters()->addModifier(EventModifierFactory::get().build(mod));
144 runParameters()->eventModifiersSequence().rbegin()->get()->initialise(*runParameters());
145 }
146
147 //----- output modules definition
148 if (const auto& out = pars.get<ParametersList>("output"); !out.keys(true).empty())
149 runParameters()->addEventExporter(EventExporterFactory::get().build(out));
150 return *this;
151 }
152 };
153 } // namespace card
154} // namespace cepgen
#define REGISTER_CARD_HANDLER(name, obj)
Add a cards handler definition to the list of handled parsers.
#define CG_FATAL(mod)
Definition Exception.h:61
#define CG_INFO(mod)
Definition Message.h:216
void setParameters(const ParametersList &) override
Set module parameters.
const ParametersList & parameters() const override
List containing all parameters handled.
void define(const ParticleProperties &)
Add a new particle definition to the library.
Definition PDG.cpp:60
static PDG & get()
Retrieve a unique instance of this particles info collection.
Definition PDG.cpp:41
A description object for parameters collection.
bool has(const std::string &key) const
Check if a given parameter is handled in this list.
U getAs(const std::string &key, const U &def=default_arg< U >::get()) const
Get a recast parameter value.
ParametersList & feed(const std::string &)
Feed a control string to the list of parameters.
bool empty() const
Is the list empty?
T get(const std::string &key, const T &def=default_arg< T >::get()) const
Get a parameter value.
ParametersList & set(const std::string &, const T &)
Set a parameter value Set a recast parameter value.
const ParametersList & fill(const std::string &key, T &value) const
Fill a variable with the key content if exists.
void setPrintEvery(size_t print_every)
Set the events display frequency.
void setNumThreads(size_t nt)
Set the number of threads for the events generation.
void setMaxGen(size_t max_gen)
Set the maximal number of events to generate.
proc::Process & process()
Process object for cross-section computation/events generation.
void addModifier(std::unique_ptr< EventModifier >)
Add a new event modification algorithm to the sequence.
void setProcess(std::unique_ptr< proc::Process >)
Set a process configuration.
bool hasProcess() const
Is this parameters collection holding any physics process?
EventModifiersSequence & eventModifiersSequence()
List of event modification algos List of event modification algos.
ParametersList & integrator()
Integrator specific user-defined parameters.
void addEventExporter(std::unique_ptr< EventExporter >)
Set a new output module definition.
void setTimeKeeper(utils::TimeKeeper *)
Initialise the timekeeper instance.
Generation & generation()
Event generation parameters.
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.
CommandLineHandler(const ParametersList &params)
Cast command line arguments into a configuration word.
static ParametersDescription description()
CommandLineHandler & parseCommands(const std::vector< std::string > &commands) override
Read configuration from command strings.
CommandLineHandler & parseFile(const std::string &filename) override
Read configuration from steering card.
Base steering card module.
Definition Handler.h:31
const RunParameters * runParameters() const
Parsed runtime parameters.
Definition Handler.h:45
static ParametersDescription description()
Definition Handler.cpp:36
const Kinematics & kinematics() const
Constant reference to the process kinematics.
Definition Process.h:50
Level
Logging threshold for the output stream.
Definition Logger.h:44
void addExceptionRule(const std::string &rule)
Add a new rule to display exceptions/messages.
Definition Logger.cpp:37
static Logger & get(std::ostream *=nullptr)
Retrieve the running instance of the logger.
Definition Logger.cpp:31
void setLevel(Level level)
Set the logging threshold.
Definition Logger.h:52
void setExtended(bool ext=true)
Set the extended information flag.
Definition Logger.h:56
Collection of clocks to benchmark execution blocks.
Definition TimeKeeper.h:35
std::string readFile(const std::string &filename)
Read the content of a file into a string buffer.
std::vector< std::string > split(const std::string &str, char delim, bool trim)
Split a string according to a separation character.
Definition String.cpp:233
Common namespace for this Monte Carlo generator.
A collection of physics constants associated to a single particle.