MaCh3  2.2.3
Reference Guide
PCAHandler.h
Go to the documentation of this file.
1 #pragma once
2 
3 // MaCh3 Includes
4 #include "Manager/Manager.h"
6 
7 #ifdef DEBUG
8  #define DEBUG_PCA 1
9 #endif
10 
11 #ifdef DEBUG_PCA
12 //KS: When debugging we produce some fancy plots, but we don't need it during normal work flow
13 #include "TCanvas.h"
14 #include "TROOT.h"
15 #include "TStyle.h"
16 #include "TColor.h"
17 #include "TLine.h"
18 #include "TText.h"
19 #include "TLegend.h"
20 
21 #if DEBUG_PCA == 2
22 #include "Eigen/Eigenvalues"
23 #endif
24 
25 #endif
26 
60 class PCAHandler{
61  public:
63  PCAHandler();
64 
66  virtual ~PCAHandler();
67 
69  void Print();
71  int GetNumberPCAedParameters() const { return NumParPCA; }
72 
76  void SetupPointers(std::vector<double>* fCurr_Val,
77  std::vector<double>* fProp_Val);
78 
85  void ConstructPCA(TMatrixDSym * CovMatrix, const int firstPCAd, const int lastPCAd,
86  const double eigen_thresh, const int _fNumPar);
87 
89  void TransferToPCA();
91  void TransferToParam();
92 
94  void ThrowParameters(const std::vector<std::unique_ptr<TRandom3>>& random_number,
95  double** throwMatrixCholDecomp,
96  double* randParams,
97  double* corr_throw,
98  const std::vector<double>& fPreFitValue,
99  const std::vector<double>& fLowBound,
100  const std::vector<double>& fUpBound,
101  int _fNumPar);
102 
104  void AcceptStep() _noexcept_;
106  void CorrelateSteps(const std::vector<double>& IndivStepScale,
107  const double GlobalStepScale,
108  const double* _restrict_ randParams,
109  const double* _restrict_ corr_throw) _noexcept_;
110 
113  void SetInitialParameters(std::vector<double>& IndStepScale);
115  void ThrowParProp(const double mag, const double* _restrict_ randParams);
117  void ThrowParCurr(const double mag, const double* _restrict_ randParams);
120  void SetBranches(TTree &tree, bool SaveProposal, const std::vector<std::string>& Names);
123  void ToggleFixAllParameters();
126  void ToggleFixParameter(const int i, const std::vector<std::string>& Names);
127 
131  void SetParametersPCA(const std::vector<double> &pars) {
132  if (int(pars.size()) != NumParPCA) {
133  MACH3LOG_ERROR("Parameter arrays of incompatible size! Not changing parameters! has size {} but was expecting {}", pars.size(), NumParPCA);
134  throw MaCh3Exception(__FILE__ , __LINE__ );
135  }
136  int parsSize = int(pars.size());
137  for (int i = 0; i < parsSize; i++) {
138  _fParPropPCA(i) = pars[i];
139  }
140  //KS: Transfer to normal base
141  TransferToParam();
142  }
143 
147  bool IsParameterFixedPCA(const int i) const {
148  if (_fErrorPCA[i] < 0) { return true; }
149  else { return false; }
150  }
151 
154  const TMatrixD GetEigenVectors() const {
155  return eigen_vectors;
156  }
157 
162  void SetParPropPCA(const int i, const double value) {
163  _fParPropPCA(i) = value;
164  // And then transfer back to the parameter basis
165  TransferToParam();
166  }
171  void SetParCurrPCA(const int i, const double value) {
172  _fParCurrPCA(i) = value;
173  // And then transfer back to the parameter basis
174  TransferToParam();
175  }
176 
180  double GetParPropPCA(const int i) const {
181  return _fParPropPCA(i);
182  }
183 
187  double GetPreFitValuePCA(const int i) const {
188  return _fPreFitValuePCA[i];
189  }
190 
194  double GetParCurrPCA(const int i) const {
195  return _fParCurrPCA(i);
196  }
197 
200  const TMatrixD GetTransferMatrix() const {
201  return TransferMat;
202  }
203 
206  const TVectorD GetEigenValues() const {
207  return eigen_values;
208  }
211  const std::vector<double> GetEigenValuesMaster() const {
212  return eigen_values_master;
213  }
214 
218  bool IsParameterDecomposed(const int i) const {
219  if(isDecomposedPCA[i] >= 0) return false;
220  else return true;
221  }
222 
223  #ifdef DEBUG_PCA
225  void DebugPCA(const double sum, TMatrixD temp, TMatrixDSym submat, int NumPar);
226  #endif
227 
228  private:
230  void SanitisePCA(TMatrixDSym* CovMatrix);
231 
233  std::vector<double> _fPreFitValuePCA;
235  TVectorD _fParPropPCA;
237  TVectorD _fParCurrPCA;
239  std::vector<double> _fErrorPCA;
241  std::vector<int> isDecomposedPCA;
244 
246  TMatrixD TransferMat;
248  TMatrixD TransferMatT;
250  TVectorD eigen_values;
252  TMatrixD eigen_vectors;
254  std::vector<double> eigen_values_master;
255 
264 
266  std::vector<double>* _pCurrVal;
268  std::vector<double>* _pPropVal;
269 };
270 
#define _noexcept_
KS: noexcept can help with performance but is terrible for debugging, this is meant to help easy way ...
Definition: Core.h:86
#define _restrict_
KS: Using restrict limits the effects of pointer aliasing, aiding optimizations. While reading I foun...
Definition: Core.h:93
#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
std::vector< int > isDecomposedPCA
If param is decomposed this will return -1, if not this will return enumerator to param in normal bas...
Definition: PCAHandler.h:241
void ThrowParameters(const std::vector< std::unique_ptr< TRandom3 >> &random_number, double **throwMatrixCholDecomp, double *randParams, double *corr_throw, const std::vector< double > &fPreFitValue, const std::vector< double > &fLowBound, const std::vector< double > &fUpBound, int _fNumPar)
Throw the parameters according to the covariance matrix. This shouldn't be used in MCMC code ase it c...
Definition: PCAHandler.cpp:352
TMatrixD TransferMat
Matrix used to converting from PCA base to normal base.
Definition: PCAHandler.h:246
std::vector< double > eigen_values_master
Eigen values which have dimension equal to _fNumParPCA, and can be used in CorrelateSteps.
Definition: PCAHandler.h:254
TMatrixD eigen_vectors
Eigen vectors only of params which are being decomposed.
Definition: PCAHandler.h:252
TVectorD _fParPropPCA
CW: Current parameter value in PCA base.
Definition: PCAHandler.h:235
void ThrowParProp(const double mag, const double *_restrict_ randParams)
Throw the proposed parameter by mag sigma.
Definition: PCAHandler.cpp:280
int NumParPCA
Number of parameters in PCA base.
Definition: PCAHandler.h:243
int GetNumberPCAedParameters() const
Retrieve number of parameters in PCA base.
Definition: PCAHandler.h:71
std::vector< double > * _pCurrVal
Pointer to current value of the parameter.
Definition: PCAHandler.h:266
TMatrixD TransferMatT
Matrix used to converting from normal base to PCA base.
Definition: PCAHandler.h:248
void AcceptStep() _noexcept_
Accepted this step.
Definition: PCAHandler.cpp:183
void CorrelateSteps(const std::vector< double > &IndivStepScale, const double GlobalStepScale, const double *_restrict_ randParams, const double *_restrict_ corr_throw) _noexcept_
Use Cholesky throw matrix for better step proposal.
Definition: PCAHandler.cpp:198
void SetInitialParameters(std::vector< double > &IndStepScale)
KS: Transfer the starting parameters to the PCA basis, you don't want to start with zero....
Definition: PCAHandler.cpp:249
void TransferToPCA()
Transfer param values from normal base to PCA base.
Definition: PCAHandler.cpp:234
virtual ~PCAHandler()
Destructor.
Definition: PCAHandler.cpp:12
void SetupPointers(std::vector< double > *fCurr_Val, std::vector< double > *fProp_Val)
KS: Setup pointers to current and proposed parameter value which we need to convert them to PCA base ...
Definition: PCAHandler.cpp:18
std::vector< double > * _pPropVal
Pointer to proposed value of the parameter.
Definition: PCAHandler.h:268
TVectorD eigen_values
Eigen value only of particles which are being decomposed.
Definition: PCAHandler.h:250
void TransferToParam()
Transfer param values from PCA base to normal base.
Definition: PCAHandler.cpp:264
int LastPCAdpar
Index of the last param that is being decomposed.
Definition: PCAHandler.h:261
int nKeptPCApars
Total number that remained after applying PCA Threshold.
Definition: PCAHandler.h:257
std::vector< double > _fErrorPCA
Tells if parameter is fixed in PCA base or not.
Definition: PCAHandler.h:239
PCAHandler()
Constructor.
Definition: PCAHandler.cpp:5
void Print()
KS: Print info about PCA parameters.
Definition: PCAHandler.cpp:303
std::vector< double > _fPreFitValuePCA
Prefit value for PCA params.
Definition: PCAHandler.h:233
void ThrowParCurr(const double mag, const double *_restrict_ randParams)
Helper function to throw the current parameter by mag sigma.
Definition: PCAHandler.cpp:292
TVectorD _fParCurrPCA
CW: Proposed parameter value in PCA base.
Definition: PCAHandler.h:237
void SanitisePCA(TMatrixDSym *CovMatrix)
KS: Make sure decomposed matrix isn't correlated with undecomposed.
Definition: PCAHandler.cpp:154
void ConstructPCA(TMatrixDSym *CovMatrix, const int firstPCAd, const int lastPCAd, const double eigen_thresh, const int _fNumPar)
CW: Calculate eigen values, prepare transition matrices and remove param based on defined threshold.
Definition: PCAHandler.cpp:26
int FirstPCAdpar
Index of the first param that is being decomposed.
Definition: PCAHandler.h:259
double eigen_threshold
CW: Threshold based on which we remove parameters in eigen base.
Definition: PCAHandler.h:263
const TMatrixD GetTransferMatrix() const
Get transfer matrix allowing to go from PCA base to normal base.
Definition: PCAHandler.h:200
const std::vector< double > GetEigenValuesMaster() const
Get eigen value of only decomposed parameters, if you want for all parameters use GetEigenValues.
Definition: PCAHandler.h:211
double GetParCurrPCA(const int i) const
Get current parameter value using PCA.
Definition: PCAHandler.h:194
bool IsParameterFixedPCA(const int i) const
Is parameter fixed in PCA base or not.
Definition: PCAHandler.h:147
double GetParPropPCA(const int i) const
Get current parameter value using PCA.
Definition: PCAHandler.h:180
double GetPreFitValuePCA(const int i) const
Get current parameter value using PCA.
Definition: PCAHandler.h:187
bool IsParameterDecomposed(const int i) const
Check if parameter in PCA base is decomposed or not.
Definition: PCAHandler.h:218
const TVectorD GetEigenValues() const
Get eigen values for all parameters, if you want for decomposed only parameters use GetEigenValuesMas...
Definition: PCAHandler.h:206
const TMatrixD GetEigenVectors() const
Get eigen vectors of covariance matrix, only works with PCA.
Definition: PCAHandler.h:154
void ToggleFixAllParameters()
fix parameters at prior values
Definition: PCAHandler.cpp:327
void SetBranches(TTree &tree, bool SaveProposal, const std::vector< std::string > &Names)
set branches for output file
Definition: PCAHandler.cpp:312
void SetParPropPCA(const int i, const double value)
Set proposed value for parameter in PCA base.
Definition: PCAHandler.h:162
void ToggleFixParameter(const int i, const std::vector< std::string > &Names)
fix parameters at prior values
Definition: PCAHandler.cpp:335
void SetParCurrPCA(const int i, const double value)
Set current value for parameter in PCA base.
Definition: PCAHandler.h:171
void SetParametersPCA(const std::vector< double > &pars)
Set values for PCA parameters in PCA base.
Definition: PCAHandler.h:131