MaCh3  2.5.1
Reference Guide
Public Member Functions | Private Member Functions | Private Attributes | List of all members
PCAHandler Class Reference

Class responsible for handling Principal Component Analysis (PCA) of covariance matrix. More...

#include <Parameters/PCAHandler.h>

Collaboration diagram for PCAHandler:
[legend]

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< M3::float_t > *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< M3::float_t > * _pPropVal
 Pointer to proposed value of the parameter. More...
 

Detailed Description

Class responsible for handling Principal Component Analysis (PCA) of covariance matrix.

Author
Clarence Wret

Introduction

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.

Dimensionality Reduction

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.

Partial Decomposition

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 38 of file PCAHandler.h.

Constructor & Destructor Documentation

◆ PCAHandler()

PCAHandler::PCAHandler ( )

Constructor.

Definition at line 4 of file PCAHandler.cpp.

4  {
5 // ********************************************
6  _pCurrVal = nullptr;
7  _pPropVal = nullptr;
8 }
std::vector< M3::float_t > * _pPropVal
Pointer to proposed value of the parameter.
Definition: PCAHandler.h:237
std::vector< double > * _pCurrVal
Pointer to current value of the parameter.
Definition: PCAHandler.h:235

◆ ~PCAHandler()

PCAHandler::~PCAHandler ( )
virtual

Destructor.

Definition at line 11 of file PCAHandler.cpp.

11  {
12 // ********************************************
13 }

Member Function Documentation

◆ AcceptStep()

void PCAHandler::AcceptStep ( )

Accepted this step.

Definition at line 180 of file PCAHandler.cpp.

180  {
181 // ********************************************
182  // Update the book-keeping for the output
183  #ifdef MULTITHREAD
184  #pragma omp parallel for
185  #endif
186  for (int i = 0; i < NumParPCA; ++i) {
187  _fParCurrPCA(i) = _fParPropPCA(i);
188  }
189  // Then update the parameter basis
190  TransferToParam();
191 }
TVectorD _fParPropPCA
CW: Current parameter value in PCA base.
Definition: PCAHandler.h:204
int NumParPCA
Number of parameters in PCA base.
Definition: PCAHandler.h:212
void TransferToParam()
Transfer param values from PCA base to normal base.
Definition: PCAHandler.cpp:256
TVectorD _fParCurrPCA
CW: Proposed parameter value in PCA base.
Definition: PCAHandler.h:206

◆ ConstructPCA()

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.

Parameters
CovMatrixSymmetric covariance matrix used for eigen decomposition.
firstPCAdIndex of the first PCA component to include.
lastPCAdIndex of the last PCA component to include.
eigen_threshThreshold for eigenvalues below which parameters are discarded.
_fNumParTotal number of parameters in the original (non-PCA) basis.

Definition at line 24 of file PCAHandler.cpp.

25  {
26 // ********************************************
27  FirstPCAdpar = firstPCAd;
28  LastPCAdpar = lastPCAd;
29  eigen_threshold = eigen_thresh;
30 
31  // Check that covariance matrix exists
32  if (CovMatrix == nullptr) {
33  MACH3LOG_ERROR("Covariance matrix for has not yet been set");
34  MACH3LOG_ERROR("Can not construct PCA until it is set");
35  throw MaCh3Exception(__FILE__ , __LINE__ );
36  }
37 
38  if(FirstPCAdpar > CovMatrix->GetNrows()-1 || LastPCAdpar>CovMatrix->GetNrows()-1) {
39  MACH3LOG_ERROR("FirstPCAdpar and LastPCAdpar are higher than the number of parameters");
40  MACH3LOG_ERROR("first: {} last: {}, params: {}", FirstPCAdpar, LastPCAdpar, CovMatrix->GetNrows()-1);
41  throw MaCh3Exception(__FILE__ , __LINE__ );
42  }
43  if(FirstPCAdpar < 0 || LastPCAdpar < 0){
44  MACH3LOG_ERROR("FirstPCAdpar and LastPCAdpar are less than 0 but not default -999");
45  MACH3LOG_ERROR("first: {} last: {}", FirstPCAdpar, LastPCAdpar);
46  throw MaCh3Exception(__FILE__ , __LINE__ );
47  }
48  MACH3LOG_INFO("PCAing parameters {} through {} inclusive", FirstPCAdpar, LastPCAdpar);
49  int NumUnPCAdPars = CovMatrix->GetNrows()-(LastPCAdpar-FirstPCAdpar+1);
50 
51  // KS: Make sure we are not doing anything silly with PCA
52  SanitisePCA(CovMatrix);
53 
54  TMatrixDSym submat(CovMatrix->GetSub(FirstPCAdpar,LastPCAdpar,FirstPCAdpar,LastPCAdpar));
55 
56  //CW: Calculate how many eigen values this threshold corresponds to
57  TMatrixDSymEigen eigen(submat);
58  eigen_values.ResizeTo(eigen.GetEigenValues());
59  eigen_vectors.ResizeTo(eigen.GetEigenVectors());
60  eigen_values = eigen.GetEigenValues();
61  eigen_vectors = eigen.GetEigenVectors();
62  double sum = 0;
63  // Loop over eigen values and sum them up
64  for (int i = 0; i < eigen_values.GetNrows(); ++i) {
65  sum += eigen_values(i);
66  }
67  nKeptPCApars = eigen_values.GetNrows();
68  //CW: Now go through again and see how many eigen values correspond to threshold
69  for (int i = 0; i < eigen_values.GetNrows(); ++i) {
70  // Get the relative size of the eigen value
71  double sig = eigen_values(i)/sum;
72  // Check against the threshold
73  if (sig < eigen_threshold) {
74  nKeptPCApars = i;
75  break;
76  }
77  }
78  NumParPCA = NumUnPCAdPars+nKeptPCApars;
79  MACH3LOG_INFO("Threshold of {} on eigen values relative sum of eigen value ({}) generates {} eigen vectors, plus we have {} unpcad pars, for a total of {}", eigen_threshold, sum, nKeptPCApars, NumUnPCAdPars, NumParPCA);
80  //DB Create array of correct size so eigen_values can be used in CorrelateSteps
81  eigen_values_master = std::vector<double>(NumParPCA, 1.0);
83 
84  // Now construct the transfer matrices
85  //These matrices will be as big as number of unPCAd pars plus number of eigenvalues kept
86  TransferMat.ResizeTo(CovMatrix->GetNrows(), NumParPCA);
87  TransferMatT.ResizeTo(CovMatrix->GetNrows(), NumParPCA);
88 
89  // Get a subset of the eigen vector matrix
90  TMatrixD temp(eigen_vectors.GetSub(0, eigen_vectors.GetNrows()-1, 0, nKeptPCApars-1));
91 
92  //Make transfer matrix which is two blocks of identity with a block of the PCA transfer matrix in between
93  TMatrixD temp2;
94  temp2.ResizeTo(CovMatrix->GetNrows(), NumParPCA);
95 
96  //First set the whole thing to 0
97  temp2.Zero();
98 
99  //Set the first identity block for non-PCAed params before PCA block, before PCA XRows == YRows
100  for(int iRow = 0; iRow < FirstPCAdpar; iRow++) {
101  temp2(iRow, iRow) = 1.0;
102  }
103 
104  //Set the transfer matrix block for the PCAd pars
105  temp2.SetSub(FirstPCAdpar,FirstPCAdpar,temp);
106 
107  //Set the second identity block starting after PCA block, remember XRows != YRows.
108  // XRows -> normal base, YRows, PCA base
109  for(int iRow = 0;iRow < (CovMatrix->GetNrows()-1)-LastPCAdpar; iRow++) {
110  const int OrigRow = LastPCAdpar + 1 + iRow;
111  const int PCARow = FirstPCAdpar + nKeptPCApars + iRow;
112  temp2(OrigRow, PCARow) = 1.;
113  }
114 
115  TransferMat = temp2;
116  // Copy the contents
118  // And then transpose
119  TransferMatT.T();
120 
121  // Make the PCA parameter arrays
122  _fParCurrPCA.ResizeTo(NumParPCA);
123  _fParPropPCA.ResizeTo(NumParPCA);
124  _fPreFitValuePCA.resize(NumParPCA);
125 
126  //KS: make easy map so we could easily find un-decomposed parameters
127  _fErrorPCA.assign(NumParPCA, 1);
128  isDecomposedPCA.assign(NumParPCA, -1);
129  // First non PCA-ed block, since this is before PCA-ed block we don't need any mapping
130  for (int i = 0; i < FirstPCAdpar; ++i) {
131  isDecomposedPCA[i] = i;
132  }
133 
134  // Second non-PCA-ed block, keep in mind this is in PCA base so we cannot use LastPCAdpar
135  // we must shift them back into the original parameter index range.
136  const int shift = _fNumPar - NumParPCA;
137  for (int i = FirstPCAdpar + nKeptPCApars; i < NumParPCA; ++i) {
138  isDecomposedPCA[i] = i + shift;
139  }
140 
141  #ifdef MACH3_DEBUG
142  //KS: Let's dump all useful matrices to properly validate PCA
144  CovMatrix->GetNrows(), FirstPCAdpar, LastPCAdpar,
146  #endif
147 }
#define MACH3LOG_ERROR
Definition: MaCh3Logger.h:37
#define MACH3LOG_INFO
Definition: MaCh3Logger.h:35
Custom exception class used throughout MaCh3.
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:210
TMatrixD TransferMat
Matrix used to converting from PCA base to normal base.
Definition: PCAHandler.h:215
std::vector< double > eigen_values_master
Eigen values which have dimension equal to _fNumParPCA, and can be used in CorrelateSteps.
Definition: PCAHandler.h:223
TMatrixD eigen_vectors
Eigen vectors only of params which are being decomposed.
Definition: PCAHandler.h:221
TMatrixD TransferMatT
Matrix used to converting from normal base to PCA base.
Definition: PCAHandler.h:217
TVectorD eigen_values
Eigen value only of particles which are being decomposed.
Definition: PCAHandler.h:219
int LastPCAdpar
Index of the last param that is being decomposed.
Definition: PCAHandler.h:230
int nKeptPCApars
Total number that remained after applying PCA Threshold.
Definition: PCAHandler.h:226
std::vector< double > _fErrorPCA
Tells if parameter is fixed in PCA base or not.
Definition: PCAHandler.h:208
std::vector< double > _fPreFitValuePCA
Prefit value for PCA params.
Definition: PCAHandler.h:202
void SanitisePCA(TMatrixDSym *CovMatrix)
KS: Make sure decomposed matrix isn't correlated with undecomposed.
Definition: PCAHandler.cpp:151
int FirstPCAdpar
Index of the first param that is being decomposed.
Definition: PCAHandler.h:228
double eigen_threshold
CW: Threshold based on which we remove parameters in eigen base.
Definition: PCAHandler.h:232
void DebugPCA(const double sum, const TMatrixD &temp, const TMatrixD &eigen_vectors, const TVectorD &eigen_values, const TMatrixD &TransferMat, const TMatrixD &TransferMatT, const int NumPar, const int FirstPCAdpar, const int LastPCAdpar, const int nKeptPCApars, const double eigen_threshold)
KS: Let's dump all useful matrices to properly validate PCA.

◆ CorrelateSteps()

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 195 of file PCAHandler.cpp.

198  {
199 // ************************************************
200  // Throw around the current step
201  #ifdef MULTITHREAD
202  #pragma omp parallel for
203  #endif
204  for (int i = 0; i < NumParPCA; ++i)
205  {
206  if (_fErrorPCA[i] > 0.)
207  {
208  double IndStepScale = 1.;
209  //KS: If undecomposed parameter apply individual step scale and Cholesky for better acceptance rate
210  if(isDecomposedPCA[i] >= 0)
211  {
212  IndStepScale *= IndivStepScale[isDecomposedPCA[i]];
213  IndStepScale *= corr_throw[isDecomposedPCA[i]];
214  }
215  //If decomposed apply only random number
216  else
217  {
218  IndStepScale *= randParams[i];
219  //KS: All PCA-ed parameters have the same step scale
220  IndStepScale *= IndivStepScale[FirstPCAdpar];
221  }
222  _fParPropPCA(i) = _fParCurrPCA(i)+GlobalStepScale*IndStepScale*eigen_values_master[i];
223  }
224  }
225  // Then update the parameter basis
226  TransferToParam();
227 }

◆ GetEigenValues()

const TVectorD PCAHandler::GetEigenValues ( ) const
inline

Get eigen values for all parameters, if you want for decomposed only parameters use GetEigenValuesMaster.

Definition at line 180 of file PCAHandler.h.

180  {
181  return eigen_values;
182  }

◆ GetEigenValuesMaster()

const std::vector<double> PCAHandler::GetEigenValuesMaster ( ) const
inline

Get eigen value of only decomposed parameters, if you want for all parameters use GetEigenValues.

Definition at line 185 of file PCAHandler.h.

185  {
186  return eigen_values_master;
187  }

◆ GetEigenVectors()

const TMatrixD PCAHandler::GetEigenVectors ( ) const
inline

Get eigen vectors of covariance matrix, only works with PCA.

Definition at line 128 of file PCAHandler.h.

128  {
129  return eigen_vectors;
130  }

◆ GetNumberPCAedParameters()

int PCAHandler::GetNumberPCAedParameters ( ) const
inline

Retrieve number of parameters in PCA base.

Definition at line 49 of file PCAHandler.h.

49 { return NumParPCA; }

◆ GetParCurrPCA()

double PCAHandler::GetParCurrPCA ( const int  i) const
inline

Get current parameter value using PCA.

Parameters
iParameter index

Definition at line 168 of file PCAHandler.h.

168  {
169  return _fParCurrPCA(i);
170  }

◆ GetParPropPCA()

double PCAHandler::GetParPropPCA ( const int  i) const
inline

Get current parameter value using PCA.

Parameters
iParameter index

Definition at line 154 of file PCAHandler.h.

154  {
155  return _fParPropPCA(i);
156  }

◆ GetPreFitValuePCA()

double PCAHandler::GetPreFitValuePCA ( const int  i) const
inline

Get current parameter value using PCA.

Parameters
iParameter index

Definition at line 161 of file PCAHandler.h.

161  {
162  return _fPreFitValuePCA[i];
163  }

◆ GetTransferMatrix()

const TMatrixD PCAHandler::GetTransferMatrix ( ) const
inline

Get transfer matrix allowing to go from PCA base to normal base.

Definition at line 174 of file PCAHandler.h.

174  {
175  return TransferMat;
176  }

◆ IsParameterDecomposed()

bool PCAHandler::IsParameterDecomposed ( const int  i) const
inline

Check if parameter in PCA base is decomposed or not.

Parameters
iParameter index

Definition at line 192 of file PCAHandler.h.

192  {
193  if(isDecomposedPCA[i] >= 0) return false;
194  else return true;
195  }

◆ IsParameterFixedPCA()

bool PCAHandler::IsParameterFixedPCA ( const int  i) const
inline

Is parameter fixed in PCA base or not.

Parameters
iParameter index

Definition at line 121 of file PCAHandler.h.

121  {
122  if (_fErrorPCA[i] < 0) { return true; }
123  else { return false; }
124  }

◆ Print()

void PCAHandler::Print ( ) const

KS: Print info about PCA parameters.

Definition at line 274 of file PCAHandler.cpp.

274  {
275 // ********************************************
276  MACH3LOG_INFO("PCA:");
277  for (int i = 0; i < NumParPCA; ++i) {
278  MACH3LOG_INFO("PCA {:<2} Current: {:<10.2f} Proposed: {:<10.2f}", i, _fParCurrPCA(i), _fParPropPCA(i));
279  }
280 }

◆ SanitisePCA()

void PCAHandler::SanitisePCA ( TMatrixDSym *  CovMatrix)
private

KS: Make sure decomposed matrix isn't correlated with undecomposed.

Definition at line 151 of file PCAHandler.cpp.

151  {
152 // ********************************************
153  constexpr double correlation_threshold = 1e-6;
154 
155  bool found_significant_correlation = false;
156 
157  int N = CovMatrix->GetNrows();
158  for (int i = FirstPCAdpar; i <= LastPCAdpar; ++i) {
159  for (int j = 0; j < N; ++j) {
160  // Skip if j is inside the decomposed range (we only want cross-correlations)
161  if (j >= FirstPCAdpar && j <= LastPCAdpar) continue;
162 
163  double corr_val = (*CovMatrix)(i, j);
164  if (std::fabs(corr_val) > correlation_threshold) {
165  found_significant_correlation = true;
166  MACH3LOG_ERROR("Significant correlation detected between decomposed parameter '{}' "
167  "and undecomposed parameter '{}': {:.6e}", i, j, corr_val);
168  }
169  }
170  }
171 
172  if (found_significant_correlation) {
173  MACH3LOG_ERROR("There are correlations between undecomposed and decomposed part of matrices, this will not work");
174  throw MaCh3Exception(__FILE__ , __LINE__);
175  }
176 }

◆ SetBranches()

void PCAHandler::SetBranches ( TTree &  tree,
const bool  SaveProposal,
const std::vector< std::string > &  Names 
)

set branches for output file

Definition at line 283 of file PCAHandler.cpp.

283  {
284 // ********************************************
285  for (int i = 0; i < NumParPCA; ++i) {
286  tree.Branch(Form("%s_PCA", Names[i].c_str()), &_fParCurrPCA.GetMatrixArray()[i], Form("%s_PCA/D", Names[i].c_str()));
287  }
288 
289  if(SaveProposal)
290  {
291  for (int i = 0; i < NumParPCA; ++i) {
292  tree.Branch(Form("%s_PCA_Prop", Names[i].c_str()), &_fParPropPCA.GetMatrixArray()[i], Form("%s_PCA_Prop/D", Names[i].c_str()));
293  }
294  }
295 }

◆ SetInitialParameters()

void PCAHandler::SetInitialParameters ( )

KS: Transfer the starting parameters to the PCA basis, you don't want to start with zero..

Parameters
IndStepScalePer parameter step scale, we set so each PCA param uses same step scale [eigen value takes care of size]

Definition at line 246 of file PCAHandler.cpp.

246  {
247 // ********************************************
248  TransferToPCA();
249  for (int i = 0; i < NumParPCA; ++i) {
251  }
252 }
void TransferToPCA()
Transfer param values from normal base to PCA base.
Definition: PCAHandler.cpp:231

◆ SetParametersPCA()

void PCAHandler::SetParametersPCA ( const std::vector< double > &  pars)
inline

Set values for PCA parameters in PCA base.

Parameters
parsvector with new values of PCA params

Definition at line 105 of file PCAHandler.h.

105  {
106  if (int(pars.size()) != NumParPCA) {
107  MACH3LOG_ERROR("Parameter arrays of incompatible size! Not changing parameters! has size {} but was expecting {}", pars.size(), NumParPCA);
108  throw MaCh3Exception(__FILE__ , __LINE__ );
109  }
110  int parsSize = int(pars.size());
111  for (int i = 0; i < parsSize; i++) {
112  _fParPropPCA(i) = pars[i];
113  }
114  //KS: Transfer to normal base
115  TransferToParam();
116  }

◆ SetParCurrPCA()

void PCAHandler::SetParCurrPCA ( const int  i,
const double  value 
)
inline

Set current value for parameter in PCA base.

Parameters
iParameter index
valuenew value

Definition at line 145 of file PCAHandler.h.

145  {
146  _fParCurrPCA(i) = value;
147  // And then transfer back to the parameter basis
148  TransferToParam();
149  }

◆ SetParPropPCA()

void PCAHandler::SetParPropPCA ( const int  i,
const double  value 
)
inline

Set proposed value for parameter in PCA base.

Parameters
iParameter index
valuenew value

Definition at line 136 of file PCAHandler.h.

136  {
137  _fParPropPCA(i) = value;
138  // And then transfer back to the parameter basis
139  TransferToParam();
140  }

◆ SetupPointers()

void PCAHandler::SetupPointers ( std::vector< double > *  fCurr_Val,
std::vector< M3::float_t > *  fProp_Val 
)

KS: Setup pointers to current and proposed parameter value which we need to convert them to PCA base each step.

Parameters
fCurr_Valpointer to current position of parameter
fProp_Valpointer to proposed position of parameter

Definition at line 16 of file PCAHandler.cpp.

17  {
18 // ********************************************
19  _pCurrVal = fCurr_Val;
20  _pPropVal = fProp_Val;
21 }

◆ ThrowParameters()

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;.

Todo:
KS: We don't check if param is out of bounds. This is more problematic for PCA params.

Definition at line 322 of file PCAHandler.cpp.

329  {
330 // ********************************************
331  #pragma GCC diagnostic push
332  #pragma GCC diagnostic ignored "-Wuseless-cast"
333  //KS: Do not multithread!
334  for (int i = 0; i < NumParPCA; ++i) {
335  // Check if parameter is fixed first: if so don't randomly throw
336  if (IsParameterFixedPCA(i)) continue;
337 
338  if(!IsParameterDecomposed(i))
339  {
340  (*_pPropVal)[i] = static_cast<M3::float_t>(fPreFitValue[i] + corr_throw[i]);
341  int throws = 0;
342  // Try again if we the initial parameter proposal falls outside of the range of the parameter
343  while ((*_pPropVal)[i] > fUpBound[i] || (*_pPropVal)[i] < fLowBound[i]) {
344  randParams[i] = random_number[M3::GetThreadIndex()]->Gaus(0, 1);
345  const double corr_throw_single = M3::MatrixVectorMultiSingle(throwMatrixCholDecomp, randParams, _fNumPar, i);
346  (*_pPropVal)[i] = static_cast<M3::float_t>(fPreFitValue[i] + corr_throw_single);
347  if (throws > 10000)
348  {
349  //KS: Since we are multithreading there is danger that those messages
350  //will be all over the place, small price to pay for faster code
351  MACH3LOG_WARN("Tried {} times to throw parameter {} but failed", throws, i);
352  MACH3LOG_WARN("Setting _fPropVal: {} to {}", (*_pPropVal)[i], fPreFitValue[i]);
353  MACH3LOG_WARN("I live at {}:{}", __FILE__, __LINE__);
354  (*_pPropVal)[i] = static_cast<M3::float_t>(fPreFitValue[i]);
355  }
356  throws++;
357  }
358  (*_pCurrVal)[i] = (*_pPropVal)[i];
359 
360  } else {
361  // KS: We have to multiply by number of parameters in PCA base
364  }
365  } // end of parameter loop
366 
368  for (int i = 0; i < _fNumPar; ++i) {
369  auto low = static_cast<M3::float_t>(fLowBound[i]);
370  auto up = static_cast<M3::float_t>(fUpBound[i]);
371  (*_pPropVal)[i] = std::max(up, std::min((*_pPropVal)[i], low));
372  (*_pCurrVal)[i] = (*_pPropVal)[i];
373  }
374  #pragma GCC diagnostic pop
375 }
#define MACH3LOG_WARN
Definition: MaCh3Logger.h:36
bool IsParameterFixedPCA(const int i) const
Is parameter fixed in PCA base or not.
Definition: PCAHandler.h:121
double GetParPropPCA(const int i) const
Get current parameter value using PCA.
Definition: PCAHandler.h:154
void SetParPropPCA(const int i, const double value)
Set proposed value for parameter in PCA base.
Definition: PCAHandler.h:136
double GetPreFitValuePCA(const int i) const
Get current parameter value using PCA.
Definition: PCAHandler.h:161
bool IsParameterDecomposed(const int i) const
Check if parameter in PCA base is decomposed or not.
Definition: PCAHandler.h:192
void SetParCurrPCA(const int i, const double value)
Set current value for parameter in PCA base.
Definition: PCAHandler.h:145
double float_t
Definition: Core.h:37
int GetThreadIndex()
thread index inside parallel loop
Definition: Monitor.h:86
double MatrixVectorMultiSingle(double **_restrict_ matrix, const double *_restrict_ vector, const int Length, const int i)
KS: Custom function to perform multiplication of matrix and single element which is thread safe.

◆ ToggleFixAllParameters()

void PCAHandler::ToggleFixAllParameters ( const std::vector< std::string > &  Names)

fix parameters at prior values

Definition at line 298 of file PCAHandler.cpp.

298  {
299 // ********************************************
300  for (int i = 0; i < NumParPCA; i++) ToggleFixParameter(i, Names);
301 }
void ToggleFixParameter(const int i, const std::vector< std::string > &Names)
fix parameters at prior values
Definition: PCAHandler.cpp:304

◆ ToggleFixParameter()

void PCAHandler::ToggleFixParameter ( const int  i,
const std::vector< std::string > &  Names 
)

fix parameters at prior values

Definition at line 304 of file PCAHandler.cpp.

304  {
305 // ********************************************
306  int isDecom = -1;
307  for (int im = 0; im < NumParPCA; ++im) {
308  if(isDecomposedPCA[im] == i) {isDecom = im;}
309  }
310  if(isDecom < 0) {
311  MACH3LOG_ERROR("Parameter {} is PCA decomposed can't fix this", Names[i]);
312  //throw MaCh3Exception(__FILE__ , __LINE__ );
313  } else {
314  _fErrorPCA[isDecom] *= -1.0;
315  if(IsParameterFixedPCA(i)) MACH3LOG_INFO("Setting un-decomposed {}(parameter {}/{} in PCA base) to fixed at {}", Names[i], i, isDecom, (*_pCurrVal)[i]);
316  else MACH3LOG_INFO("Setting un-decomposed {}(parameter {}/{} in PCA base) free", Names[i], i, isDecom);
317  }
318 }

◆ TransferToParam()

void PCAHandler::TransferToParam ( )

Transfer param values from PCA base to normal base.

Definition at line 256 of file PCAHandler.cpp.

256  {
257 // ********************************************
258  #pragma GCC diagnostic push
259  #pragma GCC diagnostic ignored "-Wuseless-cast"
260  // Make the temporary vectors
261  TVectorD fParProp_vec = TransferMat*_fParPropPCA;
262  TVectorD fParCurr_vec = TransferMat*_fParCurrPCA;
263  #ifdef MULTITHREAD
264  #pragma omp parallel for
265  #endif
266  for(int i = 0; i < static_cast<int>(_pCurrVal->size()); ++i) {
267  (*_pPropVal)[i] = static_cast<M3::float_t>(fParProp_vec(i));
268  (*_pCurrVal)[i] = fParCurr_vec(i);
269  }
270  #pragma GCC diagnostic pop
271 }

◆ TransferToPCA()

void PCAHandler::TransferToPCA ( )

Transfer param values from normal base to PCA base.

Definition at line 231 of file PCAHandler.cpp.

231  {
232 // ********************************************
233  // Make the temporary vectors
234  TVectorD fParCurr_vec(static_cast<Int_t>(_pCurrVal->size()));
235  TVectorD fParProp_vec(static_cast<Int_t>(_pCurrVal->size()));
236  for(int i = 0; i < static_cast<int>(_pCurrVal->size()); ++i) {
237  fParCurr_vec(i) = (*_pCurrVal)[i];
238  fParProp_vec(i) = (*_pPropVal)[i];
239  }
240 
241  _fParCurrPCA = TransferMatT*fParCurr_vec;
242  _fParPropPCA = TransferMatT*fParProp_vec;
243 }

Member Data Documentation

◆ _fErrorPCA

std::vector<double> PCAHandler::_fErrorPCA
private

Tells if parameter is fixed in PCA base or not.

Definition at line 208 of file PCAHandler.h.

◆ _fParCurrPCA

TVectorD PCAHandler::_fParCurrPCA
private

CW: Proposed parameter value in PCA base.

Definition at line 206 of file PCAHandler.h.

◆ _fParPropPCA

TVectorD PCAHandler::_fParPropPCA
private

CW: Current parameter value in PCA base.

Definition at line 204 of file PCAHandler.h.

◆ _fPreFitValuePCA

std::vector<double> PCAHandler::_fPreFitValuePCA
private

Prefit value for PCA params.

Definition at line 202 of file PCAHandler.h.

◆ _pCurrVal

std::vector<double>* PCAHandler::_pCurrVal
private

Pointer to current value of the parameter.

Definition at line 235 of file PCAHandler.h.

◆ _pPropVal

std::vector<M3::float_t>* PCAHandler::_pPropVal
private

Pointer to proposed value of the parameter.

Definition at line 237 of file PCAHandler.h.

◆ eigen_threshold

double PCAHandler::eigen_threshold
private

CW: Threshold based on which we remove parameters in eigen base.

Definition at line 232 of file PCAHandler.h.

◆ eigen_values

TVectorD PCAHandler::eigen_values
private

Eigen value only of particles which are being decomposed.

Definition at line 219 of file PCAHandler.h.

◆ eigen_values_master

std::vector<double> PCAHandler::eigen_values_master
private

Eigen values which have dimension equal to _fNumParPCA, and can be used in CorrelateSteps.

Definition at line 223 of file PCAHandler.h.

◆ eigen_vectors

TMatrixD PCAHandler::eigen_vectors
private

Eigen vectors only of params which are being decomposed.

Definition at line 221 of file PCAHandler.h.

◆ FirstPCAdpar

int PCAHandler::FirstPCAdpar
private

Index of the first param that is being decomposed.

Definition at line 228 of file PCAHandler.h.

◆ isDecomposedPCA

std::vector<int> PCAHandler::isDecomposedPCA
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 210 of file PCAHandler.h.

◆ LastPCAdpar

int PCAHandler::LastPCAdpar
private

Index of the last param that is being decomposed.

Definition at line 230 of file PCAHandler.h.

◆ nKeptPCApars

int PCAHandler::nKeptPCApars
private

Total number that remained after applying PCA Threshold.

Definition at line 226 of file PCAHandler.h.

◆ NumParPCA

int PCAHandler::NumParPCA
private

Number of parameters in PCA base.

Definition at line 212 of file PCAHandler.h.

◆ TransferMat

TMatrixD PCAHandler::TransferMat
private

Matrix used to converting from PCA base to normal base.

Definition at line 215 of file PCAHandler.h.

◆ TransferMatT

TMatrixD PCAHandler::TransferMatT
private

Matrix used to converting from normal base to PCA base.

Definition at line 217 of file PCAHandler.h.


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