MaCh3  2.2.3
Reference Guide
ParameterHandlerGeneric.h
Go to the documentation of this file.
1 #pragma once
2 
3 // MaCh3 includes
6 
13  public:
20  ParameterHandlerGeneric(const std::vector<std::string>& FileNames, std::string name = "xsec_cov", double threshold = -1, int FirstPCAdpar = -999, int LastPCAdpar = -999);
23 
24  // General Getter functions not split by detector
28  inline std::vector<std::string> GetParSampleID(const int i) const { return _fSampleNames[i];};
32  inline std::string GetParamTypeString(const int i) const { return SystType_ToString(_fParamType[i]); }
36  inline SystType GetParamType(const int i) const {return _fParamType[i];}
37 
40  inline SplineInterpolation GetParSplineInterpolation(const int i) const {return SplineParams.at(i)._SplineInterpolationType;}
43  const std::vector<SplineInterpolation> GetSplineInterpolationFromSampleName(const std::string& SampleName);
47  std::string GetParSplineName(const int i) const {return _fSplineNames[i];}
48 
51  const std::vector<int> GetGlobalSystIndexFromSampleName(const std::string& SampleName, const SystType Type);
55  inline double GetParSplineKnotUpperBound(const int i) const {return SplineParams.at(i)._SplineKnotUpBound;}
59  inline double GetParSplineKnotLowerBound(const int i) const {return SplineParams.at(i)._SplineKnotLowBound;}
60 
65  int GetNumParamsFromSampleName(const std::string& SampleName, const SystType Type);
70  const std::vector<std::string> GetParsNamesFromSampleName(const std::string& SampleName, const SystType Type);
75  const std::vector<int> GetParsIndexFromSampleName(const std::string& SampleName, const SystType Type);
76 
79  const std::vector<std::string> GetSplineParsNamesFromSampleName(const std::string& SampleName);
82  const std::vector<std::string> GetSplineFileParsNamesFromSampleName(const std::string& SampleName);
83 
86  const std::vector< std::vector<int> > GetSplineModeVecFromSampleName(const std::string& SampleName);
91  const std::vector<int> GetSystIndexFromSampleName(const std::string& SampleName, const SystType Type) const;
94  const std::vector<NormParameter> GetNormParsFromSampleName(const std::string& SampleName) const;
97  const std::vector<FunctionalParameter> GetFunctionalParametersFromSampleName(const std::string& SampleName) const;
100  const std::vector<SplineParameter> GetSplineParsFromSampleName(const std::string& SampleName) const;
101 
107  bool IsParFromGroup(const int i, const std::string& Group) const;
108 
111  int GetNumParFromGroup(const std::string& Group) const;
114  std::vector<std::string> GetUniqueParameterGroups();
115 
121  void SetGroupOnlyParameters(const std::string& Group, const std::vector<double>& Pars = {});
125  void SetGroupOnlyParameters(const std::vector<std::string>& Groups);
126 
130  void DumpMatrixToFile(const std::string& Name);
131 
134  std::vector<const double*> GetOscParsFromSampleName(const std::string& SampleName);
135 
136  protected:
138  void Print();
140  void PrintGlobablInfo();
142  void PrintNormParams();
144  void PrintSplineParams();
146  void PrintFunctionalParams();
148  void PrintOscillationParams();
150  void PrintParameterGroups();
151 
154 
167  template <typename FilterFunc, typename ActionFunc>
168  void IterateOverParams(const std::string& SampleName, FilterFunc filter, ActionFunc action);
169 
173  void InitParams();
174 
177  inline void InitParametersTypeFromConfig();
178 
182  inline NormParameter GetNormParameter(const YAML::Node& param, const int Index);
183 
187  inline OscillationParameter GetOscillationParameters(const YAML::Node& param, const int Index);
188 
192  inline FunctionalParameter GetFunctionalParameters(const YAML::Node& param, const int Index);
196  inline SplineParameter GetSplineParameter(const YAML::Node& param, const int Index);
201  inline void GetBaseParameter(const YAML::Node& param, const int Index, TypeParameterBase& Parameter);
202 
209  template<typename ParamT>
210  std::vector<ParamT> GetTypeParamsFromSampleName(const std::map<int, int>& indexMap, const std::vector<ParamT>& params, const std::string& SampleName) const;
211 
213  std::vector<SystType> _fParamType;
214 
216  std::vector<std::string> _fSplineNames;
217 
219  std::vector<std::string> _ParameterGroup;
220 
222  std::vector<std::map<int, int>> _fSystToGlobalSystIndexMap;
223 
225  std::vector<SplineParameter> SplineParams;
226 
228  std::vector<NormParameter> NormParams;
229 
231  std::vector<FunctionalParameter> FuncParams;
232 
234  std::vector<OscillationParameter> OscParams;
235 };
SplineInterpolation
Make an enum of the spline interpolation type.
std::string SystType_ToString(const SystType i)
Convert a Syst type type to a string.
SystType
std::vector< std::string > FileNames
Definition: ProcessMCMC.cpp:27
Base class responsible for handling of systematic error parameters. Capable of using PCA or using ada...
std::vector< std::vector< std::string > > _fSampleNames
Tells to which samples object param should be applied.
Class responsible for handling of systematic error parameters with different types defined in the con...
void PrintNormParams()
Prints normalization parameters.
FunctionalParameter GetFunctionalParameters(const YAML::Node &param, const int Index)
Get Func params.
std::vector< SplineParameter > SplineParams
Vector containing info for normalisation systematics.
void PrintParameterGroups()
Prints groups of parameters.
SplineParameter GetSplineParameter(const YAML::Node &param, const int Index)
Get Spline params.
void Print()
Print information about the whole object once it is set.
ParameterHandlerGeneric(const std::vector< std::string > &FileNames, std::string name="xsec_cov", double threshold=-1, int FirstPCAdpar=-999, int LastPCAdpar=-999)
Constructor.
std::vector< std::string > _fSplineNames
Name of spline in TTree (TBranch),.
void PrintOscillationParams()
Prints oscillation parameters.
void PrintSplineParams()
Prints spline parameters.
void InitParams()
Initializes the systematic parameters from the configuration file. This function loads parameters lik...
std::vector< NormParameter > NormParams
Vector containing info for normalisation systematics.
std::vector< std::map< int, int > > _fSystToGlobalSystIndexMap
Map between number of given parameter type with global parameter numbering. For example 2nd norm para...
std::vector< SystType > _fParamType
Type of parameter like norm, spline etc.
void PrintFunctionalParams()
Prints functional parameters.
void CheckCorrectInitialisation()
KS: Check if matrix is correctly initialised.
void PrintGlobablInfo()
Prints general information about the ParameterHandler object.
NormParameter GetNormParameter(const YAML::Node &param, const int Index)
Get Norm params.
void InitParametersTypeFromConfig()
Parses the YAML configuration to set up cross-section parameters. The YAML file defines the types of ...
std::vector< ParamT > GetTypeParamsFromSampleName(const std::map< int, int > &indexMap, const std::vector< ParamT > &params, const std::string &SampleName) const
Retrieve parameters that apply to a given sample name.
std::vector< FunctionalParameter > FuncParams
Vector containing info for functional systematics.
std::vector< OscillationParameter > OscParams
Vector containing info for functional systematics.
void DumpMatrixToFile(const std::string &Name)
Dump Matrix to ROOT file, useful when we need to pass matrix info to another fitting group.
SplineInterpolation GetParSplineInterpolation(const int i) const
Get interpolation type for a given parameter.
OscillationParameter GetOscillationParameters(const YAML::Node &param, const int Index)
Get Osc params.
void IterateOverParams(const std::string &SampleName, FilterFunc filter, ActionFunc action)
Iterates over parameters and applies a filter and action function.
void GetBaseParameter(const YAML::Node &param, const int Index, TypeParameterBase &Parameter)
Fill base parameters.
std::vector< std::string > _ParameterGroup
KS: Allow to group parameters for example to affect only cross-section or only flux etc.
const std::vector< NormParameter > GetNormParsFromSampleName(const std::string &SampleName) const
DB Get norm/func parameters depending on given SampleName.
SystType GetParamType(const int i) const
Returns enum describing our param type.
const std::vector< int > GetSystIndexFromSampleName(const std::string &SampleName, const SystType Type) const
Grab the index of the syst relative to global numbering.
const std::vector< int > GetGlobalSystIndexFromSampleName(const std::string &SampleName, const SystType Type)
DB Get spline parameters depending on given SampleName.
std::vector< std::string > GetParSampleID(const int i) const
ETA - just return the int of the SampleName, this can be removed to do a string comp at some point.
double GetParSplineKnotUpperBound(const int i) const
EM: value at which we cap spline knot weight.
const std::vector< std::string > GetSplineFileParsNamesFromSampleName(const std::string &SampleName)
DB Get spline parameters depending on given SampleName.
const std::vector< int > GetParsIndexFromSampleName(const std::string &SampleName, const SystType Type)
DB Grab the parameter indices for the relevant SampleName.
const std::vector< SplineParameter > GetSplineParsFromSampleName(const std::string &SampleName) const
KS: Grab the Spline parameters for the relevant SampleName.
std::string GetParSplineName(const int i) const
Get the name of the spline associated with the spline at index i.
const std::vector< FunctionalParameter > GetFunctionalParametersFromSampleName(const std::string &SampleName) const
HH Get functional parameters for the relevant SampleName.
int GetNumParamsFromSampleName(const std::string &SampleName, const SystType Type)
DB Grab the number of parameters for the relevant SampleName.
const std::vector< std::string > GetParsNamesFromSampleName(const std::string &SampleName, const SystType Type)
DB Grab the parameter names for the relevant SampleName.
std::vector< std::string > GetUniqueParameterGroups()
KS: Get names of all unique parameter groups.
const std::vector< std::string > GetSplineParsNamesFromSampleName(const std::string &SampleName)
DB Get spline parameters depending on given SampleName.
std::string GetParamTypeString(const int i) const
ETA - just return a string of "spline", "norm" or "functional".
bool IsParFromGroup(const int i, const std::string &Group) const
Checks if parameter belongs to a given group.
const std::vector< SplineInterpolation > GetSplineInterpolationFromSampleName(const std::string &SampleName)
Get the interpolation types for splines affecting a particular SampleName.
double GetParSplineKnotLowerBound(const int i) const
EM: value at which we cap spline knot weight.
std::vector< const double * > GetOscParsFromSampleName(const std::string &SampleName)
Get pointers to Osc params from Sample name.
const std::vector< std::vector< int > > GetSplineModeVecFromSampleName(const std::string &SampleName)
DB Grab the Spline Modes for the relevant SampleName.
int GetNumParFromGroup(const std::string &Group) const
KS: Check how many parameters are associated with given group.
void SetGroupOnlyParameters(const std::string &Group, const std::vector< double > &Pars={})
KS Function to set to prior parameters of a given group or values from vector.
HH - Functional parameters Carrier for whether you want to apply a systematic to an event or not.
ETA - Normalisations for cross-section parameters Carrier for whether you want to apply a systematic ...
KS: Struct holding info about oscillation Systematics.
KS: Struct holding info about Spline Systematics.
Base class storing info for parameters types, helping unify codebase.