MaCh3  2.5.1
Reference Guide
samples.h
Go to the documentation of this file.
1 #pragma once
2 
5 
6 #include <pybind11/pybind11.h>
7 #include <pybind11/stl.h>
8 #include <pybind11/numpy.h>
9 
11 #include "TH1.h"
12 #include "TH2.h"
13 
14 #include "histutils.h"
15 
16 namespace py = pybind11;
17 
18 
21 public:
22  /* Inherit the constructors */
24 
25  /* Trampoline (need one for each virtual function) */
26  std::string GetName() const override {
27  PYBIND11_OVERRIDE_PURE(
28  std::string, /* Return type */
29  SampleHandlerInterface, /* Parent class */
30  GetName, /* Name of function in C++ (must match Python name) */
31  );
32  }
33 
34  /* Trampoline (need one for each virtual function) */
35  std::string GetSampleTitle(const int iSample) const override {
36  PYBIND11_OVERRIDE_PURE(
37  std::string, /* Return type */
38  SampleHandlerInterface, /* Parent class */
39  GetSampleTitle, /* Name of function in C++ (must match Python name) */
40  iSample /* Argument(s) */
41  );
42  }
43 
44  /* Trampoline (need one for each virtual function) */
45  int GetNOscChannels(const int iSample) const override {
46  PYBIND11_OVERRIDE_PURE(
47  int, /* Return type */
48  SampleHandlerInterface, /* Parent class */
49  GetNOscChannels, /* Name of function in C++ (must match Python name) */
50  iSample /* Argument(s) */
51  );
52  }
53 
54  /* Trampoline (need one for each virtual function) */
55  void Reweight() override {
56  PYBIND11_OVERRIDE_PURE_NAME(
57  void, /* Return type */
58  SampleHandlerInterface, /* Parent class */
59  "reweight",
60  Reweight /* Name of function in C++ (must match Python name) */
61  );
62  }
63 
64 
65  /* Trampoline (need one for each virtual function) */
66  double GetSampleLikelihood(const int iSample) const override {
67  PYBIND11_OVERRIDE_PURE_NAME(
68  double, /* Return type */
69  SampleHandlerInterface, /* Parent class */
70  "get_sample_likelihood",
71  GetSampleLikelihood, /* Name of function in C++ (must match Python name) */
72  iSample /* Argument(s) */
73  );
74  }
75 
76  /* Trampoline (need one for each virtual function) */
77  void CleanMemoryBeforeFit() override {
78  PYBIND11_OVERRIDE_PURE_NAME(
79  void, /* Return type */
80  SampleHandlerInterface, /* Parent class */
81  "clean_memory_before_fit",
82  CleanMemoryBeforeFit /* Name of function in C++ (must match Python name) */
83  );
84  }
85 
86  /* Trampoline (need one for each virtual function) */
87  void PrintRates(const bool DataOnly = false) override {
88  PYBIND11_OVERRIDE_PURE(
89  void, /* Return type */
90  SampleHandlerInterface, /* Parent class */
91  PrintRates, /* Name of function in C++ (must match Python name) */
92  DataOnly /* Argument(s) */
93  );
94  }
95 
96  /* Trampoline (need one for each virtual function) */
97  std::string GetKinVarName(const int iSample, const int Dimension) const override {
98  PYBIND11_OVERRIDE_PURE(
99  std::string, /* Return type */
100  SampleHandlerInterface, /* Parent class */
101  GetKinVarName, /* Name of function in C++ (must match Python name) */
102  iSample, /* Argument(s) */
103  Dimension /* Argument(s) */
104  );
105  }
106 
107  /* Trampoline (need one for each virtual function) */
108  std::vector<double> ReturnKinematicParameterBinning(const int Sample, const std::string &KinematicParameter) const override {
109  PYBIND11_OVERRIDE_PURE(
110  std::vector<double>, /* Return type */
111  SampleHandlerInterface, /* Parent class */
112  GetKinVarName, /* Name of function in C++ (must match Python name) */
113  Sample, /* Argument(s) */
114  KinematicParameter /* Argument(s) */
115  );
116  }
117 
118  TH1* GetDataHist(const int Sample) override {
119  PYBIND11_OVERRIDE_PURE(
120  TH1*, /* Return type */
121  SampleHandlerInterface, /* Parent class */
122  GetDataHist, /* Name of function in C++ (must match Python name) */
123  Sample /* Argument(s) */
124  );
125  }
126 
127  TH1* GetMCHist(const int Sample) override {
128  PYBIND11_OVERRIDE_PURE(
129  TH1*, /* Return type */
130  SampleHandlerInterface, /* Parent class */
131  GetMCHist, /* Name of function in C++ (must match Python name) */
132  Sample /* Argument(s) */
133  );
134  }
135 
136  TH1* GetW2Hist(const int Sample) override {
137  PYBIND11_OVERRIDE_PURE(
138  TH1*, /* Return type */
139  SampleHandlerInterface, /* Parent class */
140  GetW2Hist, /* Name of function in C++ (must match Python name) */
141  Sample /* Argument(s) */
142  );
143  }
144 
145  double GetLikelihood() const override {
146  PYBIND11_OVERRIDE_PURE_NAME(
147  double, /* Return type */
148  SampleHandlerInterface, /* Parent class */
149  "get_likelihood", /* Python name*/
150  GetLikelihood /* Name of function in C++ (must match Python name) */
151  /* Argument(s) */
152  );
153  }
154 
155  std::unique_ptr<TH1> Get1DVarHistByModeAndChannel(const int iSample,
156  const std::string& ProjectionVar_Str,
157  const int kModeToFill = -1,
158  const int kChannelToFill = -1,
159  const int WeightStyle = 0) override {
160  (void) iSample;
161  (void) ProjectionVar_Str;
162  (void) kModeToFill;
163  (void) kChannelToFill;
164  (void) WeightStyle;
165  return nullptr;
166  }
167 
168  std::unique_ptr<TH2> Get2DVarHistByModeAndChannel(const int iSample,
169  const std::string& ProjectionVar_StrX,
170  const std::string& ProjectionVar_StrY,
171  const int kModeToFill = -1,
172  const int kChannelToFill = -1,
173  const int WeightStyle = 0) override {
174  (void) iSample;
175  (void) ProjectionVar_StrX;
176  (void) ProjectionVar_StrY;
177  (void) kModeToFill;
178  (void) kChannelToFill;
179  (void) WeightStyle;
180  return nullptr;
181  }
182 
183  std::unique_ptr<TH1> Get1DVarHist(const int iSample,
184  const std::string &ProjectionVar,
185  const std::vector<KinematicCut> &EventSelectionVec = {},
186  const int WeightStyle = 0,
187  const std::vector<KinematicCut> &SubEventSelectionVec = {}) override {
188  (void) iSample;
189  (void) ProjectionVar;
190  (void) EventSelectionVec;
191  (void) WeightStyle;
192  (void) SubEventSelectionVec;
193  return nullptr;
194  }
195 
196  std::unique_ptr<TH2> Get2DVarHist(const int iSample,
197  const std::string& ProjectionVarX,
198  const std::string& ProjectionVarY,
199  const std::vector<KinematicCut>& EventSelectionVec = {},
200  const int WeightStyle = 0,
201  const std::vector<KinematicCut>& SubEventSelectionVec = {}) override {
202  (void) iSample;
203  (void) ProjectionVarX;
204  (void) ProjectionVarY;
205  (void) EventSelectionVec;
206  (void) WeightStyle;
207  (void) SubEventSelectionVec;
208  return nullptr;
209  }
210 
211  int GetNDim(const int Sample) const override {
212  PYBIND11_OVERRIDE_PURE(
213  int, /* Return type */
214  SampleHandlerInterface, /* Parent class */
215  GetNDim, /* Name of function in C++ */
216  Sample
217  );
218  }
219 
220  std::string GetFlavourName(const int iSample,
221  const int iChannel) const override {
222  PYBIND11_OVERRIDE_PURE(
223  std::string, /* Return type */
224  SampleHandlerInterface, /* Parent class */
225  GetFlavourName, /* Name of function in C++ */
226  iSample,
227  iChannel
228  );
229  }
230 };
231 
232 
235 public:
236  /* Inherit the constructors */
238 
239  /* Trampoline (need one for each virtual function) */
240  void AddAdditionalWeightPointers() override {
241  PYBIND11_OVERRIDE_PURE_NAME(
242  void, /* Return type */
243  SampleHandlerBase, /* Parent class */
244  "add_additional_weight_pointers", /*python name*/
245  AddAdditionalWeightPointers, /* Name of function in C++ */
246  /* Argument(s) */
247  );
248  }
249 
250  /* Trampoline (need one for each virtual function) */
251  void CleanMemoryBeforeFit() override {
252  PYBIND11_OVERRIDE_PURE(
253  void, /* Return type */
254  SampleHandlerBase, /* Parent class */
255  CleanMemoryBeforeFit /* Name of function in C++ (must match Python name) */
256  );
257  }
258 
259  /* Trampoline (need one for each virtual function) */
260  void SetupSplines() override {
261  PYBIND11_OVERRIDE_PURE_NAME(
262  void, /* Return type */
263  SampleHandlerBase, /* Parent class */
264  "setup_splines", /*python name*/
265  SetupSplines, /* Name of function in C++ */
266  /* Argument(s) */
267  );
268  }
269 
270  /* Trampoline (need one for each virtual function) */
271  void Init() override {
272  PYBIND11_OVERRIDE_PURE_NAME(
273  void, /* Return type */
274  SampleHandlerBase, /* Parent class */
275  "init", /*python name*/
276  Init, /* Name of function in C++ */
277  /* Argument(s) */
278  );
279  }
280 
281  /* Trampoline (need one for each virtual function) */
282  int SetupExperimentMC() override {
283  PYBIND11_OVERRIDE_PURE_NAME(
284  int, /* Return type */
285  SampleHandlerBase, /* Parent class */
286  "setup_experiment_MC", /*python name*/
287  SetupExperimentMC, /* Name of function in C++ */
288  );
289  }
290 
291  /* Trampoline (need one for each virtual function) */
292  void InititialiseData() override {
293  PYBIND11_OVERRIDE_PURE_NAME(
294  void, /* Return type */
295  SampleHandlerBase, /* Parent class */
296  "inititialis_data", /*python name*/
297  InititialiseData, /* Name of function in C++ */
298  );
299  }
300 
301  /* Trampoline (need one for each virtual function) */
302  void SetupMC() override {
303  PYBIND11_OVERRIDE_PURE_NAME(
304  void, /* Return type */
305  SampleHandlerBase, /* Parent class */
306  "setup_MC", /*python name*/
307  SetupMC, /* Name of function in C++ */
308  );
309  }
310 
311  double ReturnKinematicParameter(int, int) const override {
312  PYBIND11_OVERRIDE_PURE_NAME(
313  double, /* Return type */
314  SampleHandlerBase, /* Parent class */
315  "get_event_kinematic_value", /* python name*/
316  ReturnKinematicParameter, /* Name of function in C++ (must match Python name) */
317  py::arg("variable"),
318  py::arg("event") /* Argument(s) */
319  );
320  }
321 
322  const double *GetPointerToKinematicParameter(int, int) const override {
323  PYBIND11_OVERRIDE_PURE_NAME(
324  const double *, /* Return type */
325  SampleHandlerBase, /* Parent class */
326  "get_event_kinematic_value_reference", /* python name*/
327  GetPointerToKinematicParameter, /* Name of function in C++ (must match Python name) */
328  py::arg("variable"), /* Argument(s) */
329  py::arg("event") /* Argument(s) */
330  );
331  }
332 
334  PYBIND11_OVERRIDE_PURE_NAME(
335  void, /* Return type */
336  SampleHandlerBase, /* Parent class */
337  "register_functional_parameters",
338  RegisterFunctionalParameters /* Name of function in C++ (must match Python name) */
339  );
340  }
341 };
342 
343 void initSamplesModule(py::module &m_samples){
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
void SetupMC() override
Function which translates experiment struct into core struct.
Definition: samples.h:302
void InititialiseData() override
Function responsible for loading data from file or loading from file.
Definition: samples.h:292
void CleanMemoryBeforeFit() override
Allow to clean not used memory before fit starts.
Definition: samples.h:251
void Init() override
Initialise any variables that your experiment specific SampleHandler needs.
Definition: samples.h:271
int SetupExperimentMC() override
Experiment specific setup, returns the number of events which were loaded.
Definition: samples.h:282
const double * GetPointerToKinematicParameter(int, int) const override
Definition: samples.h:322
void SetupSplines() override
initialise your splineXX object and then use InitialiseSplineObject to conviently setup everything up
Definition: samples.h:260
double ReturnKinematicParameter(int, int) const override
Definition: samples.h:311
void RegisterFunctionalParameters() override
HH - a experiment-specific function where the maps to actual functions are set up.
Definition: samples.h:333
void AddAdditionalWeightPointers() override
DB Function to determine which weights apply to which types of samples.
Definition: samples.h:240
EW: As SampleHandlerBase is an abstract base class we have to do some gymnastics to get it to get it ...
Definition: samples.h:20
void Reweight() override
main routine modifying MC prediction based on proposed parameter values
Definition: samples.h:55
std::string GetKinVarName(const int iSample, const int Dimension) const override
Return Kinematic Variable name for specified sample and dimension for example "Reconstructed_Neutrino...
Definition: samples.h:97
std::unique_ptr< TH2 > Get2DVarHist(const int iSample, const std::string &ProjectionVarX, const std::string &ProjectionVarY, const std::vector< KinematicCut > &EventSelectionVec={}, const int WeightStyle=0, const std::vector< KinematicCut > &SubEventSelectionVec={}) override
Build a 2D projection of MC events into specified variables.
Definition: samples.h:196
std::unique_ptr< TH1 > Get1DVarHistByModeAndChannel(const int iSample, const std::string &ProjectionVar_Str, const int kModeToFill=-1, const int kChannelToFill=-1, const int WeightStyle=0) override
Build a 1D histogram for a given variable, optionally filtered by mode and channel.
Definition: samples.h:155
std::vector< double > ReturnKinematicParameterBinning(const int Sample, const std::string &KinematicParameter) const override
Return the binning used to draw a kinematic parameter.
Definition: samples.h:108
TH1 * GetDataHist(const int Sample) override
Get Data histogram.
Definition: samples.h:118
void CleanMemoryBeforeFit() override
Allow to clean not used memory before fit starts.
Definition: samples.h:77
std::string GetSampleTitle(const int iSample) const override
Get fancy title for specified samples.
Definition: samples.h:35
double GetSampleLikelihood(const int iSample) const override
Get likelihood (-logL) for a single sample.
Definition: samples.h:66
int GetNOscChannels(const int iSample) const override
Get number of oscillation channels for a single sample.
Definition: samples.h:45
std::string GetName() const override
Get name for Sample Handler.
Definition: samples.h:26
TH1 * GetMCHist(const int Sample) override
Get MC histogram.
Definition: samples.h:127
std::string GetFlavourName(const int iSample, const int iChannel) const override
Get the flavour name for a given sample and oscillation channel.
Definition: samples.h:220
std::unique_ptr< TH1 > Get1DVarHist(const int iSample, const std::string &ProjectionVar, const std::vector< KinematicCut > &EventSelectionVec={}, const int WeightStyle=0, const std::vector< KinematicCut > &SubEventSelectionVec={}) override
Return 1D projection of MC into given 1D variable (doesn't have to be variable used in the fit)
Definition: samples.h:183
int GetNDim(const int Sample) const override
DB Get what dimensionality binning for given sample has.
Definition: samples.h:211
std::unique_ptr< TH2 > Get2DVarHistByModeAndChannel(const int iSample, const std::string &ProjectionVar_StrX, const std::string &ProjectionVar_StrY, const int kModeToFill=-1, const int kChannelToFill=-1, const int WeightStyle=0) override
Build a 2D histogram for given variables, optionally filtered by mode and channel.
Definition: samples.h:168
void PrintRates(const bool DataOnly=false) override
Helper function to print rates for the samples with LLH.
Definition: samples.h:87
double GetLikelihood() const override
Return likelihood (-logL) for all samples.
Definition: samples.h:145
TH1 * GetW2Hist(const int Sample) override
Get W2 histogram.
Definition: samples.h:136
Class responsible for handling implementation of samples used in analysis, reweighting and returning ...
SampleHandlerBase(std::string ConfigFileName, ParameterHandlerGeneric *xsec_cov, const std::shared_ptr< OscillationHandler > &OscillatorObj_=nullptr)
Constructor.
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
SampleHandlerInterface()
The main constructor.
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.
HW: Histogram utilities for converting ROOT histograms to numpy arrays for use in Python.
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.
void initSamplesModule(py::module &m_samples)
Definition: samples.h:343
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.