MaCh3  2.4.2
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 
32 
33  // ETA - maybe need to add checks to index on the setters? i.e. if( i > _fPropVal.size()){throw;}
36  void SetCovMatrix(TMatrixDSym *cov);
38  void SetName(const std::string& name) { matrixName = name; }
42  void SetParName(const int i, const std::string& name) { _fNames.at(i) = name; }
44  void SetSingleParameter(const int parNo, const double parVal);
48  void SetPar(const int i, const double val);
52  void SetParCurrProp(const int i, const double val);
56  void SetParProp(const int i, const double val) {
57  _fPropVal[i] = val;
58  if (pca) PCAObj->TransferToPCA();
59  }
62  void SetParameters(const std::vector<double>& pars = {});
66  void SetFlatPrior(const int i, const bool eL);
67 
71  void SetRandomThrow(const int i, const double rand) { randParams[i] = rand;}
74  double GetRandomThrow(const int i) const { return randParams[i];}
75 
79  void SetBranches(TTree &tree, const bool SaveProposal = false);
84  void SetStepScale(const double scale, const bool verbose=true);
88  void SetIndivStepScale(const int ParameterIndex, const double StepScale){ _fIndivStepScale.at(ParameterIndex) = StepScale; }
91  void SetIndivStepScale(const std::vector<double>& stepscale);
93  inline void SetPrintLength(const unsigned int PriLen) { PrintLength = PriLen; }
94 
97 
99  void ThrowParameters();
101  void RandomConfiguration();
102 
104  int CheckBounds() const _noexcept_;
118  double CalcLikelihood() const _noexcept_;
120  virtual double GetLikelihood();
121 
123  TMatrixDSym *GetCovMatrix() const { return covMatrix; }
125  TMatrixDSym *GetInvCovMatrix() const { return invCovMatrix; }
127  double GetInvCovMatrix(const int i, const int j) const { return InvertCovMatrix[i][j]; }
128 
131  double GetCorrThrows(const int i) const { return corr_throw[i]; }
132 
135  bool GetFlatPrior(const int i) const { return _fFlatPrior[i]; }
136 
138  std::string GetName() const { return matrixName; }
141  std::string GetParName(const int i) const {return _fNames[i];}
142 
144  int GetParIndex(const std::string& name) const;
145 
148  std::string GetParFancyName(const int i) const {return _fFancyNames[i];}
150  std::string GetInputFile() const { return inputFile; }
151 
154  inline double GetDiagonalError(const int i) const { return std::sqrt((*covMatrix)(i,i)); }
157  inline double GetError(const int i) const {return _fError[i];}
158 
160  void ResetIndivStepScale();
161 
164 
167  void InitialiseAdaption(const YAML::Node& adapt_manager);
169  void SaveAdaptiveToFile(const std::string& outFileName, const std::string& systematicName) {
170  AdaptiveHandler->SaveAdaptiveToFile(outFileName, systematicName); }
171 
173  bool GetDoAdaption() const {return use_adaptive;}
175  void SetThrowMatrix(TMatrixDSym *cov);
176  void SetSubThrowMatrix(int first_index, int last_index, TMatrixDSym const &subcov);
178  void UpdateThrowMatrix(TMatrixDSym *cov);
180  inline void SetNumberOfSteps(const int nsteps) {
181  AdaptiveHandler->SetTotalSteps(nsteps);
182  if(AdaptiveHandler->AdaptionUpdate()) ResetIndivStepScale();
183  }
184 
186  inline TMatrixDSym *GetThrowMatrix() const {return throwMatrix;}
188  double GetThrowMatrix(const int i, const int j) const { return throwMatrixCholDecomp[i][j];}
189 
194  TH2D* GetCorrelationMatrix();
195 
205  inline const double* RetPointer(const int iParam) {return &(_fPropVal.data()[iParam]);}
206 
209  inline const std::vector<double> &GetParPropVec() {return _fPropVal;}
210 
212  int GetNumParams() const {return _fNumPar;}
214  std::vector<double> GetPreFitValues() const {return _fPreFitValue;}
216  std::vector<double> GetProposed() const;
219  double GetParProp(const int i) const { return _fPropVal[i]; }
222  inline double GetParCurr(const int i) const { return _fCurrVal[i]; }
224  inline const std::vector<double> &GetParCurrVec() const { return _fCurrVal; }
225 
228  inline double GetParInit(const int i) const { return _fPreFitValue[i]; }
231  inline double GetUpperBound(const int i) const { return _fUpBound[i];}
234  inline double GetLowerBound(const int i) const { return _fLowBound[i]; }
237  inline double GetIndivStepScale(const int ParameterIndex) const {return _fIndivStepScale.at(ParameterIndex); }
239  inline double GetGlobalStepScale() const {return _fGlobalStepScale; }
240 
242  inline int GetNParameters() const {
243  if (pca) return PCAObj->GetNumberPCAedParameters();
244  else return _fNumPar;
245  }
246 
248  void PrintNominal() const;
250  void PrintNominalCurrProp() const;
255  void PrintIndivStepScale() const;
256 
258  virtual void ProposeStep();
260  void Randomize() _noexcept_;
262  void CorrelateSteps() _noexcept_;
266 
268  void AcceptStep() _noexcept_;
269 
271  void SetFixAllParameters();
274  void SetFixParameter(const int i);
277  void SetFixParameter(const std::string& name);
278 
280  void SetFreeAllParameters();
283  void SetFreeParameter(const int i);
286  void SetFreeParameter(const std::string& name);
287 
289  void ToggleFixAllParameters();
292  void ToggleFixParameter(const int i);
295  void ToggleFixParameter(const std::string& name);
298  bool IsParameterFixed(const int i) const {
299  if (_fError[i] < 0) { return true; }
300  else { return false; }
301  }
304  bool IsParameterFixed(const std::string& name) const;
305 
311  void ConstructPCA(const double eigen_threshold, int FirstPCAdpar, int LastPCAdpar);
312 
314  inline bool IsPCA() const { return pca; }
315 
317  YAML::Node GetConfig() const { return _fYAMLDoc; }
318 
321  if (!use_adaptive) {
322  MACH3LOG_ERROR("Am not running in Adaptive mode");
323  throw MaCh3Exception(__FILE__ , __LINE__ );
324  }
325  return AdaptiveHandler.get();
326  }
327 
329  void SetTune(const std::string& TuneName);
330 
332  inline PCAHandler* GetPCAHandler() const {
333  if (!pca) {
334  MACH3LOG_ERROR("Am not running in PCA mode");
335  throw MaCh3Exception(__FILE__ , __LINE__ );
336  }
337  return PCAObj.get();
338  }
339 
353  void MatchMaCh3OutputBranches(TTree *PosteriorFile,
354  std::vector<double>& BranchValues,
355  std::vector<std::string>& BranchNames,
356  const std::vector<std::string>& FancyNames = {});
357 
358 protected:
360  void Init(const std::string& name, const std::string& file);
363  void Init(const std::vector<std::string>& YAMLFile);
366  void ReserveMemory(const int size);
367 
370  void MakePosDef(TMatrixDSym *cov = nullptr);
371 
373  void MakeClosestPosDef(TMatrixDSym *cov);
374 
379  void SetThrowMatrixFromFile(const std::string& matrix_file_name, const std::string& matrix_name, const std::string& means_name);
380 
384  bool AppliesToSample(const int SystIndex, const std::string& SampleName) const;
385 
389  void FlipParameterValue(const int index, const double FlipPoint);
390 
393  void CircularParBounds(const int i, const double LowBound, const double UpBound);
394 
396  void EnableSpecialProposal(const YAML::Node& param, const int Index);
397 
400  void SpecialStepProposal();
401 
404 
406  const std::string inputFile;
407 
409  std::string matrixName;
411  TMatrixDSym *covMatrix;
413  TMatrixDSym *invCovMatrix;
415  std::vector<std::vector<double>> InvertCovMatrix;
416 
418  std::vector<std::unique_ptr<TRandom3>> random_number;
419 
421  double* randParams;
423  double* corr_throw;
426 
429 
431  std::vector<std::string> _fNames;
433  std::vector<std::string> _fFancyNames;
435  YAML::Node _fYAMLDoc;
437  int _fNumPar;
439  std::vector<double> _fPreFitValue;
441  std::vector<double> _fCurrVal;
443  std::vector<double> _fPropVal;
445  std::vector<double> _fError;
447  std::vector<double> _fLowBound;
449  std::vector<double> _fUpBound;
451  std::vector<double> _fIndivStepScale;
453  std::vector<bool> _fFlatPrior;
455  std::vector<std::vector<std::string>> _fSampleNames;
456 
458  std::vector<double> _fIndivStepScaleInitial;
459 
462 
464  std::vector<bool> param_skip_adapt_flags;
465 
467  TMatrixDSym* throwMatrix;
470 
472  bool pca;
475 
477  std::unique_ptr<PCAHandler> PCAObj;
479  std::unique_ptr<adaptive_mcmc::AdaptiveMCMCHandler> AdaptiveHandler;
481  std::unique_ptr<ParameterTunes> Tunes;
482 
484  std::vector<int> FlipParameterIndex;
486  std::vector<double> FlipParameterPoint;
488  std::vector<int> CircularBoundsIndex;
490  std::vector<std::pair<double,double>> CircularBoundsValues;
491 };
#define _noexcept_
KS: noexcept can help with performance but is terrible for debugging, this is meant to help easy way ...
Definition: Core.h:96
#define MACH3LOG_ERROR
Definition: MaCh3Logger.h:37
std::vector< TString > BranchNames
Definition: RHat.cpp:54
Custom exception class used throughout MaCh3.
Class responsible for handling Principal Component Analysis (PCA) of covariance matrix.
Definition: PCAHandler.h:58
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.
int GetNumParams() const
Get total number of parameters.
std::vector< bool > param_skip_adapt_flags
Flags telling if parameter should be skipped during adaption.
void ToggleFixParameter(const int i)
Toggle fixing parameter at prior values.
bool use_adaptive
Are we using AMCMC?
void SetFixAllParameters()
Set all parameters to be fixed at prior values.
virtual void ProposeStep()
Generate a new proposed state.
double * randParams
Random number taken from gaussian around prior error used for corr_throw.
void SetName(const std::string &name)
Set matrix name.
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.
std::string GetInputFile() const
Get name of input file.
void SetStepScale(const double scale, const bool verbose=true)
Set global step scale for covariance object.
virtual ~ParameterHandlerBase()
Destructor.
void SetRandomThrow(const int i, const double rand)
Set random value useful for debugging/CI.
TMatrixDSym * invCovMatrix
The inverse covariance matrix.
void SetCovMatrix(TMatrixDSym *cov)
Set covariance matrix.
const double * RetPointer(const int iParam)
DB Pointer return to param position.
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.
void SetFlatPrior(const int i, const bool eL)
Set if parameter should have flat prior or not.
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
void SetParProp(const int i, const double val)
Set proposed parameter value.
bool AppliesToSample(const int SystIndex, const std::string &SampleName) const
Check if parameter is affecting given sample name.
void AcceptStep() _noexcept_
Accepted this step.
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
std::unique_ptr< ParameterTunes > Tunes
Struct containing information about adaption.
void InitialiseAdaption(const YAML::Node &adapt_manager)
Initialise adaptive MCMC.
int GetNParameters() const
Get number of params which will be different depending if using Eigen decomposition or not.
TMatrixDSym * throwMatrix
Matrix which we use for step proposal before Cholesky decomposition (not actually used for step propo...
double _fGlobalStepScale
Global step scale applied to all params in this class.
int GetParIndex(const std::string &name) const
Get index based on name.
void ReserveMemory(const int size)
Initialise vectors with parameters information.
std::vector< std::string > _fFancyNames
Fancy name for example rather than param_0 it is MAQE, useful for human reading.
bool doSpecialStepProposal
Check if any of special step proposal were enabled.
TMatrixDSym * GetInvCovMatrix() const
Return inverted covariance matrix.
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...
double _fGlobalStepScaleInitial
Backup of _fGlobalStepScale for parameters which are skipped during adaption.
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.
TMatrixDSym * GetCovMatrix() const
Return covariance matrix.
double GetThrowMatrix(const int i, const int j) const
Get matrix used for step proposal.
std::vector< double > FlipParameterPoint
Central points around which parameters are flipped.
void UpdateAdaptiveCovariance()
Method to update adaptive MCMC .
void SetThrowMatrix(TMatrixDSym *cov)
Use new throw matrix, used in adaptive MCMC.
std::string GetParFancyName(const int i) const
Get fancy name of the Parameter.
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...
PCAHandler * GetPCAHandler() const
Get pointer for PCAHandler.
void SetFreeAllParameters()
Set all parameters to be treated as free.
void SetFixParameter(const int i)
Set parameter to be fixed at prior value.
void PrintNominal() const
Print prior value for every parameter.
void SetParameters(const std::vector< double > &pars={})
Set parameter values using vector, it has to have same size as covariance class.
const std::vector< double > & GetParCurrVec() const
Get vector of current parameter values.
YAML::Node _fYAMLDoc
Stores config describing systematics.
void SetPrintLength(const unsigned int PriLen)
KS: In case someone really want to change this.
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...
double GetRandomThrow(const int i) const
Get random value useful for debugging/CI.
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.
void UpdateThrowMatrix(TMatrixDSym *cov)
Replaces old throw matrix with new one.
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 SetBranches(TTree &tree, const bool SaveProposal=false)
set branches for output file
void SaveUpdatedMatrixConfig()
KS: After step scale, prefit etc. value were modified save this modified config.
double GetLowerBound(const int i) const
Get lower parameter bound in which it is physically valid.
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 > _fIndivStepScaleInitial
Backup of _fIndivStepScale for parameters which are skipped during adaption.
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.
void SetParName(const int i, const std::string &name)
change parameter name
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.
double GetGlobalStepScale() const
Get global step scale for covariance object.
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 param_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 GetFlatPrior(const int i) const
Get if param has flat prior or not.
bool IsPCA() const
is PCA, can use to query e.g. LLH scans
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.
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 SetFreeParameter(const int i)
Set parameter to be treated as free.
std::vector< double > _fUpBound
Upper physical bound, parameter will not be able to go beyond it.
void SetSingleParameter(const int parNo, const double parVal)
Set value of single param to a given value.
double GetParInit(const int i) const
Get prior parameter 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 SetIndivStepScaleForSkippedAdaptParams()
Set individual step scale for parameters which are skipped during adaption to initial values.
const std::string inputFile
The input root file we read in.
double GetDiagonalError(const int i) const
Get diagonal error for ith parameter.
void SetSubThrowMatrix(int first_index, int last_index, TMatrixDSym const &subcov)
void SetTune(const std::string &TuneName)
KS: Set proposed parameter values vector to be base on tune values, for example set proposed values t...
YAML::Node GetConfig() const
Getter to return a copy of the YAML node.
double GetParProp(const int i) const
Get proposed parameter value.
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.
double GetCorrThrows(const int i) const
Return correlated throws.
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.
double GetParCurr(const int i) const
Get current parameter value.
int CheckBounds() const _noexcept_
Check if parameters were proposed outside physical boundary.
void ToggleFixAllParameters()
Toggle fixing parameters at prior values.
TMatrixDSym * GetThrowMatrix() const
Get matrix used for step proposal.
adaptive_mcmc::AdaptiveMCMCHandler * GetAdaptiveHandler() const
Get pointer for AdaptiveHandler.
Contains information about adaptive covariance matrix .