MaCh3  2.4.2
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 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 ...