MaCh3  2.4.2
Reference Guide
parameters.cpp
Go to the documentation of this file.
1 // pybind includes
2 #include <pybind11/pybind11.h>
3 #include <pybind11/stl.h>
4 #include <pybind11/numpy.h>
5 // MaCh3 includes
8 
9 namespace py = pybind11;
10 
13 public:
14  /* Inherit the constructors */
16 
17  /* Trampoline (need one for each virtual function) */
18  double GetLikelihood() override {
19  PYBIND11_OVERRIDE_NAME(
20  double, /* Return type */
21  ParameterHandlerBase, /* Parent class */
22  "get_likelihood", /* Name in python*/
23  GetLikelihood /* Name of function in C++ (must match Python name) */
24  );
25  }
26 
27  void ProposeStep() override {
28  PYBIND11_OVERRIDE_NAME(
29  void, /* Return type */
30  ParameterHandlerBase, /* Parent class */
31  "propose_step", /* Name in python*/
32  ProposeStep /* Name of function in C++ (must match Python name) */
33  );
34  }
35 };
36 
37 
38 void initParameters(py::module &m) {
39  auto m_parameters = m.def_submodule("parameters");
40  m_parameters.doc() =
41  "This is a Python binding of MaCh3s C++ parameters library.";
42 
43 
44  // Bind the systematic type enum that lets us set different types of systematics
45  py::enum_<SystType>(m_parameters, "SystematicType")
46  .value("Normalisation", SystType::kNorm)
47  .value("Spline", SystType::kSpline)
48  .value("Functional", SystType::kFunc)
49  .value("N_Systematic_Types", SystType::kSystTypes);
50 
51 
52  py::class_<ParameterHandlerBase, PyParameterHandlerBase /* <--- trampoline*/>(m_parameters, "ParameterHandlerBase")
53  .def(
54  py::init<const std::vector<std::string>&, const char *, double, int, int>(),
55  "Construct a parameters object from a set of yaml files that define the systematic parameters \n\
56  :param yaml_files: The name of the yaml file to initialise from. \n\
57  :param name: the name of this ParameterHandler object. \n\
58  :param threshold: threshold PCA threshold from 0 to 1. Default is -1 and means no PCA. \n\
59  :param first_PCA_par: FirstPCAdpar First PCA parameter that will be decomposed. \n\
60  :param last_PCA_par: LastPCAdpar First PCA parameter that will be decomposed.",
61  py::arg("yaml_files"),
62  py::arg("name"),
63  py::arg("threshold") = -1.0,
64  py::arg("firs_PCA_par") = -999,
65  py::arg("last_PCA_par") = -999
66  )
67 
68  .def(
69  "calculate_likelihood",
71  "Calculate penalty term based on inverted covariance matrix."
72  )
73 
74  .def(
75  "get_internal_par_name",
76  [](ParameterHandlerBase &self, int index)
77  {
78  // do this to disambiguate between the std::string and const char* version of this fn
79  std::string ret;
80  ret = self.GetParName(index);
81  return ret;
82  },
83  "Get the internally used name of this parameter. \n\
84  :param index: The global index of the parameter",
85  py::arg("index")
86  )
87 
88  .def(
89  "get_fancy_par_name",
90  [](ParameterHandlerBase &self, int index)
91  {
92  // do this to disambiguate between the std::string and const char* version of this fn
93  std::string ret;
94  ret = self.GetParFancyName(index);
95  return ret;
96  },
97  "Get the name of this parameter. \n\
98  :param index: The global index of the parameter",
99  py::arg("index")
100  )
101 
102  .def(
103  "get_n_pars",
105  "Get the number of parameters that this ParameterHandler object knows about."
106  )
107 
108  .def(
109  "propose_step",
111  "Propose a step based on the covariances. Also feel free to overwrite if you want something more funky."
112  )
113 
114  .def(
115  "get_proposal_array",
116  [](ParameterHandlerBase &self)
117  {
118  return py::memoryview::from_buffer<double>(
119  self.GetParPropVec().data(), // the data pointer
120  {self.GetNParameters()}, // shape
121  {sizeof(double)} // shape
122  );
123  },
124  "Bind a python array to the parameter proposal values for this ParameterHandler object. \n\
125  This allows you to set e.g. a numpy array to 'track' the parameter proposal values. You could either use this to directly set the proposals, or to just read the values proposed by e.g. throw_par_prop() \n\
126  :warning: This should be set *AFTER* all of the parameters have been read in from the config file as it resizes the array to fit the number of parameters. \n\
127  :param array: This is the array that will be set. Size and contents don't matter as it will be changed to fit the parameters. "
128  )
129 
130 
131  ; // End of ParameterHandlerBase binding
132 
133 
134  py::class_<ParameterHandlerGeneric, ParameterHandlerBase /* <--- trampoline*/>(m_parameters, "ParameterHandlerGeneric")
135  .def(
136  py::init<const std::vector<std::string>&, const char *, double, int, int>(),
137  "Construct a systematic ParameterHandler object from a set of yaml files that define the systematic parameters \n\
138  :param yaml_files: The name of the yaml file to initialise from. \n\
139  :param name: the name of this ParameterHandler object. \n\
140  :param threshold: threshold PCA threshold from 0 to 1. Default is -1 and means no PCA. \n\
141  :param first_PCA_par: FirstPCAdpar First PCA parameter that will be decomposed. \n\
142  :param last_PCA_par: LastPCAdpar First PCA parameter that will be decomposed.",
143  py::arg("yaml_files"),
144  py::arg("name") = "xsec_cov",
145  py::arg("threshold") = -1.0,
146  py::arg("firs_PCA_par") = -999,
147  py::arg("last_PCA_par") = -999
148  )
149 
150  .def(
151  "get_par_type",
153  "Get what type of systematic this parameters is (see :py:class:`pyMaCh3.m_parameters.SystematicType` for possible types). \n\
154  :param index: The global index of the parameter",
155  py::arg("index")
156  )
157 
158  .def(
159  "get_par_spline_type",
161  "Get what type of spline this parameter is set to use (assuming that it is a spline type parameter). \n\
162  :param index: The index of the spline parameter",
163  py::arg("index")
164  )
165 
166  .def(
167  "get_par_spline_name",
169  "Get the name of the spline associated with a spline parameter. This is generally what it is called in input spline files and can in principle be different to the parameters name. \n\
170  :param index: The index of the spline parameter",
171  py::arg("index")
172  )
173 
174  ; // End of ParameterHandlerGeneric binding
175 }
@ kNorm
For normalisation parameters.
@ kSpline
For splined parameters (1D)
@ kSystTypes
This only enumerates.
@ kFunc
For functional parameters.
Base class responsible for handling of systematic error parameters. Capable of using PCA or using ada...
virtual void ProposeStep()
Generate a new proposed state.
int GetNParameters() const
Get number of params which will be different depending if using Eigen decomposition or not.
double CalcLikelihood() const _noexcept_
Calc penalty term based on inverted covariance matrix.
ParameterHandlerBase(const std::vector< std::string > &YAMLFile, std::string name, double threshold=-1, int FirstPCAdpar=-999, int LastPCAdpar=-999)
ETA - constructor for a YAML file.
Class responsible for handling of systematic error parameters with different types defined in the con...
SystType GetParamType(const int i) const
Returns enum describing our param type.
std::string GetParSplineName(const int i) const
Get the name of the spline associated with the spline at index i.
SplineInterpolation GetParSplineInterpolation(const int i) const
Get interpolation type for a given parameter.
EW: As ParameterHandlerBase is an abstract base class we have to do some gymnastics to get it to get ...
Definition: parameters.cpp:12
double GetLikelihood() override
Return CalcLikelihood if some params were thrown out of boundary return LARGE_LOGL
Definition: parameters.cpp:18
void ProposeStep() override
Generate a new proposed state.
Definition: parameters.cpp:27
void initParameters(py::module &m)
Definition: parameters.cpp:38