MaCh3  2.5.1
Reference Guide
SampleHandlerInterface.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 
35 {
36  public:
40  virtual ~SampleHandlerInterface();
41 
43  virtual M3::int_t GetNSamples(){ return nSamples; };
46  virtual std::string GetSampleTitle(const int iSample) const = 0;
48  virtual std::string GetName() const = 0;
51  virtual double GetSampleLikelihood(const int iSample) const = 0;
53  virtual void CleanMemoryBeforeFit() = 0;
56  virtual void SaveAdditionalInfo(TDirectory* Dir) {(void) Dir;};
58  MaCh3Modes* GetMaCh3Modes() const { return Modes.get(); }
60  virtual void Reweight()=0;
62  virtual double GetLikelihood() const = 0;
63 
66  virtual void PrintRates(const bool DataOnly = false) = 0;
67 
69  unsigned int GetNEvents() const {return nEvents;}
72  virtual int GetNOscChannels(const int iSample) const = 0;
73 
77  virtual std::string GetKinVarName(const int iSample, const int Dimension) const = 0;
78 
81  virtual const TH1* GetDataHist(const int Sample) = 0;
84  virtual const TH1* GetMCHist(const int Sample) = 0;
87  virtual const TH1* GetW2Hist(const int Sample) = 0;
88 
91  virtual int GetNDim(const int Sample) const = 0;
95  virtual std::string GetFlavourName(const int iSample, const int iChannel) const = 0;
96 
100  virtual std::vector<double> ReturnKinematicParameterBinning(const int Sample, const std::string &KinematicParameter) const = 0;
101 
108  virtual std::unique_ptr<TH1> Get1DVarHistByModeAndChannel(const int iSample, const std::string& ProjectionVar_Str,
109  const int kModeToFill = -1, const int kChannelToFill = -1,
110  const int WeightStyle = 0) = 0;
118  virtual std::unique_ptr<TH2> Get2DVarHistByModeAndChannel(const int iSample, const std::string& ProjectionVar_StrX,
119  const std::string& ProjectionVar_StrY, int kModeToFill = -1,
120  const int kChannelToFill = -1, const int WeightStyle = 0) = 0;
127  virtual std::unique_ptr<TH1> Get1DVarHist(const int iSample, const std::string &ProjectionVar,
128  const std::vector<KinematicCut> &EventSelectionVec = {}, int WeightStyle = 0,
129  const std::vector<KinematicCut> &SubEventSelectionVec = {}) = 0;
137  virtual std::unique_ptr<TH2> Get2DVarHist(const int iSample, const std::string& ProjectionVarX,
138  const std::string& ProjectionVarY,
139  const std::vector< KinematicCut >& EventSelectionVec = {},
140  const int WeightStyle = 0, const std::vector< KinematicCut >& SubEventSelectionVec = {}) = 0;
141 
145  double GetPoissonLLH(const double data, const double mc) const;
146 
231  double GetTestStatLLH(const double data, const double mc, const double w2) const;
232 
235  void SetTestStatistic(TestStatistic testStat){ fTestStatistic = testStat; }
238 
239 protected:
241  void QuietPlease();
243  void NowTalk();
244 
246  template <typename T>
247  bool MatchCondition(const std::vector<T>& allowedValues, const T& value) {
248  if (allowedValues.empty()) {
249  return true; // Apply to all if no specific values are specified
250  }
251  return std::find(allowedValues.begin(), allowedValues.end(), value) != allowedValues.end();
252  }
253 
256 
258  std::streambuf *buf;
260  std::streambuf *errbuf;
261 
264 
266  unsigned int nEvents;
267 
269  std::unique_ptr<MaCh3Modes> Modes;
270 };
#define _MaCh3_Safe_Include_Start_
KS: Avoiding warning checking for headers.
Definition: Core.h:126
#define _MaCh3_Safe_Include_End_
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:142
Class responsible for handling implementation of samples used in analysis, reweighting and returning ...
virtual 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)=0
Build a 1D histogram for a given variable, optionally filtered by mode and channel.
virtual std::string GetName() const =0
Get name for Sample Handler.
virtual void SaveAdditionalInfo(TDirectory *Dir)
Store additional info in a chain.
virtual std::unique_ptr< TH1 > Get1DVarHist(const int iSample, const std::string &ProjectionVar, const std::vector< KinematicCut > &EventSelectionVec={}, int WeightStyle=0, const std::vector< KinematicCut > &SubEventSelectionVec={})=0
Return 1D projection of MC into given 1D variable (doesn't have to be variable used in the fit)
virtual const TH1 * GetMCHist(const int Sample)=0
Get MC histogram.
TestStatistic fTestStatistic
Test statistic tells what kind of likelihood sample is using.
virtual const TH1 * GetW2Hist(const int Sample)=0
Get W2 histogram.
TestStatistic GetTestStatistic() const
Get the test statistic used when calculating the binned likelihoods.
std::streambuf * errbuf
Keep the cerr buffer.
virtual double GetLikelihood() const =0
Return likelihood (-logL) for all samples.
bool MatchCondition(const std::vector< T > &allowedValues, const T &value)
check if event is affected by following conditions, for example pdg, or modes etc
std::streambuf * buf
Keep the cout buffer.
virtual std::string GetFlavourName(const int iSample, const int iChannel) const =0
Get the flavour name for a given sample and oscillation channel.
double GetPoissonLLH(const double data, const double mc) const
Calculate test statistic for a single bin using Poisson.
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 const TH1 * GetDataHist(const int Sample)=0
Get Data histogram.
virtual std::vector< double > ReturnKinematicParameterBinning(const int Sample, const std::string &KinematicParameter) const =0
Return the binning used to draw a kinematic parameter.
MaCh3Modes * GetMaCh3Modes() const
Return pointer to MaCh3 modes.
virtual int GetNOscChannels(const int iSample) const =0
Get number of oscillation channels for a single sample.
M3::int_t nSamples
Contains how many samples we've got.
unsigned int nEvents
Number of MC events are there.
virtual M3::int_t GetNSamples()
returns total number of samples
std::unique_ptr< MaCh3Modes > Modes
Holds information about used Generator and MaCh3 modes.
virtual std::unique_ptr< TH2 > Get2DVarHistByModeAndChannel(const int iSample, const std::string &ProjectionVar_StrX, const std::string &ProjectionVar_StrY, int kModeToFill=-1, const int kChannelToFill=-1, const int WeightStyle=0)=0
Build a 2D histogram for given variables, optionally filtered by mode and channel.
unsigned int GetNEvents() const
Return total number of events.
virtual double GetSampleLikelihood(const int iSample) const =0
Get likelihood (-logL) for a single sample.
virtual void Reweight()=0
main routine modifying MC prediction based on proposed parameter values
virtual ~SampleHandlerInterface()
destructor
void NowTalk()
CW: Redirect std::cout to silence some experiment specific libraries.
virtual void PrintRates(const bool DataOnly=false)=0
Helper function to print rates for the samples with LLH.
virtual std::string GetSampleTitle(const int iSample) const =0
Get fancy title for specified samples.
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...
SampleHandlerInterface()
The main constructor.
void QuietPlease()
CW: Redirect std::cout to silence some experiment specific libraries.
void SetTestStatistic(TestStatistic testStat)
Set the test statistic to be used when calculating the binned likelihoods.
virtual 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={})=0
Build a 2D projection of MC events into specified variables.
virtual int GetNDim(const int Sample) const =0
DB Get what dimensionality binning for given sample has.
virtual void CleanMemoryBeforeFit()=0
Allow to clean not used memory before fit starts.
int int_t
Definition: Core.h:38