MaCh3  2.5.1
Reference Guide
Classes | Functions
samples.h File Reference
#include <pybind11/pybind11.h>
#include <pybind11/stl.h>
#include <pybind11/numpy.h>
#include "Samples/SampleHandlerBase.h"
#include "TH1.h"
#include "TH2.h"
#include "histutils.h"
Include dependency graph for samples.h:
This graph shows which files directly or indirectly include this file:

Go to the source code of this file.

Classes

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

Functions

void initSamplesModule (py::module &m_samples)
 

Detailed Description

Author
Ewan Miller

Definition in file samples.h.

Function Documentation

◆ initSamplesModule()

void initSamplesModule ( py::module &  m_samples)

Definition at line 343 of file samples.h.

343  {
344 
345  // Bind the systematic type enum that lets us set different types of systematics
346  py::enum_<TestStatistic>(m_samples, "TestStatistic")
347  .value("Poisson", TestStatistic::kPoisson)
348  .value("Barlow_Beeston", TestStatistic::kBarlowBeeston)
349  .value("Ice_Cube", TestStatistic::kIceCube)
350  .value("Pearson", TestStatistic::kPearson)
351  .value("Dembinski_Abdelmottele", TestStatistic::kDembinskiAbdelmotteleb)
352  .value("N_Test_Statistics", TestStatistic::kNTestStatistics);
353 
354  py::class_<SampleHandlerInterface, PySampleHandlerInterface /* <--- trampoline*/>(m_samples, "SampleHandlerInterface")
355  .def(py::init())
356 
357  .def(
358  "reweight",
360  "reweight the MC events in this sample. You will need to override this.")
361 
362  .def(
363  "get_n_samples",
365  "Get the total number of samples"
366  )
367 
368  .def(
369  "get_n_dim",
371  py::arg("sample"),
372  "Get the dimension of a given sample"
373  )
374 
375  .def(
376  "get_mc_hist",
377  [](SampleHandlerBase &self, const int sample)
378  {
379  auto hist_original = M3::Clone(self.GetMCHist(sample));
380 
381  auto edges = HistToNumpy(hist_original);
382  return edges;
383  },
384 
385  py::return_value_policy::reference_internal,
386  py::arg("sample"),
387  "Get MC histogram as numpy arrays.\n"
388  "For 1D: Returns (contents, edges)\n"
389  "For 2D: Returns (contents, edgesX, edgesY)\n"
390  "where contents is shape (nbinsY, nbinsX) for 2D")
391 
392  .def(
393  "get_likelihood",
395  "Get the sample likelihood at the current point in your model space. You will need to override this.")
396 
397  .def(
398  "set_test_stat",
400  "Set the test statistic that should be used when calculating likelihoods. \n\
401  :param test_stat: The new test statistic to use",
402  py::arg("test_stat"))
403 
404  .def(
405  "get_bin_LLH",
406  py::overload_cast<double, double, double>(&SampleHandlerInterface::GetTestStatLLH, py::const_),
407  "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.SampleHandlerInterface.set_test_stat` \n\
408  :param data: The data content of the bin. \n\
409  :param mc: The mc content of the bin \n\
410  :param w2: The Sum(w_{i}^2) (sum of weights squared) in the bin, which is sigma^2_{MC stats}",
411  py::arg("data"),
412  py::arg("mc"),
413  py::arg("w2"))
414 
415  ; // End of SampleHandlerInterface binding
416 
417  py::class_<SampleHandlerBase, PySampleHandlerBase /* <--- trampoline*/,
418  SampleHandlerInterface>(m_samples, "SampleHandlerBase")
419  .def(
420  py::init<std::string, ParameterHandlerGeneric *>(),
421  "This should never be called directly as SampleHandlerBase is an abstract base class. \n\
422  However when creating a derived class, in the __init__() method, you should call the parent constructor i.e. this one by doing:: \n\
423  \n\
424  \tsuper(<your derived SampleHandler class>, self).__init__(*args) \n\
425  \n ",
426  py::arg("mc_version"), py::arg("xsec_cov"))
427 
428  .def(
429  "add_data",
430  py::overload_cast<const int, const std::vector<double>&>(&SampleHandlerBase::AddData),
431  py::arg("sample"),
432  py::arg("data_array"),
433  "Set the data for your sample handler (assumes the binning is the same as your MC!)"
434  )
435 
436  // ================
437  // Useful getters
438  // ===============
439  .def(
440  "get_sample_title",
442  py::arg("sample"),
443  "Get the title for a given sample"
444  )
445  .def(
446  "get_data_array",
447  py::overload_cast<const int>(&SampleHandlerBase::GetDataArray, py::const_),
448  py::arg("sample"),
449  "Returns the contents of the MC histogram as a flat list"
450  )
451 
452  .def(
453  "get_mc_array",
454  py::overload_cast<const int>(&SampleHandlerBase::GetMCArray, py::const_),
455  py::arg("sample"),
456  "Returns the contents of the MC histogram as a flat list"
457  )
458 
459  .def(
460  "get_w2_array",
461  py::overload_cast<const int>(&SampleHandlerBase::GetW2Array, py::const_),
462  py::arg("sample"),
463  "Returns the contents of the W2 histogram as a flat list"
464  )
465 
466  .def(
467  "get_data_hist",
468  [](SampleHandlerBase &self, const int sample)
469  {
470  auto hist_original = M3::Clone(self.GetDataHist(sample));
471 
472  auto edges = HistToNumpy(hist_original);
473  return edges;
474 
475  },
476  py::arg("Dimension"),
477  "Get Data histogram as numpy arrays.\n"
478  "For 1D: Returns (contents, edges)\n"
479  "For 2D: Returns (contents, edgesX, edgesY)\n"
480  "where contents is shape (nbinsY, nbinsX) for 2D")
481 
482  .def(
483  "get_w2_hist",
484  [](SampleHandlerBase &self, const int sample)
485  {
486  auto hist_original = M3::Clone(self.GetW2Hist(sample));
487 
488  auto edges = HistToNumpy(hist_original);
489  return edges;
490  },
491 
492 
493  py::arg("sample"),
494  "Get W2 histogram as numpy arrays.\n"
495  "For 1D: Returns (contents, edges)\n"
496  "For 2D: Returns (contents, edgesX, edgesY)\n"
497  "where contents is shape (nbinsY, nbinsX) for 2D")
498 
499  .def("get_var_hist",
500  [](SampleHandlerBase &self, const int iSample,
501  const std::string& ProjectionVarX,
502  const std::string& ProjectionVarY="",
503  const std::vector<KinematicCut> &EventSelectionVec = {},
504  int WeightStyle = 0,
505  const std::vector< KinematicCut >& SubEventSelectionVec = {})
506  {
507 
508  py::array_t<M3::float_t> edgesY, edgesX, contents;
509  std::unique_ptr<TH1> hist;
510  int dim;
511  if(ProjectionVarY==""){
512  hist = self.Get1DVarHist(iSample,
513  ProjectionVarX,
514  EventSelectionVec,
515  WeightStyle,
516  SubEventSelectionVec);
517  } else{
518  hist = self.Get2DVarHist(iSample,
519  ProjectionVarX,
520  ProjectionVarY,
521  EventSelectionVec,
522  WeightStyle,
523  SubEventSelectionVec);
524  }
525  return HistToNumpy(hist);
526  }
527  )
528  ; // End of SampleHandler Base
529 
530 
531  py::class_<KinematicCut>(m_samples, "KinematicCut")
532  .def(py::init<>(), "Simple wrapper around Kinematic cuts")
533  .def_readwrite("param_name", &KinematicCut::ParamToCutOnIt, "Parameter to cut on")
534  .def_readwrite("lower_bound", &KinematicCut::LowerBound, "Lower bound")
535  .def_readwrite("upper_bound", &KinematicCut::UpperBound, "Upper Bound");
536 
537 
538  /* Not sure if this will be needed in future versions of MaCh3 so leaving commented for now
539  py::class_<fdmc_base>(m_samples, "MCstruct")
540  .def(py::init())
541 
542  // Because a lot of the variables in fdmc_base use c style arrays,
543  // we need to provide some setter functions to be able to set them using more
544  // "pythony" objects, e.g. lists and numpy arrays
545  .def(
546  "set_event_variable_values",
547  [](fdmc_base &self, int dim, py::array_t<M3::float_t, py::array::c_style> &array)
548  {
549  py::buffer_info bufInfo = array.request();
550 
551  if ( dim > 2 )
552  throw MaCh3Exception(__FILE__, __LINE__, "Currently only dimensions of 1 or 2 are supported sorry :(");
553 
554  if ( bufInfo.ndim != 1 )
555  throw MaCh3Exception(__FILE__, __LINE__, "Number of dimensions in parameter array must be one if setting only one of the event variable arrays!");
556 
557  if( dim ==1 )
558  self.x_var = array.data();
559 
560  else if ( dim == 2)
561  self.y_var = array.data();
562  }
563  )
564  ;
565  */
566 }
@ 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 .
As SampleHandlerBase is an abstract base class we have to do some gymnastics to get it to get it into...
Definition: samples.h:234
EW: As SampleHandlerBase is an abstract base class we have to do some gymnastics to get it to get it ...
Definition: samples.h:20
Class responsible for handling implementation of samples used in analysis, reweighting and returning ...
void AddData(const int Sample, TH1 *Data)
DB: Add data for a given sample from a ROOT histogram.
auto GetDataArray() const
Return array storing data entries for every bin.
std::string GetSampleTitle(const int Sample) const final
Get fancy title for specified samples.
auto GetMCArray() const
Return array storing MC entries for every bin.
auto GetW2Array() const
Return array storing W2 entries for every bin.
Class responsible for handling implementation of samples used in analysis, reweighting and returning ...
virtual double GetLikelihood() const =0
Return likelihood (-logL) for all samples.
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 M3::int_t GetNSamples()
returns total number of samples
virtual void Reweight()=0
main routine modifying MC prediction based on proposed parameter values
void SetTestStatistic(TestStatistic testStat)
Set the test statistic to be used when calculating the binned likelihoods.
virtual int GetNDim(const int Sample) const =0
DB Get what dimensionality binning for given sample has.
py::tuple HistToNumpy(std::unique_ptr< TH1 > &hist)
Definition: histutils.h:56
std::unique_ptr< ObjectType > Clone(const ObjectType *obj, const std::string &name="")
KS: Creates a copy of a ROOT-like object and wraps it in a smart pointer.
double UpperBound
Upper bound on which we apply cut.
double LowerBound
Lower bound on which we apply cut.
int ParamToCutOnIt
Index or enum value identifying the kinematic variable to cut on.