MaCh3  2.2.3
Reference Guide
ParameterHandlerBase.h
Go to the documentation of this file.
1 #pragma once
2 
3 // MaCh3 includes
4 #include "Manager/Manager.h"
9 
16  public:
23  ParameterHandlerBase(const std::vector<std::string>& YAMLFile, std::string name, double threshold = -1, int FirstPCAdpar = -999, int LastPCAdpar = -999);
27  ParameterHandlerBase(std::string name, std::string file, double threshold = -1, int FirstPCAdpar = -999, int LastPCAdpar = -999);
28 
30  virtual ~ParameterHandlerBase();
31 
34 
37 
38  // ETA - maybe need to add checks to index on the setters? i.e. if( i > _fPropVal.size()){throw;}
42  void SetCovMatrix(TMatrixDSym *cov);
45  void SetName(const std::string& name) { matrixName = name; }
50  void SetParName(const int i, const std::string& name) { _fNames.at(i) = name; }
53  void SetSingleParameter(const int parNo, const double parVal);
58  void SetPar(const int i, const double val);
63  void SetParCurrProp(const int i, const double val);
68  void SetParProp(const int i, const double val) {
69  _fPropVal[i] = val;
70  if (pca) PCAObj->TransferToPCA();
71  }
75  void SetParameters(const std::vector<double>& pars = {});
80  void SetFlatPrior(const int i, const bool eL);
81 
86  void SetRandomThrow(const int i, const double rand) { randParams[i] = rand;}
90  double GetRandomThrow(const int i) const { return randParams[i];}
91 
96  void SetBranches(TTree &tree, const bool SaveProposal = false);
102  void SetStepScale(const double scale, const bool verbose=true);
107  void SetIndivStepScale(const int ParameterIndex, const double StepScale){ _fIndivStepScale.at(ParameterIndex) = StepScale; }
111  void SetIndivStepScale(const std::vector<double>& stepscale);
114  inline void SetPrintLength(const unsigned int PriLen) { PrintLength = PriLen; }
115 
118 
120  void ThrowParProp(const double mag = 1.);
121 
123  void ThrowParCurr(const double mag = 1.);
125  void ThrowParameters();
127  void RandomConfiguration();
128 
130  int CheckBounds() const _noexcept_;
132  double CalcLikelihood() const _noexcept_;
135  virtual double GetLikelihood();
136 
139  TMatrixDSym *GetCovMatrix() const { return covMatrix; }
142  TMatrixDSym *GetInvCovMatrix() const { return invCovMatrix; }
145  double GetInvCovMatrix(const int i, const int j) const { return InvertCovMatrix[i][j]; }
146 
150  double GetCorrThrows(const int i) const { return corr_throw[i]; }
151 
155  inline bool GetFlatPrior(const int i) const { return _fFlatPrior[i]; }
156 
159  std::string GetName() const { return matrixName; }
163  std::string GetParName(const int i) const {return _fNames[i];}
164 
167  int GetParIndex(const std::string& name) const;
168 
172  std::string GetParFancyName(const int i) const {return _fFancyNames[i];}
175  std::string GetInputFile() const { return inputFile; }
176 
180  inline double GetDiagonalError(const int i) const { return std::sqrt((*covMatrix)(i,i)); }
184  inline double GetError(const int i) const {return _fError[i];}
185 
187  void ResetIndivStepScale();
188 
191  void InitialiseAdaption(const YAML::Node& adapt_manager);
193  void SaveAdaptiveToFile(const std::string& outFileName, const std::string& systematicName) {
194  AdaptiveHandler->SaveAdaptiveToFile(outFileName, systematicName); }
195 
198  bool GetDoAdaption() const {return use_adaptive;}
201  void SetThrowMatrix(TMatrixDSym *cov);
202  void UpdateThrowMatrix(TMatrixDSym *cov);
205  inline void SetNumberOfSteps(const int nsteps) {
206  AdaptiveHandler->SetTotalSteps(nsteps);
207  if(AdaptiveHandler->AdaptionUpdate()) ResetIndivStepScale();
208  }
209 
212  inline TMatrixDSym *GetThrowMatrix() const {return throwMatrix;}
215  double GetThrowMatrix(const int i, const int j) const { return throwMatrixCholDecomp[i][j];}
216 
222  TH2D* GetCorrelationMatrix();
223 
233  inline const double* RetPointer(const int iParam) {return &(_fPropVal.data()[iParam]);}
234 
237  inline const std::vector<double> &GetParPropVec() {return _fPropVal;}
238 
239  //Some Getters
242  inline int GetNumParams() const {return _fNumPar;}
245  std::vector<double> GetPreFitValues() const {return _fPreFitValue;}
248  std::vector<double> GetProposed() const;
252  inline double GetParProp(const int i) const { return _fPropVal[i]; }
256  inline double GetParCurr(const int i) const { return _fCurrVal[i]; }
259  inline const std::vector<double> &GetParCurrVec() const { return _fCurrVal; }
260 
264  inline double GetParInit(const int i) const { return _fPreFitValue[i]; }
268  inline double GetUpperBound(const int i) const { return _fUpBound[i];}
272  inline double GetLowerBound(const int i) const { return _fLowBound[i]; }
276  inline double GetIndivStepScale(const int ParameterIndex) const {return _fIndivStepScale.at(ParameterIndex); }
279  inline double GetGlobalStepScale() const {return _fGlobalStepScale; }
280 
283  inline int GetNParameters() const {
284  if (pca) return PCAObj->GetNumberPCAedParameters();
285  else return _fNumPar;
286  }
287 
289  void PrintNominal() const;
291  void PrintNominalCurrProp() const;
296  void PrintIndivStepScale() const;
297 
299  virtual void ProposeStep();
301  void Randomize() _noexcept_;
303  void CorrelateSteps() _noexcept_;
307 
309  void AcceptStep() _noexcept_;
310 
312  void ToggleFixAllParameters();
315  void ToggleFixParameter(const int i);
318  void ToggleFixParameter(const std::string& name);
321  bool IsParameterFixed(const int i) const {
322  if (_fError[i] < 0) { return true; }
323  else { return false; }
324  }
327  bool IsParameterFixed(const std::string& name) const;
328 
334  void ConstructPCA(const double eigen_threshold, int FirstPCAdpar, int LastPCAdpar);
335 
337  inline bool IsPCA() const { return pca; }
338 
341  YAML::Node GetConfig() const { return _fYAMLDoc; }
342 
346  if (!use_adaptive) {
347  MACH3LOG_ERROR("Am not running in Adaptive mode");
348  throw MaCh3Exception(__FILE__ , __LINE__ );
349  }
350  return AdaptiveHandler.get();
351  }
352 
355  void SetTune(const std::string& TuneName);
356 
358  inline PCAHandler* GetPCAHandler() const {
359  if (!pca) {
360  MACH3LOG_ERROR("Am not running in PCA mode");
361  throw MaCh3Exception(__FILE__ , __LINE__ );
362  }
363  return PCAObj.get();
364  }
365 
373  void MatchMaCh3OutputBranches(TTree *PosteriorFile,
374  std::vector<double>& BranchValues,
375  std::vector<std::string>& BranchNames);
376 
377 protected:
379  void Init(const std::string& name, const std::string& file);
382  void Init(const std::vector<std::string>& YAMLFile);
385  void ReserveMemory(const int size);
386 
389  void MakePosDef(TMatrixDSym *cov = nullptr);
390 
392  void MakeClosestPosDef(TMatrixDSym *cov);
393 
398  void SetThrowMatrixFromFile(const std::string& matrix_file_name, const std::string& matrix_name, const std::string& means_name);
399 
403  bool AppliesToSample(const int SystIndex, const std::string& SampleName) const;
404 
408  void FlipParameterValue(const int index, const double FlipPoint);
409 
412  void CircularParBounds(const int i, const double LowBound, const double UpBound);
413 
415  void EnableSpecialProposal(const YAML::Node& param, const int Index);
416 
418  void SpecialStepProposal();
419 
422 
424  const std::string inputFile;
425 
427  std::string matrixName;
429  TMatrixDSym *covMatrix;
431  TMatrixDSym *invCovMatrix;
433  std::vector<std::vector<double>> InvertCovMatrix;
434 
436  std::vector<std::unique_ptr<TRandom3>> random_number;
437 
439  double* randParams;
441  double* corr_throw;
444 
447 
449  std::vector<std::string> _fNames;
451  std::vector<std::string> _fFancyNames;
453  YAML::Node _fYAMLDoc;
455  int _fNumPar;
457  std::vector<double> _fPreFitValue;
459  std::vector<double> _fCurrVal;
461  std::vector<double> _fPropVal;
463  std::vector<double> _fError;
465  std::vector<double> _fLowBound;
467  std::vector<double> _fUpBound;
469  std::vector<double> _fIndivStepScale;
471  std::vector<bool> _fFlatPrior;
473  std::vector<std::vector<std::string>> _fSampleNames;
474 
476  TMatrixDSym* throwMatrix;
479 
481  bool pca;
484 
486  std::unique_ptr<PCAHandler> PCAObj;
488  std::unique_ptr<adaptive_mcmc::AdaptiveMCMCHandler> AdaptiveHandler;
490  std::unique_ptr<ParameterTunes> Tunes;
491 
493  std::vector<int> FlipParameterIndex;
495  std::vector<double> FlipParameterPoint;
497  std::vector<int> CircularBoundsIndex;
499  std::vector<std::pair<double,double>> CircularBoundsValues;
500 };
#define _noexcept_
KS: noexcept can help with performance but is terrible for debugging, this is meant to help easy way ...
Definition: Core.h:86
int size
std::vector< TString > BranchNames
#define MACH3LOG_ERROR
Definition: MaCh3Logger.h:27
Custom exception class for MaCh3 errors.
Class responsible for handling Principal Component Analysis (PCA) of covariance matrix.
Definition: PCAHandler.h:60
Base class responsible for handling of systematic error parameters. Capable of using PCA or using ada...
std::vector< int > CircularBoundsIndex
Indices of parameters with circular bounds.
void ToggleFixParameter(const int i)
fix parameter at prior values
bool use_adaptive
Are we using AMCMC?
virtual void ProposeStep()
Generate a new proposed state.
double * randParams
Random number taken from gaussian around prior error used for corr_throw.
double ** throwMatrixCholDecomp
Throw matrix that is being used in the fit, much faster as TMatrixDSym cache miss.
std::unique_ptr< adaptive_mcmc::AdaptiveMCMCHandler > AdaptiveHandler
Struct containing information about adaption.
virtual ~ParameterHandlerBase()
Destructor.
TMatrixDSym * invCovMatrix
The inverse covariance matrix.
const double * RetPointer(const int iParam)
DB Pointer return to param position.
std::vector< int > FlipParameterIndex
Indices of parameters with flip symmetry.
void SetThrowMatrixFromFile(const std::string &matrix_file_name, const std::string &matrix_name, const std::string &means_name)
sets throw matrix from a file
bool AppliesToSample(const int SystIndex, const std::string &SampleName) const
Check if parameter is affecting given sample name.
void AcceptStep() _noexcept_
Accepted this step.
void ThrowParProp(const double mag=1.)
Throw the proposed parameter by mag sigma. Should really just have the user specify this throw by hav...
std::unique_ptr< ParameterTunes > Tunes
Struct containing information about adaption.
void InitialiseAdaption(const YAML::Node &adapt_manager)
Initialise adaptive MCMC.
TMatrixDSym * throwMatrix
Matrix which we use for step proposal before Cholesky decomposition (not actually used for step propo...
void ThrowParCurr(const double mag=1.)
Helper function to throw the current parameter by mag sigma. Can study bias in MCMC with this; put di...
double _fGlobalStepScale
Global step scale applied to all params in this class.
void ReserveMemory(const int size)
Initialise vectors with parameters information.
std::vector< std::string > _fFancyNames
Fancy name for example rather than xsec_0 it is MAQE, useful for human reading.
bool doSpecialStepProposal
Check if any of special step proposal were enabled.
std::vector< bool > _fFlatPrior
Whether to apply flat prior or not.
void ThrowParameters()
Throw the parameters according to the covariance matrix. This shouldn't be used in MCMC code ase it c...
void Init(const std::string &name, const std::string &file)
Initialisation of the class using matrix from root file.
std::string matrixName
Name of cov matrix.
void MakeClosestPosDef(TMatrixDSym *cov)
HW: Finds closest possible positive definite matrix in Frobenius Norm ||.||_frob Where ||X||_frob=sqr...
void RandomConfiguration()
Randomly throw the parameters in their 1 sigma range.
void SaveAdaptiveToFile(const std::string &outFileName, const std::string &systematicName)
Save adaptive throw matrix to file.
std::vector< double > FlipParameterPoint
Central points around which parameters are flipped.
void UpdateAdaptiveCovariance()
Method to update adaptive MCMC .
std::vector< double > _fError
Prior error on the parameter.
PCAHandler * GetPCAHandler() const
Get pointer for PCAHandler.
void PrintNominal() const
Print prior value for every parameter.
YAML::Node _fYAMLDoc
Stores config describing systematics.
bool IsParameterFixed(const int i) const
Is parameter fixed or not.
const std::vector< double > & GetParPropVec()
Get a reference to the proposed parameter values Can be useful if you want to track these without hav...
void UpdateThrowMatrix(TMatrixDSym *cov)
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.
void SaveUpdatedMatrixConfig()
KS: After step scale, prefit etc. value were modified save this modified config.
void ResetIndivStepScale()
Adaptive Step Tuning Stuff.
std::unique_ptr< PCAHandler > PCAObj
Struct containing information about PCA.
int _fNumPar
Number of systematic parameters.
double CalcLikelihood() const _noexcept_
Calc penalty term based on inverted covariance matrix.
void FlipParameterValue(const int index, const double FlipPoint)
KS: Flip parameter around given value, for example mass ordering around 0.
void ConstructPCA(const double eigen_threshold, int FirstPCAdpar, int LastPCAdpar)
CW: Calculate eigen values, prepare transition matrices and remove param based on defined threshold.
std::vector< double > _fLowBound
Lowest physical bound, parameter will not be able to go beyond it.
std::vector< double > _fCurrVal
Current value of the parameter.
void PrintNominalCurrProp() const
Print prior, current and proposed value for each parameter.
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::vector< double > > InvertCovMatrix
KS: Same as above but much faster as TMatrixDSym cache miss.
TMatrixDSym * covMatrix
The covariance matrix.
std::vector< double > _fPreFitValue
Parameter value dictated by the prior model. Based on it penalty term is calculated.
void Randomize() _noexcept_
"Randomize" the parameters in the covariance class for the proposed step. Used the proposal kernel an...
void PrintIndivStepScale() const
Print step scale for each parameter.
std::vector< std::string > _fNames
ETA _fNames is set automatically in the covariance class to be something like xsec_i,...
std::vector< std::pair< double, double > > CircularBoundsValues
Circular bounds for each parameter (lower, upper)
void EnableSpecialProposal(const YAML::Node &param, const int Index)
Enable special proposal.
double * corr_throw
Result of multiplication of Cholesky matrix and randParams.
bool IsPCA() const
is PCA, can use to query e.g. LLH scans
void MatchMaCh3OutputBranches(TTree *PosteriorFile, std::vector< double > &BranchValues, std::vector< std::string > &BranchNames)
Matches branches in a TTree to parameters in a systematic handler.
std::vector< double > _fUpBound
Upper physical bound, parameter will not be able to go beyond it.
const std::string inputFile
The input root file we read in.
bool pca
perform PCA or not
std::vector< double > _fIndivStepScale
Individual step scale used by MCMC algorithm.
std::vector< std::unique_ptr< TRandom3 > > random_number
KS: Set Random numbers for each thread so each thread has different seed.
void MakePosDef(TMatrixDSym *cov=nullptr)
Make matrix positive definite by adding small values to diagonal, necessary for inverting matrix.
void SpecialStepProposal()
Perform Special Step Proposal.
int PrintLength
KS: This is used when printing parameters, sometimes we have super long parameters name,...
void CorrelateSteps() _noexcept_
Use Cholesky throw matrix for better step proposal.
std::vector< double > _fPropVal
Proposed value of the parameter.
std::vector< std::vector< std::string > > _fSampleNames
Tells to which samples object param should be applied.
int CheckBounds() const _noexcept_
Check if parameters were proposed outside physical boundary.
void ToggleFixAllParameters()
fix parameters at prior values
Contains information about adaptive covariance matrix .
int GetNumParams() const
Get total number of parameters.
std::string GetInputFile() const
Get name of input file.
std::string GetParName(const int i) const
Get name of parameter.
double GetInvCovMatrix(const int i, const int j) const
Return inverted covariance matrix.
bool GetDoAdaption() const
Do we adapt or not.
std::vector< double > GetProposed() const
Get vector of all proposed parameter values.
double GetUpperBound(const int i) const
Get upper parameter bound in which it is physically valid.
virtual double GetLikelihood()
Return CalcLikelihood if some params were thrown out of boundary return LARGE_LOGL
int GetNParameters() const
Get number of params which will be different depending if using Eigen decomposition or not.
int GetParIndex(const std::string &name) const
Get index based on name.
TMatrixDSym * GetInvCovMatrix() const
Return inverted covariance matrix.
TMatrixDSym * GetCovMatrix() const
Return covariance matrix.
double GetThrowMatrix(const int i, const int j) const
Get matrix used for step proposal.
std::string GetParFancyName(const int i) const
Get fancy name of the Parameter.
TH2D * GetCorrelationMatrix()
KS: Convert covariance matrix to correlation matrix and return TH2D which can be used for fancy plott...
const std::vector< double > & GetParCurrVec() const
Get vector of current parameter values.
double GetRandomThrow(const int i) const
Get random value useful for debugging/CI.
double GetLowerBound(const int i) const
Get lower parameter bound in which it is physically valid.
double GetGlobalStepScale() const
Get global step scale for covariance object.
bool GetFlatPrior(const int i) const
Get if param has flat prior or not.
std::vector< double > GetPreFitValues() const
Get the pre-fit values of the parameters.
double GetIndivStepScale(const int ParameterIndex) const
Get individual step scale for selected parameter.
double GetError(const int i) const
Get the error for the ith parameter.
std::string GetName() const
Get name of covariance.
double GetParInit(const int i) const
Get prior parameter value.
double GetDiagonalError(const int i) const
Get diagonal error for ith parameter.
YAML::Node GetConfig() const
Getter to return a copy of the YAML node.
double GetParProp(const int i) const
Get proposed parameter value.
double GetCorrThrows(const int i) const
Return correlated throws.
double GetParCurr(const int i) const
Get current parameter value.
TMatrixDSym * GetThrowMatrix() const
Get matrix used for step proposal.
adaptive_mcmc::AdaptiveMCMCHandler * GetAdaptiveHandler() const
Get pointer for AdaptiveHandler.
void SetName(const std::string &name)
Set matrix name.
void SetStepScale(const double scale, const bool verbose=true)
Set global step scale for covariance object.
void SetRandomThrow(const int i, const double rand)
Set random value useful for debugging/CI.
void SetCovMatrix(TMatrixDSym *cov)
Set covariance matrix.
void SetFlatPrior(const int i, const bool eL)
Set if parameter should have flat prior or not.
void SetParProp(const int i, const double val)
Set proposed parameter value.
void SetThrowMatrix(TMatrixDSym *cov)
Use new throw matrix, used in adaptive MCMC.
void SetParameters(const std::vector< double > &pars={})
Set parameter values using vector, it has to have same size as covariance class.
void SetPrintLength(const unsigned int PriLen)
KS: In case someone really want to change this.
void SetBranches(TTree &tree, const bool SaveProposal=false)
set branches for output file
void SetParName(const int i, const std::string &name)
change parameter name
void SetNumberOfSteps(const int nsteps)
Set number of MCMC step, when running adaptive MCMC it is updated with given frequency....
void SetParCurrProp(const int i, const double val)
Set current parameter value.
void SetSingleParameter(const int parNo, const double parVal)
Set value of single param to a given value.
void SetPar(const int i, const double val)
Set all the covariance matrix parameters to a user-defined value.
void SetIndivStepScale(const int ParameterIndex, const double StepScale)
DB Function to set fIndivStepScale from a vector (Can be used from execs and inside covariance constr...
void SetTune(const std::string &TuneName)
KS: Set proposed parameter values vector to be base on tune values, for example set proposed values t...