MaCh3  2.2.3
Reference Guide
Classes | Functions
samples.cpp File Reference
#include <pybind11/pybind11.h>
#include <pybind11/stl.h>
#include <pybind11/numpy.h>
#include "Samples/SampleHandlerFD.h"
Include dependency graph for samples.cpp:

Go to the source code of this file.

Classes

class  PySampleHandlerBase
 EW: As SampleHandlerBase is an abstract base class we have to do some gymnastics to get it to get it into python. More...
 
class  PySampleHandlerFD
 As SampleHandlerFD is an abstract base class we have to do some gymnastics to get it to get it into python. More...
 

Functions

void initSamples (py::module &m)
 

Function Documentation

◆ initSamples()

void initSamples ( py::module &  m)

Definition at line 196 of file samples.cpp.

196  {
197 
198  auto m_samples = m.def_submodule("samples");
199  m_samples.doc() =
200  "This is a Python binding of MaCh3s C++ based samples library.";
201 
202  // Bind the systematic type enum that lets us set different types of systematics
203  py::enum_<TestStatistic>(m_samples, "TestStatistic")
204  .value("Poisson", TestStatistic::kPoisson)
205  .value("Barlow_Beeston", TestStatistic::kBarlowBeeston)
206  .value("Ice_Cube", TestStatistic::kIceCube)
207  .value("Pearson", TestStatistic::kPearson)
208  .value("Dembinski_Abdelmottele", TestStatistic::kDembinskiAbdelmotteleb)
209  .value("N_Test_Statistics", TestStatistic::kNTestStatistics);
210 
211  py::class_<SampleHandlerBase, PySampleHandlerBase /* <--- trampoline*/>(m_samples, "SampleHandlerBase")
212  .def(py::init())
213 
214  .def(
215  "reweight",
217  "reweight the MC events in this sample. You will need to override this."
218  )
219 
220  .def(
221  "get_likelihood",
223  "Get the sample likelihood at the current point in your model space. You will need to override this."
224  )
225 
226  .def(
227  "set_test_stat",
229  "Set the test statistic that should be used when calculating likelihoods. \n\
230  :param test_stat: The new test statistic to use",
231  py::arg("test_stat")
232  )
233 
234  .def(
235  "get_bin_LLH",
236  py::overload_cast<double, double, double>(&SampleHandlerBase::GetTestStatLLH, py::const_),
237  "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\
238  :param data: The data content of the bin. \n\
239  :param mc: The mc content of the bin \n\
240  :param w2: The Sum(w_{i}^2) (sum of weights squared) in the bin, which is sigma^2_{MC stats}",
241  py::arg("data"),
242  py::arg("mc"),
243  py::arg("w2")
244  )
245  ; // End of SampleHandlerBase binding
246 
247  py::class_<SampleHandlerFD, PySampleHandlerFD /* <--- trampoline*/, SampleHandlerBase>(m_samples, "SampleHandlerFD")
248  .def(
249  py::init<std::string, ParameterHandlerGeneric*>(),
250  "This should never be called directly as SampleHandlerFD is an abstract base class. \n\
251  However when creating a derived class, in the __init__() method, you should call the parent constructor i.e. this one by doing:: \n\
252  \n\
253  \tsuper(<your derived SampleHandler class>, self).__init__(*args) \n\
254  \n ",
255  py::arg("mc_version"),
256  py::arg("xsec_cov")
257  )
258  ;
259 
260  /* Not sure if this will be needed in future versions of MaCh3 so leaving commented for now
261  py::class_<fdmc_base>(m_samples, "MCstruct")
262  .def(py::init())
263 
264  // Because a lot of the variables in fdmc_base use c style arrays,
265  // we need to provide some setter functions to be able to set them using more
266  // "pythony" objects, e.g. lists and numpy arrays
267  .def(
268  "set_event_variable_values",
269  [](fdmc_base &self, int dim, py::array_t<double, py::array::c_style> &array)
270  {
271  py::buffer_info bufInfo = array.request();
272 
273  if ( dim > 2 )
274  throw MaCh3Exception(__FILE__, __LINE__, "Currently only dimensions of 1 or 2 are supported sorry :(");
275 
276  if ( bufInfo.ndim != 1 )
277  throw MaCh3Exception(__FILE__, __LINE__, "Number of dimensions in parameter array must be one if setting only one of the event variable arrays!");
278 
279  if( dim ==1 )
280  self.x_var = array.data();
281 
282  else if ( dim == 2)
283  self.y_var = array.data();
284  }
285  )
286  ;
287  */
288 }
@ 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:56
Class responsible for handling implementation of samples used in analysis, reweighting and returning ...
virtual void Reweight()=0
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.
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()=0