MaCh3  2.2.3
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 ()
 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 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 double * > 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 ThrowParProp (const double mag=1.)
 Throw the proposed parameter by mag sigma. Should really just have the user specify this throw by having argument double. More...
 
void ThrowParCurr (const double mag=1.)
 Helper function to throw the current parameter by mag sigma. Can study bias in MCMC with this; put different starting parameters. 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 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 UpdateThrowMatrix (TMatrixDSym *cov)
 
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 double * RetPointer (const int iParam)
 DB Pointer return to param position. More...
 
const std::vector< double > & 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...
 
double 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 [12]. More...
 
void AcceptStep () _noexcept_
 Accepted this step. More...
 
void ToggleFixAllParameters ()
 fix parameters at prior values More...
 
void ToggleFixParameter (const int i)
 fix parameter at prior values More...
 
void ToggleFixParameter (const std::string &name)
 Fix 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...
 
adaptive_mcmc::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)
 Matches branches in a TTree to parameters in a systematic handler. More...
 

Protected Member Functions

void Print ()
 Print information about the whole object once it is set. More...
 
void PrintGlobablInfo ()
 Prints general information about the ParameterHandler object. More...
 
void PrintNormParams ()
 Prints normalization parameters. More...
 
void PrintSplineParams ()
 Prints spline parameters. More...
 
void PrintFunctionalParams ()
 Prints functional parameters. More...
 
void PrintOscillationParams ()
 Prints oscillation parameters. More...
 
void PrintParameterGroups ()
 Prints groups of parameters. More...
 
void CheckCorrectInitialisation ()
 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 > _fSplineNames
 Name of spline in TTree (TBranch),. 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 xsec_i, this is currently to make things compatible with the Diagnostic tools. More...
 
std::vector< std::string > _fFancyNames
 Fancy name for example rather than xsec_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< double > _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...
 
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< adaptive_mcmc::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.

See also
For more details, visit the Wiki.
Author
Dan Barrow
Ed Atkin
Kamil Skwarczynski

Definition at line 12 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:24
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 xsec_i,...
int PrintLength
KS: This is used when printing parameters, sometimes we have super long parameters name,...
void Print()
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 108 of file ParameterHandlerGeneric.cpp.

108  {
109 // ********************************************
110  MACH3LOG_DEBUG("Deleting ParameterHandler");
111 }

Member Function Documentation

◆ CheckCorrectInitialisation()

void ParameterHandlerGeneric::CheckCorrectInitialisation ( )
protected

KS: Check if matrix is correctly initialised.

Definition at line 658 of file ParameterHandlerGeneric.cpp.

658  {
659 // ********************************************
660  // KS: Lambda Function which simply checks if there are no duplicates in std::vector
661  auto CheckForDuplicates = [](const std::vector<std::string>& names, const std::string& nameType) {
662  std::unordered_map<std::string, size_t> seenStrings;
663  for (size_t i = 0; i < names.size(); ++i) {
664  const auto& name = names[i];
665  if (seenStrings.find(name) != seenStrings.end()) {
666  size_t firstIndex = seenStrings[name];
667  MACH3LOG_CRITICAL("There are two systematics with the same {} '{}', first at index {}, and again at index {}", nameType, name, firstIndex, i);
668  throw MaCh3Exception(__FILE__, __LINE__);
669  }
670  seenStrings[name] = i;
671  }
672  };
673 
674  // KS: Checks if there are no duplicates in fancy names etc, this can happen if we merge configs etc
675  CheckForDuplicates(_fFancyNames, "_fFancyNames");
676  CheckForDuplicates(_fSplineNames, "_fSplineNames");
677 }
#define MACH3LOG_CRITICAL
Definition: MaCh3Logger.h:28
Custom exception class for MaCh3 errors.
std::vector< std::string > _fFancyNames
Fancy name for example rather than xsec_0 it is MAQE, useful for human reading.
std::vector< std::string > _fSplineNames
Name of spline in TTree (TBranch),.

◆ 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 759 of file ParameterHandlerGeneric.cpp.

759  {
760 // ********************************************
761  TFile* outputFile = new TFile(Name.c_str(), "RECREATE");
762 
763  TObjArray* xsec_param_names = new TObjArray();
764  TObjArray* xsec_spline_interpolation = new TObjArray();
765  TObjArray* xsec_spline_names = new TObjArray();
766 
767  TVectorD* xsec_param_prior = new TVectorD(_fNumPar);
768  TVectorD* xsec_flat_prior = new TVectorD(_fNumPar);
769  TVectorD* xsec_stepscale = new TVectorD(_fNumPar);
770  TVectorD* xsec_param_lb = new TVectorD(_fNumPar);
771  TVectorD* xsec_param_ub = new TVectorD(_fNumPar);
772 
773  TVectorD* xsec_param_knot_weight_lb = new TVectorD(_fNumPar);
774  TVectorD* xsec_param_knot_weight_ub = new TVectorD(_fNumPar);
775  TVectorD* xsec_error = new TVectorD(_fNumPar);
776 
777  for(int i = 0; i < _fNumPar; ++i)
778  {
779  TObjString* nameObj = new TObjString(_fFancyNames[i].c_str());
780  xsec_param_names->AddLast(nameObj);
781 
782  TObjString* splineType = new TObjString("TSpline3");
783  xsec_spline_interpolation->AddLast(splineType);
784 
785  TObjString* splineName = new TObjString("");
786  xsec_spline_names->AddLast(splineName);
787 
788  (*xsec_param_prior)[i] = _fPreFitValue[i];
789  (*xsec_flat_prior)[i] = _fFlatPrior[i];
790  (*xsec_stepscale)[i] = _fIndivStepScale[i];
791  (*xsec_error)[i] = _fError[i];
792 
793  (*xsec_param_lb)[i] = _fLowBound[i];
794  (*xsec_param_ub)[i] = _fUpBound[i];
795 
796  //Default values
797  (*xsec_param_knot_weight_lb)[i] = -9999;
798  (*xsec_param_knot_weight_ub)[i] = +9999;
799  }
800 
801  for (auto &pair : _fSystToGlobalSystIndexMap[SystType::kSpline]) {
802  auto &SplineIndex = pair.first;
803  auto &SystIndex = pair.second;
804 
805  (*xsec_param_knot_weight_lb)[SystIndex] = SplineParams.at(SplineIndex)._SplineKnotLowBound;
806  (*xsec_param_knot_weight_ub)[SystIndex] = SplineParams.at(SplineIndex)._SplineKnotUpBound;
807 
808  TObjString* splineType = new TObjString(SplineInterpolation_ToString(SplineParams.at(SplineIndex)._SplineInterpolationType).c_str());
809  xsec_spline_interpolation->AddAt(splineType, SystIndex);
810 
811  TObjString* splineName = new TObjString(_fSplineNames[SplineIndex].c_str());
812  xsec_spline_names->AddAt(splineName, SystIndex);
813  }
814  xsec_param_names->Write("xsec_param_names", TObject::kSingleKey);
815  delete xsec_param_names;
816  xsec_spline_interpolation->Write("xsec_spline_interpolation", TObject::kSingleKey);
817  delete xsec_spline_interpolation;
818  xsec_spline_names->Write("xsec_spline_names", TObject::kSingleKey);
819  delete xsec_spline_names;
820 
821  xsec_param_prior->Write("xsec_param_prior");
822  delete xsec_param_prior;
823  xsec_flat_prior->Write("xsec_flat_prior");
824  delete xsec_flat_prior;
825  xsec_stepscale->Write("xsec_stepscale");
826  delete xsec_stepscale;
827  xsec_param_lb->Write("xsec_param_lb");
828  delete xsec_param_lb;
829  xsec_param_ub->Write("xsec_param_ub");
830  delete xsec_param_ub;
831 
832  xsec_param_knot_weight_lb->Write("xsec_param_knot_weight_lb");
833  delete xsec_param_knot_weight_lb;
834  xsec_param_knot_weight_ub->Write("xsec_param_knot_weight_ub");
835  delete xsec_param_knot_weight_ub;
836  xsec_error->Write("xsec_error");
837  delete xsec_error;
838 
839  covMatrix->Write("xsec_cov");
840  TH2D* CorrMatrix = GetCorrelationMatrix();
841  CorrMatrix->Write("hcov");
842  delete CorrMatrix;
843 
844  outputFile->Close();
845  delete outputFile;
846 
847  MACH3LOG_INFO("Finished dumping ParameterHandler object");
848 }
#define MACH3LOG_INFO
Definition: MaCh3Logger.h:25
std::string SplineInterpolation_ToString(const SplineInterpolation i)
Convert a LLH type to a string.
@ kSpline
For splined parameters (1D)
std::vector< bool > _fFlatPrior
Whether to apply flat prior or not.
std::vector< double > _fError
Prior error on the parameter.
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.
std::vector< SplineParameter > SplineParams
Vector containing info for normalisation systematics.
std::vector< std::map< int, int > > _fSystToGlobalSystIndexMap
Map between number of given parameter type with global parameter numbering. For example 2nd norm para...
TH2D * GetCorrelationMatrix()
KS: Convert covariance matrix to correlation matrix and return TH2D which can be used for fancy plott...

◆ 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 218 of file ParameterHandlerGeneric.cpp.

218  {
219 // ********************************************
220  // KS: For now we don't use so avoid compilation error
221  (void) param;
222 
223  Parameter.name = GetParFancyName(Index);
224 
225  // Set the global parameter index of the normalisation parameter
226  Parameter.index = Index;
227 }
std::string GetParFancyName(const int i) const
Get fancy name of the Parameter.
int index
Parameter number of this normalisation in current systematic model.
std::string name
Name of parameters.

◆ 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 292 of file ParameterHandlerGeneric.cpp.

292  {
293 // ********************************************
294  FunctionalParameter func;
295  GetBaseParameter(param, Index, func);
296 
297  func.pdgs = GetFromManager<std::vector<int>>(param["NeutrinoFlavour"], std::vector<int>(), __FILE__ , __LINE__);
298  func.targets = GetFromManager<std::vector<int>>(param["TargetNuclei"], std::vector<int>(), __FILE__ , __LINE__);
299  func.modes = GetFromManager<std::vector<int>>(param["Mode"], std::vector<int>(), __FILE__ , __LINE__);
300  func.preoscpdgs = GetFromManager<std::vector<int>>(param["NeutrinoFlavourUnosc"], std::vector<int>(), __FILE__ , __LINE__);
301 
302  // HH - Copied from GetXsecNorm
303  int NumKinematicCuts = 0;
304  if(param["KinematicCuts"]){
305 
306  NumKinematicCuts = int(param["KinematicCuts"].size());
307 
308  std::vector<std::string> TempKinematicStrings;
309  std::vector<std::vector<std::vector<double>>> TempKinematicBounds;
310  //First element of TempKinematicBounds is always -999, and size is then 3
311  for(int KinVar_i = 0 ; KinVar_i < NumKinematicCuts ; ++KinVar_i){
312  //ETA: This is a bit messy, Kinematic cuts is a list of maps
313  for (YAML::const_iterator it = param["KinematicCuts"][KinVar_i].begin();it!=param["KinematicCuts"][KinVar_i].end();++it) {
314  TempKinematicStrings.push_back(it->first.as<std::string>());
315  TempKinematicBounds.push_back(Get2DBounds(it->second));
316  }
317  if(TempKinematicStrings.size() == 0) {
318  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++))");
319  MACH3LOG_ERROR("For Param {}", func.name);
320  throw MaCh3Exception(__FILE__, __LINE__);
321  }
322  }//KinVar_i
323  func.KinematicVarStr = TempKinematicStrings;
324  func.Selection = TempKinematicBounds;
325  }
326  func.valuePtr = RetPointer(Index);
327  return func;
328 }
int size
#define MACH3LOG_ERROR
Definition: MaCh3Logger.h:27
#define Get2DBounds(filename)
Definition: YamlHelper.h:563
const double * 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.
std::vector< std::string > KinematicVarStr
const double * valuePtr
Parameter value pointer.
std::vector< int > targets
Targets which parameter applies to.
std::vector< int > pdgs
PDG which parameter applies to.
std::vector< std::vector< std::vector< double > > > Selection
std::vector< int > preoscpdgs
Preosc PDG which parameter applies to.

◆ GetNormParameter()

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

Get Norm params.

Parameters
paramYaml node describing param
IndexGlobal parameter index

ETA size 0 to mean apply to all Ultimately all this information ends up in the @NormParams vector

Definition at line 162 of file ParameterHandlerGeneric.cpp.

162  {
163 // ********************************************
164  NormParameter norm;
165 
166  GetBaseParameter(param, Index, norm);
167 
170  norm.modes = GetFromManager<std::vector<int>>(param["Mode"], {}, __FILE__ , __LINE__);
171  norm.pdgs = GetFromManager<std::vector<int>>(param["NeutrinoFlavour"], {}, __FILE__ , __LINE__);
172  norm.preoscpdgs = GetFromManager<std::vector<int>>(param["NeutrinoFlavourUnosc"], {}, __FILE__ , __LINE__);
173  norm.targets = GetFromManager<std::vector<int>>(param["TargetNuclei"], {}, __FILE__ , __LINE__);
174 
175  if(_fLowBound[Index] < 0.) {
176  MACH3LOG_ERROR("Normalisation Parameter {} ({}), has lower parameters bound which can go below 0 and is equal {}",
177  GetParFancyName(Index), Index, _fLowBound[Index]);
178  MACH3LOG_ERROR("Normalisation parameters can't go bellow 0 as this is unphysical");
179  throw MaCh3Exception(__FILE__, __LINE__);
180  }
181  int NumKinematicCuts = 0;
182  if(param["KinematicCuts"]) {
183  NumKinematicCuts = int(param["KinematicCuts"].size());
184 
185  std::vector<std::string> TempKinematicStrings;
186  std::vector<std::vector<std::vector<double>>> TempKinematicBounds;
187  //First element of TempKinematicBounds is always -999, and size is then 3
188  for(int KinVar_i = 0 ; KinVar_i < NumKinematicCuts ; ++KinVar_i) {
189  //ETA: This is a bit messy, Kinematic cuts is a list of maps
190  for (YAML::const_iterator it = param["KinematicCuts"][KinVar_i].begin();it!=param["KinematicCuts"][KinVar_i].end();++it) {
191  TempKinematicStrings.push_back(it->first.as<std::string>());
192  TempKinematicBounds.push_back(Get2DBounds(it->second));
193  }
194  if(TempKinematicStrings.size() == 0) {
195  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++))");
196  MACH3LOG_ERROR("For Param {}", norm.name);
197  throw MaCh3Exception(__FILE__, __LINE__);
198  }
199  }//KinVar_i
200  norm.KinematicVarStr = TempKinematicStrings;
201  norm.Selection = TempKinematicBounds;
202  }
203 
204  //Next ones are kinematic bounds on where normalisation parameter should apply
205  //We set a bool to see if any bounds exist so we can short-circuit checking all of them every step
206  bool HasKinBounds = false;
207 
208  if(norm.KinematicVarStr.size() > 0) HasKinBounds = true;
209 
210  norm.hasKinBounds = HasKinBounds;
211  //End of kinematic bound checking
212 
213  return norm;
214 }
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.
bool hasKinBounds
Does this parameter have kinematic bounds.
std::vector< std::vector< std::vector< double > > > Selection
std::vector< int > modes
Mode which parameter applies to.
std::vector< int > pdgs
PDG which parameter applies to.
std::vector< std::string > KinematicVarStr
std::vector< int > targets
Targets which parameter applies to.

◆ 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 332 of file ParameterHandlerGeneric.cpp.

332  {
333 // ********************************************
334  OscillationParameter OscParamInfo;
335  GetBaseParameter(param, Index, OscParamInfo);
336 
337  return OscParamInfo;
338 }
KS: Struct holding info about oscillation Systematics.

◆ 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 40 of file ParameterHandlerGeneric.h.

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

◆ 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 263 of file ParameterHandlerGeneric.cpp.

263  {
264 // ********************************************
265  SplineParameter Spline;
266 
267  GetBaseParameter(param, Index, Spline);
268  //Now get the Spline interpolation type
269  if (param["SplineInformation"]["InterpolationType"]){
270  for(int InterpType = 0; InterpType < kSplineInterpolations ; ++InterpType){
271  if(param["SplineInformation"]["InterpolationType"].as<std::string>() == SplineInterpolation_ToString(SplineInterpolation(InterpType)))
272  Spline._SplineInterpolationType = SplineInterpolation(InterpType);
273  }
274  } else { //KS: By default use TSpline3
276  }
277  Spline._SplineKnotUpBound = GetFromManager<double>(param["SplineInformation"]["SplineKnotUpBound"], M3::DefSplineKnotUpBound, __FILE__ , __LINE__);
278  Spline._SplineKnotLowBound = GetFromManager<double>(param["SplineInformation"]["SplineKnotLowBound"], M3::DefSplineKnotLowBound, __FILE__ , __LINE__);
279 
281  MACH3LOG_WARN("Spline knot capping enabled with bounds [{}, {}]. For reliable fits, consider modifying the input generation instead.",
282  Spline._SplineKnotLowBound, Spline._SplineKnotUpBound);
283  }
284  //If there is no mode information given then this will be an empty vector
285  Spline._fSplineModes = GetFromManager(param["SplineInformation"]["Mode"], std::vector<int>(), __FILE__ , __LINE__);
286 
287  return Spline;
288 }
#define MACH3LOG_WARN
Definition: MaCh3Logger.h:26
SplineInterpolation
Make an enum of the spline interpolation type.
@ kTSpline3
Default TSpline3 interpolation.
@ kSplineInterpolations
This only enumerates.
Type GetFromManager(const YAML::Node &node, 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:299
constexpr static const double DefSplineKnotUpBound
Default value for spline knot capping, default mean not capping is being applied.
Definition: Core.h:79
constexpr static const double DefSplineKnotLowBound
Default value for spline knot capping, default mean not capping is being applied.
Definition: Core.h:81
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)

◆ 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 363 of file ParameterHandlerGeneric.cpp.

363  {
364 // ********************************************
365  std::vector<ParamT> returnVec;
366  for (const auto& pair : indexMap) {
367  const auto& localIndex = pair.first;
368  const auto& globalIndex = pair.second;
369  if (AppliesToSample(globalIndex, SampleName)) {
370  returnVec.push_back(params[localIndex]);
371  }
372  }
373  return returnVec;
374 }
bool AppliesToSample(const int SystIndex, const std::string &SampleName) const
Check if parameter is affecting given sample name.

◆ InitParametersTypeFromConfig()

void ParameterHandlerGeneric::InitParametersTypeFromConfig ( )
inlineprotected

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.find(SystType_ToString(SystType::kSpline)) != std::string::npos)
52  {
53  //Set param type
55  // Fill Spline info
56  SplineParams.push_back(GetSplineParameter(param["Systematic"], i));
57 
58  if (param["Systematic"]["SplineInformation"]["SplineName"]) {
59  _fSplineNames.push_back(param["Systematic"]["SplineInformation"]["SplineName"].as<std::string>());
60  }
61 
62  //Insert the mapping from the spline index i.e. the length of _fSplineNames etc
63  //to the Systematic index i.e. the counter for things like _fSampleID
64  _fSystToGlobalSystIndexMap[SystType::kSpline].insert(std::make_pair(ParamCounter[SystType::kSpline], i));
65  ParamCounter[SystType::kSpline]++;
66  } else if(param["Systematic"]["Type"].as<std::string>() == SystType_ToString(SystType::kNorm)) {
68  NormParams.push_back(GetNormParameter(param["Systematic"], i));
69  _fSystToGlobalSystIndexMap[SystType::kNorm].insert(std::make_pair(ParamCounter[SystType::kNorm], i));
70  ParamCounter[SystType::kNorm]++;
71  } else if(param["Systematic"]["Type"].as<std::string>() == SystType_ToString(SystType::kFunc)){
73  FuncParams.push_back(GetFunctionalParameters(param["Systematic"], i));
74  _fSystToGlobalSystIndexMap[SystType::kFunc].insert(std::make_pair(ParamCounter[SystType::kFunc], i));
75  ParamCounter[SystType::kFunc]++;
76  } else if(param["Systematic"]["Type"].as<std::string>() == SystType_ToString(SystType::kOsc)){
78  OscParams.push_back(GetOscillationParameters(param["Systematic"], i));
79  _fSystToGlobalSystIndexMap[SystType::kOsc].insert(std::make_pair(ParamCounter[SystType::kOsc], i));
80  ParamCounter[SystType::kOsc]++;
81  } else{
82  MACH3LOG_ERROR("Given unrecognised systematic type: {}", param["Systematic"]["Type"].as<std::string>());
83  std::string expectedTypes = "Expecting ";
84  for (int s = 0; s < SystType::kSystTypes; ++s) {
85  if (s > 0) expectedTypes += ", ";
86  expectedTypes += SystType_ToString(static_cast<SystType>(s)) + "\"";
87  }
88  expectedTypes += ".";
89  MACH3LOG_ERROR(expectedTypes);
90  throw MaCh3Exception(__FILE__, __LINE__);
91  }
92  i++;
93  } //end loop over params
94 
95  //Add a sanity check,
96  if(_fSplineNames.size() != ParamCounter[SystType::kSpline]){
97  MACH3LOG_ERROR("_fSplineNames is of size {} but found {} spline parameters", _fSplineNames.size(), ParamCounter[SystType::kSpline]);
98  throw MaCh3Exception(__FILE__, __LINE__);
99  }
100  //KS We resized them above to all params to fight memory fragmentation, now let's resize to fit only allocated memory to save RAM
101  NormParams.shrink_to_fit();
102  SplineParams.shrink_to_fit();
103  FuncParams.shrink_to_fit();
104  OscParams.shrink_to_fit();
105 }
std::string SystType_ToString(const SystType i)
Convert a Syst type type to a string.
SystType
@ kNorm
For normalisation parameters.
@ kSystTypes
This only enumerates.
@ kOsc
For oscillation parameters.
@ kFunc
For functional parameters.
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.
std::vector< NormParameter > NormParams
Vector containing info for normalisation systematics.
std::vector< SystType > _fParamType
Type of parameter like norm, spline etc.
NormParameter GetNormParameter(const YAML::Node &param, const int Index)
Get Norm params.
std::vector< FunctionalParameter > FuncParams
Vector containing info for functional systematics.
std::vector< OscillationParameter > OscParams
Vector containing info for functional systematics.
OscillationParameter GetOscillationParameters(const YAML::Node &param, const int Index)
Get Osc params.
std::vector< std::string > _ParameterGroup
KS: Allow to group parameters for example to affect only cross-section or only flux etc.

◆ 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 424 of file ParameterHandlerGeneric.cpp.

424  {
425 // ********************************************
426  for (int i = 0; i < _fNumPar; ++i) {
427  //ETA - set the name to be xsec_% as this is what ProcessorMCMC expects
428  _fNames[i] = "xsec_"+std::to_string(i);
429 
430  // KS: Plenty
431  if(_fParamType[i] == kOsc){
432  _fNames[i] = _fFancyNames[i];
433 
434  if(_ParameterGroup[i] != "Osc"){
435  MACH3LOG_ERROR("Parameter {}, is of type Oscillation but doesn't belong to Osc group", _fFancyNames[i]);
436  MACH3LOG_ERROR("It belongs to {} group", _ParameterGroup[i]);
437  throw MaCh3Exception(__FILE__ , __LINE__ );
438  }
439  }
440  // Set ParameterHandler parameters (Curr = current, Prop = proposed, Sigma = step)
441  _fCurrVal[i] = _fPreFitValue[i];
442  _fPropVal[i] = _fCurrVal[i];
443  }
444  Randomize();
445  //KS: Transfer the starting parameters to the PCA basis, you don't want to start with zero..
446  if (pca) {
447  PCAObj->SetInitialParameters(_fIndivStepScale);
448  }
449 }
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
std::vector< double > _fPropVal
Proposed value of the parameter.

◆ 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 414 of file ParameterHandlerGeneric.cpp.

414  {
415 // ********************************************
416  for (int i = 0; i < _fNumPar; ++i) {
417  if ((AppliesToSample(i, SampleName)) && filter(i)) { // Common filter logic
418  action(i); // Specific action for each function
419  }
420  }
421 }

◆ Print()

void ParameterHandlerGeneric::Print ( )
protected

Print information about the whole object once it is set.

Definition at line 453 of file ParameterHandlerGeneric.cpp.

453  {
454 // ********************************************
455  MACH3LOG_INFO("#################################################");
456  MACH3LOG_INFO("Printing ParameterHandlerGeneric:");
457 
459 
460  PrintNormParams();
461 
463 
465 
467 
469 
470  MACH3LOG_INFO("Finished");
471  MACH3LOG_INFO("#################################################");
472 
474 } // End
void PrintNormParams()
Prints normalization parameters.
void PrintParameterGroups()
Prints groups of parameters.
void PrintOscillationParams()
Prints oscillation parameters.
void PrintSplineParams()
Prints spline parameters.
void PrintFunctionalParams()
Prints functional parameters.
void CheckCorrectInitialisation()
KS: Check if matrix is correctly initialised.
void PrintGlobablInfo()
Prints general information about the ParameterHandler object.

◆ PrintFunctionalParams()

void ParameterHandlerGeneric::PrintFunctionalParams ( )
protected

Prints functional parameters.

Definition at line 592 of file ParameterHandlerGeneric.cpp.

592  {
593 // ********************************************
594  MACH3LOG_INFO("Functional parameters: {}", _fSystToGlobalSystIndexMap[SystType::kFunc].size());
596  MACH3LOG_INFO("┌────┬──────────┬────────────────────────────────────────┐");
597  MACH3LOG_INFO("│{0:4}│{1:10}│{2:40}│", "#", "Global #", "Name");
598  MACH3LOG_INFO("├────┼──────────┼────────────────────────────────────────┤");
599  for (auto &pair : _fSystToGlobalSystIndexMap[SystType::kFunc]) {
600  auto &FuncIndex = pair.first;
601  auto &GlobalIndex = pair.second;
602  MACH3LOG_INFO("│{0:4}│{1:<10}│{2:40}│", std::to_string(FuncIndex), GlobalIndex, GetParFancyName(GlobalIndex));
603  }
604  MACH3LOG_INFO("└────┴──────────┴────────────────────────────────────────┘");
605 }

◆ PrintGlobablInfo()

void ParameterHandlerGeneric::PrintGlobablInfo ( )
protected

Prints general information about the ParameterHandler object.

Definition at line 477 of file ParameterHandlerGeneric.cpp.

477  {
478 // ********************************************
479  MACH3LOG_INFO("============================================================================================================================================================");
480  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");
481  MACH3LOG_INFO("------------------------------------------------------------------------------------------------------------------------------------------------------------");
482  for (int i = 0; i < GetNumParams(); i++) {
483  std::string ErrString = fmt::format("{:.2f}", _fError[i]);
484  std::string SampleNameString = "";
485  for (const auto& SampleName : _fSampleNames[i]) {
486  if (!SampleNameString.empty()) {
487  SampleNameString += ", ";
488  }
489  SampleNameString += SampleName;
490  }
491  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]));
492  }
493  MACH3LOG_INFO("============================================================================================================================================================");
494 }
std::vector< std::vector< std::string > > _fSampleNames
Tells to which samples object param should be applied.
int GetNumParams() const
Get total number of parameters.

◆ PrintNormParams()

void ParameterHandlerGeneric::PrintNormParams ( )
protected

Prints normalization parameters.

Definition at line 497 of file ParameterHandlerGeneric.cpp.

497  {
498 // ********************************************
499  // Output the normalisation parameters as a sanity check!
500  MACH3LOG_INFO("Normalisation parameters: {}", NormParams.size());
502 
503  bool have_parameter_with_kin_bounds = false;
504 
505  //KS: Consider making some class producing table..
506  MACH3LOG_INFO("┌────┬──────────┬────────────────────────────────────────┬────────────────────┬────────────────────┬────────────────────┐");
507  MACH3LOG_INFO("│{0:4}│{1:10}│{2:40}│{3:20}│{4:20}│{5:20}│", "#", "Global #", "Name", "Int. mode", "Target", "pdg");
508  MACH3LOG_INFO("├────┼──────────┼────────────────────────────────────────┼────────────────────┼────────────────────┼────────────────────┤");
509 
510  for (unsigned int i = 0; i < NormParams.size(); ++i)
511  {
512  std::string intModeString;
513  for (unsigned int j = 0; j < NormParams[i].modes.size(); j++) {
514  intModeString += std::to_string(NormParams[i].modes[j]);
515  intModeString += " ";
516  }
517  if (NormParams[i].modes.empty()) intModeString += "all";
518 
519  std::string targetString;
520  for (unsigned int j = 0; j < NormParams[i].targets.size(); j++) {
521  targetString += std::to_string(NormParams[i].targets[j]);
522  targetString += " ";
523  }
524  if (NormParams[i].targets.empty()) targetString += "all";
525 
526  std::string pdgString;
527  for (unsigned int j = 0; j < NormParams[i].pdgs.size(); j++) {
528  pdgString += std::to_string(NormParams[i].pdgs[j]);
529  pdgString += " ";
530  }
531  if (NormParams[i].pdgs.empty()) pdgString += "all";
532 
533  MACH3LOG_INFO("│{: <4}│{: <10}│{: <40}│{: <20}│{: <20}│{: <20}│", i, NormParams[i].index, NormParams[i].name, intModeString, targetString, pdgString);
534 
535  if(NormParams[i].hasKinBounds) have_parameter_with_kin_bounds = true;
536  }
537  MACH3LOG_INFO("└────┴──────────┴────────────────────────────────────────┴────────────────────┴────────────────────┴────────────────────┘");
538 
539  if(have_parameter_with_kin_bounds) {
540  MACH3LOG_INFO("Normalisation parameters KinematicCuts information");
541  MACH3LOG_INFO("┌────┬──────────┬────────────────────────────────────────┬────────────────────┬────────────────────────────────────────┐");
542  MACH3LOG_INFO("│{0:4}│{1:10}│{2:40}│{3:20}│{4:40}│", "#", "Global #", "Name", "KinematicCut", "Value");
543  MACH3LOG_INFO("├────┼──────────┼────────────────────────────────────────┼────────────────────┼────────────────────────────────────────┤");
544  for (unsigned int i = 0; i < NormParams.size(); ++i)
545  {
546  //skip parameters with no KinematicCuts
547  if(!NormParams[i].hasKinBounds) continue;
548 
549  const long unsigned int ncuts = NormParams[i].KinematicVarStr.size();
550  for(long unsigned int icut = 0; icut < ncuts; icut++) {
551  std::string kinematicCutValueString;
552  for(const auto & value : NormParams[i].Selection[icut]) {
553  for (const auto& v : value) {
554  kinematicCutValueString += fmt::format("{:.2f} ", v);
555  }
556  }
557  if(icut == 0)
558  MACH3LOG_INFO("│{: <4}│{: <10}│{: <40}│{: <20}│{: <40}│", i, NormParams[i].index, NormParams[i].name, NormParams[i].KinematicVarStr[icut], kinematicCutValueString);
559  else
560  MACH3LOG_INFO("│{: <4}│{: <10}│{: <40}│{: <20}│{: <40}│", "", "", "", NormParams[i].KinematicVarStr[icut], kinematicCutValueString);
561  }//icut
562  }//i
563  MACH3LOG_INFO("└────┴──────────┴────────────────────────────────────────┴────────────────────┴────────────────────────────────────────┘");
564  }
565  else
566  MACH3LOG_INFO("No normalisation parameters have KinematicCuts defined");
567 }

◆ PrintOscillationParams()

void ParameterHandlerGeneric::PrintOscillationParams ( )
protected

Prints oscillation parameters.

Definition at line 608 of file ParameterHandlerGeneric.cpp.

608  {
609 // ********************************************
610  MACH3LOG_INFO("Oscillation parameters: {}", _fSystToGlobalSystIndexMap[SystType::kOsc].size());
612  MACH3LOG_INFO("┌────┬──────────┬────────────────────────────────────────┐");
613  MACH3LOG_INFO("│{0:4}│{1:10}│{2:40}│", "#", "Global #", "Name");
614  MACH3LOG_INFO("├────┼──────────┼────────────────────────────────────────┤");
615  for (auto &pair : _fSystToGlobalSystIndexMap[SystType::kOsc]) {
616  auto &OscIndex = pair.first;
617  auto &GlobalIndex = pair.second;
618  MACH3LOG_INFO("│{0:4}│{1:<10}│{2:40}│", std::to_string(OscIndex), GlobalIndex, GetParFancyName(GlobalIndex));
619  }
620  MACH3LOG_INFO("└────┴──────────┴────────────────────────────────────────┘");
621 }

◆ PrintParameterGroups()

void ParameterHandlerGeneric::PrintParameterGroups ( )
protected

Prints groups of parameters.

Definition at line 624 of file ParameterHandlerGeneric.cpp.

624  {
625 // ********************************************
626  // KS: Create a map to store the counts of unique strings, in principle this could be in header file
627  std::unordered_map<std::string, int> paramCounts;
628 
629  std::for_each(_ParameterGroup.begin(), _ParameterGroup.end(),
630  [&paramCounts](const std::string& param) {
631  paramCounts[param]++;
632  });
633 
634  MACH3LOG_INFO("Printing parameter groups");
635  // Output the counts
636  for (const auto& pair : paramCounts) {
637  MACH3LOG_INFO("Found {}: {} params", pair.second, pair.first);
638  }
639 }

◆ PrintSplineParams()

void ParameterHandlerGeneric::PrintSplineParams ( )
protected

Prints spline parameters.

Definition at line 570 of file ParameterHandlerGeneric.cpp.

570  {
571 // ********************************************
572  MACH3LOG_INFO("Spline parameters: {}", _fSystToGlobalSystIndexMap[SystType::kSpline].size());
574  MACH3LOG_INFO("=====================================================================================================================================================================");
575  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", "|");
576  MACH3LOG_INFO("---------------------------------------------------------------------------------------------------------------------------------------------------------------------");
577  for (auto &pair : _fSystToGlobalSystIndexMap[SystType::kSpline]) {
578  auto &SplineIndex = pair.first;
579  auto &GlobalIndex = pair.second;
580 
581  MACH3LOG_INFO("{:<4} {:<2} {:<40} {:<2} {:<40} {:<2} {:<20} {:<2} {:<20} {:<2} {:<20} {:<2}",
582  SplineIndex, "|", GetParFancyName(GlobalIndex), "|",
583  _fSplineNames[SplineIndex], "|",
585  GetParSplineKnotLowerBound(SplineIndex), "|",
586  GetParSplineKnotUpperBound(SplineIndex), "|");
587  }
588  MACH3LOG_INFO("=====================================================================================================================================================================");
589 }
SplineInterpolation GetParSplineInterpolation(const int i) const
Get interpolation type for a given parameter.
double GetParSplineKnotUpperBound(const int i) const
EM: value at which we cap spline knot weight.
double GetParSplineKnotLowerBound(const int i) const
EM: value at which we cap spline knot weight.

Member Data Documentation

◆ _fParamType

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

Type of parameter like norm, spline etc.

Definition at line 213 of file ParameterHandlerGeneric.h.

◆ _fSplineNames

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

Name of spline in TTree (TBranch),.

Definition at line 216 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 222 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 219 of file ParameterHandlerGeneric.h.

◆ FuncParams

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

Vector containing info for functional systematics.

Definition at line 231 of file ParameterHandlerGeneric.h.

◆ NormParams

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

Vector containing info for normalisation systematics.

Definition at line 228 of file ParameterHandlerGeneric.h.

◆ OscParams

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

Vector containing info for functional systematics.

Definition at line 234 of file ParameterHandlerGeneric.h.

◆ SplineParams

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

Vector containing info for normalisation systematics.

Definition at line 225 of file ParameterHandlerGeneric.h.


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