MaCh3  2.5.0
Reference Guide
BinnedSplineHandler.h
Go to the documentation of this file.
1 #pragma once
2 
3 //MaCh3 includes
4 #include "Splines/SplineBase.h"
5 #include "Manager/MaCh3Modes.h"
6 
13  public:
18  virtual ~BinnedSplineHandler();
19 
23  void Evaluate() final;
24 
26  void AddSample(const std::string& SampleName,
27  const std::string& SampleTitle,
28  const std::vector<std::string>& OscChanFileNames,
29  const std::vector<std::string>& SplineVarNames);
31  void TransferToMonolith();
33  void cleanUpMemory();
34 
37  virtual void FillSampleArray(const std::string& SampleTitle, const std::vector<std::string>& OscChanFileNames);
40  std::vector< std::vector<int> > StripDuplicatedModes(const std::vector< std::vector<int> >& InputVector) const;
42  std::vector< std::vector<int> > GetEventSplines(const std::string& SampleTitle, int iOscChan, int EventMode, double Var1Val, double Var2Val, double Var3Val);
44  void SynchroniseMemTransfer() const final {return;}
46  std::vector<TAxis*> FindSplineBinning(const std::string& FileName, const std::string& SampleTitle);
47 
48  int CountNumberOfLoadedSplines(bool NonFlat=false, int Verbosity=0);
49  std::string getDimLabel(const int BinningOpt, const unsigned int Axis) const;
52  int GetSampleIndex(const std::string& SampleTitle) const;
54  bool isValidSplineIndex(const std::string& SampleTitle, int iSyst, int iOscChan, int iMode, int iVar1, int iVar2, int iVar3) const;
55 
56  void BuildSampleIndexingArray(const std::string& SampleTitle);
57  void PrepForReweight();
58  void getSplineCoeff_SepMany(int splineindex, M3::float_t *& xArray, M3::float_t *&manyArray);
59  void PrintBinning(TAxis* Axis) const;
61  void PrintSampleDetails(const std::string& SampleTitle) const;
62  void PrintArrayDetails(const std::string& SampleTitle) const;
63 
65  const M3::float_t* RetPointer(const int sample, const int oscchan, const int syst, const int mode,
66  const int var1bin, const int var2bin, const int var3bin) const {
67  int index = indexvec[sample][oscchan][syst][mode][var1bin][var2bin][var3bin];
68  return &weightvec_Monolith[index];
69  }
71  void PrepareSplineFile(std::string FileName) final;
74  void LoadSplineFile(std::string FileName) final;
75 
76  protected:
78  void CalcSplineWeights() final;
81 
82  //And now the actual member variables
83  std::vector<std::string> SampleNames;
84  std::vector<std::string> SampleTitles;
85  std::vector<int> Dimensions;
86  std::vector<std::vector<std::string>> DimensionLabels;
87  std::vector<int> nSplineParams;
88  std::vector<int> nOscChans;
89 
90  std::vector< std::vector< std::vector<TAxis*> > > SplineBinning;
91  std::vector< std::vector<std::string> > SplineFileParPrefixNames;
95  std::vector< std::vector< std::vector<int> > > SplineModeVecs;
99  std::vector< std::vector<int> > GlobalSystIndex;
102  std::vector< std::vector<SplineInterpolation> > SplineInterpolationTypes;
103 
105  std::vector<std::string> UniqueSystNames;
107  std::vector<int> UniqueSystIndices;
108 
110  std::vector< std::vector< std::vector< std::vector< std::vector< std::vector< std::vector< int > > > > > > > indexvec;
111  std::vector<int > coeffindexvec;
113  std::vector<int> uniquecoeffindices;
114 
117 
121 
123  bool *isflatarray;
128 
132  std::vector<int> uniquesplinevec_Monolith;
133 
137  virtual std::vector<std::string> GetTokensFromSplineName(const std::string& FullSplineName) = 0;
138 
139  private:
141  void InvestigateMissingSplines() const;
142 
145  void LoadSettingsDir(std::unique_ptr<TFile>& SplineFile);
148  void LoadMonolithDir(std::unique_ptr<TFile>& SplineFile);
151  void LoadIndexDir(std::unique_ptr<TFile>& SplineFile);
152 
155  void PrepareSettingsDir(std::unique_ptr<TFile>& SplineFile) const;
158  void PrepareMonolithDir(std::unique_ptr<TFile>& SplineFile) const;
161  void PrepareIndexDir(std::unique_ptr<TFile>& SplineFile) const;
164  void PrepareOtherInfoDir(std::unique_ptr<TFile>& SplineFile) const;
165 };
SplineInterpolation
Make an enum of the spline interpolation type.
Bin-by-bin class calculating response for spline parameters.
void getSplineCoeff_SepMany(int splineindex, M3::float_t *&xArray, M3::float_t *&manyArray)
virtual void FillSampleArray(const std::string &SampleTitle, const std::vector< std::string > &OscChanFileNames)
Loads and processes splines from ROOT files for a given sample.
std::string getDimLabel(const int BinningOpt, const unsigned int Axis) const
bool * isflatarray
Need to keep track of which splines are flat and which aren't.
void LoadIndexDir(std::unique_ptr< TFile > &SplineFile)
KS: Load preprocessed Index.
int CountNumberOfLoadedSplines(bool NonFlat=false, int Verbosity=0)
std::vector< std::vector< std::vector< std::vector< std::vector< std::vector< std::vector< int > > > > > > > indexvec
Variables related to determined which modes have splines and which piggy-back of other modes.
ParameterHandlerGeneric * ParHandler
Pointer to covariance from which we get information about spline params.
void CalcSplineWeights() final
CPU based code which eval weight for each spline.
std::vector< std::vector< std::string > > SplineFileParPrefixNames
void LoadSettingsDir(std::unique_ptr< TFile > &SplineFile)
KS: Load preprocessed Settings.
std::vector< int > coeffindexvec
void cleanUpMemory()
Remove setup variables not needed for spline evaluations.
bool isValidSplineIndex(const std::string &SampleTitle, int iSyst, int iOscChan, int iMode, int iVar1, int iVar2, int iVar3) const
Ensure we have spline for a given bin.
void PrintBinning(TAxis *Axis) const
void PrepareSplineFile(std::string FileName) final
KS: Prepare spline file that can be used for fast loading.
void InvestigateMissingSplines() const
This function will find missing splines in file.
std::vector< int > nSplineParams
std::vector< std::vector< std::vector< int > > > SplineModeVecs
std::vector< int > nOscChans
M3::float_t * xcoeff_arr
x coefficients for each spline
std::vector< std::string > SampleTitles
void PrintSampleDetails(const std::string &SampleTitle) const
Print info like Sample ID of spline params etc.
std::vector< std::string > UniqueSystNames
name of each spline parameter
std::vector< std::vector< SplineInterpolation > > SplineInterpolationTypes
spline interpolation types for each sample. These vectors are from a call to GetSplineInterpolationFr...
std::vector< std::vector< int > > GetEventSplines(const std::string &SampleTitle, int iOscChan, int EventMode, double Var1Val, double Var2Val, double Var3Val)
Return the splines which affect a given event.
void Evaluate() final
CW: This Eval should be used when using two separate x,{y,a,b,c,d} arrays to store the weights; proba...
std::vector< int > Dimensions
std::vector< int > uniquesplinevec_Monolith
Maps single spline object with single parameter.
void PrepareIndexDir(std::unique_ptr< TFile > &SplineFile) const
KS: Prepare Index Info within SplineFile.
void BuildSampleIndexingArray(const std::string &SampleTitle)
std::vector< M3::float_t > weightvec_Monolith
Stores weight from spline evaluation for each single spline.
const M3::float_t * RetPointer(const int sample, const int oscchan, const int syst, const int mode, const int var1bin, const int var2bin, const int var3bin) const
get pointer to spline weight based on bin variables
void AddSample(const std::string &SampleName, const std::string &SampleTitle, const std::vector< std::string > &OscChanFileNames, const std::vector< std::string > &SplineVarNames)
add oscillation channel to spline monolith
std::vector< TAxis * > FindSplineBinning(const std::string &FileName, const std::string &SampleTitle)
Grab histograms with spline binning.
std::vector< std::vector< std::string > > DimensionLabels
std::vector< std::vector< int > > StripDuplicatedModes(const std::vector< std::vector< int > > &InputVector) const
Check if there are any repeated modes. This is used to reduce the number of modes in case many intera...
void LoadSplineFile(std::string FileName) final
KS: Load preprocessed spline file.
void PrepareMonolithDir(std::unique_ptr< TFile > &SplineFile) const
KS: Prepare Monolith Info within SplineFile.
BinnedSplineHandler(ParameterHandlerGeneric *xsec_, MaCh3Modes *Modes_)
Constructor.
void LoadMonolithDir(std::unique_ptr< TFile > &SplineFile)
KS: Load preprocessed Monolith.
std::vector< std::vector< std::vector< TAxis * > > > SplineBinning
std::vector< TSpline3_red * > splinevec_Monolith
holds each spline object before stripping into coefficient monolith
std::vector< int > uniquecoeffindices
Unique coefficient indices.
void PrepareOtherInfoDir(std::unique_ptr< TFile > &SplineFile) const
KS: Prepare Other Info within SplineFile.
void PrepareSettingsDir(std::unique_ptr< TFile > &SplineFile) const
KS: Prepare Settings Info within SplineFile.
MaCh3Modes * Modes
pointer to MaCh3 Mode from which we get spline suffix
void PrintArrayDetails(const std::string &SampleTitle) const
std::vector< std::vector< int > > GlobalSystIndex
This holds the global spline index and is used to grab the current parameter value to evaluate spline...
virtual ~BinnedSplineHandler()
Destructor.
virtual std::vector< std::string > GetTokensFromSplineName(const std::string &FullSplineName)=0
std::vector< std::string > SampleNames
void TransferToMonolith()
flatten multidimensional spline array into proper monolith
std::vector< int > UniqueSystIndices
Global index of each spline param, it allows us to match spline ordering with global.
int GetSampleIndex(const std::string &SampleTitle) const
Get index of sample based on name.
void SynchroniseMemTransfer() const final
KS: After calculations are done on GPU we copy memory to CPU. This operation is asynchronous meaning ...
M3::float_t * manycoeff_arr
ybcd coefficients for each spline
KS: Class describing MaCh3 modes used in the analysis, it is being initialised from config.
Definition: MaCh3Modes.h:135
Class responsible for handling of systematic error parameters with different types defined in the con...
Base class for calculating weight from spline.
Definition: SplineBase.h:27
CW: Reduced TSpline3 class.
Main namespace for MaCh3 software.
double float_t
Definition: Core.h:37