MaCh3  2.5.0
Reference Guide
Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
ParameterHandlerGeneric Class Reference

Class responsible for handling of systematic error parameters with different types defined in the config. Like spline, normalisation parameters etc. More...

#include <Parameters/ParameterHandlerGeneric.h>

Inheritance diagram for ParameterHandlerGeneric:
[legend]
Collaboration diagram for ParameterHandlerGeneric:
[legend]

Public Member Functions

 ParameterHandlerGeneric (const std::vector< std::string > &FileNames, std::string name="xsec_cov", double threshold=-1, int FirstPCAdpar=-999, int LastPCAdpar=-999)
 Constructor. More...
 
 ~ParameterHandlerGeneric ()
 Destructor. More...
 
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. More...
 
std::string GetParamTypeString (const int i) const
 ETA - just return a string of "spline", "norm" or "functional". More...
 
SystType GetParamType (const int i) const
 Returns enum describing our param type. More...
 
SplineInterpolation GetParSplineInterpolation (const int i) const
 Get interpolation type for a given parameter. More...
 
const std::vector< SplineInterpolationGetSplineInterpolationFromSampleName (const std::string &SampleName)
 Get the interpolation types for splines affecting a particular SampleName. More...
 
std::string GetParSplineName (const int i) const
 Get the name of the spline associated with the spline at index i. More...
 
const std::vector< int > GetGlobalSystIndexFromSampleName (const std::string &SampleName, const SystType Type)
 DB Get spline parameters depending on given SampleName. More...
 
double GetParSplineKnotUpperBound (const int i) const
 EM: value at which we cap spline knot weight. More...
 
double GetParSplineKnotLowerBound (const int i) const
 EM: value at which we cap spline knot weight. More...
 
int GetNumParamsFromSampleName (const std::string &SampleName, const SystType Type)
 DB Grab the number of parameters for the relevant SampleName. More...
 
const std::vector< std::string > GetParsNamesFromSampleName (const std::string &SampleName, const SystType Type)
 DB Grab the parameter names for the relevant SampleName. More...
 
const std::vector< int > GetParsIndexFromSampleName (const std::string &SampleName, const SystType Type)
 DB Grab the parameter indices for the relevant SampleName. More...
 
const std::vector< std::string > GetSplineParsNamesFromSampleName (const std::string &SampleName)
 DB Get spline parameters depending on given SampleName. More...
 
const std::vector< std::string > GetSplineFileParsNamesFromSampleName (const std::string &SampleName)
 DB Get spline parameters depending on given SampleName. More...
 
const std::vector< std::vector< int > > GetSplineModeVecFromSampleName (const std::string &SampleName)
 DB Grab the Spline Modes for the relevant SampleName. More...
 
const std::vector< int > GetSystIndexFromSampleName (const std::string &SampleName, const SystType Type) const
 Grab the index of the syst relative to global numbering. More...
 
const std::vector< NormParameterGetNormParsFromSampleName (const std::string &SampleName) const
 DB Get norm/func parameters depending on given SampleName. More...
 
const std::vector< FunctionalParameterGetFunctionalParametersFromSampleName (const std::string &SampleName) const
 HH Get functional parameters for the relevant SampleName. More...
 
const std::vector< SplineParameterGetSplineParsFromSampleName (const std::string &SampleName) const
 KS: Grab the Spline parameters for the relevant SampleName. More...
 
bool IsParFromGroup (const int i, const std::string &Group) const
 Checks if parameter belongs to a given group. More...
 
int GetNumParFromGroup (const std::string &Group) const
 KS: Check how many parameters are associated with given group. More...
 
std::vector< std::string > GetUniqueParameterGroups () const
 KS: Get names of all unique parameter groups. More...
 
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. More...
 
void SetGroupOnlyParameters (const std::vector< std::string > &Groups)
 KS Function to set to prior parameters of a given groups or values from vector. More...
 
void SetFixGroupOnlyParameters (const std::string &Group)
 TN Method to set parameters within a group to be fixed to their prior values. More...
 
void SetFixGroupOnlyParameters (const std::vector< std::string > &Groups)
 TN Method to set parameters of certain groups to be fixed to their prior values. More...
 
void SetFreeGroupOnlyParameters (const std::string &Group)
 TN Method to set parameters within a group to be treated as free. More...
 
void SetFreeGroupOnlyParameters (const std::vector< std::string > &Groups)
 TN Method to set parameters of certain groups to be treated as free. More...
 
void ToggleFixGroupOnlyParameters (const std::string &Group)
 TN Method to toggle fix/free parameters within a group. More...
 
void ToggleFixGroupOnlyParameters (const std::vector< std::string > &Groups)
 TN Method to toggle fix/free parameters within given groups. More...
 
void DumpMatrixToFile (const std::string &Name)
 Dump Matrix to ROOT file, useful when we need to pass matrix info to another fitting group. More...
 
std::vector< const M3::float_t * > GetOscParsFromSampleName (const std::string &SampleName)
 Get pointers to Osc params from Sample name. More...
 
- Public Member Functions inherited from ParameterHandlerBase
 ParameterHandlerBase (const std::vector< std::string > &YAMLFile, std::string name, double threshold=-1, int FirstPCAdpar=-999, int LastPCAdpar=-999)
 ETA - constructor for a YAML file. More...
 
 ParameterHandlerBase (std::string name, std::string file, double threshold=-1, int FirstPCAdpar=-999, int LastPCAdpar=-999)
 "Usual" constructors from root file More...
 
virtual ~ParameterHandlerBase ()
 Destructor. More...
 
void SetCovMatrix (TMatrixDSym *cov)
 Set covariance matrix. More...
 
void SetName (const std::string &name)
 Set matrix name. More...
 
void SetParName (const int i, const std::string &name)
 change parameter name More...
 
void SetSingleParameter (const int parNo, const double parVal)
 Set value of single param to a given value. More...
 
void SetPar (const int i, const double val)
 Set all the covariance matrix parameters to a user-defined value. More...
 
void SetParCurrProp (const int i, const double val)
 Set current parameter value. More...
 
void SetParProp (const int i, const double val)
 Set proposed parameter value. More...
 
void SetParameters (const std::vector< double > &pars={})
 Set parameter values using vector, it has to have same size as covariance class. More...
 
void SetFlatPrior (const int i, const bool eL)
 Set if parameter should have flat prior or not. More...
 
void SetRandomThrow (const int i, const double rand)
 Set random value useful for debugging/CI. More...
 
double GetRandomThrow (const int i) const
 Get random value useful for debugging/CI. More...
 
void SetBranches (TTree &tree, const bool SaveProposal=false)
 set branches for output file More...
 
void SetStepScale (const double scale, const bool verbose=true)
 Set global step scale for covariance object. More...
 
void SetIndivStepScale (const int ParameterIndex, const double StepScale)
 DB Function to set fIndivStepScale from a vector (Can be used from execs and inside covariance constructors) More...
 
void SetIndivStepScale (const std::vector< double > &stepscale)
 DB Function to set fIndivStepScale from a vector (Can be used from execs and inside covariance constructors) More...
 
void SetPrintLength (const unsigned int PriLen)
 KS: In case someone really want to change this. More...
 
void SaveUpdatedMatrixConfig ()
 KS: After step scale, prefit etc. value were modified save this modified config. More...
 
void ThrowParameters ()
 Throw the parameters according to the covariance matrix. This shouldn't be used in MCMC code ase it can break Detailed Balance;. More...
 
void RandomConfiguration ()
 Randomly throw the parameters in their 1 sigma range. More...
 
int CheckBounds () const _noexcept_
 Check if parameters were proposed outside physical boundary. More...
 
double CalcLikelihood () const _noexcept_
 Calc penalty term based on inverted covariance matrix. More...
 
virtual double GetLikelihood ()
 Return CalcLikelihood if some params were thrown out of boundary return LARGE_LOGL More...
 
TMatrixDSym * GetCovMatrix () const
 Return covariance matrix. More...
 
TMatrixDSym * GetInvCovMatrix () const
 Return inverted covariance matrix. More...
 
double GetInvCovMatrix (const int i, const int j) const
 Return inverted covariance matrix. More...
 
double GetCorrThrows (const int i) const
 Return correlated throws. More...
 
bool GetFlatPrior (const int i) const
 Get if param has flat prior or not. More...
 
std::string GetName () const
 Get name of covariance. More...
 
std::string GetParName (const int i) const
 Get name of parameter. More...
 
int GetParIndex (const std::string &name) const
 Get index based on name. More...
 
std::string GetParFancyName (const int i) const
 Get fancy name of the Parameter. More...
 
std::string GetInputFile () const
 Get name of input file. More...
 
double GetDiagonalError (const int i) const
 Get diagonal error for ith parameter. More...
 
double GetError (const int i) const
 Get the error for the ith parameter. More...
 
void ResetIndivStepScale ()
 Adaptive Step Tuning Stuff. More...
 
void SetIndivStepScaleForSkippedAdaptParams ()
 Set individual step scale for parameters which are skipped during adaption to initial values. More...
 
void InitialiseAdaption (const YAML::Node &adapt_manager)
 Initialise adaptive MCMC. More...
 
void SaveAdaptiveToFile (const std::string &outFileName, const std::string &systematicName)
 Save adaptive throw matrix to file. More...
 
bool GetDoAdaption () const
 Do we adapt or not. More...
 
void SetThrowMatrix (TMatrixDSym *cov)
 Use new throw matrix, used in adaptive MCMC. More...
 
void SetSubThrowMatrix (int first_index, int last_index, TMatrixDSym const &subcov)
 
void UpdateThrowMatrix (TMatrixDSym *cov)
 Replaces old throw matrix with new one. More...
 
void SetNumberOfSteps (const int nsteps)
 Set number of MCMC step, when running adaptive MCMC it is updated with given frequency. We need number of steps to determine frequency. More...
 
TMatrixDSym * GetThrowMatrix () const
 Get matrix used for step proposal. More...
 
double GetThrowMatrix (const int i, const int j) const
 Get matrix used for step proposal. More...
 
TH2D * GetCorrelationMatrix ()
 KS: Convert covariance matrix to correlation matrix and return TH2D which can be used for fancy plotting. More...
 
const M3::float_tRetPointer (const int iParam)
 DB Pointer return to param position. More...
 
const std::vector< M3::float_t > & GetParPropVec ()
 Get a reference to the proposed parameter values Can be useful if you want to track these without having to copy values using getProposed() More...
 
int GetNumParams () const
 Get total number of parameters. More...
 
std::vector< double > GetPreFitValues () const
 Get the pre-fit values of the parameters. More...
 
std::vector< double > GetProposed () const
 Get vector of all proposed parameter values. More...
 
M3::float_t GetParProp (const int i) const
 Get proposed parameter value. More...
 
double GetParCurr (const int i) const
 Get current parameter value. More...
 
const std::vector< double > & GetParCurrVec () const
 Get vector of current parameter values. More...
 
double GetParInit (const int i) const
 Get prior parameter value. More...
 
double GetUpperBound (const int i) const
 Get upper parameter bound in which it is physically valid. More...
 
double GetLowerBound (const int i) const
 Get lower parameter bound in which it is physically valid. More...
 
double GetIndivStepScale (const int ParameterIndex) const
 Get individual step scale for selected parameter. More...
 
double GetGlobalStepScale () const
 Get global step scale for covariance object. More...
 
int GetNParameters () const
 Get number of params which will be different depending if using Eigen decomposition or not. More...
 
void PrintNominal () const
 Print prior value for every parameter. More...
 
void PrintNominalCurrProp () const
 Print prior, current and proposed value for each parameter. More...
 
void PrintParameters () const
 
void PrintIndivStepScale () const
 Print step scale for each parameter. More...
 
virtual void ProposeStep ()
 Generate a new proposed state. More...
 
void Randomize () _noexcept_
 "Randomize" the parameters in the covariance class for the proposed step. Used the proposal kernel and the current parameter value to set proposed step More...
 
void CorrelateSteps () _noexcept_
 Use Cholesky throw matrix for better step proposal. More...
 
void UpdateAdaptiveCovariance ()
 Method to update adaptive MCMC [15]. More...
 
void AcceptStep () _noexcept_
 Accepted this step. More...
 
void SetFixAllParameters ()
 Set all parameters to be fixed at prior values. More...
 
void SetFixParameter (const int i)
 Set parameter to be fixed at prior value. More...
 
void SetFixParameter (const std::string &name)
 Set parameter to be fixed at prior value. More...
 
void SetFreeAllParameters ()
 Set all parameters to be treated as free. More...
 
void SetFreeParameter (const int i)
 Set parameter to be treated as free. More...
 
void SetFreeParameter (const std::string &name)
 Set parameter to be treated as free. More...
 
void ToggleFixAllParameters ()
 Toggle fixing parameters at prior values. More...
 
void ToggleFixParameter (const int i)
 Toggle fixing parameter at prior values. More...
 
void ToggleFixParameter (const std::string &name)
 Toggle fixing parameter at prior values. More...
 
bool IsParameterFixed (const int i) const
 Is parameter fixed or not. More...
 
bool IsParameterFixed (const std::string &name) const
 Is parameter fixed or not. More...
 
void ConstructPCA (const double eigen_threshold, int FirstPCAdpar, int LastPCAdpar)
 CW: Calculate eigen values, prepare transition matrices and remove param based on defined threshold. More...
 
bool IsPCA () const
 is PCA, can use to query e.g. LLH scans More...
 
YAML::Node GetConfig () const
 Getter to return a copy of the YAML node. More...
 
AdaptiveMCMCHandlerGetAdaptiveHandler () const
 Get pointer for AdaptiveHandler. More...
 
void SetTune (const std::string &TuneName)
 KS: Set proposed parameter values vector to be base on tune values, for example set proposed values to be of generated or maybe PostND. More...
 
PCAHandlerGetPCAHandler () const
 Get pointer for PCAHandler. More...
 
void MatchMaCh3OutputBranches (TTree *PosteriorFile, std::vector< double > &BranchValues, std::vector< std::string > &BranchNames, const std::vector< std::string > &FancyNames={})
 Matches branches in a TTree to parameters in a systematic handler. More...
 

Protected Member Functions

void Print () const
 Print information about the whole object once it is set. More...
 
void PrintGlobablInfo () const
 Prints general information about the ParameterHandler object. More...
 
void PrintNormParams () const
 Prints normalization parameters. More...
 
void PrintSplineParams () const
 Prints spline parameters. More...
 
void PrintFunctionalParams () const
 Prints functional parameters. More...
 
void PrintOscillationParams () const
 Prints oscillation parameters. More...
 
void PrintParameterGroups () const
 Prints groups of parameters. More...
 
void CheckCorrectInitialisation () const
 KS: Check if matrix is correctly initialised. More...
 
template<typename FilterFunc , typename ActionFunc >
void IterateOverParams (const std::string &SampleName, FilterFunc filter, ActionFunc action)
 Iterates over parameters and applies a filter and action function. More...
 
void InitParams ()
 Initializes the systematic parameters from the configuration file. This function loads parameters like normalizations and splines from the provided YAML file. More...
 
void InitParametersTypeFromConfig ()
 Parses the YAML configuration to set up cross-section parameters. The YAML file defines the types of systematic errors, interpolation types, and bounds for splines. More...
 
NormParameter GetNormParameter (const YAML::Node &param, const int Index)
 Get Norm params. More...
 
OscillationParameter GetOscillationParameters (const YAML::Node &param, const int Index)
 Get Osc params. More...
 
FunctionalParameter GetFunctionalParameters (const YAML::Node &param, const int Index)
 Get Func params. More...
 
SplineParameter GetSplineParameter (const YAML::Node &param, const int Index)
 Get Spline params. More...
 
void GetBaseParameter (const YAML::Node &param, const int Index, TypeParameterBase &Parameter)
 Fill base parameters. More...
 
template<typename ParamT >
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. More...
 
- Protected Member Functions inherited from ParameterHandlerBase
void Init (const std::string &name, const std::string &file)
 Initialisation of the class using matrix from root file. More...
 
void Init (const std::vector< std::string > &YAMLFile)
 Initialisation of the class using config. More...
 
void ReserveMemory (const int size)
 Initialise vectors with parameters information. More...
 
void MakePosDef (TMatrixDSym *cov=nullptr)
 Make matrix positive definite by adding small values to diagonal, necessary for inverting matrix. More...
 
void MakeClosestPosDef (TMatrixDSym *cov)
 HW: Finds closest possible positive definite matrix in Frobenius Norm ||.||_frob Where ||X||_frob=sqrt[sum_ij(x_ij^2)] (basically just turns an n,n matrix into vector in n^2 space then does Euclidean norm) More...
 
void SetThrowMatrixFromFile (const std::string &matrix_file_name, const std::string &matrix_name, const std::string &means_name)
 sets throw matrix from a file More...
 
bool AppliesToSample (const int SystIndex, const std::string &SampleName) const
 Check if parameter is affecting given sample name. More...
 
void FlipParameterValue (const int index, const double FlipPoint)
 KS: Flip parameter around given value, for example mass ordering around 0. More...
 
void CircularParBounds (const int i, const double LowBound, const double UpBound)
 HW :: This method is a tad hacky but modular arithmetic gives me a headache. More...
 
void EnableSpecialProposal (const YAML::Node &param, const int Index)
 Enable special proposal. More...
 
void SpecialStepProposal ()
 Perform Special Step Proposal. More...
 

Protected Attributes

std::vector< SystType_fParamType
 Type of parameter like norm, spline etc. More...
 
std::vector< std::string > _ParameterGroup
 KS: Allow to group parameters for example to affect only cross-section or only flux etc. More...
 
std::vector< std::map< int, int > > _fSystToGlobalSystIndexMap
 Map between number of given parameter type with global parameter numbering. For example 2nd norm param may be 10-th global param. More...
 
std::vector< SplineParameterSplineParams
 Vector containing info for normalisation systematics. More...
 
std::vector< NormParameterNormParams
 Vector containing info for normalisation systematics. More...
 
std::vector< FunctionalParameterFuncParams
 Vector containing info for functional systematics. More...
 
std::vector< OscillationParameterOscParams
 Vector containing info for functional systematics. More...
 
- Protected Attributes inherited from ParameterHandlerBase
bool doSpecialStepProposal
 Check if any of special step proposal were enabled. More...
 
const std::string inputFile
 The input root file we read in. More...
 
std::string matrixName
 Name of cov matrix. More...
 
TMatrixDSym * covMatrix
 The covariance matrix. More...
 
TMatrixDSym * invCovMatrix
 The inverse covariance matrix. More...
 
std::vector< std::vector< double > > InvertCovMatrix
 KS: Same as above but much faster as TMatrixDSym cache miss. More...
 
std::vector< std::unique_ptr< TRandom3 > > random_number
 KS: Set Random numbers for each thread so each thread has different seed. More...
 
double * randParams
 Random number taken from gaussian around prior error used for corr_throw. More...
 
double * corr_throw
 Result of multiplication of Cholesky matrix and randParams. More...
 
double _fGlobalStepScale
 Global step scale applied to all params in this class. More...
 
int PrintLength
 KS: This is used when printing parameters, sometimes we have super long parameters name, we want to flexibly adjust couts. More...
 
std::vector< std::string > _fNames
 ETA _fNames is set automatically in the covariance class to be something like param_i, this is currently to make things compatible with the Diagnostic tools. More...
 
std::vector< std::string > _fFancyNames
 Fancy name for example rather than param_0 it is MAQE, useful for human reading. More...
 
YAML::Node _fYAMLDoc
 Stores config describing systematics. More...
 
int _fNumPar
 Number of systematic parameters. More...
 
std::vector< double > _fPreFitValue
 Parameter value dictated by the prior model. Based on it penalty term is calculated. More...
 
std::vector< double > _fCurrVal
 Current value of the parameter. More...
 
std::vector< M3::float_t_fPropVal
 Proposed value of the parameter. More...
 
std::vector< double > _fError
 Prior error on the parameter. More...
 
std::vector< double > _fLowBound
 Lowest physical bound, parameter will not be able to go beyond it. More...
 
std::vector< double > _fUpBound
 Upper physical bound, parameter will not be able to go beyond it. More...
 
std::vector< double > _fIndivStepScale
 Individual step scale used by MCMC algorithm. More...
 
std::vector< bool > _fFlatPrior
 Whether to apply flat prior or not. More...
 
std::vector< std::vector< std::string > > _fSampleNames
 Tells to which samples object param should be applied. More...
 
std::vector< double > _fIndivStepScaleInitial
 Backup of _fIndivStepScale for parameters which are skipped during adaption. More...
 
double _fGlobalStepScaleInitial
 Backup of _fGlobalStepScale for parameters which are skipped during adaption. More...
 
std::vector< bool > param_skip_adapt_flags
 Flags telling if parameter should be skipped during adaption. More...
 
TMatrixDSym * throwMatrix
 Matrix which we use for step proposal before Cholesky decomposition (not actually used for step proposal) More...
 
double ** throwMatrixCholDecomp
 Throw matrix that is being used in the fit, much faster as TMatrixDSym cache miss. More...
 
bool pca
 perform PCA or not More...
 
bool use_adaptive
 Are we using AMCMC? More...
 
std::unique_ptr< PCAHandlerPCAObj
 Struct containing information about PCA. More...
 
std::unique_ptr< AdaptiveMCMCHandlerAdaptiveHandler
 Struct containing information about adaption. More...
 
std::unique_ptr< ParameterTunesTunes
 Struct containing information about adaption. More...
 
std::vector< int > FlipParameterIndex
 Indices of parameters with flip symmetry. More...
 
std::vector< double > FlipParameterPoint
 Central points around which parameters are flipped. More...
 
std::vector< int > CircularBoundsIndex
 Indices of parameters with circular bounds. More...
 
std::vector< std::pair< double, double > > CircularBoundsValues
 Circular bounds for each parameter (lower, upper) More...
 

Detailed Description

Class responsible for handling of systematic error parameters with different types defined in the config. Like spline, normalisation parameters etc.

Author
Dan Barrow
Ed Atkin
Kamil Skwarczynski

Definition at line 13 of file ParameterHandlerGeneric.h.

Constructor & Destructor Documentation

◆ ParameterHandlerGeneric()

ParameterHandlerGeneric::ParameterHandlerGeneric ( const std::vector< std::string > &  FileNames,
std::string  name = "xsec_cov",
double  threshold = -1,
int  FirstPCAdpar = -999,
int  LastPCAdpar = -999 
)

Constructor.

Parameters
FileNamesA vector of strings representing the YAML files used for initialisation of matrix
nameMatrix name
thresholdPCA threshold from 0 to 1. Default is -1 and means no PCA
FirstPCAdparFirst PCA parameter that will be decomposed.
LastPCAdparFirst PCA parameter that will be decomposed.

Definition at line 7 of file ParameterHandlerGeneric.cpp.

8  : ParameterHandlerBase(YAMLFile, name, threshold, FirstPCA, LastPCA){
9 // ********************************************
11 
12  //ETA - again this really doesn't need to be hear...
13  for (int i = 0; i < _fNumPar; i++)
14  {
15  // Sort out the print length
16  if(int(_fNames[i].length()) > PrintLength) PrintLength = int(_fNames[i].length());
17  } // end the for loop
18 
19  MACH3LOG_DEBUG("Constructing instance of ParameterHandler");
20  InitParams();
21  // Print
22  Print();
23 }
#define MACH3LOG_DEBUG
Definition: MaCh3Logger.h:34
int _fNumPar
Number of systematic parameters.
ParameterHandlerBase(const std::vector< std::string > &YAMLFile, std::string name, double threshold=-1, int FirstPCAdpar=-999, int LastPCAdpar=-999)
ETA - constructor for a YAML file.
std::vector< std::string > _fNames
ETA _fNames is set automatically in the covariance class to be something like param_i,...
int PrintLength
KS: This is used when printing parameters, sometimes we have super long parameters name,...
void Print() const
Print information about the whole object once it is set.
void InitParams()
Initializes the systematic parameters from the configuration file. This function loads parameters lik...
void InitParametersTypeFromConfig()
Parses the YAML configuration to set up cross-section parameters. The YAML file defines the types of ...

◆ ~ParameterHandlerGeneric()

ParameterHandlerGeneric::~ParameterHandlerGeneric ( )

Destructor.

Definition at line 99 of file ParameterHandlerGeneric.cpp.

99  {
100 // ********************************************
101  MACH3LOG_DEBUG("Deleting ParameterHandler");
102 }

Member Function Documentation

◆ CheckCorrectInitialisation()

void ParameterHandlerGeneric::CheckCorrectInitialisation ( ) const
protected

KS: Check if matrix is correctly initialised.

Definition at line 625 of file ParameterHandlerGeneric.cpp.

625  {
626 // ********************************************
627  // KS: Lambda Function which simply checks if there are no duplicates in std::vector
628  auto CheckForDuplicates = [](const std::vector<std::string>& names, const std::string& nameType, bool Warning) {
629  std::unordered_map<std::string, size_t> seenStrings;
630  for (size_t i = 0; i < names.size(); ++i) {
631  const auto& name = names[i];
632  if (seenStrings.find(name) != seenStrings.end()) {
633  size_t firstIndex = seenStrings[name];
634  if(Warning){
635  MACH3LOG_WARN("There are two systematics with the same {} '{}', first at index {}, and again at index {}", nameType, name, firstIndex, i);
636  MACH3LOG_WARN("Is this desired?");
637  return;
638  } else {
639  MACH3LOG_CRITICAL("There are two systematics with the same {} '{}', first at index {}, and again at index {}", nameType, name, firstIndex, i);
640  throw MaCh3Exception(__FILE__, __LINE__);
641  }
642  }
643  seenStrings[name] = i;
644  }
645  };
646  std::vector<std::string> SplineNameTemp(SplineParams.size());
647  for(size_t it = 0; it < SplineParams.size(); it++){
648  SplineNameTemp[it] = SplineParams[it]._fSplineNames;
649  }
650  // KS: Checks if there are no duplicates in fancy names etc, this can happen if we merge configs etc
651  CheckForDuplicates(_fFancyNames, "_fFancyNames", false);
652  CheckForDuplicates(SplineNameTemp, "_fSplineNames", true);
653 }
#define MACH3LOG_CRITICAL
Definition: MaCh3Logger.h:38
#define MACH3LOG_WARN
Definition: MaCh3Logger.h:36
Custom exception class used throughout MaCh3.
std::vector< std::string > _fFancyNames
Fancy name for example rather than param_0 it is MAQE, useful for human reading.
std::vector< SplineParameter > SplineParams
Vector containing info for normalisation systematics.

◆ DumpMatrixToFile()

void ParameterHandlerGeneric::DumpMatrixToFile ( const std::string &  Name)

Dump Matrix to ROOT file, useful when we need to pass matrix info to another fitting group.

Parameters
NameName of TFile to which we save stuff
Warning
This is mostly used for backward compatibility

Definition at line 786 of file ParameterHandlerGeneric.cpp.

786  {
787 // ********************************************
788  // KS: Ideally we remove it eventually...
789  TH2D* CorrMatrix = GetCorrelationMatrix();
792  _fError,
793  _fLowBound,
794  _fUpBound,
796  _fFancyNames,
797  _fFlatPrior,
798  SplineParams,
799  covMatrix,
800  CorrMatrix,
801  Name);
802  delete CorrMatrix;
803  MACH3LOG_INFO("Finished dumping ParameterHandler object");
804 }
#define MACH3LOG_INFO
Definition: MaCh3Logger.h:35
std::vector< bool > _fFlatPrior
Whether to apply flat prior or not.
std::vector< double > _fError
Prior error on the parameter.
TH2D * GetCorrelationMatrix()
KS: Convert covariance matrix to correlation matrix and return TH2D which can be used for fancy plott...
std::vector< double > _fLowBound
Lowest physical bound, parameter will not be able to go beyond it.
TMatrixDSym * covMatrix
The covariance matrix.
std::vector< double > _fPreFitValue
Parameter value dictated by the prior model. Based on it penalty term is calculated.
std::vector< double > _fUpBound
Upper physical bound, parameter will not be able to go beyond it.
std::vector< double > _fIndivStepScale
Individual step scale used by MCMC algorithm.
void DumpParamHandlerToFile(const int _fNumPar, const std::vector< double > &_fPreFitValue, const std::vector< double > &_fError, const std::vector< double > &_fLowBound, const std::vector< double > &_fUpBound, const std::vector< double > &_fIndivStepScale, const std::vector< std::string > &_fFancyNames, const std::vector< bool > &_fFlatPrior, const std::vector< SplineParameter > &SplineParams, TMatrixDSym *covMatrix, TH2D *CorrMatrix, const std::string &Name)
Dump Matrix to ROOT file, useful when we need to pass matrix info to another fitting group.

◆ GetBaseParameter()

void ParameterHandlerGeneric::GetBaseParameter ( const YAML::Node &  param,
const int  Index,
TypeParameterBase Parameter 
)
inlineprotected

Fill base parameters.

Parameters
paramYaml node describing param
IndexGlobal parameter index
ParameterObject storing info

Definition at line 178 of file ParameterHandlerGeneric.cpp.

178  {
179 // ********************************************
180  Parameter.name = GetParFancyName(Index);
181 
182  // Set the global parameter index of the normalisation parameter
183  Parameter.index = Index;
184 
185  int NumKinematicCuts = 0;
186  if(param["KinematicCuts"]) {
187  NumKinematicCuts = int(param["KinematicCuts"].size());
188 
189  std::vector<std::string> TempKinematicStrings;
190  std::vector<std::vector<std::vector<double>>> TempKinematicBounds;
191  //First element of TempKinematicBounds is always -999, and size is then 3
192  for(int KinVar_i = 0 ; KinVar_i < NumKinematicCuts ; ++KinVar_i) {
193  //ETA: This is a bit messy, Kinematic cuts is a list of maps
194  for (YAML::const_iterator it = param["KinematicCuts"][KinVar_i].begin();it!=param["KinematicCuts"][KinVar_i].end();++it) {
195  TempKinematicStrings.push_back(it->first.as<std::string>());
196  TempKinematicBounds.push_back(Get2DBounds(it->second));
197  }
198  if(TempKinematicStrings.size() == 0) {
199  MACH3LOG_ERROR("Received a KinematicCuts node but couldn't read the contents (it's a list of single-element dictionaries (python) = map of pairs (C++))");
200  MACH3LOG_ERROR("For Param {}", Parameter.name);
201  throw MaCh3Exception(__FILE__, __LINE__);
202  }
203  }//KinVar_i
204  Parameter.KinematicVarStr = TempKinematicStrings;
205  Parameter.Selection = TempKinematicBounds;
206  }
207 
208  //Next ones are kinematic bounds on where normalisation parameter should apply
209  //We set a bool to see if any bounds exist so we can short-circuit checking all of them every step
210  bool HasKinBounds = false;
211 
212  if(Parameter.KinematicVarStr.size() > 0) HasKinBounds = true;
213 
214  Parameter.hasKinBounds = HasKinBounds;
215  //End of kinematic bound checking
216 }
#define MACH3LOG_ERROR
Definition: MaCh3Logger.h:37
#define Get2DBounds(filename)
Definition: YamlHelper.h:591
std::string GetParFancyName(const int i) const
Get fancy name of the Parameter.
int index
Parameter number of this normalisation in current systematic model.
bool hasKinBounds
Does this parameter have kinematic bounds.
std::string name
Name of parameters.
std::vector< std::string > KinematicVarStr
std::vector< std::vector< std::vector< double > > > Selection

◆ GetFunctionalParameters()

FunctionalParameter ParameterHandlerGeneric::GetFunctionalParameters ( const YAML::Node &  param,
const int  Index 
)
inlineprotected

Get Func params.

Parameters
paramYaml node describing param
IndexGlobal parameter index

Definition at line 283 of file ParameterHandlerGeneric.cpp.

283  {
284 // ********************************************
285  FunctionalParameter func;
286  GetBaseParameter(param, Index, func);
287 
288  func.modes = GetFromManager<std::vector<int>>(param["Mode"], std::vector<int>(), __FILE__ , __LINE__);
289 
290  func.valuePtr = RetPointer(Index);
291  return func;
292 }
const M3::float_t * RetPointer(const int iParam)
DB Pointer return to param position.
void GetBaseParameter(const YAML::Node &param, const int Index, TypeParameterBase &Parameter)
Fill base parameters.
HH - Functional parameters Carrier for whether you want to apply a systematic to an event or not.
std::vector< int > modes
Mode which parameter applies to.
const M3::float_t * valuePtr
Parameter value pointer.

◆ GetFunctionalParametersFromSampleName()

const std::vector< FunctionalParameter > ParameterHandlerGeneric::GetFunctionalParametersFromSampleName ( const std::string &  SampleName) const

HH Get functional parameters for the relevant SampleName.

Definition at line 306 of file ParameterHandlerGeneric.cpp.

306  {
307 // ********************************************
309 }
@ kFunc
For functional parameters.
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< 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.

◆ GetGlobalSystIndexFromSampleName()

const std::vector< int > ParameterHandlerGeneric::GetGlobalSystIndexFromSampleName ( const std::string &  SampleName,
const SystType  Type 
)

DB Get spline parameters depending on given SampleName.

Definition at line 221 of file ParameterHandlerGeneric.cpp.

221  {
222 // ********************************************
223  std::vector<int> returnVec;
224  for (auto &pair : _fSystToGlobalSystIndexMap[Type]) {
225  auto &SystIndex = pair.second;
226  if (AppliesToSample(SystIndex, SampleName)) { //If parameter applies to required SampleID
227  returnVec.push_back(SystIndex);
228  }
229  }
230  return returnVec;
231 }
bool AppliesToSample(const int SystIndex, const std::string &SampleName) const
Check if parameter is affecting given sample name.

◆ GetNormParameter()

NormParameter ParameterHandlerGeneric::GetNormParameter ( const YAML::Node &  param,
const int  Index 
)
inlineprotected

Get Norm params.

Parameters
paramYaml node describing param
IndexGlobal parameter index

Definition at line 153 of file ParameterHandlerGeneric.cpp.

153  {
154 // ********************************************
155  NormParameter norm;
156 
157  GetBaseParameter(param, Index, norm);
158 
159  // ETA size 0 to mean apply to all
160  // Ultimately all this information ends up in the NormParams vector
161  norm.modes = GetFromManager<std::vector<int>>(param["Mode"], {}, __FILE__ , __LINE__);
162  norm.pdgs = GetFromManager<std::vector<int>>(param["NeutrinoFlavour"], {}, __FILE__ , __LINE__);
163  norm.preoscpdgs = GetFromManager<std::vector<int>>(param["NeutrinoFlavourUnosc"], {}, __FILE__ , __LINE__);
164  norm.targets = GetFromManager<std::vector<int>>(param["TargetNuclei"], {}, __FILE__ , __LINE__);
165 
166  if(_fLowBound[Index] < 0.) {
167  MACH3LOG_ERROR("Normalisation Parameter {} ({}), has lower parameters bound which can go below 0 and is equal {}",
168  GetParFancyName(Index), Index, _fLowBound[Index]);
169  MACH3LOG_ERROR("Normalisation parameters can't go bellow 0 as this is unphysical");
170  throw MaCh3Exception(__FILE__, __LINE__);
171  }
172 
173  return norm;
174 }
ETA - Normalisations for cross-section parameters Carrier for whether you want to apply a systematic ...
std::vector< int > preoscpdgs
Preosc PDG which parameter applies to.
std::vector< int > modes
Mode which parameter applies to.
std::vector< int > pdgs
PDG which parameter applies to.
std::vector< int > targets
Targets which parameter applies to.

◆ GetNormParsFromSampleName()

const std::vector< NormParameter > ParameterHandlerGeneric::GetNormParsFromSampleName ( const std::string &  SampleName) const

DB Get norm/func parameters depending on given SampleName.

Definition at line 313 of file ParameterHandlerGeneric.cpp.

313  {
314 // ********************************************
316 }
@ kNorm
For normalisation parameters.
std::vector< NormParameter > NormParams
Vector containing info for normalisation systematics.

◆ GetNumParamsFromSampleName()

int ParameterHandlerGeneric::GetNumParamsFromSampleName ( const std::string &  SampleName,
const SystType  Type 
)

DB Grab the number of parameters for the relevant SampleName.

Parameters
SampleNameproperty of SampleHandler class based on which we select whether to apply uncertainties or not
TypeType of syst, for example kNorm, kSpline etc

Definition at line 342 of file ParameterHandlerGeneric.cpp.

342  {
343 // ********************************************
344  int returnVal = 0;
345  IterateOverParams(SampleName,
346  [&](int i) { return GetParamType(i) == Type; }, // Filter condition
347  [&](int) { returnVal += 1; } // Action to perform if filter passes
348  );
349  return returnVal;
350 }
SystType GetParamType(const int i) const
Returns enum describing our param type.
void IterateOverParams(const std::string &SampleName, FilterFunc filter, ActionFunc action)
Iterates over parameters and applies a filter and action function.

◆ GetNumParFromGroup()

int ParameterHandlerGeneric::GetNumParFromGroup ( const std::string &  Group) const

KS: Check how many parameters are associated with given group.

Definition at line 761 of file ParameterHandlerGeneric.cpp.

761  {
762 // ********************************************
763  int Counter = 0;
764  for (int i = 0; i < _fNumPar; i++) {
765  if(IsParFromGroup(i, Group)) Counter++;
766  }
767  return Counter;
768 }
bool IsParFromGroup(const int i, const std::string &Group) const
Checks if parameter belongs to a given group.

◆ GetOscillationParameters()

OscillationParameter ParameterHandlerGeneric::GetOscillationParameters ( const YAML::Node &  param,
const int  Index 
)
inlineprotected

Get Osc params.

Parameters
paramYaml node describing param
IndexGlobal parameter index

Definition at line 296 of file ParameterHandlerGeneric.cpp.

296  {
297 // ********************************************
298  OscillationParameter OscParamInfo;
299  GetBaseParameter(param, Index, OscParamInfo);
300 
301  return OscParamInfo;
302 }
KS: Struct holding info about oscillation Systematics.

◆ GetOscParsFromSampleName()

std::vector< const M3::float_t * > ParameterHandlerGeneric::GetOscParsFromSampleName ( const std::string &  SampleName)

Get pointers to Osc params from Sample name.

Definition at line 772 of file ParameterHandlerGeneric.cpp.

772  {
773 // ********************************************
774  std::vector<const M3::float_t*> returnVec;
775  for (const auto& pair : _fSystToGlobalSystIndexMap[SystType::kOsc]) {
776  const auto& globalIndex = pair.second;
777  if (AppliesToSample(globalIndex, SampleName)) {
778  returnVec.push_back(RetPointer(globalIndex));
779  }
780  }
781  return returnVec;
782 }
@ kOsc
For oscillation parameters.

◆ GetParamType()

SystType ParameterHandlerGeneric::GetParamType ( const int  i) const
inline

Returns enum describing our param type.

Parameters
iparameter index

Definition at line 34 of file ParameterHandlerGeneric.h.

34 {return _fParamType[i];}
std::vector< SystType > _fParamType
Type of parameter like norm, spline etc.

◆ GetParamTypeString()

std::string ParameterHandlerGeneric::GetParamTypeString ( const int  i) const
inline

ETA - just return a string of "spline", "norm" or "functional".

Parameters
iparameter index

Definition at line 31 of file ParameterHandlerGeneric.h.

31 { return SystType_ToString(_fParamType[i]); }
std::string SystType_ToString(const SystType i)
Convert a Syst type type to a string.

◆ GetParSampleID()

std::vector<std::string> ParameterHandlerGeneric::GetParSampleID ( const int  i) const
inline

ETA - just return the int of the SampleName, this can be removed to do a string comp at some point.

Parameters
iparameter index

Definition at line 28 of file ParameterHandlerGeneric.h.

28 { return _fSampleNames[i];};
std::vector< std::vector< std::string > > _fSampleNames
Tells to which samples object param should be applied.

◆ GetParsIndexFromSampleName()

const std::vector< int > ParameterHandlerGeneric::GetParsIndexFromSampleName ( const std::string &  SampleName,
const SystType  Type 
)

DB Grab the parameter indices for the relevant SampleName.

Parameters
SampleNameproperty of SampleHandler class based on which we select whether to apply uncertainties or not
TypeType of syst, for example kNorm, kSpline etc

Definition at line 366 of file ParameterHandlerGeneric.cpp.

366  {
367 // ********************************************
368  std::vector<int> returnVec;
369  IterateOverParams(SampleName,
370  [&](int i) { return GetParamType(i) == Type; }, // Filter condition
371  [&](int i) { returnVec.push_back(i); } // Action to perform if filter passes
372  );
373  return returnVec;
374 }

◆ GetParsNamesFromSampleName()

const std::vector< std::string > ParameterHandlerGeneric::GetParsNamesFromSampleName ( const std::string &  SampleName,
const SystType  Type 
)

DB Grab the parameter names for the relevant SampleName.

Parameters
SampleNameproperty of SampleHandler class based on which we select whether to apply uncertainties or not
TypeType of syst, for example kNorm, kSpline etc

Definition at line 354 of file ParameterHandlerGeneric.cpp.

354  {
355 // ********************************************
356  std::vector<std::string> returnVec;
357  IterateOverParams(SampleName,
358  [&](int i) { return GetParamType(i) == Type; }, // Filter condition
359  [&](int i) { returnVec.push_back(GetParFancyName(i)); } // Action to perform if filter passes
360  );
361  return returnVec;
362 }

◆ GetParSplineInterpolation()

SplineInterpolation ParameterHandlerGeneric::GetParSplineInterpolation ( const int  i) const
inline

Get interpolation type for a given parameter.

Parameters
ispline parameter index, not confuse with global index

Definition at line 38 of file ParameterHandlerGeneric.h.

38 {return SplineParams.at(i)._SplineInterpolationType;}

◆ GetParSplineKnotLowerBound()

double ParameterHandlerGeneric::GetParSplineKnotLowerBound ( const int  i) const
inline

EM: value at which we cap spline knot weight.

Parameters
ispline parameter index, not confuse with global index

Definition at line 52 of file ParameterHandlerGeneric.h.

52 {return SplineParams.at(i)._SplineKnotLowBound;}

◆ GetParSplineKnotUpperBound()

double ParameterHandlerGeneric::GetParSplineKnotUpperBound ( const int  i) const
inline

EM: value at which we cap spline knot weight.

Parameters
ispline parameter index, not confuse with global index

Definition at line 49 of file ParameterHandlerGeneric.h.

49 {return SplineParams.at(i)._SplineKnotUpBound;}

◆ GetParSplineName()

std::string ParameterHandlerGeneric::GetParSplineName ( const int  i) const
inline

Get the name of the spline associated with the spline at index i.

Parameters
ispline parameter index, not to be confused with global index

Definition at line 43 of file ParameterHandlerGeneric.h.

43 {return SplineParams.at(i)._fSplineNames;}

◆ GetSplineFileParsNamesFromSampleName()

const std::vector<std::string> ParameterHandlerGeneric::GetSplineFileParsNamesFromSampleName ( const std::string &  SampleName)

DB Get spline parameters depending on given SampleName.

◆ GetSplineInterpolationFromSampleName()

const std::vector< SplineInterpolation > ParameterHandlerGeneric::GetSplineInterpolationFromSampleName ( const std::string &  SampleName)

Get the interpolation types for splines affecting a particular SampleName.

Definition at line 120 of file ParameterHandlerGeneric.cpp.

120  {
121 // ********************************************
122  std::vector<SplineInterpolation> returnVec;
123  for (auto &pair : _fSystToGlobalSystIndexMap[SystType::kSpline]) {
124  auto &SplineIndex = pair.first;
125  auto &SystIndex = pair.second;
126 
127  if (AppliesToSample(SystIndex, SampleName)) { //If parameter applies to required SampleID
128  returnVec.push_back(SplineParams.at(SplineIndex)._SplineInterpolationType);
129  }
130  }
131  return returnVec;
132 }
@ kSpline
For splined parameters (1D)

◆ GetSplineModeVecFromSampleName()

const std::vector< std::vector< int > > ParameterHandlerGeneric::GetSplineModeVecFromSampleName ( const std::string &  SampleName)

DB Grab the Spline Modes for the relevant SampleName.

Definition at line 136 of file ParameterHandlerGeneric.cpp.

136  {
137 // ********************************************
138  std::vector< std::vector<int> > returnVec;
139  //Need a counter or something to correctly get the index in _fSplineModes since it's not of length nPars
140  //Should probably just make a std::map<std::string, int> for param name to FD spline index
141  for (auto &pair : _fSystToGlobalSystIndexMap[SystType::kSpline]) {
142  auto &SplineIndex = pair.first;
143  auto &SystIndex = pair.second;
144  if (AppliesToSample(SystIndex, SampleName)) { //If parameter applies to required SampleID
145  returnVec.push_back(SplineParams.at(SplineIndex)._fSplineModes);
146  }
147  }
148  return returnVec;
149 }

◆ GetSplineParameter()

SplineParameter ParameterHandlerGeneric::GetSplineParameter ( const YAML::Node &  param,
const int  Index 
)
inlineprotected

Get Spline params.

Parameters
paramYaml node describing param
IndexGlobal parameter index

Definition at line 251 of file ParameterHandlerGeneric.cpp.

251  {
252 // ********************************************
253  SplineParameter Spline;
254 
255  GetBaseParameter(param, Index, Spline);
256  auto& SplinePar = param["SplineInformation"];
257  //Now get the Spline interpolation type
258  if (SplinePar["InterpolationType"]) {
259  for(int InterpType = 0; InterpType < kSplineInterpolations ; ++InterpType){
260  if(SplinePar["InterpolationType"].as<std::string>() == SplineInterpolation_ToString(SplineInterpolation(InterpType)))
261  Spline._SplineInterpolationType = SplineInterpolation(InterpType);
262  }
263  } else { //KS: By default use TSpline3
265  }
266 
267  Spline._fSplineNames = SplinePar["SplineName"].as<std::string>();
268  Spline._SplineKnotUpBound = GetFromManager<double>(SplinePar["SplineKnotUpBound"], M3::DefSplineKnotUpBound, __FILE__ , __LINE__);
269  Spline._SplineKnotLowBound = GetFromManager<double>(SplinePar["SplineKnotLowBound"], M3::DefSplineKnotLowBound, __FILE__ , __LINE__);
270 
272  MACH3LOG_WARN("Spline knot capping enabled with bounds [{}, {}]. For reliable fits, consider modifying the input generation instead.",
273  Spline._SplineKnotLowBound, Spline._SplineKnotUpBound);
274  }
275  //If there is no mode information given then this will be an empty vector
276  Spline._fSplineModes = GetFromManager(SplinePar["Mode"], std::vector<int>(), __FILE__ , __LINE__);
277 
278  return Spline;
279 }
std::string SplineInterpolation_ToString(const SplineInterpolation i)
Convert a LLH type to a string.
SplineInterpolation
Make an enum of the spline interpolation type.
@ kTSpline3
Default TSpline3 interpolation.
@ kSplineInterpolations
This only enumerates.
Type GetFromManager(const YAML::Node &node, const Type defval, const std::string File="", const int Line=1)
Get content of config file if node is not found take default value specified.
Definition: YamlHelper.h:329
constexpr static const double DefSplineKnotUpBound
Default value for spline knot capping, default mean not capping is being applied.
Definition: Core.h:86
constexpr static const double DefSplineKnotLowBound
Default value for spline knot capping, default mean not capping is being applied.
Definition: Core.h:88
KS: Struct holding info about Spline Systematics.
SplineInterpolation _SplineInterpolationType
Spline interpolation vector.
double _SplineKnotUpBound
EM: Cap spline knot higher value.
double _SplineKnotLowBound
EM: Cap spline knot lower value.
std::vector< int > _fSplineModes
Modes to which spline applies (valid only for binned splines)
std::string _fSplineNames
Name of spline in TTree (TBranch),.

◆ GetSplineParsFromSampleName()

const std::vector< SplineParameter > ParameterHandlerGeneric::GetSplineParsFromSampleName ( const std::string &  SampleName) const

KS: Grab the Spline parameters for the relevant SampleName.

Definition at line 320 of file ParameterHandlerGeneric.cpp.

320  {
321 // ********************************************
323 }

◆ GetSplineParsNamesFromSampleName()

const std::vector< std::string > ParameterHandlerGeneric::GetSplineParsNamesFromSampleName ( const std::string &  SampleName)

DB Get spline parameters depending on given SampleName.

Definition at line 106 of file ParameterHandlerGeneric.cpp.

106  {
107 // ********************************************
108  std::vector<std::string> returnVec;
109  for (auto &pair : _fSystToGlobalSystIndexMap[SystType::kSpline]) {
110  auto &SplineIndex = pair.first;
111  auto &SystIndex = pair.second;
112  if (AppliesToSample(SystIndex, SampleName)) { //If parameter applies to required Sample
113  returnVec.push_back(SplineParams[SplineIndex]._fSplineNames);
114  }
115  }
116  return returnVec;
117 }

◆ GetSystIndexFromSampleName()

const std::vector< int > ParameterHandlerGeneric::GetSystIndexFromSampleName ( const std::string &  SampleName,
const SystType  Type 
) const

Grab the index of the syst relative to global numbering.

Parameters
SampleNameproperty of SampleHandler class based on which we select whether to apply uncertainties or not
TypeType of syst, for example kNorm, kSpline etc

Definition at line 236 of file ParameterHandlerGeneric.cpp.

236  {
237 // ********************************************
238  std::vector<int> returnVec;
239  for (auto &pair : _fSystToGlobalSystIndexMap[Type]) {
240  auto &SplineIndex = pair.first;
241  auto &systIndex = pair.second;
242  if (AppliesToSample(systIndex, SampleName)) { //If parameter applies to required SampleID
243  returnVec.push_back(SplineIndex);
244  }
245  }
246  return returnVec;
247 }

◆ GetTypeParamsFromSampleName()

template<typename ParamT >
std::vector< ParamT > ParameterHandlerGeneric::GetTypeParamsFromSampleName ( const std::map< int, int > &  indexMap,
const std::vector< ParamT > &  params,
const std::string &  SampleName 
) const
protected

Retrieve parameters that apply to a given sample name.

Template Parameters
ParamTType of parameter (e.g., FunctionalParameter, NormParameter).
Parameters
indexMapMap from local to global parameter indices.
paramsVector of all parameters of the given type.
SampleNameThe name of the sample to filter applicable parameters for.
Returns
Vector of parameters of type ParamT that apply to the specified sample.

Definition at line 327 of file ParameterHandlerGeneric.cpp.

327  {
328 // ********************************************
329  std::vector<ParamT> returnVec;
330  for (const auto& pair : indexMap) {
331  const auto& localIndex = pair.first;
332  const auto& globalIndex = pair.second;
333  if (AppliesToSample(globalIndex, SampleName)) {
334  returnVec.push_back(params[localIndex]);
335  }
336  }
337  return returnVec;
338 }

◆ GetUniqueParameterGroups()

std::vector< std::string > ParameterHandlerGeneric::GetUniqueParameterGroups ( ) const

KS: Get names of all unique parameter groups.

Definition at line 609 of file ParameterHandlerGeneric.cpp.

609  {
610 // ********************************************
611  std::unordered_set<std::string> uniqueGroups;
612 
613  // Fill the set with unique values
614  for (const auto& param : _ParameterGroup) {
615  uniqueGroups.insert(param);
616  }
617 
618  // Convert to vector and return
619  std::vector<std::string> result(uniqueGroups.begin(), uniqueGroups.end());
620  return result;
621 }
std::vector< std::string > _ParameterGroup
KS: Allow to group parameters for example to affect only cross-section or only flux etc.

◆ InitParametersTypeFromConfig()

void ParameterHandlerGeneric::InitParametersTypeFromConfig ( )
protected

Parses the YAML configuration to set up cross-section parameters. The YAML file defines the types of systematic errors, interpolation types, and bounds for splines.

Definition at line 26 of file ParameterHandlerGeneric.cpp.

26  {
27 // ********************************************
29 
30  _fParamType = std::vector<SystType>(_fNumPar);
31  _ParameterGroup = std::vector<std::string>(_fNumPar);
32 
33  //KS: We know at most how params we expect so reserve memory for max possible params. Later we will shrink to size to not waste memory. Reserving means slightly faster loading and possible less memory fragmentation.
34  NormParams.reserve(_fNumPar);
35  SplineParams.reserve(_fNumPar);
36  FuncParams.reserve(_fNumPar);
37  OscParams.reserve(_fNumPar);
38 
39  int i = 0;
40  unsigned int ParamCounter[SystType::kSystTypes] = {0};
41  //ETA - read in the systematics. Would be good to add in some checks to make sure
42  //that there are the correct number of entries i.e. are the _fNumPars for Names,
43  //PreFitValues etc etc.
44  for (auto const &param : _fYAMLDoc["Systematics"])
45  {
46  _ParameterGroup[i] = Get<std::string>(param["Systematic"]["ParameterGroup"], __FILE__ , __LINE__);
47 
48  //Fill the map to get the correlations later as well
49  auto ParamType = Get<std::string>(param["Systematic"]["Type"], __FILE__ , __LINE__);
50  //Now load in variables for spline systematics only
51  if (ParamType == SystType_ToString(SystType::kSpline))
52  {
53  //Set param type
55  // Fill Spline info
56  SplineParams.push_back(GetSplineParameter(param["Systematic"], i));
57 
58  //Insert the mapping from the spline index i.e. the length of _fSplineNames etc
59  //to the Systematic index i.e. the counter for things like _fSampleID
60  _fSystToGlobalSystIndexMap[SystType::kSpline].insert(std::make_pair(ParamCounter[SystType::kSpline], i));
61  ParamCounter[SystType::kSpline]++;
62  } else if(ParamType == SystType_ToString(SystType::kNorm)) {
64  NormParams.push_back(GetNormParameter(param["Systematic"], i));
65  _fSystToGlobalSystIndexMap[SystType::kNorm].insert(std::make_pair(ParamCounter[SystType::kNorm], i));
66  ParamCounter[SystType::kNorm]++;
67  } else if(ParamType == SystType_ToString(SystType::kFunc)){
69  FuncParams.push_back(GetFunctionalParameters(param["Systematic"], i));
70  _fSystToGlobalSystIndexMap[SystType::kFunc].insert(std::make_pair(ParamCounter[SystType::kFunc], i));
71  ParamCounter[SystType::kFunc]++;
72  } else if(ParamType == SystType_ToString(SystType::kOsc)){
74  OscParams.push_back(GetOscillationParameters(param["Systematic"], i));
75  _fSystToGlobalSystIndexMap[SystType::kOsc].insert(std::make_pair(ParamCounter[SystType::kOsc], i));
76  ParamCounter[SystType::kOsc]++;
77  } else{
78  MACH3LOG_ERROR("Given unrecognised systematic type: {}", ParamType);
79  std::string expectedTypes = "Expecting ";
80  for (int s = 0; s < SystType::kSystTypes; ++s) {
81  if (s > 0) expectedTypes += ", ";
82  expectedTypes += SystType_ToString(static_cast<SystType>(s)) + "\"";
83  }
84  expectedTypes += ".";
85  MACH3LOG_ERROR(expectedTypes);
86  throw MaCh3Exception(__FILE__, __LINE__);
87  }
88  i++;
89  } //end loop over params
90 
91  //KS We resized them above to all params to fight memory fragmentation, now let's resize to fit only allocated memory to save RAM
92  NormParams.shrink_to_fit();
93  SplineParams.shrink_to_fit();
94  FuncParams.shrink_to_fit();
95  OscParams.shrink_to_fit();
96 }
SystType
@ kSystTypes
This only enumerates.
YAML::Node _fYAMLDoc
Stores config describing systematics.
FunctionalParameter GetFunctionalParameters(const YAML::Node &param, const int Index)
Get Func params.
SplineParameter GetSplineParameter(const YAML::Node &param, const int Index)
Get Spline params.
NormParameter GetNormParameter(const YAML::Node &param, const int Index)
Get Norm params.
std::vector< OscillationParameter > OscParams
Vector containing info for functional systematics.
OscillationParameter GetOscillationParameters(const YAML::Node &param, const int Index)
Get Osc params.

◆ InitParams()

void ParameterHandlerGeneric::InitParams ( )
protected

Initializes the systematic parameters from the configuration file. This function loads parameters like normalizations and splines from the provided YAML file.

Note
This is used internally during the object's initialization process.

Definition at line 388 of file ParameterHandlerGeneric.cpp.

388  {
389 // ********************************************
390  for (int i = 0; i < _fNumPar; ++i) {
391  //ETA - set the name to be param_% as this is what ProcessorMCMC expects
392  _fNames[i] = "param_"+std::to_string(i);
393 
394  // KS: Plenty of MaCh3 processing script rely on osc params having "fancy name" this is to maintain backward compatibility with them
395  if(_fParamType[i] == kOsc) {
396  _fNames[i] = _fFancyNames[i];
397 
398  if(_ParameterGroup[i] != "Osc"){
399  MACH3LOG_ERROR("Parameter {}, is of type Oscillation but doesn't belong to Osc group", _fFancyNames[i]);
400  MACH3LOG_ERROR("It belongs to {} group", _ParameterGroup[i]);
401  throw MaCh3Exception(__FILE__ , __LINE__ );
402  }
403  }
404  #pragma GCC diagnostic push
405  #pragma GCC diagnostic ignored "-Wuseless-cast"
406  // Set ParameterHandler parameters (Curr = current, Prop = proposed, Sigma = step)
407  _fCurrVal[i] = _fPreFitValue[i];
408  _fPropVal[i] = static_cast<M3::float_t>(_fCurrVal[i]);
409  #pragma GCC diagnostic pop
410  }
411  Randomize();
412  //KS: Transfer the starting parameters to the PCA basis, you don't want to start with zero..
413  if (pca) {
414  PCAObj->SetInitialParameters();
415  }
416 }
std::vector< M3::float_t > _fPropVal
Proposed value of the parameter.
std::unique_ptr< PCAHandler > PCAObj
Struct containing information about PCA.
std::vector< double > _fCurrVal
Current value of the parameter.
void Randomize() _noexcept_
"Randomize" the parameters in the covariance class for the proposed step. Used the proposal kernel an...
bool pca
perform PCA or not
double float_t
Definition: Core.h:37

◆ IsParFromGroup()

bool ParameterHandlerGeneric::IsParFromGroup ( const int  i,
const std::string &  Group 
) const

Checks if parameter belongs to a given group.

Parameters
iparameter index
Groupname of group, like Xsec or Flux
Returns
bool telling whether param is part of group

Definition at line 748 of file ParameterHandlerGeneric.cpp.

748  {
749 // ********************************************
750  std::string groupLower = Group;
751  std::string paramGroupLower = _ParameterGroup[i];
752 
753  // KS: Convert both strings to lowercase, this way comparison will be case insensitive
754  std::transform(groupLower.begin(), groupLower.end(), groupLower.begin(), ::tolower);
755  std::transform(paramGroupLower.begin(), paramGroupLower.end(), paramGroupLower.begin(), ::tolower);
756 
757  return groupLower == paramGroupLower;
758 }

◆ IterateOverParams()

template<typename FilterFunc , typename ActionFunc >
void ParameterHandlerGeneric::IterateOverParams ( const std::string &  SampleName,
FilterFunc  filter,
ActionFunc  action 
)
protected

Iterates over parameters and applies a filter and action function.

This template function provides a way to iterate over parameters associated with a specific Sample ID (SampleName). It applies a filter function to determine which parameters to process and an action function to define what to do with the selected parameters.

Template Parameters
FilterFuncThe type of the filter function used to determine which parameters to include.
ActionFuncThe type of the action function applied to each selected parameter.
Parameters
SampleNameThe Sample ID used to filter parameters.

Definition at line 378 of file ParameterHandlerGeneric.cpp.

378  {
379 // ********************************************
380  for (int i = 0; i < _fNumPar; ++i) {
381  if ((AppliesToSample(i, SampleName)) && filter(i)) { // Common filter logic
382  action(i); // Specific action for each function
383  }
384  }
385 }

◆ Print()

void ParameterHandlerGeneric::Print ( ) const
protected

Print information about the whole object once it is set.

Definition at line 420 of file ParameterHandlerGeneric.cpp.

420  {
421 // ********************************************
422  MACH3LOG_INFO("#################################################");
423  MACH3LOG_INFO("Printing ParameterHandlerGeneric:");
424 
426 
427  PrintNormParams();
428 
430 
432 
434 
436 
437  MACH3LOG_INFO("Finished");
438  MACH3LOG_INFO("#################################################");
439 
441 } // End
void PrintGlobablInfo() const
Prints general information about the ParameterHandler object.
void PrintSplineParams() const
Prints spline parameters.
void PrintOscillationParams() const
Prints oscillation parameters.
void PrintFunctionalParams() const
Prints functional parameters.
void CheckCorrectInitialisation() const
KS: Check if matrix is correctly initialised.
void PrintNormParams() const
Prints normalization parameters.
void PrintParameterGroups() const
Prints groups of parameters.

◆ PrintFunctionalParams()

void ParameterHandlerGeneric::PrintFunctionalParams ( ) const
protected

Prints functional parameters.

Definition at line 559 of file ParameterHandlerGeneric.cpp.

559  {
560 // ********************************************
561  MACH3LOG_INFO("Functional parameters: {}", _fSystToGlobalSystIndexMap[SystType::kFunc].size());
562  if(_fSystToGlobalSystIndexMap[SystType::kFunc].size() == 0) return;
563  MACH3LOG_INFO("┌────┬──────────┬────────────────────────────────────────┐");
564  MACH3LOG_INFO("│{0:4}│{1:10}│{2:40}│", "#", "Global #", "Name");
565  MACH3LOG_INFO("├────┼──────────┼────────────────────────────────────────┤");
566  for (auto &pair : _fSystToGlobalSystIndexMap[SystType::kFunc]) {
567  auto &FuncIndex = pair.first;
568  auto &GlobalIndex = pair.second;
569  MACH3LOG_INFO("│{0:4}│{1:<10}│{2:40}│", std::to_string(FuncIndex), GlobalIndex, GetParFancyName(GlobalIndex));
570  }
571  MACH3LOG_INFO("└────┴──────────┴────────────────────────────────────────┘");
572 }

◆ PrintGlobablInfo()

void ParameterHandlerGeneric::PrintGlobablInfo ( ) const
protected

Prints general information about the ParameterHandler object.

Definition at line 444 of file ParameterHandlerGeneric.cpp.

444  {
445 // ********************************************
446  MACH3LOG_INFO("============================================================================================================================================================");
447  MACH3LOG_INFO("{:<5} {:2} {:<40} {:2} {:<10} {:2} {:<10} {:2} {:<10} {:2} {:<10} {:2} {:<10} {:2} {:<20} {:2} {:<10}", "#", "|", "Name", "|", "Prior", "|", "Error", "|", "Lower", "|", "Upper", "|", "StepScale", "|", "SampleNames", "|", "Type");
448  MACH3LOG_INFO("------------------------------------------------------------------------------------------------------------------------------------------------------------");
449  for (int i = 0; i < GetNumParams(); i++) {
450  std::string ErrString = fmt::format("{:.2f}", _fError[i]);
451  std::string SampleNameString = "";
452  for (const auto& SampleName : _fSampleNames[i]) {
453  if (!SampleNameString.empty()) {
454  SampleNameString += ", ";
455  }
456  SampleNameString += SampleName;
457  }
458  MACH3LOG_INFO("{:<5} {:2} {:<40} {:2} {:<10} {:2} {:<10} {:2} {:<10} {:2} {:<10} {:2} {:<10} {:2} {:<20} {:2} {:<10}", i, "|", GetParFancyName(i), "|", _fPreFitValue[i], "|", "+/- " + ErrString, "|", _fLowBound[i], "|", _fUpBound[i], "|", _fIndivStepScale[i], "|", SampleNameString, "|", SystType_ToString(_fParamType[i]));
459  }
460  MACH3LOG_INFO("============================================================================================================================================================");
461 }
int GetNumParams() const
Get total number of parameters.

◆ PrintNormParams()

void ParameterHandlerGeneric::PrintNormParams ( ) const
protected

Prints normalization parameters.

Definition at line 464 of file ParameterHandlerGeneric.cpp.

464  {
465 // ********************************************
466  // Output the normalisation parameters as a sanity check!
467  MACH3LOG_INFO("Normalisation parameters: {}", NormParams.size());
468  if(_fSystToGlobalSystIndexMap[SystType::kNorm].size() == 0) return;
469 
470  bool have_parameter_with_kin_bounds = false;
471 
472  //KS: Consider making some class producing table..
473  MACH3LOG_INFO("┌────┬──────────┬────────────────────────────────────────┬────────────────────┬────────────────────┬────────────────────┐");
474  MACH3LOG_INFO("│{0:4}│{1:10}│{2:40}│{3:20}│{4:20}│{5:20}│", "#", "Global #", "Name", "Int. mode", "Target", "pdg");
475  MACH3LOG_INFO("├────┼──────────┼────────────────────────────────────────┼────────────────────┼────────────────────┼────────────────────┤");
476 
477  for (unsigned int i = 0; i < NormParams.size(); ++i)
478  {
479  std::string intModeString;
480  for (unsigned int j = 0; j < NormParams[i].modes.size(); j++) {
481  intModeString += std::to_string(NormParams[i].modes[j]);
482  intModeString += " ";
483  }
484  if (NormParams[i].modes.empty()) intModeString += "all";
485 
486  std::string targetString;
487  for (unsigned int j = 0; j < NormParams[i].targets.size(); j++) {
488  targetString += std::to_string(NormParams[i].targets[j]);
489  targetString += " ";
490  }
491  if (NormParams[i].targets.empty()) targetString += "all";
492 
493  std::string pdgString;
494  for (unsigned int j = 0; j < NormParams[i].pdgs.size(); j++) {
495  pdgString += std::to_string(NormParams[i].pdgs[j]);
496  pdgString += " ";
497  }
498  if (NormParams[i].pdgs.empty()) pdgString += "all";
499 
500  MACH3LOG_INFO("│{: <4}│{: <10}│{: <40}│{: <20}│{: <20}│{: <20}│", i, NormParams[i].index, NormParams[i].name, intModeString, targetString, pdgString);
501 
502  if(NormParams[i].hasKinBounds) have_parameter_with_kin_bounds = true;
503  }
504  MACH3LOG_INFO("└────┴──────────┴────────────────────────────────────────┴────────────────────┴────────────────────┴────────────────────┘");
505 
506  if(have_parameter_with_kin_bounds) {
507  MACH3LOG_INFO("Normalisation parameters KinematicCuts information");
508  MACH3LOG_INFO("┌────┬──────────┬────────────────────────────────────────┬────────────────────┬────────────────────────────────────────┐");
509  MACH3LOG_INFO("│{0:4}│{1:10}│{2:40}│{3:20}│{4:40}│", "#", "Global #", "Name", "KinematicCut", "Value");
510  MACH3LOG_INFO("├────┼──────────┼────────────────────────────────────────┼────────────────────┼────────────────────────────────────────┤");
511  for (unsigned int i = 0; i < NormParams.size(); ++i)
512  {
513  //skip parameters with no KinematicCuts
514  if(!NormParams[i].hasKinBounds) continue;
515 
516  const long unsigned int ncuts = NormParams[i].KinematicVarStr.size();
517  for(long unsigned int icut = 0; icut < ncuts; icut++) {
518  std::string kinematicCutValueString;
519  for(const auto & value : NormParams[i].Selection[icut]) {
520  for (const auto& v : value) {
521  kinematicCutValueString += fmt::format("{:.2f} ", v);
522  }
523  }
524  if(icut == 0)
525  MACH3LOG_INFO("│{: <4}│{: <10}│{: <40}│{: <20}│{: <40}│", i, NormParams[i].index, NormParams[i].name, NormParams[i].KinematicVarStr[icut], kinematicCutValueString);
526  else
527  MACH3LOG_INFO("│{: <4}│{: <10}│{: <40}│{: <20}│{: <40}│", "", "", "", NormParams[i].KinematicVarStr[icut], kinematicCutValueString);
528  }//icut
529  }//i
530  MACH3LOG_INFO("└────┴──────────┴────────────────────────────────────────┴────────────────────┴────────────────────────────────────────┘");
531  }
532  else
533  MACH3LOG_INFO("No normalisation parameters have KinematicCuts defined");
534 }

◆ PrintOscillationParams()

void ParameterHandlerGeneric::PrintOscillationParams ( ) const
protected

Prints oscillation parameters.

Definition at line 575 of file ParameterHandlerGeneric.cpp.

575  {
576 // ********************************************
577  MACH3LOG_INFO("Oscillation parameters: {}", _fSystToGlobalSystIndexMap[SystType::kOsc].size());
578  if(_fSystToGlobalSystIndexMap[SystType::kOsc].size() == 0) return;
579  MACH3LOG_INFO("┌────┬──────────┬────────────────────────────────────────┐");
580  MACH3LOG_INFO("│{0:4}│{1:10}│{2:40}│", "#", "Global #", "Name");
581  MACH3LOG_INFO("├────┼──────────┼────────────────────────────────────────┤");
582  for (auto &pair : _fSystToGlobalSystIndexMap[SystType::kOsc]) {
583  auto &OscIndex = pair.first;
584  auto &GlobalIndex = pair.second;
585  MACH3LOG_INFO("│{0:4}│{1:<10}│{2:40}│", std::to_string(OscIndex), GlobalIndex, GetParFancyName(GlobalIndex));
586  }
587  MACH3LOG_INFO("└────┴──────────┴────────────────────────────────────────┘");
588 }

◆ PrintParameterGroups()

void ParameterHandlerGeneric::PrintParameterGroups ( ) const
protected

Prints groups of parameters.

Definition at line 591 of file ParameterHandlerGeneric.cpp.

591  {
592 // ********************************************
593  // KS: Create a map to store the counts of unique strings, in principle this could be in header file
594  std::unordered_map<std::string, int> paramCounts;
595 
596  std::for_each(_ParameterGroup.begin(), _ParameterGroup.end(),
597  [&paramCounts](const std::string& param) {
598  paramCounts[param]++;
599  });
600 
601  MACH3LOG_INFO("Printing parameter groups");
602  // Output the counts
603  for (const auto& pair : paramCounts) {
604  MACH3LOG_INFO("Found {}: {} params", pair.second, pair.first);
605  }
606 }

◆ PrintSplineParams()

void ParameterHandlerGeneric::PrintSplineParams ( ) const
protected

Prints spline parameters.

Definition at line 537 of file ParameterHandlerGeneric.cpp.

537  {
538 // ********************************************
539  MACH3LOG_INFO("Spline parameters: {}", _fSystToGlobalSystIndexMap[SystType::kSpline].size());
540  if(_fSystToGlobalSystIndexMap[SystType::kSpline].size() == 0) return;
541  MACH3LOG_INFO("=====================================================================================================================================================================");
542  MACH3LOG_INFO("{:<4} {:<2} {:<40} {:<2} {:<40} {:<2} {:<20} {:<2} {:<20} {:<2} {:<20} {:<2}", "#", "|", "Name", "|", "Spline Name", "|", "Spline Interpolation", "|", "Low Knot Bound", "|", "Up Knot Bound", "|");
543  MACH3LOG_INFO("---------------------------------------------------------------------------------------------------------------------------------------------------------------------");
544  for (auto &pair : _fSystToGlobalSystIndexMap[SystType::kSpline]) {
545  auto &SplineIndex = pair.first;
546  auto &GlobalIndex = pair.second;
547 
548  MACH3LOG_INFO("{:<4} {:<2} {:<40} {:<2} {:<40} {:<2} {:<20} {:<2} {:<20} {:<2} {:<20} {:<2}",
549  SplineIndex, "|", GetParFancyName(GlobalIndex), "|",
550  SplineParams[SplineIndex]._fSplineNames, "|",
552  GetParSplineKnotLowerBound(SplineIndex), "|",
553  GetParSplineKnotUpperBound(SplineIndex), "|");
554  }
555  MACH3LOG_INFO("=====================================================================================================================================================================");
556 }
double GetParSplineKnotUpperBound(const int i) const
EM: value at which we cap spline knot weight.
SplineInterpolation GetParSplineInterpolation(const int i) const
Get interpolation type for a given parameter.
double GetParSplineKnotLowerBound(const int i) const
EM: value at which we cap spline knot weight.

◆ SetFixGroupOnlyParameters() [1/2]

void ParameterHandlerGeneric::SetFixGroupOnlyParameters ( const std::string &  Group)

TN Method to set parameters within a group to be fixed to their prior values.

Parameters
Groupname of the parameter group (Xsec, Flux, Osc, etc.)

Definition at line 716 of file ParameterHandlerGeneric.cpp.

716  {
717 // ********************************************
718  for (int i = 0; i < _fNumPar; ++i)
719  if(IsParFromGroup(i, Group)) SetFixParameter(i);
720 }
void SetFixParameter(const int i)
Set parameter to be fixed at prior value.

◆ SetFixGroupOnlyParameters() [2/2]

void ParameterHandlerGeneric::SetFixGroupOnlyParameters ( const std::vector< std::string > &  Groups)

TN Method to set parameters of certain groups to be fixed to their prior values.

Parameters
Groupsvector of group names (e.g. {"Xsec", "Flux"})

Definition at line 724 of file ParameterHandlerGeneric.cpp.

724  {
725 // ********************************************
726  for(size_t i = 0; i < Groups.size(); i++)
727  SetFixGroupOnlyParameters(Groups[i]);
728 }
void SetFixGroupOnlyParameters(const std::string &Group)
TN Method to set parameters within a group to be fixed to their prior values.

◆ SetFreeGroupOnlyParameters() [1/2]

void ParameterHandlerGeneric::SetFreeGroupOnlyParameters ( const std::string &  Group)

TN Method to set parameters within a group to be treated as free.

Parameters
Groupname of the parameter group (Xsec, Flux, Osc, etc.)

Definition at line 732 of file ParameterHandlerGeneric.cpp.

732  {
733 // ********************************************
734  for (int i = 0; i < _fNumPar; ++i)
735  if(IsParFromGroup(i, Group)) SetFreeParameter(i);
736 }
void SetFreeParameter(const int i)
Set parameter to be treated as free.

◆ SetFreeGroupOnlyParameters() [2/2]

void ParameterHandlerGeneric::SetFreeGroupOnlyParameters ( const std::vector< std::string > &  Groups)

TN Method to set parameters of certain groups to be treated as free.

Parameters
Groupsvector of group names (e.g. {"Xsec", "Flux"})

Definition at line 740 of file ParameterHandlerGeneric.cpp.

740  {
741 // ********************************************
742  for(size_t i = 0; i < Groups.size(); i++)
743  SetFreeGroupOnlyParameters(Groups[i]);
744 }
void SetFreeGroupOnlyParameters(const std::string &Group)
TN Method to set parameters within a group to be treated as free.

◆ SetGroupOnlyParameters() [1/2]

void ParameterHandlerGeneric::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.

Parameters
Groupname of group, like Xsec or Flux
ParsValues which will overwrite proposed step
Note
This mimics functionality of SetParameters

Definition at line 666 of file ParameterHandlerGeneric.cpp.

666  {
667 // ********************************************
668  #pragma GCC diagnostic push
669  #pragma GCC diagnostic ignored "-Wuseless-cast"
670  // If empty, set the proposed to prior
671  if (Pars.empty()) {
672  for (int i = 0; i < _fNumPar; i++) {
673  if(IsParFromGroup(i, Group)) _fPropVal[i] = static_cast<M3::float_t>(_fPreFitValue[i]);
674  }
675  } else{
676  const size_t ExpectedSize = static_cast<size_t>(GetNumParFromGroup(Group));
677  if (Pars.size() != ExpectedSize) {
678  MACH3LOG_ERROR("Number of param in group {} is {}, while you passed {}", Group, ExpectedSize, Pars.size());
679  throw MaCh3Exception(__FILE__, __LINE__);
680  }
681  int Counter = 0;
682  for (int i = 0; i < _fNumPar; i++) {
683  // If belongs to group set value from parsed vector, otherwise use propose value
684  if(IsParFromGroup(i, Group)){
685  _fPropVal[i] = static_cast<M3::float_t>(Pars[Counter]);
686  Counter++;
687  }
688  }
689  }
690  // And if pca make the transfer
691  if (pca) {
692  PCAObj->TransferToPCA();
693  PCAObj->TransferToParam();
694  }
695  #pragma GCC diagnostic pop
696 }
int GetNumParFromGroup(const std::string &Group) const
KS: Check how many parameters are associated with given group.

◆ SetGroupOnlyParameters() [2/2]

void ParameterHandlerGeneric::SetGroupOnlyParameters ( const std::vector< std::string > &  Groups)

KS Function to set to prior parameters of a given groups or values from vector.

Parameters
Groupsvector of group names, like Xsec or Flux

Definition at line 657 of file ParameterHandlerGeneric.cpp.

657  {
658 // ********************************************
659  for(size_t i = 0; i < Groups.size(); i++){
660  SetGroupOnlyParameters(Groups[i]);
661  }
662 }
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.

◆ ToggleFixGroupOnlyParameters() [1/2]

void ParameterHandlerGeneric::ToggleFixGroupOnlyParameters ( const std::string &  Group)

TN Method to toggle fix/free parameters within a group.

Parameters
Groupname of the parameter group (Xsec, Flux, Osc, etc.)

Definition at line 700 of file ParameterHandlerGeneric.cpp.

700  {
701 // ********************************************
702  for (int i = 0; i < _fNumPar; ++i)
703  if(IsParFromGroup(i, Group)) ToggleFixParameter(i);
704 }
void ToggleFixParameter(const int i)
Toggle fixing parameter at prior values.

◆ ToggleFixGroupOnlyParameters() [2/2]

void ParameterHandlerGeneric::ToggleFixGroupOnlyParameters ( const std::vector< std::string > &  Groups)

TN Method to toggle fix/free parameters within given groups.

Parameters
Groupsvector of group names (e.g. {"Xsec", "Flux"})

Definition at line 708 of file ParameterHandlerGeneric.cpp.

708  {
709 // ********************************************
710  for(size_t i = 0; i < Groups.size(); i++)
711  ToggleFixGroupOnlyParameters(Groups[i]);
712 }
void ToggleFixGroupOnlyParameters(const std::string &Group)
TN Method to toggle fix/free parameters within a group.

Member Data Documentation

◆ _fParamType

std::vector<SystType> ParameterHandlerGeneric::_fParamType
protected

Type of parameter like norm, spline etc.

Definition at line 211 of file ParameterHandlerGeneric.h.

◆ _fSystToGlobalSystIndexMap

std::vector<std::map<int, int> > ParameterHandlerGeneric::_fSystToGlobalSystIndexMap
protected

Map between number of given parameter type with global parameter numbering. For example 2nd norm param may be 10-th global param.

Definition at line 217 of file ParameterHandlerGeneric.h.

◆ _ParameterGroup

std::vector<std::string> ParameterHandlerGeneric::_ParameterGroup
protected

KS: Allow to group parameters for example to affect only cross-section or only flux etc.

Definition at line 214 of file ParameterHandlerGeneric.h.

◆ FuncParams

std::vector<FunctionalParameter> ParameterHandlerGeneric::FuncParams
protected

Vector containing info for functional systematics.

Definition at line 226 of file ParameterHandlerGeneric.h.

◆ NormParams

std::vector<NormParameter> ParameterHandlerGeneric::NormParams
protected

Vector containing info for normalisation systematics.

Definition at line 223 of file ParameterHandlerGeneric.h.

◆ OscParams

std::vector<OscillationParameter> ParameterHandlerGeneric::OscParams
protected

Vector containing info for functional systematics.

Definition at line 229 of file ParameterHandlerGeneric.h.

◆ SplineParams

std::vector<SplineParameter> ParameterHandlerGeneric::SplineParams
protected

Vector containing info for normalisation systematics.

Definition at line 220 of file ParameterHandlerGeneric.h.


The documentation for this class was generated from the following files: