MaCh3  2.5.1
Reference Guide
ParameterHandlerGeneric.h
Go to the documentation of this file.
1 #pragma once
2 
3 // MaCh3 includes
6 
14  public:
21  ParameterHandlerGeneric(const std::vector<std::string>& FileNames, std::string name = "xsec_cov", double threshold = -1, int FirstPCAdpar = -999, int LastPCAdpar = -999);
24 
25  // General Getter functions not split by detector
28  std::vector<std::string> GetParSampleID(const int i) const { return _fSampleNames[i];};
31  std::string GetParamTypeString(const int i) const { return SystType_ToString(_fParamType[i]); }
34  SystType GetParamType(const int i) const {return _fParamType[i];}
35 
38  SplineInterpolation GetParSplineInterpolation(const int i) const {return SplineParams.at(i)._SplineInterpolationType;}
40  const std::vector<SplineInterpolation> GetSplineInterpolationFromSampleName(const std::string& SampleName) const;
43  std::string GetParSplineName(const int i) const {return SplineParams.at(i)._fSplineNames;}
44 
46  const std::vector<int> GetGlobalSystIndexFromSampleName(const std::string& SampleName, const SystType Type) const;
49  double GetParSplineKnotUpperBound(const int i) const {return SplineParams.at(i)._SplineKnotUpBound;}
52  double GetParSplineKnotLowerBound(const int i) const {return SplineParams.at(i)._SplineKnotLowBound;}
53 
57  int GetNumParamsFromSampleName(const std::string& SampleName, const SystType Type) const;
61  const std::vector<std::string> GetParsNamesFromSampleName(const std::string& SampleName, const SystType Type) const;
65  const std::vector<int> GetParsIndexFromSampleName(const std::string& SampleName, const SystType Type) const;
66 
68  const std::vector<std::string> GetSplineParsNamesFromSampleName(const std::string& SampleName) const;
70  const std::vector<std::string> GetSplineFileParsNamesFromSampleName(const std::string& SampleName) const;
71 
73  const std::vector< std::vector<int> > GetSplineModeVecFromSampleName(const std::string& SampleName) const;
77  const std::vector<int> GetSystIndexFromSampleName(const std::string& SampleName, const SystType Type) const;
79  const std::vector<NormParameter> GetNormParsFromSampleName(const std::string& SampleName) const;
81  const std::vector<FunctionalParameter> GetFunctionalParametersFromSampleName(const std::string& SampleName) const;
83  const std::vector<SplineParameter> GetSplineParsFromSampleName(const std::string& SampleName) const;
84 
89  bool IsParFromGroup(const int i, const std::string& Group) const;
90 
92  int GetNumParFromGroup(const std::string& Group) const;
94  std::vector<std::string> GetUniqueParameterGroups() const;
95 
100  void SetGroupOnlyParameters(const std::string& Group, const std::vector<double>& Pars = {});
103  void SetGroupOnlyParameters(const std::vector<std::string>& Groups);
104 
107  void SetFixGroupOnlyParameters(const std::string& Group);
110  void SetFixGroupOnlyParameters(const std::vector<std::string>& Groups);
111 
114  void SetFreeGroupOnlyParameters(const std::string& Group);
117  void SetFreeGroupOnlyParameters(const std::vector<std::string>& Groups);
118 
122  void DumpMatrixToFile(const std::string& Name);
123 
125  std::vector<const M3::float_t*> GetOscParsFromSampleName(const std::string& SampleName) const;
126 
127  protected:
130  void InitialiseFromConfig(const std::vector<std::string>& YAMLFile);
132  void Print() const;
134  void PrintGlobablInfo() const;
136  void PrintNormParams() const;
138  void PrintSplineParams() const;
140  void PrintFunctionalParams() const;
142  void PrintOscillationParams() const;
144  void PrintParameterGroups() const;
145 
147  void CheckCorrectInitialisation() const;
149  void LoadCorrelationFromConfig(std::vector<std::map<std::string,double>>& Correlations,
150  std::map<std::string, int>& CorrNamesMap);
152  void LoadAndMergeYAML(const std::vector<std::string>& YAMLFile,
153  std::map<std::pair<int,int>, std::unique_ptr<TMatrixDSym>>& overrides);
166  template <typename FilterFunc, typename ActionFunc>
167  void IterateOverParams(const std::string& SampleName, FilterFunc filter, ActionFunc action) const;
168 
172  void InitParameters();
173 
177 
181  NormParameter GetNormParameter(const YAML::Node& param, const int Index) const;
182 
186  OscillationParameter GetOscillationParameters(const YAML::Node& param, const int Index) const;
187 
191  FunctionalParameter GetFunctionalParameters(const YAML::Node& param, const int Index) const;
195  SplineParameter GetSplineParameter(const YAML::Node& param, const int Index) const;
200  void GetBaseParameter(const YAML::Node& param, const int Index, TypeParameterBase& Parameter) const;
201 
205  bool AppliesToSample(const int SystIndex, const std::string& SampleName) const;
206 
213  template<typename ParamT>
214  std::vector<ParamT> GetTypeParamsFromSampleName(const std::map<int, int>& indexMap, const std::vector<ParamT>& params, const std::string& SampleName) const;
215 
217  std::vector<std::vector<std::string>> _fSampleNames;
218 
220  std::vector<SystType> _fParamType;
221 
223  std::vector<std::string> _ParameterGroup;
224 
226  std::vector<std::map<int, int>> _fSystToGlobalSystIndexMap;
227 
229  std::vector<SplineParameter> SplineParams;
230 
232  std::vector<NormParameter> NormParams;
233 
235  std::vector<FunctionalParameter> FuncParams;
236 
238  std::vector<OscillationParameter> OscParams;
239 };
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:28
Base class responsible for handling of systematic error parameters. Capable of using PCA or using ada...
Class responsible for handling of systematic error parameters with different types defined in the con...
void PrintGlobablInfo() const
Prints general information about the ParameterHandler object.
void Print() const
Print information about the whole object once it is set.
std::vector< SplineParameter > SplineParams
Vector containing info for normalisation systematics.
FunctionalParameter GetFunctionalParameters(const YAML::Node &param, const int Index) const
Get Func params.
const std::vector< NormParameter > GetNormParsFromSampleName(const std::string &SampleName) const
DB Get norm/func parameters depending on given SampleName.
std::vector< std::string > GetUniqueParameterGroups() const
KS: Get names of all unique parameter groups.
ParameterHandlerGeneric(const std::vector< std::string > &FileNames, std::string name="xsec_cov", double threshold=-1, int FirstPCAdpar=-999, int LastPCAdpar=-999)
Constructor.
void PrintSplineParams() const
Prints spline parameters.
std::vector< const M3::float_t * > GetOscParsFromSampleName(const std::string &SampleName) const
Get pointers to Osc params from Sample name.
void IterateOverParams(const std::string &SampleName, FilterFunc filter, ActionFunc action) const
Iterates over parameters and applies a filter and action function.
const std::vector< int > GetParsIndexFromSampleName(const std::string &SampleName, const SystType Type) const
DB Grab the parameter indices for the relevant SampleName.
SplineParameter GetSplineParameter(const YAML::Node &param, const int Index) const
Get Spline params.
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.
void GetBaseParameter(const YAML::Node &param, const int Index, TypeParameterBase &Parameter) const
Fill base parameters.
std::vector< NormParameter > NormParams
Vector containing info for normalisation systematics.
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.
std::vector< std::map< int, int > > _fSystToGlobalSystIndexMap
Map between number of given parameter type with global parameter numbering. For example 2nd norm para...
void LoadCorrelationFromConfig(std::vector< std::map< std::string, double >> &Correlations, std::map< std::string, int > &CorrNamesMap)
Load correlation from yaml config.
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.
std::vector< SystType > _fParamType
Type of parameter like norm, spline etc.
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.
void PrintOscillationParams() const
Prints oscillation parameters.
const std::vector< std::string > GetSplineFileParsNamesFromSampleName(const std::string &SampleName) const
DB Get spline parameters depending on given SampleName.
const std::vector< FunctionalParameter > GetFunctionalParametersFromSampleName(const std::string &SampleName) const
HH Get functional parameters for the relevant SampleName.
void SetFreeGroupOnlyParameters(const std::string &Group)
TN Method to set parameters within a group to be treated as free.
void InitParametersTypeFromConfig()
Parses the YAML configuration to set up cross-section parameters. The YAML file defines the types of ...
std::vector< std::vector< std::string > > _fSampleNames
Tells to which samples object param should be applied.
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.
int GetNumParamsFromSampleName(const std::string &SampleName, const SystType Type) const
DB Grab the number of parameters for the relevant SampleName.
void PrintFunctionalParams() const
Prints functional parameters.
const std::vector< std::string > GetParsNamesFromSampleName(const std::string &SampleName, const SystType Type) const
DB Grab the parameter names for the relevant SampleName.
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.
const std::vector< std::string > GetSplineParsNamesFromSampleName(const std::string &SampleName) const
DB Get spline parameters depending on given SampleName.
SplineInterpolation GetParSplineInterpolation(const int i) const
Get interpolation type for a given parameter.
void InitParameters()
Initializes the systematic parameters from the configuration file. This function loads parameters lik...
const std::vector< int > GetGlobalSystIndexFromSampleName(const std::string &SampleName, const SystType Type) const
DB Get spline parameters depending on given SampleName.
OscillationParameter GetOscillationParameters(const YAML::Node &param, const int Index) const
Get Osc params.
NormParameter GetNormParameter(const YAML::Node &param, const int Index) const
Get Norm params.
bool AppliesToSample(const int SystIndex, const std::string &SampleName) const
Check if parameter is affecting given sample name.
const std::vector< SplineInterpolation > GetSplineInterpolationFromSampleName(const std::string &SampleName) const
Get the interpolation types for splines affecting a particular 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.
void CheckCorrectInitialisation() const
KS: Check if matrix is correctly initialised.
void InitialiseFromConfig(const std::vector< std::string > &YAMLFile)
Initialisation of the class using config.
void SetFixGroupOnlyParameters(const std::string &Group)
TN Method to set parameters within a group to be fixed to their prior values.
std::vector< std::string > _ParameterGroup
KS: Allow to group parameters for example to affect only cross-section or only flux etc.
void LoadAndMergeYAML(const std::vector< std::string > &YAMLFile, std::map< std::pair< int, int >, std::unique_ptr< TMatrixDSym >> &overrides)
Initialise single yaml based on several yaml files.
double GetParSplineKnotLowerBound(const int i) const
EM: value at which we cap spline knot weight.
void PrintNormParams() const
Prints normalization parameters.
int GetNumParFromGroup(const std::string &Group) const
KS: Check how many parameters are associated with given group.
const std::vector< std::vector< int > > GetSplineModeVecFromSampleName(const std::string &SampleName) const
DB Grab the Spline Modes for the relevant SampleName.
void PrintParameterGroups() const
Prints groups of parameters.
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.