MaCh3  2.4.2
Reference Guide
SampleHandlerBase.h
Go to the documentation of this file.
1 #pragma once
2 
3 //C++ includes
4 #include <assert.h>
5 
6 //MaCh3 includes
9 #include "Manager/Manager.h"
10 #include "Manager/MaCh3Modes.h"
11 
13 //ROOT includes
14 #include "TTree.h"
15 #include "TH1D.h"
16 #include "TH2D.h"
17 #include "TMath.h"
18 #include "TFile.h"
19 #include "TROOT.h"
20 #include "TRandom3.h"
21 #include "TString.h"
23 
27 {
28  public:
32  virtual ~SampleHandlerBase();
33 
34  virtual inline M3::int_t GetNsamples(){ return nSamples; };
35  virtual std::string GetSampleTitle(const int Sample) const = 0;
36  virtual std::string GetName() const = 0;
37  virtual double GetSampleLikelihood(const int isample) const = 0;
39  virtual void CleanMemoryBeforeFit() = 0;
41  virtual void SaveAdditionalInfo(TDirectory* Dir) {(void) Dir;};
43  MaCh3Modes* GetMaCh3Modes() const { return Modes.get(); }
44 
45  virtual void Reweight()=0;
46  virtual double GetLikelihood() const = 0;
47 
50  virtual void PrintRates(const bool DataOnly = false) = 0;
51 
52  unsigned int GetNEvents() const {return nEvents;}
53  virtual int GetNOscChannels(const int iSample) const = 0;
54 
58  virtual std::string GetKinVarName(const int iSample, const int Dimension) const = 0;
59 
61  virtual TH1* GetDataHist(const int Sample) = 0;
63  virtual TH1* GetMCHist(const int Sample) = 0;
65  virtual TH1* GetW2Hist(const int Sample) = 0;
66 
68  virtual int GetNDim(const int Sample) const = 0;
69  virtual std::string GetFlavourName(const int iSample, const int iChannel) const = 0;
70 
72  virtual std::vector<double> ReturnKinematicParameterBinning(const int Sample, const std::string &KinematicParameter) const = 0;
73 
74  virtual TH1* Get1DVarHistByModeAndChannel(const int iSample, const std::string& ProjectionVar_Str,
75  int kModeToFill = -1, int kChannelToFill = -1, int WeightStyle = 0, TAxis* Axis = nullptr) = 0;
76 
77  virtual TH2* Get2DVarHistByModeAndChannel(const int iSample, const std::string& ProjectionVar_StrX,
78  const std::string& ProjectionVar_StrY, int kModeToFill = -1,
79  int kChannelToFill = -1, int WeightStyle = 0,
80  TAxis* AxisX = nullptr, TAxis* AxisY = nullptr) = 0;
81  virtual TH1 *Get1DVarHist(const int iSample, const std::string &ProjectionVar,
82  const std::vector<KinematicCut> &EventSelectionVec = {}, int WeightStyle = 0,
83  TAxis *Axis = nullptr, const std::vector<KinematicCut> &SubEventSelectionVec = {}) = 0;
84  virtual TH2* Get2DVarHist(const int iSample, const std::string& ProjectionVarX, const std::string& ProjectionVarY,
85  const std::vector< KinematicCut >& EventSelectionVec = {},
86  int WeightStyle = 0, TAxis* AxisX = nullptr, TAxis* AxisY = nullptr,
87  const std::vector< KinematicCut >& SubEventSelectionVec = {}) = 0;
88 
92  double GetPoissonLLH(const double data, const double mc) const;
93 
178  double GetTestStatLLH(const double data, const double mc, const double w2) const;
179 
182  void SetTestStatistic(TestStatistic testStat){ fTestStatistic = testStat; }
185 
186 protected:
188  void QuietPlease();
190  void NowTalk();
191 
193  template <typename T>
194  bool MatchCondition(const std::vector<T>& allowedValues, const T& value) {
195  if (allowedValues.empty()) {
196  return true; // Apply to all if no specific values are specified
197  }
198  return std::find(allowedValues.begin(), allowedValues.end(), value) != allowedValues.end();
199  }
200 
203 
205  std::streambuf *buf;
207  std::streambuf *errbuf;
208 
211 
213  unsigned int nEvents;
214 
216  std::unique_ptr<MaCh3Modes> Modes;
217 };
#define _MaCh3_Safe_Include_Start_
KS: Avoiding warning checking for headers.
Definition: Core.h:126
#define _MaCh3_Safe_Include_End_
KS: Restore warning checking after including external headers.
Definition: Core.h:140
TestStatistic
Make an enum of the test statistic that we're using.
KS: Class describing MaCh3 modes used in the analysis, it is being initialised from config.
Definition: MaCh3Modes.h:135
Class responsible for handling implementation of samples used in analysis, reweighting and returning ...
virtual ~SampleHandlerBase()
destructor
void SetTestStatistic(TestStatistic testStat)
Set the test statistic to be used when calculating the binned likelihoods.
virtual std::string GetKinVarName(const int iSample, const int Dimension) const =0
Return Kinematic Variable name for specified sample and dimension for example "Reconstructed_Neutrino...
void NowTalk()
CW: Redirect std::cout to silence some experiment specific libraries.
std::unique_ptr< MaCh3Modes > Modes
Holds information about used Generator and MaCh3 modes.
virtual void Reweight()=0
std::streambuf * buf
Keep the cout buffer.
MaCh3Modes * GetMaCh3Modes() const
Return pointer to MaCh3 modes.
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 TH1 * GetW2Hist(const int Sample)=0
Get W2 histogram.
virtual std::string GetFlavourName(const int iSample, const int iChannel) const =0
virtual double GetSampleLikelihood(const int isample) const =0
virtual TH2 * Get2DVarHist(const int iSample, const std::string &ProjectionVarX, const std::string &ProjectionVarY, const std::vector< KinematicCut > &EventSelectionVec={}, int WeightStyle=0, TAxis *AxisX=nullptr, TAxis *AxisY=nullptr, const std::vector< KinematicCut > &SubEventSelectionVec={})=0
TestStatistic GetTestStatistic() const
Get the test statistic used when calculating the binned likelihoods.
SampleHandlerBase()
The main constructor.
std::streambuf * errbuf
Keep the cerr buffer.
unsigned int GetNEvents() const
virtual TH1 * Get1DVarHistByModeAndChannel(const int iSample, const std::string &ProjectionVar_Str, int kModeToFill=-1, int kChannelToFill=-1, int WeightStyle=0, TAxis *Axis=nullptr)=0
virtual TH1 * GetMCHist(const int Sample)=0
Get MC histogram.
virtual TH1 * Get1DVarHist(const int iSample, const std::string &ProjectionVar, const std::vector< KinematicCut > &EventSelectionVec={}, int WeightStyle=0, TAxis *Axis=nullptr, const std::vector< KinematicCut > &SubEventSelectionVec={})=0
TestStatistic fTestStatistic
Test statistic tells what kind of likelihood sample is using.
virtual int GetNOscChannels(const int iSample) const =0
M3::int_t nSamples
Contains how many samples we've got.
virtual std::vector< double > ReturnKinematicParameterBinning(const int Sample, const std::string &KinematicParameter) const =0
Return the binning used to draw a kinematic parameter.
unsigned int nEvents
Number of MC events are there.
virtual void CleanMemoryBeforeFit()=0
Allow to clean not used memory before fit starts.
virtual void SaveAdditionalInfo(TDirectory *Dir)
Store additional info in a chan.
virtual M3::int_t GetNsamples()
virtual double GetLikelihood() const =0
bool MatchCondition(const std::vector< T > &allowedValues, const T &value)
check if event is affected by following conditions, for example pdg, or modes etc
void QuietPlease()
CW: Redirect std::cout to silence some experiment specific libraries.
virtual std::string GetSampleTitle(const int Sample) const =0
virtual std::string GetName() const =0
virtual TH2 * Get2DVarHistByModeAndChannel(const int iSample, const std::string &ProjectionVar_StrX, const std::string &ProjectionVar_StrY, int kModeToFill=-1, int kChannelToFill=-1, int WeightStyle=0, TAxis *AxisX=nullptr, TAxis *AxisY=nullptr)=0
virtual int GetNDim(const int Sample) const =0
DB Function to differentiate 1D or 2D binning.
virtual TH1 * GetDataHist(const int Sample)=0
Get Data histogram.
double GetPoissonLLH(const double data, const double mc) const
Calculate test statistic for a single bin using Poisson.
virtual void PrintRates(const bool DataOnly=false)=0
Helper function to print rates for the samples with LLH.
int int_t
Definition: Core.h:38