![]() |
MaCh3
2.4.2
Reference Guide
|
Class responsible for handling Principal Component Analysis (PCA) of covariance matrix. More...
#include <Parameters/PCAHandler.h>
Public Member Functions | |
| PCAHandler () | |
| Constructor. More... | |
| virtual | ~PCAHandler () |
| Destructor. More... | |
| void | Print () const |
| KS: Print info about PCA parameters. More... | |
| int | GetNumberPCAedParameters () const |
| Retrieve number of parameters in PCA base. More... | |
| 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 each step. More... | |
| 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. More... | |
| void | TransferToPCA () |
| Transfer param values from normal base to PCA base. More... | |
| void | TransferToParam () |
| Transfer param values from PCA base to normal base. More... | |
| 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 can break Detailed Balance;. More... | |
| void | AcceptStep () _noexcept_ |
| Accepted this step. More... | |
| 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. More... | |
| void | SetInitialParameters () |
| KS: Transfer the starting parameters to the PCA basis, you don't want to start with zero.. More... | |
| void | SetBranches (TTree &tree, const bool SaveProposal, const std::vector< std::string > &Names) |
| set branches for output file More... | |
| void | ToggleFixAllParameters (const std::vector< std::string > &Names) |
| fix parameters at prior values More... | |
| void | ToggleFixParameter (const int i, const std::vector< std::string > &Names) |
| fix parameters at prior values More... | |
| void | SetParametersPCA (const std::vector< double > &pars) |
| Set values for PCA parameters in PCA base. More... | |
| bool | IsParameterFixedPCA (const int i) const |
| Is parameter fixed in PCA base or not. More... | |
| const TMatrixD | GetEigenVectors () const |
| Get eigen vectors of covariance matrix, only works with PCA. More... | |
| void | SetParPropPCA (const int i, const double value) |
| Set proposed value for parameter in PCA base. More... | |
| void | SetParCurrPCA (const int i, const double value) |
| Set current value for parameter in PCA base. More... | |
| double | GetParPropPCA (const int i) const |
| Get current parameter value using PCA. More... | |
| double | GetPreFitValuePCA (const int i) const |
| Get current parameter value using PCA. More... | |
| double | GetParCurrPCA (const int i) const |
| Get current parameter value using PCA. More... | |
| const TMatrixD | GetTransferMatrix () const |
| Get transfer matrix allowing to go from PCA base to normal base. More... | |
| const TVectorD | GetEigenValues () const |
| Get eigen values for all parameters, if you want for decomposed only parameters use GetEigenValuesMaster. More... | |
| const std::vector< double > | GetEigenValuesMaster () const |
| Get eigen value of only decomposed parameters, if you want for all parameters use GetEigenValues. More... | |
| bool | IsParameterDecomposed (const int i) const |
| Check if parameter in PCA base is decomposed or not. More... | |
Private Member Functions | |
| void | SanitisePCA (TMatrixDSym *CovMatrix) |
| KS: Make sure decomposed matrix isn't correlated with undecomposed. More... | |
Private Attributes | |
| std::vector< double > | _fPreFitValuePCA |
| Prefit value for PCA params. More... | |
| TVectorD | _fParPropPCA |
| CW: Current parameter value in PCA base. More... | |
| TVectorD | _fParCurrPCA |
| CW: Proposed parameter value in PCA base. More... | |
| std::vector< double > | _fErrorPCA |
| Tells if parameter is fixed in PCA base or not. More... | |
| std::vector< int > | isDecomposedPCA |
| If param is decomposed this will return -1, if not this will return enumerator to param in normal base. This way we can map stuff like step scale etc between normal base and undecomposed param in eigen base. More... | |
| int | NumParPCA |
| Number of parameters in PCA base. More... | |
| TMatrixD | TransferMat |
| Matrix used to converting from PCA base to normal base. More... | |
| TMatrixD | TransferMatT |
| Matrix used to converting from normal base to PCA base. More... | |
| TVectorD | eigen_values |
| Eigen value only of particles which are being decomposed. More... | |
| TMatrixD | eigen_vectors |
| Eigen vectors only of params which are being decomposed. More... | |
| std::vector< double > | eigen_values_master |
| Eigen values which have dimension equal to _fNumParPCA, and can be used in CorrelateSteps. More... | |
| int | nKeptPCApars |
| Total number that remained after applying PCA Threshold. More... | |
| int | FirstPCAdpar |
| Index of the first param that is being decomposed. More... | |
| int | LastPCAdpar |
| Index of the last param that is being decomposed. More... | |
| double | eigen_threshold |
| CW: Threshold based on which we remove parameters in eigen base. More... | |
| std::vector< double > * | _pCurrVal |
| Pointer to current value of the parameter. More... | |
| std::vector< double > * | _pPropVal |
| Pointer to proposed value of the parameter. More... | |
Class responsible for handling Principal Component Analysis (PCA) of covariance matrix.
This class applies Principal Component Analysis (PCA) to covariance matrices using eigen-decomposition. The main benefit is that PCA rotates the parameter space into a basis where parameters are uncorrelated. In practice, this helps MCMC fits move more efficiently, since correlated directions in the original space often slow down convergence.
PCA makes it possible to drop directions in parameter space with very small eigenvalues. These correspond to combinations of parameters that are poorly constrained or dominated by statistical noise. By applying a threshold (e.g. discarding eigenvalues smaller than 1e-4 of the largest one), the effective number of parameters can be reduced.
This is particularly useful in large MCMC fits, where many parameters may contribute little information. Removing them reduces dimensionality while keeping the dominant structure of the parameter space intact, often leading to faster and more stable fits.
PCA does not need to be applied to the full covariance matrix. This class supports eigen-decomposition of only a submatrix. This is useful when only a subset of parameters is strongly correlated and benefits from PCA, while the rest remain in the original parameter basis.
Definition at line 58 of file PCAHandler.h.
| PCAHandler::PCAHandler | ( | ) |
Constructor.
Definition at line 5 of file PCAHandler.cpp.
|
virtual |
Destructor.
Definition at line 12 of file PCAHandler.cpp.
| void PCAHandler::AcceptStep | ( | ) |
Accepted this step.
Definition at line 179 of file PCAHandler.cpp.
| void PCAHandler::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.
| CovMatrix | Symmetric covariance matrix used for eigen decomposition. |
| firstPCAd | Index of the first PCA component to include. |
| lastPCAd | Index of the last PCA component to include. |
| eigen_thresh | Threshold for eigenvalues below which parameters are discarded. |
| _fNumPar | Total number of parameters in the original (non-PCA) basis. |
Definition at line 25 of file PCAHandler.cpp.
| void PCAHandler::CorrelateSteps | ( | const std::vector< double > & | IndivStepScale, |
| const double | GlobalStepScale, | ||
| const double *_restrict_ | randParams, | ||
| const double *_restrict_ | corr_throw | ||
| ) |
Use Cholesky throw matrix for better step proposal.
Definition at line 194 of file PCAHandler.cpp.
|
inline |
Get eigen values for all parameters, if you want for decomposed only parameters use GetEigenValuesMaster.
Definition at line 200 of file PCAHandler.h.
|
inline |
Get eigen value of only decomposed parameters, if you want for all parameters use GetEigenValues.
Definition at line 205 of file PCAHandler.h.
|
inline |
Get eigen vectors of covariance matrix, only works with PCA.
Definition at line 148 of file PCAHandler.h.
|
inline |
|
inline |
Get current parameter value using PCA.
| i | Parameter index |
Definition at line 188 of file PCAHandler.h.
|
inline |
Get current parameter value using PCA.
| i | Parameter index |
Definition at line 174 of file PCAHandler.h.
|
inline |
Get current parameter value using PCA.
| i | Parameter index |
Definition at line 181 of file PCAHandler.h.
|
inline |
Get transfer matrix allowing to go from PCA base to normal base.
Definition at line 194 of file PCAHandler.h.
|
inline |
Check if parameter in PCA base is decomposed or not.
| i | Parameter index |
Definition at line 212 of file PCAHandler.h.
|
inline |
Is parameter fixed in PCA base or not.
| i | Parameter index |
Definition at line 141 of file PCAHandler.h.
| void PCAHandler::Print | ( | ) | const |
KS: Print info about PCA parameters.
Definition at line 270 of file PCAHandler.cpp.
|
private |
KS: Make sure decomposed matrix isn't correlated with undecomposed.
Definition at line 150 of file PCAHandler.cpp.
| void PCAHandler::SetBranches | ( | TTree & | tree, |
| const bool | SaveProposal, | ||
| const std::vector< std::string > & | Names | ||
| ) |
set branches for output file
Definition at line 279 of file PCAHandler.cpp.
| void PCAHandler::SetInitialParameters | ( | ) |
KS: Transfer the starting parameters to the PCA basis, you don't want to start with zero..
| IndStepScale | Per parameter step scale, we set so each PCA param uses same step scale [eigen value takes care of size] |
Definition at line 245 of file PCAHandler.cpp.
|
inline |
Set values for PCA parameters in PCA base.
| pars | vector with new values of PCA params |
Definition at line 125 of file PCAHandler.h.
|
inline |
Set current value for parameter in PCA base.
| i | Parameter index |
| value | new value |
Definition at line 165 of file PCAHandler.h.
|
inline |
Set proposed value for parameter in PCA base.
| i | Parameter index |
| value | new value |
Definition at line 156 of file PCAHandler.h.
| void PCAHandler::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 each step.
| fCurr_Val | pointer to current position of parameter |
| fProp_Val | pointer to proposed position of parameter |
Definition at line 17 of file PCAHandler.cpp.
| void PCAHandler::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 can break Detailed Balance;.
Definition at line 318 of file PCAHandler.cpp.
| void PCAHandler::ToggleFixAllParameters | ( | const std::vector< std::string > & | Names | ) |
fix parameters at prior values
Definition at line 294 of file PCAHandler.cpp.
| void PCAHandler::ToggleFixParameter | ( | const int | i, |
| const std::vector< std::string > & | Names | ||
| ) |
fix parameters at prior values
Definition at line 300 of file PCAHandler.cpp.
| void PCAHandler::TransferToParam | ( | ) |
| void PCAHandler::TransferToPCA | ( | ) |
|
private |
Tells if parameter is fixed in PCA base or not.
Definition at line 233 of file PCAHandler.h.
|
private |
CW: Proposed parameter value in PCA base.
Definition at line 231 of file PCAHandler.h.
|
private |
CW: Current parameter value in PCA base.
Definition at line 229 of file PCAHandler.h.
|
private |
Prefit value for PCA params.
Definition at line 227 of file PCAHandler.h.
|
private |
Pointer to current value of the parameter.
Definition at line 260 of file PCAHandler.h.
|
private |
Pointer to proposed value of the parameter.
Definition at line 262 of file PCAHandler.h.
|
private |
CW: Threshold based on which we remove parameters in eigen base.
Definition at line 257 of file PCAHandler.h.
|
private |
Eigen value only of particles which are being decomposed.
Definition at line 244 of file PCAHandler.h.
|
private |
Eigen values which have dimension equal to _fNumParPCA, and can be used in CorrelateSteps.
Definition at line 248 of file PCAHandler.h.
|
private |
Eigen vectors only of params which are being decomposed.
Definition at line 246 of file PCAHandler.h.
|
private |
Index of the first param that is being decomposed.
Definition at line 253 of file PCAHandler.h.
|
private |
If param is decomposed this will return -1, if not this will return enumerator to param in normal base. This way we can map stuff like step scale etc between normal base and undecomposed param in eigen base.
Definition at line 235 of file PCAHandler.h.
|
private |
Index of the last param that is being decomposed.
Definition at line 255 of file PCAHandler.h.
|
private |
Total number that remained after applying PCA Threshold.
Definition at line 251 of file PCAHandler.h.
|
private |
Number of parameters in PCA base.
Definition at line 237 of file PCAHandler.h.
|
private |
Matrix used to converting from PCA base to normal base.
Definition at line 240 of file PCAHandler.h.
|
private |
Matrix used to converting from normal base to PCA base.
Definition at line 242 of file PCAHandler.h.