MaCh3  2.4.2
Reference Guide
Functions
pyMaCh3.cpp File Reference
#include <pybind11/pybind11.h>
#include <pybind11/stl.h>
Include dependency graph for pyMaCh3.cpp:

Go to the source code of this file.

Functions

void initPlotting (py::module &)
 
void initFitters (py::module &)
 
void initSamples (py::module &)
 
void initManager (py::module &)
 
void initParameters (py::module &)
 
void initSplines (py::module &)
 
 PYBIND11_MODULE (_pyMaCh3, m)
 

Detailed Description

Author
Ewan Miller

Definition in file pyMaCh3.cpp.

Function Documentation

◆ initFitters()

void initFitters ( py::module &  m)

Definition at line 47 of file fitters.cpp.

47  {
48  auto m_fitters = m.def_submodule("fitters");
49  m_fitters.doc() =
50  "This is a Python binding of MaCh3s C++ fitters library.";
51 
52 
53  py::class_<FitterBase, PyFitterBase /* <--- trampoline*/>(m_fitters, "FitterBase")
54  .def(py::init<Manager* const>())
55 
56  .def(
57  "run",
59  "The implementation of the fitter, you should override this with your own desired fitting algorithm"
60  )
61 
62  .def(
63  "get_name",
65  " The name of the algorithm, you should override this with something like:: \n"
66  "\n"
67  " return 'mySuperCoolAlgoName' \n"
68  )
69 
70  .def(
71  "run_LLH_scan",
73  "Perform a 1D likelihood scan"
74  )
75 
76  .def(
77  "get_step_scale_from_LLH_scan",
79  "LLH scan is good first estimate of step scale, this will get the rough estimates for the step scales based on running an LLH scan"
80  )
81 
82  .def(
83  "run_2d_LLH_scan",
85  " Perform a 2D likelihood scan. \n"
86  " :param warning: This operation may take a significant amount of time, especially for complex models."
87  )
88 
89  .def(
90  "run_sigma_var",
92  " Perform a 2D and 1D sigma var for all samples. \n"
93  " :param warning: Code uses TH2Poly"
94  )
95 
96  .def(
97  "drag_race",
99  " Calculates the required time for each sample or covariance object in a drag race simulation. Inspired by Dan's feature \n"
100  " :param NLaps: number of laps, every part of Fitter will be tested with given number of laps and you will get total and average time",
101  py::arg("NLaps") = 100
102  )
103 
104  // stuff for registering other objects with the fitter
105 
106  .def(
107  "add_sample_handler",
109  " This function adds a sample handler object to the analysis framework. The sample handler object will be utilized in fitting procedures or likelihood scans. \n"
110  " :param sample: A sample handler object derived from SampleHandlerBase. ",
111  py::arg("sample")
112  )
113 
114  .def(
115  "add_syst_object",
117  " This function adds a Covariance object to the analysis framework. The Covariance object will be utilized in fitting procedures or likelihood scans. \n"
118  " :param cov: A Covariance object derived from ParameterHandlerBase. ",
119  py::arg("cov")
120  )
121 
122  ; // End of FitterBase class binding
123 
124  py::class_<MR2T2, FitterBase>(m_fitters, "mcmc")
125  .def(py::init<Manager *const>())
126 
127  .def(
128  "set_chain_length",
130  "Set how long chain should be.",
131  py::arg("length")); // end of MCMC class binding
132 
133  py::class_<DelayedMR2T2, FitterBase>(m_fitters, "DelayedMCMC")
134  .def(py::init<Manager *const>())
135 
136  .def(
137  "set_chain_length",
139  "Set how long chain should be.",
140  py::arg("length")); // end of MCMC class binding
141 
142  py::class_<LikelihoodFit, PyLikelihoodFit /* <--- trampoline*/, FitterBase>(m_fitters, "LikelihoodFit")
143  .def(py::init<Manager* const>())
144 
145  .def(
146  "caluclate_chi2",
147  [](LikelihoodFit &self, const std::vector<double> &parameterVals)
148  {
149  return self.CalcChi2(parameterVals.data());
150  },
151  "Get the Chi2 calculation over all included samples and syst objects for the specified parameter_values \n\
152  :param parameter_valuse: The location to evaluate the chi2 at.",
153  py::arg("parameter_values")
154  )
155 
156  .def(
157  "get_n_params",
159  "Get The total number of parameters across all known covariance objects associated with this LikelihoodFit object."
160  )
161  ; // end of LikelihoodFit class binding
162 
163  py::class_<MinuitFit, LikelihoodFit>(m_fitters, "MinuitFit")
164  .def(py::init<Manager* const>())
165 
166  ; // end of MinuitFit class binding
167 
168  py::class_<PSO, LikelihoodFit>(m_fitters, "PSO")
169  .def(py::init<Manager* const>())
170 
171  .def(
172  "init",
173  &PSO::init,
174  "Initialise the fitter"
175  )
176 
177  ; // end of PSO class binding
178 }
Base class for implementing fitting algorithms.
Definition: FitterBase.h:26
void RunLLHScan()
Perform a 1D likelihood scan.
Definition: FitterBase.cpp:622
void AddSystObj(ParameterHandlerBase *cov)
This function adds a Covariance object to the analysis framework. The Covariance object will be utili...
Definition: FitterBase.cpp:298
std::string GetName() const
Get name of class.
Definition: FitterBase.h:70
void AddSampleHandler(SampleHandlerBase *sample)
This function adds a sample PDF object to the analysis framework. The sample PDF object will be utili...
Definition: FitterBase.cpp:262
virtual void RunMCMC()=0
The specific fitting algorithm implemented in this function depends on the derived class....
void GetStepScaleBasedOnLLHScan()
LLH scan is good first estimate of step scale.
Definition: FitterBase.cpp:887
void RunSigmaVar()
Perform a 1D/2D sigma var for all samples.
void DragRace(const int NLaps=100)
Calculates the required time for each sample or covariance object in a drag race simulation....
Definition: FitterBase.cpp:461
void Run2DLLHScan()
Perform a 2D likelihood scan.
Definition: FitterBase.cpp:936
Implementation of base Likelihood Fit class, it is mostly responsible for likelihood calculation whil...
Definition: LikelihoodFit.h:6
int GetNPars()
Get total number of params, this sums over all covariance objects.
Definition: LikelihoodFit.h:16
void setChainLength(unsigned int L)
Set how long chain should be.
Definition: MCMCBase.h:26
void init()
Definition: PSO.cpp:51
EW: As FitterBase is an abstract base class we have to do some gymnastics to get it to get it into py...
Definition: fitters.cpp:14
EW: As LikelihoodFit is an abstract base class we have to do some gymnastics to get it to get it into...
Definition: fitters.cpp:31

◆ initManager()

void initManager ( py::module &  m)

Definition at line 11 of file manager.cpp.

11  {
12  auto m_manager = m.def_submodule("manager");
13  m_manager.doc() =
14  "This is a Python binding of MaCh3s C++ based manager library.";
15 
16  // Bind some of the cpp-yaml library
17  // shamelessly stolen from stackoverflow: https://stackoverflow.com/questions/62347521/using-pybind11-to-wrap-yaml-cpp-iterator
18  py::enum_<YAML::NodeType::value>(m_manager, "NodeType")
19  .value("Undefined", YAML::NodeType::Undefined)
20  .value("Null", YAML::NodeType::Null)
21  .value("Scalar", YAML::NodeType::Scalar)
22  .value("Sequence", YAML::NodeType::Sequence)
23  .value("Map", YAML::NodeType::Map);
24 
25  py::class_<YAML::Node>(m_manager, "YamlNode")
26  .def(py::init<const std::string &>())
27 
28  .def("data",
29  [](const YAML::Node node){
30  if ( node.Type() != YAML::NodeType::Scalar )
31  {
32  throw MaCh3Exception(__FILE__, __LINE__, "Attempting to access the data of non-scalar yaml node. This is undefined.");
33  }
34  return node.Scalar();
35  },
36  "Access the data stored in the node. This is only valid if the node is a 'scalar' type, i.e. it is a leaf of the yaml tree structure.")
37 
38  .def("__getitem__",
39  [](const YAML::Node node, const std::string& key){
40  return node[key];
41  })
42 
43  .def("__getitem__",
44  [](const YAML::Node node, const int& key){
45  if ( node.Type() != YAML::NodeType::Sequence)
46  {
47  throw MaCh3Exception(__FILE__, __LINE__, "Trying to access a non sequence yaml node with integer index");
48  }
49  return node[key];
50  })
51 
52  .def("__iter__",
53  [](const YAML::Node &node) {
54  return py::make_iterator(node.begin(), node.end());},
55  py::keep_alive<0, 1>())
56 
57  .def("__str__",
58  [](const YAML::Node& node) {
59  YAML::Emitter out;
60  out << node;
61  return std::string(out.c_str());
62  })
63 
64  .def("type", &YAML::Node::Type)
65 
66  .def("__len__", &YAML::Node::size)
67  ;
68 
69  py::class_<YAML::detail::iterator_value, YAML::Node>(m_manager, "_YamlDetailIteratorValue")
70  .def(py::init<>())
71  .def("first", [](YAML::detail::iterator_value& val) { return val.first;})
72  .def("second", [](YAML::detail::iterator_value& val) { return val.second;})
73  ;
74 
75  m.def("load_file", &YAML::LoadFile, "");
76 
77  py::class_<Manager>(m_manager, "Manager")
78  .def(
79  py::init<std::string const &>(),
80  "create a Manager object with a specified *config_file*",
81  py::arg("config_file")
82  )
83 
84  .def(
85  "print",
87  "Print currently used config."
88  )
89 
90  .def(
91  "raw",
92  &Manager::raw,
93  "Get the raw yaml config."
94  )
95 
96  .def(
97  "get_test_stat",
99  "Get the test statistic that was specified in the config file."
100  )
101  ;
102 }
int GetMCStatLLH() const
Get likelihood type defined in the config.
Definition: Manager.cpp:98
YAML::Node const & raw() const
Return config.
Definition: Manager.h:41
void Print() const
Print currently used config.
Definition: Manager.cpp:90

◆ initParameters()

void initParameters ( py::module &  m)

Definition at line 38 of file parameters.cpp.

38  {
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.
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

◆ initPlotting()

void initPlotting ( py::module &  m)

Definition at line 15 of file plotting.cpp.

15  {
16  auto m_plotting = m.def_submodule("plotting");
17  m_plotting.doc() = "This is a Python binding of MaCh3s C++ based plotting library.";
18 
19  py::class_<MaCh3Plotting::PlottingManager>(m_plotting, "PlottingManager")
20  .def(
21  py::init<const std::string &>(),
22  "construct a PlottingManager using the specified config file",
23  py::arg("config_file_name")
24  )
25 
26  .def(
27  py::init(),
28  "default constructor, will initialise the PlottingManager with the default plotting config"
29  )
30 
31  .def(
32  "initialise",
33  &MaCh3Plotting::PlottingManager::initialise,
34  "initalise this PlottingManager"
35  )
36 
37  .def(
38  "usage",
40  "Print a usage message for the current executable"
41  )
42 
43  .def(
44  "parse_inputs",
45  &MaCh3Plotting::PlottingManager::parseInputsVec,
46  "Parse command line variables",
47  py::arg("arguments")
48  )
49 
50  .def(
51  "set_exec",
52  &MaCh3Plotting::PlottingManager::setExec,
53  "Set the name of the current executable, which will be used when getting executable specific options from the plotting config file",
54  py::arg("exec_name")
55  )
56 
57  .def(
58  "get_file_name",
59  &MaCh3Plotting::PlottingManager::getFileName,
60  "Get the path to a particular file",
61  py::arg("input_file_id")
62  )
63 
64  .def(
65  "get_file_label",
66  &MaCh3Plotting::PlottingManager::getFileLabel,
67  "Get the specified label of a particular input file",
68  py::arg("input_file_id")
69  )
70 
71  .def(
72  "get_draw_options",
73  &MaCh3Plotting::PlottingManager::getDrawOptions,
74  "Get any additional root drawing options specified by the user"
75  )
76 
77  .def(
78  "get_output_name",
79  py::overload_cast<const std::string &>(&MaCh3Plotting::PlottingManager::getOutputName),
80  "Get the output name specified by the user, can specify an additional *suffix* to append to the file name but before the file extension",
81  py::arg("suffix") = ""
82  )
83 
84  .def(
85  "get_file_names",
86  &MaCh3Plotting::PlottingManager::getFileNames,
87  "Get the list of all file names"
88  )
89 
90  .def(
91  "get_file_labels",
92  &MaCh3Plotting::PlottingManager::getFileLabels,
93  "Get the list of all file labels"
94  )
95 
96  .def(
97  "get_n_files",
98  &MaCh3Plotting::PlottingManager::getNFiles,
99  "Get the number of specified files"
100  )
101 
102  .def(
103  "get_split_by_sample",
104  &MaCh3Plotting::PlottingManager::getSplitBySample,
105  "Get whether or not the user has set the 'split by sample' (-s) option"
106  )
107 
108  .def(
109  "get_plot_ratios",
110  &MaCh3Plotting::PlottingManager::getPlotRatios,
111  "Get whether or not the user specified the 'plot ratios' (-r) option"
112  )
113 
114  .def(
115  "get_draw_grid",
116  &MaCh3Plotting::PlottingManager::getDrawGrid,
117  "Get wheter or not the user has specified the 'draw grid' (-g) option"
118  )
119 
120  .def(
121  "style",
122  &MaCh3Plotting::PlottingManager::style, py::return_value_policy::reference,
123  "Get the StyleManager associated with this PlottingManager"
124  )
125 
126  .def(
127  "input",
128  &MaCh3Plotting::PlottingManager::input, py::return_value_policy::reference,
129  "Get the InputManager associated with this PlottingManager"
130  )
131 
132  // EM: I can't figure out how to add the getOption methods as theres not really a way to deduce the return type from yaml so leaving them out for now :/
133  // I think one solution would be to extend the PlottingManager class inside of python (add a pyPlottingManager (or something like that) that derives
134  // from this one and add the functions to that) and then fold that python code into the module somehow. But that is currently beyond my pybinding abilities
135  ;
136 
137  py::class_<MaCh3Plotting::InputManager>(m_plotting, "InputManager")
138  .def(
139  "print",
140  &MaCh3Plotting::InputManager::print,
141  "Print a summary of everything this manager knows"
142  )
143 
144  .def(
145  "get_llh_scan",
146  &MaCh3Plotting::InputManager::getLLHScan,
147  "Get the LLH scan for a particular parameter from a particular file",
148  py::arg("input_file_id"),
149  py::arg("param_name"),
150  py::arg("LLH_type") = "total"
151  )
152 
153  .def(
154  "get_llh_scan_by_sample",
155  &MaCh3Plotting::InputManager::getSampleSpecificLLHScan,
156  "Get the LLH scan for a particular parameter from a particular file for a particular sample",
157  py::arg("input_file_id"),
158  py::arg("param"),
159  py::arg("sample")
160  )
161 
162  .def(
163  "get_enabled_llh",
164  &MaCh3Plotting::InputManager::getEnabledLLH,
165  "Get whether a particular file has LLH scans for a particular parameter",
166  py::arg("input_file_id"),
167  py::arg("param"),
168  py::arg("LLH_type") = "total"
169  )
170 
171  .def(
172  "get_enabled_llh_by_sample",
173  &MaCh3Plotting::InputManager::getEnabledLLHBySample,
174  "Get whether a particular file has LLH scans for a particular parameter for a particular sample",
175  py::arg("input_file_id"),
176  py::arg("param"),
177  py::arg("sample")
178  )
179 
180  .def(
181  "get_post_fit_error",
182  &MaCh3Plotting::InputManager::getPostFitError,
183  "Get the post fit error for a parameter from a particular file",
184  py::arg("input_file_id"),
185  py::arg("param"),
186  py::arg("error_type") = ""
187  )
188 
189  .def(
190  "get_post_fit_value",
191  &MaCh3Plotting::InputManager::getPostFitValue,
192  "Get the post fit value for a parameter from a particular file",
193  py::arg("input_file_id"),
194  py::arg("param"),
195  py::arg("error_type") = ""
196  )
197 
198  .def(
199  "get_known_parameters",
200  &MaCh3Plotting::InputManager::getKnownParameters,
201  "Get all the parameters that this manager knows about. Useful for iterating over"
202  )
203 
204  .def(
205  "get_known_samples",
206  &MaCh3Plotting::InputManager::getKnownSamples,
207  "Get all the samples that this manager knows about. Useful for iterating over"
208  )
209 
210  .def(
211  "get_tagged_parameters",
212  &MaCh3Plotting::InputManager::getTaggedParameters,
213  "Get all the parameters whose tags match some specified list",
214  py::arg("tags"),
215  py::arg("check_type") = "all"
216  )
217 
218  .def(
219  "get_tagged_samples",
220  &MaCh3Plotting::InputManager::getTaggedSamples,
221  "Get all the samples whose tags match some specified list",
222  py::arg("tags"),
223  py::arg("check_type") = "all"
224  )
225 
226  .def(
227  "get_n_input_files",
228  &MaCh3Plotting::InputManager::getNInputFiles,
229  "Get the number of input files registered with this manager"
230  )
231 
232  .def(
233  "get_known_llh_parameters",
234  &MaCh3Plotting::InputManager::getKnownLLHParameters,
235  "Get all the parameters that a file has LLH scans for",
236  py::arg("file_id")
237  )
238 
239  .def(
240  "get_known_llh_samples",
241  &MaCh3Plotting::InputManager::getKnownLLHSamples,
242  "Get all the samples that a file has individual LLH scans for",
243  py::arg("file_id")
244  )
245 
246  .def(
247  "get_known_post_fit_parameters",
248  &MaCh3Plotting::InputManager::getKnownPostFitParameters,
249  "Get all the parameters that a file has post fit values and errors for",
250  py::arg("file_id")
251  )
252 
253  .def(
254  "get_known_MCMC_parameters",
255  &MaCh3Plotting::InputManager::getKnownMCMCParameters,
256  "Get all the parameters that a file has MCMC chain entries for",
257  py::arg("file_id")
258  )
259 
260  .def(
261  "get_known_1d_posterior_parameters",
262  &MaCh3Plotting::InputManager::getKnown1dPosteriorParameters,
263  "Get all the parameters that a file has processed 1d posteriors for",
264  py::arg("file_id")
265  )
266 
267  .def(
268  "get_MCMC_entry",
269  &MaCh3Plotting::InputManager::getMCMCentry,
270  "Load up a particular step in the MCMC chain for a particular input file",
271  py::arg("file_id"),
272  py::arg("step")
273  )
274 
275  .def(
276  "get_MCMC_value",
277  &MaCh3Plotting::InputManager::getMCMCvalue,
278  "Get the value of a particular parameter for the current entry (set by set_MCMC_entry) in the chain for a particular file",
279  py::arg("file_id"),
280  py::arg("param")
281  )
282 
283  .def(
284  "get_n_MCMC_entries",
285  &MaCh3Plotting::InputManager::getnMCMCentries,
286  "Get the number of entries in the MCMC chain in a particular file"
287  )
288 
289  .def(
290  "get_1d_posterior",
291  &MaCh3Plotting::InputManager::get1dPosterior,
292  "Get the 1d posterior for a particular parameter from a particular file",
293  py::arg("file_id"),
294  py::arg("param")
295  )
296  ;
297 
298  py::class_<MaCh3Plotting::StyleManager>(m_plotting, "StyleManager")
299  .def(
300  "prettify_parameter_name",
301  &MaCh3Plotting::StyleManager::prettifyParamName,
302  "Convert internally used parameter name to a nice pretty name that can be used in plots",
303  py::arg("param")
304  )
305 
306  .def(
307  "prettify_sample_name",
308  &MaCh3Plotting::StyleManager::prettifyParamName,
309  "Convert internally used sample name to a nice pretty name that can be used in plots",
310  py::arg("sample")
311  )
312  ;
313 }
void usage()

◆ initSamples()

void initSamples ( py::module &  m)

Definition at line 392 of file samples.cpp.

392  {
393  auto m_samples = m.def_submodule("samples");
394  m_samples.doc() =
395  "This is a Python binding of MaCh3s C++ based samples library.";
396 
397  // Bind the systematic type enum that lets us set different types of systematics
398  py::enum_<TestStatistic>(m_samples, "TestStatistic")
399  .value("Poisson", TestStatistic::kPoisson)
400  .value("Barlow_Beeston", TestStatistic::kBarlowBeeston)
401  .value("Ice_Cube", TestStatistic::kIceCube)
402  .value("Pearson", TestStatistic::kPearson)
403  .value("Dembinski_Abdelmottele", TestStatistic::kDembinskiAbdelmotteleb)
404  .value("N_Test_Statistics", TestStatistic::kNTestStatistics);
405 
406  py::class_<SampleHandlerBase, PySampleHandlerBase /* <--- trampoline*/>(m_samples, "SampleHandlerBase")
407  .def(py::init())
408 
409  .def(
410  "reweight",
412  "reweight the MC events in this sample. You will need to override this."
413  )
414 
415  .def(
416  "get_likelihood",
418  "Get the sample likelihood at the current point in your model space. You will need to override this."
419  )
420 
421  .def(
422  "set_test_stat",
424  "Set the test statistic that should be used when calculating likelihoods. \n\
425  :param test_stat: The new test statistic to use",
426  py::arg("test_stat")
427  )
428 
429  .def(
430  "get_bin_LLH",
431  py::overload_cast<double, double, double>(&SampleHandlerBase::GetTestStatLLH, py::const_),
432  "Get the LLH for a bin by comparing the data and MC. The result depends on having previously set the test statistic using :py:meth:`pyMaCh3.samples.SampleHandlerBase.set_test_stat` \n\
433  :param data: The data content of the bin. \n\
434  :param mc: The mc content of the bin \n\
435  :param w2: The Sum(w_{i}^2) (sum of weights squared) in the bin, which is sigma^2_{MC stats}",
436  py::arg("data"),
437  py::arg("mc"),
438  py::arg("w2")
439  )
440  ; // End of SampleHandlerBase binding
441 
442  py::class_<SampleHandlerFD, PySampleHandlerFD /* <--- trampoline*/, SampleHandlerBase>(m_samples, "SampleHandlerFD")
443  .def(
444  py::init<std::string, ParameterHandlerGeneric*>(),
445  "This should never be called directly as SampleHandlerFD is an abstract base class. \n\
446  However when creating a derived class, in the __init__() method, you should call the parent constructor i.e. this one by doing:: \n\
447  \n\
448  \tsuper(<your derived SampleHandler class>, self).__init__(*args) \n\
449  \n ",
450  py::arg("mc_version"),
451  py::arg("xsec_cov")
452  )
453  ;
454 
455  /* Not sure if this will be needed in future versions of MaCh3 so leaving commented for now
456  py::class_<fdmc_base>(m_samples, "MCstruct")
457  .def(py::init())
458 
459  // Because a lot of the variables in fdmc_base use c style arrays,
460  // we need to provide some setter functions to be able to set them using more
461  // "pythony" objects, e.g. lists and numpy arrays
462  .def(
463  "set_event_variable_values",
464  [](fdmc_base &self, int dim, py::array_t<double, py::array::c_style> &array)
465  {
466  py::buffer_info bufInfo = array.request();
467 
468  if ( dim > 2 )
469  throw MaCh3Exception(__FILE__, __LINE__, "Currently only dimensions of 1 or 2 are supported sorry :(");
470 
471  if ( bufInfo.ndim != 1 )
472  throw MaCh3Exception(__FILE__, __LINE__, "Number of dimensions in parameter array must be one if setting only one of the event variable arrays!");
473 
474  if( dim ==1 )
475  self.x_var = array.data();
476 
477  else if ( dim == 2)
478  self.y_var = array.data();
479  }
480  )
481  ;
482  */
483 }
@ kNTestStatistics
Number of test statistics.
@ kPearson
Standard Pearson likelihood .
@ kBarlowBeeston
Barlow-Beeston () following Conway approximation ()
@ kIceCube
Based on .
@ kDembinskiAbdelmotteleb
Based on .
@ kPoisson
Standard Poisson likelihood .
EW: As SampleHandlerBase is an abstract base class we have to do some gymnastics to get it to get it ...
Definition: samples.cpp:10
As SampleHandlerFD is an abstract base class we have to do some gymnastics to get it to get it into p...
Definition: samples.cpp:252
Class responsible for handling implementation of samples used in analysis, reweighting and returning ...
void SetTestStatistic(TestStatistic testStat)
Set the test statistic to be used when calculating the binned likelihoods.
virtual void Reweight()=0
double GetTestStatLLH(const double data, const double mc, const double w2) const
Calculate test statistic for a single bin. Calculation depends on setting of fTestStatistic....
virtual double GetLikelihood() const =0
Class responsible for handling implementation of samples used in analysis, reweighting and returning ...

◆ initSplines()

void initSplines ( py::module &  m)

Definition at line 72 of file splines.cpp.

72  {
73  auto m_splines = m.def_submodule("splines");
74  m_splines.doc() =
75  "This is a Python binding of MaCh3s C++ based spline library.";
76 
77  // Bind the interpolation type enum that lets us set different interpolation types for our splines
78  py::enum_<SplineInterpolation>(m_splines, "InterpolationType")
79  .value("Linear", SplineInterpolation::kLinear, "Linear interpolation between the knots")
80  .value("Linear_Func", SplineInterpolation::kLinearFunc, "Same as 'Linear'")
81  .value("Cubic_TSpline3", SplineInterpolation::kTSpline3, "Use same coefficients as `ROOT's TSpline3 <https://root.cern.ch/doc/master/classTSpline3.html>`_ implementation")
82  .value("Cubic_Monotonic", SplineInterpolation::kMonotonic, "Coefficients are calculated such that the segments between knots are forced to be monotonic. The implementation we use is based on `this method <https://www.jstor.org/stable/2156610>`_ by Fritsch and Carlson.")
83  .value("Cubic_Akima", SplineInterpolation::kAkima, "The second derivative is not required to be continuous at the knots. This means that these splines are useful if the second derivative is rapidly varying. The implementation we used is based on `this paper <http://www.leg.ufpr.br/lib/exe/fetch.php/wiki:internas:biblioteca:akima.pdf>`_ by Akima.")
84  .value("N_Interpolation_Types", SplineInterpolation::kSplineInterpolations, "This is only to be used when iterating and is not a valid interpolation type.");
85 
86 
87  py::class_<SplineBase, PySplineBase /* <--- trampoline*/>(m_splines, "SplineBase");
88 
89  py::class_<TResponseFunction_red>(m_splines, "_ResponseFunctionBase")
90  .doc() = "Base class of the response function, this binding only exists for consistency with the inheritance structure of the c++ code. Just pretend it doesn't exist and don't worry about it...";
91 
92  // Bind the TSpline3_red class. Decided to go with a clearer name of ResponseFunction for the python binding
93  // and make the interface a bit more python-y. Additionally remove passing root stuff so we don't need to deal
94  // with root python binding and can just pass it native python objects.
95  py::class_<TSpline3_red, TResponseFunction_red, std::unique_ptr<TSpline3_red, py::nodelete>>(m_splines, "ResponseFunction")
96  .def(
97  // define a more python friendly constructor that massages the inputs and passes them
98  // through to the c++ constructor
99  py::init
100  (
101  // Just take in some vectors, then build a TSpline3 and pass this to the constructor
102  [](std::vector<double> &xVals, std::vector<double> &yVals, SplineInterpolation interpType)
103  {
104  if ( xVals.size() != yVals.size() )
105  {
106  throw MaCh3Exception(__FILE__, __LINE__, "Different number of x values and y values!");
107  }
108 
109  int length = int(xVals.size());
110 
111  if (length == 1)
112  {
113  M3::float_t xKnot = M3::float_t(xVals[0]);
114  M3::float_t yKnot = M3::float_t(yVals[0]);
115 
116  std::vector<M3::float_t *> pars;
117  pars.resize(3);
118  pars[0] = new M3::float_t(0.0);
119  pars[1] = new M3::float_t(0.0);
120  pars[2] = new M3::float_t(0.0);
121  delete pars[0];
122  delete pars[1];
123  delete pars[2];
124 
125  return new TSpline3_red(&xKnot, &yKnot, 1, pars.data());
126  }
127 
128  TSpline3 *splineTmp = new TSpline3( "spline_tmp", xVals.data(), yVals.data(), length );
129  return new TSpline3_red(splineTmp, interpType);
130  }
131  )
132  )
133 
134  .def(
135  "find_segment",
137  "Find the segment that a particular *value* lies in. \n"
138  ":param value: The value to test",
139  py::arg("value")
140  )
141 
142  .def(
143  "evaluate",
145  "Evaluate the response function at a particular *value*. \n"
146  ":param value: The value to evaluate at.",
147  py::arg("value")
148  )
149  ; // End of binding for ResponseFunction
150 
151  py::class_<SMonolith, SplineBase>(m_splines, "EventSplineMonolith")
152  .def(
153  py::init(
154  [](std::vector<std::vector<TResponseFunction_red*>> &responseFns, const bool saveFlatTree)
155  {
156  std::vector<RespFuncType> respFnTypes;
157  for(uint i = 0; i < responseFns[0].size(); i++)
158  {
159  // ** WARNING **
160  // Right now I'm only pushing back TSpline3_reds as that's all that's supported right now
161  // In the future there might be more
162  // I think what would be best to do would be to store the interpolation type somehow in the ResponseFunction objects
163  // then just read them here and pass through to the constructor
164  respFnTypes.push_back(RespFuncType::kTSpline3_red);
165  }
166  return new SMonolith(responseFns, respFnTypes, saveFlatTree);
167  }
168  ),
169  "Create an EventSplineMonolith \n"
170  ":param master_splines: These are the 'knot' values to make splines from. This should be an P x E 2D list where P is the number of parameters and E is the number of events. \n"
171  ":param save_flat_tree: Whether we want to save monolith into speedy flat tree",
172  py::arg("master_splines"),
173  py::arg("save_flat_tree") = false
174  )
175 
176  .def(
177  py::init<std::string>(),
178  "Constructor where you pass path to preprocessed root FileName which is generated by creating an EventSplineMonolith with the `save_flat_tree` flag set to True. \n"
179  ":param file_name: The name of the file to read from.",
180  py::arg("file_name")
181  )
182 
183  .def(
184  "evaluate",
186  "Evaluate the splines at their current values."
187  )
188 
189  .def(
190  "sync_mem_transfer",
192  "This is important when running on GPU. After calculations are done on GPU we copy memory to CPU. This operation is asynchronous meaning while memory is being copied some operations are being carried. Memory must be copied before actual reweight. This function make sure all has been copied."
193  )
194 
195  .def(
196  "get_event_weight",
198  py::return_value_policy::reference,
199  "Get the weight of a particular event. \n"
200  ":param event: The index of the event whose weight you would like.",
201  py::arg("event")
202  )
203 
204  .def(
205  "set_param_value_array",
206  // Wrap up the setSplinePointers method so that we can take in a numpy array and get
207  // pointers to it's sweet sweet data and use those pointers in the splineMonolith
208  [](SMonolith &self, py::array_t<double, py::array::c_style> &array)
209  {
210  py::buffer_info bufInfo = array.request();
211 
212  if ( bufInfo.ndim != 1)
213  {
214  throw MaCh3Exception(__FILE__, __LINE__, "Number of dimensions in parameter array must be one!");
215  }
216 
217  if ( bufInfo.shape[0] != self.GetNParams() )
218  {
219  throw MaCh3Exception(__FILE__, __LINE__, "Number of entries in parameter array must equal the number of parameters!");
220  }
221 
222  std::vector<const double *> paramVec;
223  paramVec.resize(self.GetNParams());
224 
225  for( int idx = 0; idx < self.GetNParams(); idx++ )
226  {
227  // booooo pointer arithmetic
228  paramVec[idx] = array.data() + idx;
229  }
230 
231  self.setSplinePointers(paramVec);
232  },
233  "Set the array that the monolith should use to read parameter values from. \n"
234  "Usage of this might vary a bit from what you're used to in python. \n"
235  "Rather than just setting the values here, what you're really doing is setting pointers in the underlying c++ code. \n"
236  "What that means is that you pass an array to this function like:: \n"
237  "\n event_spline_monolith_instance.set_param_value_array(array) \n\n"
238  "Then when you set values in that array as normal, they will also be updated inside of the event_spline_monolith_instance.",
239  py::arg("array")
240 
241  )
242 
243  .doc() = "This 'monolith' deals with event by event weighting using splines."
244 
245  ; // End of binding for EventSplineMonolith
246 }
SplineInterpolation
Make an enum of the spline interpolation type.
@ kTSpline3
Default TSpline3 interpolation.
@ kMonotonic
EM: DOES NOT make the entire spline monotonic, only the segments.
@ kSplineInterpolations
This only enumerates.
@ kLinear
Linear interpolation between knots.
@ kLinearFunc
Liner interpolation using TF1 not spline.
@ kAkima
EM: Akima spline iis allowed to be discontinuous in 2nd derivative and coefficients in any segment.
@ kTSpline3_red
Uses TSpline3_red for interpolation.
Custom exception class used throughout MaCh3.
EW: As SplineBase is an abstract base class we have to do some gymnastics to get it to get it into py...
Definition: splines.cpp:19
Even-by-event class calculating response for spline parameters. It is possible to use GPU acceleratio...
void SynchroniseMemTransfer() const override
KS: After calculations are done on GPU we copy memory to CPU. This operation is asynchronous meaning ...
void Evaluate() override
CW: This Eval should be used when using two separate x,{y,a,b,c,d} arrays to store the weights; proba...
const float * retPointer(const int event) const
KS: Get pointer to total weight to make fit faster wrooom!
Base class for calculating weight from spline.
Definition: SplineBase.h:25
CW: Reduced TSpline3 class.
double Eval(double var) override
CW: Evaluate the weight from a variation.
int FindX(double x)
Find the segment relevant to this variation in x.
double float_t
Definition: Core.h:37

◆ PYBIND11_MODULE()

PYBIND11_MODULE ( _pyMaCh3  ,
 
)

Definition at line 16 of file pyMaCh3.cpp.

16  {
17  initPlotting(m);
18  initFitters(m);
19  initSamples(m);
20  initManager(m);
21  initParameters(m);
22  initSplines(m);
23 }
void initFitters(py::module &)
Definition: fitters.cpp:47
void initParameters(py::module &)
Definition: parameters.cpp:38
void initSplines(py::module &)
Definition: splines.cpp:72
void initManager(py::module &)
Definition: manager.cpp:11
void initSamples(py::module &)
Definition: samples.cpp:392
void initPlotting(py::module &)
Definition: plotting.cpp:15