cepgen is hosted by Hepforge, IPPP Durham
CepGen 1.2.5
Central exclusive processes event generator
Loading...
Searching...
No Matches
ProcessVariablesAnalyser.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
21#include "CepGen/Utils/Drawer.h"
25
26namespace cepgen {
27 namespace utils {
29 : SteeredObject(params), proc_(proc) {
30 for (const auto& var : proc_.mapped_variables_)
31 if (auto hist = steer<ParametersList>(var.name); !hist.empty())
32 hists_.insert(std::make_pair(var.name, new Hist1D(hist.set<std::string>("name", var.name))));
33 else
34 hists_.insert(std::make_pair(var.name,
36 .set<std::string>("name", var.name)
37 .set<int>("nbinsX", 50)
38 .set<Limits>("xrange", var.limits))));
39 }
40
41 void ProcessVariablesAnalyser::feed(double weight) {
42 for (const auto& var : proc_.mapped_variables_)
43 if (hists_.count(var.name))
44 hists_.at(var.name)->fill(var.value, weight);
45 }
46
48 auto drawer = DrawerFactory::get().build(steer<ParametersList>("drawer"));
49 for (const auto& var : hists_)
50 drawer->draw(*var.second);
51 }
52
54 auto desc = ParametersDescription();
55 ParametersDescription hist_desc;
56 hist_desc.add<std::vector<double> >("xbins", {}).setDescription("x-axis bins definition");
57 hist_desc.add<int>("nbinsX", 25).setDescription("Bins multiplicity for x-axis");
58 hist_desc.add<Limits>("xrange", Limits{0., 1.}).setDescription("Minimum-maximum range for x-axis");
59 desc.addParametersDescriptionVector("histVariables", hist_desc, {}).setDescription("Histogram definition");
60 desc.add<ParametersDescription>("drawer", DrawerFactory::get().describeParameters("root"));
61 return desc;
62 }
63 } // namespace utils
64} // namespace cepgen
Validity interval for a variable.
Definition Limits.h:28
A description object for parameters collection.
ParametersDescription & add(const std::string &name, const T &def)
Add the description to a new parameter.
ParametersList & set(const std::string &, const T &)
Set a parameter value Set a recast parameter value.
Base user-steerable object.
Class template to define any process to compute using this MC integrator/events generator.
Definition Process.h:34
1D histogram container
Definition Histogram.h:72
ProcessVariablesAnalyser(const proc::Process &, const ParametersList &)
Common namespace for this Monte Carlo generator.