MaCh3  2.5.0
Reference Guide
PCAHandler.cpp
Go to the documentation of this file.
2 
3 #ifdef MACH3_DEBUG
5 //KS: When debugging we produce some fancy plots, but we don't need it during normal work flow
6 #include "TCanvas.h"
7 #include "TROOT.h"
8 #include "TStyle.h"
9 #include "TColor.h"
10 #include "TLine.h"
11 #include "TText.h"
12 #include "TLegend.h"
14 #endif
15 
16 // ********************************************
18 // ********************************************
19  _pCurrVal = nullptr;
20  _pPropVal = nullptr;
21 }
22 
23 // ********************************************
25 // ********************************************
26 }
27 
28 // ********************************************
29 void PCAHandler::SetupPointers(std::vector<double>* fCurr_Val,
30  std::vector<M3::float_t>* fProp_Val) {
31 // ********************************************
32  _pCurrVal = fCurr_Val;
33  _pPropVal = fProp_Val;
34 }
35 
36 // ********************************************
37 void PCAHandler::ConstructPCA(TMatrixDSym* CovMatrix, const int firstPCAd, const int lastPCAd,
38  const double eigen_thresh, const int _fNumPar) {
39 // ********************************************
40  FirstPCAdpar = firstPCAd;
41  LastPCAdpar = lastPCAd;
42  eigen_threshold = eigen_thresh;
43 
44  // Check that covariance matrix exists
45  if (CovMatrix == nullptr) {
46  MACH3LOG_ERROR("Covariance matrix for has not yet been set");
47  MACH3LOG_ERROR("Can not construct PCA until it is set");
48  throw MaCh3Exception(__FILE__ , __LINE__ );
49  }
50 
51  if(FirstPCAdpar > CovMatrix->GetNrows()-1 || LastPCAdpar>CovMatrix->GetNrows()-1) {
52  MACH3LOG_ERROR("FirstPCAdpar and LastPCAdpar are higher than the number of parameters");
53  MACH3LOG_ERROR("first: {} last: {}, params: {}", FirstPCAdpar, LastPCAdpar, CovMatrix->GetNrows()-1);
54  throw MaCh3Exception(__FILE__ , __LINE__ );
55  }
56  if(FirstPCAdpar < 0 || LastPCAdpar < 0){
57  MACH3LOG_ERROR("FirstPCAdpar and LastPCAdpar are less than 0 but not default -999");
58  MACH3LOG_ERROR("first: {} last: {}", FirstPCAdpar, LastPCAdpar);
59  throw MaCh3Exception(__FILE__ , __LINE__ );
60  }
61  MACH3LOG_INFO("PCAing parameters {} through {} inclusive", FirstPCAdpar, LastPCAdpar);
62  int NumUnPCAdPars = CovMatrix->GetNrows()-(LastPCAdpar-FirstPCAdpar+1);
63 
64  // KS: Make sure we are not doing anything silly with PCA
65  SanitisePCA(CovMatrix);
66 
67  TMatrixDSym submat(CovMatrix->GetSub(FirstPCAdpar,LastPCAdpar,FirstPCAdpar,LastPCAdpar));
68 
69  //CW: Calculate how many eigen values this threshold corresponds to
70  TMatrixDSymEigen eigen(submat);
71  eigen_values.ResizeTo(eigen.GetEigenValues());
72  eigen_vectors.ResizeTo(eigen.GetEigenVectors());
73  eigen_values = eigen.GetEigenValues();
74  eigen_vectors = eigen.GetEigenVectors();
75  double sum = 0;
76  // Loop over eigen values and sum them up
77  for (int i = 0; i < eigen_values.GetNrows(); ++i) {
78  sum += eigen_values(i);
79  }
80  nKeptPCApars = eigen_values.GetNrows();
81  //CW: Now go through again and see how many eigen values correspond to threshold
82  for (int i = 0; i < eigen_values.GetNrows(); ++i) {
83  // Get the relative size of the eigen value
84  double sig = eigen_values(i)/sum;
85  // Check against the threshold
86  if (sig < eigen_threshold) {
87  nKeptPCApars = i;
88  break;
89  }
90  }
91  NumParPCA = NumUnPCAdPars+nKeptPCApars;
92  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);
93  //DB Create array of correct size so eigen_values can be used in CorrelateSteps
94  eigen_values_master = std::vector<double>(NumParPCA, 1.0);
96 
97  // Now construct the transfer matrices
98  //These matrices will be as big as number of unPCAd pars plus number of eigenvalues kept
99  TransferMat.ResizeTo(CovMatrix->GetNrows(), NumParPCA);
100  TransferMatT.ResizeTo(CovMatrix->GetNrows(), NumParPCA);
101 
102  // Get a subset of the eigen vector matrix
103  TMatrixD temp(eigen_vectors.GetSub(0, eigen_vectors.GetNrows()-1, 0, nKeptPCApars-1));
104 
105  //Make transfer matrix which is two blocks of identity with a block of the PCA transfer matrix in between
106  TMatrixD temp2;
107  temp2.ResizeTo(CovMatrix->GetNrows(), NumParPCA);
108 
109  //First set the whole thing to 0
110  temp2.Zero();
111 
112  //Set the first identity block for non-PCAed params before PCA block, before PCA XRows == YRows
113  for(int iRow = 0; iRow < FirstPCAdpar; iRow++) {
114  temp2(iRow, iRow) = 1.0;
115  }
116 
117  //Set the transfer matrix block for the PCAd pars
118  temp2.SetSub(FirstPCAdpar,FirstPCAdpar,temp);
119 
120  //Set the second identity block starting after PCA block, remember XRows != YRows.
121  // XRows -> normal base, YRows, PCA base
122  for(int iRow = 0;iRow < (CovMatrix->GetNrows()-1)-LastPCAdpar; iRow++) {
123  const int OrigRow = LastPCAdpar + 1 + iRow;
124  const int PCARow = FirstPCAdpar + nKeptPCApars + iRow;
125  temp2(OrigRow, PCARow) = 1.;
126  }
127 
128  TransferMat = temp2;
129  // Copy the contents
131  // And then transpose
132  TransferMatT.T();
133 
134  // Make the PCA parameter arrays
135  _fParCurrPCA.ResizeTo(NumParPCA);
136  _fParPropPCA.ResizeTo(NumParPCA);
137  _fPreFitValuePCA.resize(NumParPCA);
138 
139  //KS: make easy map so we could easily find un-decomposed parameters
140  _fErrorPCA.assign(NumParPCA, 1);
141  isDecomposedPCA.assign(NumParPCA, -1);
142  // First non PCA-ed block, since this is before PCA-ed block we don't need any mapping
143  for (int i = 0; i < FirstPCAdpar; ++i) {
144  isDecomposedPCA[i] = i;
145  }
146 
147  // Second non-PCA-ed block, keep in mind this is in PCA base so we cannot use LastPCAdpar
148  // we must shift them back into the original parameter index range.
149  const int shift = _fNumPar - NumParPCA;
150  for (int i = FirstPCAdpar + nKeptPCApars; i < NumParPCA; ++i) {
151  isDecomposedPCA[i] = i + shift;
152  }
153 
154  #ifdef MACH3_DEBUG_PCA
155  //KS: Let's dump all useful matrices to properly validate PCA
156  DebugPCA(sum, temp, CovMatrix->GetNrows());
157  #endif
158 }
159 
160 // ********************************************
161 // Make sure decomposed matrix isn't correlated with undecomposed
162 void PCAHandler::SanitisePCA(TMatrixDSym* CovMatrix) {
163 // ********************************************
164  constexpr double correlation_threshold = 1e-6;
165 
166  bool found_significant_correlation = false;
167 
168  int N = CovMatrix->GetNrows();
169  for (int i = FirstPCAdpar; i <= LastPCAdpar; ++i) {
170  for (int j = 0; j < N; ++j) {
171  // Skip if j is inside the decomposed range (we only want cross-correlations)
172  if (j >= FirstPCAdpar && j <= LastPCAdpar) continue;
173 
174  double corr_val = (*CovMatrix)(i, j);
175  if (std::fabs(corr_val) > correlation_threshold) {
176  found_significant_correlation = true;
177  MACH3LOG_ERROR("Significant correlation detected between decomposed parameter '{}' "
178  "and undecomposed parameter '{}': {:.6e}", i, j, corr_val);
179  }
180  }
181  }
182 
183  if (found_significant_correlation) {
184  MACH3LOG_ERROR("There are correlations between undecomposed and decomposed part of matrices, this will not work");
185  throw MaCh3Exception(__FILE__ , __LINE__);
186  }
187 }
188 
189 // ********************************************
190 // Update so that current step becomes the previously proposed step
192 // ********************************************
193  // Update the book-keeping for the output
194  #ifdef MULTITHREAD
195  #pragma omp parallel for
196  #endif
197  for (int i = 0; i < NumParPCA; ++i) {
198  _fParCurrPCA(i) = _fParPropPCA(i);
199  }
200  // Then update the parameter basis
201  TransferToParam();
202 }
203 
204 // ************************************************
205 // Correlate the steps by setting the proposed step of a parameter to its current value + some correlated throw
206 void PCAHandler::CorrelateSteps(const std::vector<double>& IndivStepScale,
207  const double GlobalStepScale,
208  const double* _restrict_ randParams,
209  const double* _restrict_ corr_throw) _noexcept_ {
210 // ************************************************
211  // Throw around the current step
212  #ifdef MULTITHREAD
213  #pragma omp parallel for
214  #endif
215  for (int i = 0; i < NumParPCA; ++i)
216  {
217  if (_fErrorPCA[i] > 0.)
218  {
219  double IndStepScale = 1.;
220  //KS: If undecomposed parameter apply individual step scale and Cholesky for better acceptance rate
221  if(isDecomposedPCA[i] >= 0)
222  {
223  IndStepScale *= IndivStepScale[isDecomposedPCA[i]];
224  IndStepScale *= corr_throw[isDecomposedPCA[i]];
225  }
226  //If decomposed apply only random number
227  else
228  {
229  IndStepScale *= randParams[i];
230  //KS: All PCA-ed parameters have the same step scale
231  IndStepScale *= IndivStepScale[FirstPCAdpar];
232  }
233  _fParPropPCA(i) = _fParCurrPCA(i)+GlobalStepScale*IndStepScale*eigen_values_master[i];
234  }
235  }
236  // Then update the parameter basis
237  TransferToParam();
238 }
239 
240 // ********************************************
241 // Transfer a parameter variation in the parameter basis to the eigen basis
243 // ********************************************
244  // Make the temporary vectors
245  TVectorD fParCurr_vec(static_cast<Int_t>(_pCurrVal->size()));
246  TVectorD fParProp_vec(static_cast<Int_t>(_pCurrVal->size()));
247  for(int i = 0; i < static_cast<int>(_pCurrVal->size()); ++i) {
248  fParCurr_vec(i) = (*_pCurrVal)[i];
249  fParProp_vec(i) = (*_pPropVal)[i];
250  }
251 
252  _fParCurrPCA = TransferMatT*fParCurr_vec;
253  _fParPropPCA = TransferMatT*fParProp_vec;
254 }
255 
256 // ********************************************
258 // ********************************************
259  TransferToPCA();
260  for (int i = 0; i < NumParPCA; ++i) {
262  }
263 }
264 
265 // ********************************************
266 // Transfer a parameter variation in the eigen basis to the parameter basis
268 // ********************************************
269  #pragma GCC diagnostic push
270  #pragma GCC diagnostic ignored "-Wuseless-cast"
271  // Make the temporary vectors
272  TVectorD fParProp_vec = TransferMat*_fParPropPCA;
273  TVectorD fParCurr_vec = TransferMat*_fParCurrPCA;
274  #ifdef MULTITHREAD
275  #pragma omp parallel for
276  #endif
277  for(int i = 0; i < static_cast<int>(_pCurrVal->size()); ++i) {
278  (*_pPropVal)[i] = static_cast<M3::float_t>(fParProp_vec(i));
279  (*_pCurrVal)[i] = fParCurr_vec(i);
280  }
281  #pragma GCC diagnostic pop
282 }
283 
284 // ********************************************
285 void PCAHandler::Print() const {
286 // ********************************************
287  MACH3LOG_INFO("PCA:");
288  for (int i = 0; i < NumParPCA; ++i) {
289  MACH3LOG_INFO("PCA {:<2} Current: {:<10.2f} Proposed: {:<10.2f}", i, _fParCurrPCA(i), _fParPropPCA(i));
290  }
291 }
292 
293 // ********************************************
294 void PCAHandler::SetBranches(TTree &tree, const bool SaveProposal, const std::vector<std::string>& Names) {
295 // ********************************************
296  for (int i = 0; i < NumParPCA; ++i) {
297  tree.Branch(Form("%s_PCA", Names[i].c_str()), &_fParCurrPCA.GetMatrixArray()[i], Form("%s_PCA/D", Names[i].c_str()));
298  }
299 
300  if(SaveProposal)
301  {
302  for (int i = 0; i < NumParPCA; ++i) {
303  tree.Branch(Form("%s_PCA_Prop", Names[i].c_str()), &_fParPropPCA.GetMatrixArray()[i], Form("%s_PCA_Prop/D", Names[i].c_str()));
304  }
305  }
306 }
307 
308 // ********************************************
309 void PCAHandler::ToggleFixAllParameters(const std::vector<std::string>& Names) {
310 // ********************************************
311  for (int i = 0; i < NumParPCA; i++) ToggleFixParameter(i, Names);
312 }
313 
314 // ********************************************
315 void PCAHandler::ToggleFixParameter(const int i, const std::vector<std::string>& Names) {
316 // ********************************************
317  int isDecom = -1;
318  for (int im = 0; im < NumParPCA; ++im) {
319  if(isDecomposedPCA[im] == i) {isDecom = im;}
320  }
321  if(isDecom < 0) {
322  MACH3LOG_ERROR("Parameter {} is PCA decomposed can't fix this", Names[i]);
323  //throw MaCh3Exception(__FILE__ , __LINE__ );
324  } else {
325  _fErrorPCA[isDecom] *= -1.0;
326  if(IsParameterFixedPCA(i)) MACH3LOG_INFO("Setting un-decomposed {}(parameter {}/{} in PCA base) to fixed at {}", Names[i], i, isDecom, (*_pCurrVal)[i]);
327  else MACH3LOG_INFO("Setting un-decomposed {}(parameter {}/{} in PCA base) free", Names[i], i, isDecom);
328  }
329 }
330 
331 
332 // ********************************************
333 void PCAHandler::ThrowParameters(const std::vector<std::unique_ptr<TRandom3>>& random_number,
334  double** throwMatrixCholDecomp,
335  double* randParams,
336  double* corr_throw,
337  const std::vector<double>& fPreFitValue,
338  const std::vector<double>& fLowBound,
339  const std::vector<double>& fUpBound,
340  const int _fNumPar) {
341 // ********************************************
342  #pragma GCC diagnostic push
343  #pragma GCC diagnostic ignored "-Wuseless-cast"
344  //KS: Do not multithread!
345  for (int i = 0; i < NumParPCA; ++i) {
346  // Check if parameter is fixed first: if so don't randomly throw
347  if (IsParameterFixedPCA(i)) continue;
348 
349  if(!IsParameterDecomposed(i))
350  {
351  (*_pPropVal)[i] = static_cast<M3::float_t>(fPreFitValue[i] + corr_throw[i]);
352  int throws = 0;
353  // Try again if we the initial parameter proposal falls outside of the range of the parameter
354  while ((*_pPropVal)[i] > fUpBound[i] || (*_pPropVal)[i] < fLowBound[i]) {
355  randParams[i] = random_number[M3::GetThreadIndex()]->Gaus(0, 1);
356  const double corr_throw_single = M3::MatrixVectorMultiSingle(throwMatrixCholDecomp, randParams, _fNumPar, i);
357  (*_pPropVal)[i] = static_cast<M3::float_t>(fPreFitValue[i] + corr_throw_single);
358  if (throws > 10000)
359  {
360  //KS: Since we are multithreading there is danger that those messages
361  //will be all over the place, small price to pay for faster code
362  MACH3LOG_WARN("Tried {} times to throw parameter {} but failed", throws, i);
363  MACH3LOG_WARN("Setting _fPropVal: {} to {}", (*_pPropVal)[i], fPreFitValue[i]);
364  MACH3LOG_WARN("I live at {}:{}", __FILE__, __LINE__);
365  (*_pPropVal)[i] = static_cast<M3::float_t>(fPreFitValue[i]);
366  }
367  throws++;
368  }
369  (*_pCurrVal)[i] = (*_pPropVal)[i];
370 
371  } else {
372  // KS: We have to multiply by number of parameters in PCA base
375  }
376  } // end of parameter loop
377 
379  for (int i = 0; i < _fNumPar; ++i) {
380  auto low = static_cast<M3::float_t>(fLowBound[i]);
381  auto up = static_cast<M3::float_t>(fUpBound[i]);
382  (*_pPropVal)[i] = std::max(up, std::min((*_pPropVal)[i], low));
383  (*_pCurrVal)[i] = (*_pPropVal)[i];
384  }
385  #pragma GCC diagnostic pop
386 }
387 
388 #ifdef MACH3_DEBUG
389 #pragma GCC diagnostic ignored "-Wfloat-conversion"
390 // ********************************************
391 //KS: Let's dump all useful matrices to properly validate PCA
392 void PCAHandler::DebugPCA(const double sum, TMatrixD temp, int NumPar) {
393 // ********************************************
394  int originalErrorWarning = gErrorIgnoreLevel;
395  gErrorIgnoreLevel = kFatal;
396 
397  for (int i = 0; i < NumParPCA; ++i) {
398  MACH3LOG_DEBUG("Param {} isDecomposedPCA={}", i, isDecomposedPCA[i]);
399  }
400 
401  TFile *PCA_Debug = new TFile("Debug_PCA.root", "RECREATE");
402  PCA_Debug->cd();
403 
404  bool PlotText = true;
405  //KS: If we have more than 200 plot becomes unreadable :(
406  if(NumPar > 200) PlotText = false;
407 
408  auto heigen_values = std::make_unique<TH1D>("eigen_values", "Eigen Values", eigen_values.GetNrows(), 0.0, eigen_values.GetNrows());
409  heigen_values->SetDirectory(nullptr);
410  auto heigen_cumulative = std::make_unique<TH1D>("heigen_cumulative", "heigen_cumulative", eigen_values.GetNrows(), 0.0, eigen_values.GetNrows());
411  heigen_cumulative->SetDirectory(nullptr);
412  auto heigen_frac = std::make_unique<TH1D>("heigen_fractional", "heigen_fractional", eigen_values.GetNrows(), 0.0, eigen_values.GetNrows());
413  heigen_frac->SetDirectory(nullptr);
414  heigen_values->GetXaxis()->SetTitle("Eigen Vector");
415  heigen_values->GetYaxis()->SetTitle("Eigen Value");
416 
417  double Cumulative = 0;
418  for(int i = 0; i < eigen_values.GetNrows(); i++)
419  {
420  heigen_values->SetBinContent(i+1, (eigen_values)(i));
421  heigen_cumulative->SetBinContent(i+1, (eigen_values)(i)/sum + Cumulative);
422  heigen_frac->SetBinContent(i+1, (eigen_values)(i)/sum);
423  Cumulative += (eigen_values)(i)/sum;
424  }
425  heigen_values->Write("heigen_values");
426  eigen_values.Write("eigen_values");
427  heigen_cumulative->Write("heigen_values_cumulative");
428  heigen_frac->Write("heigen_values_frac");
429 
430  TH2D* heigen_vectors = new TH2D(eigen_vectors);
431  heigen_vectors->GetXaxis()->SetTitle("Parameter in Normal Base");
432  heigen_vectors->GetYaxis()->SetTitle("Parameter in Decomposed Base");
433  heigen_vectors->Write("heigen_vectors");
434  eigen_vectors.Write("eigen_vectors");
435 
436  TH2D* SubsetPCA = new TH2D(temp);
437  SubsetPCA->GetXaxis()->SetTitle("Parameter in Normal Base");
438  SubsetPCA->GetYaxis()->SetTitle("Parameter in Decomposed Base");
439 
440  SubsetPCA->Write("hSubsetPCA");
441  temp.Write("SubsetPCA");
442  TH2D* hTransferMat = new TH2D(TransferMat);
443  hTransferMat->GetXaxis()->SetTitle("Parameter in Normal Base");
444  hTransferMat->GetYaxis()->SetTitle("Parameter in Decomposed Base");
445  TH2D* hTransferMatT = new TH2D(TransferMatT);
446 
447  hTransferMatT->GetXaxis()->SetTitle("Parameter in Decomposed Base");
448  hTransferMatT->GetYaxis()->SetTitle("Parameter in Normal Base");
449 
450  hTransferMat->Write("hTransferMat");
451  TransferMat.Write("TransferMat");
452  hTransferMatT->Write("hTransferMatT");
453  TransferMatT.Write("TransferMatT");
454 
455  auto c1 = std::make_unique<TCanvas>("c1", " ", 0, 0, 1024, 1024);
456  c1->SetBottomMargin(0.1);
457  c1->SetTopMargin(0.05);
458  c1->SetRightMargin(0.05);
459  c1->SetLeftMargin(0.12);
460  c1->SetGrid();
461 
462  gStyle->SetPaintTextFormat("4.1f");
463  gStyle->SetOptFit(0);
464  gStyle->SetOptStat(0);
465  // Make pretty correlation colors (red to blue)
466  constexpr int NRGBs = 5;
467  TColor::InitializeColors();
468  Double_t stops[NRGBs] = { 0.00, 0.25, 0.50, 0.75, 1.00 };
469  Double_t red[NRGBs] = { 0.00, 0.25, 1.00, 1.00, 0.50 };
470  Double_t green[NRGBs] = { 0.00, 0.25, 1.00, 0.25, 0.00 };
471  Double_t blue[NRGBs] = { 0.50, 1.00, 1.00, 0.25, 0.00 };
472  TColor::CreateGradientColorTable(5, stops, red, green, blue, 255);
473  gStyle->SetNumberContours(255);
474 
475  double maxz = 0;
476  double minz = 0;
477 
478  c1->Print("Debug_PCA.pdf[");
479  auto EigenLine = std::make_unique<TLine>(nKeptPCApars, 0, nKeptPCApars, heigen_cumulative->GetMaximum());
480  EigenLine->SetLineColor(kPink);
481  EigenLine->SetLineWidth(2);
482  EigenLine->SetLineStyle(kSolid);
483 
484  auto text = std::make_unique<TText>(0.5, 0.5, Form("Threshold = %g", eigen_threshold));
485  text->SetTextFont (43);
486  text->SetTextSize (40);
487 
488  heigen_values->SetLineColor(kRed);
489  heigen_values->SetLineWidth(2);
490  heigen_cumulative->SetLineColor(kGreen);
491  heigen_cumulative->SetLineWidth(2);
492  heigen_frac->SetLineColor(kBlue);
493  heigen_frac->SetLineWidth(2);
494 
495  c1->SetLogy();
496  heigen_values->SetMaximum(heigen_cumulative->GetMaximum()+heigen_cumulative->GetMaximum()*0.4);
497  heigen_values->Draw();
498  heigen_frac->Draw("SAME");
499  heigen_cumulative->Draw("SAME");
500  EigenLine->Draw("Same");
501  text->DrawTextNDC(0.42, 0.84,Form("Threshold = %g", eigen_threshold));
502 
503  auto leg = std::make_unique<TLegend>(0.2, 0.2, 0.6, 0.5);
504  leg->SetTextSize(0.04);
505  leg->AddEntry(heigen_values.get(), "Absolute", "l");
506  leg->AddEntry(heigen_frac.get(), "Fractional", "l");
507  leg->AddEntry(heigen_cumulative.get(), "Cumulative", "l");
508 
509  leg->SetLineColor(0);
510  leg->SetLineStyle(0);
511  leg->SetFillColor(0);
512  leg->SetFillStyle(0);
513  leg->Draw("Same");
514 
515  c1->Print("Debug_PCA.pdf");
516  c1->SetRightMargin(0.15);
517  c1->SetLogy(0);
518 
519  heigen_vectors->SetMarkerSize(0.2);
520  minz = heigen_vectors->GetMinimum();
521  if (fabs(0-maxz)>fabs(0-minz)) heigen_vectors->GetZaxis()->SetRangeUser(0-fabs(0-maxz),0+fabs(0-maxz));
522  else heigen_vectors->GetZaxis()->SetRangeUser(0-fabs(0-minz),0+fabs(0-minz));
523  if(PlotText) heigen_vectors->Draw("COLZ TEXT");
524  else heigen_vectors->Draw("COLZ");
525 
526  auto Eigen_Line = std::make_unique<TLine>(0, nKeptPCApars, LastPCAdpar - FirstPCAdpar, nKeptPCApars);
527  Eigen_Line->SetLineColor(kGreen);
528  Eigen_Line->SetLineWidth(2);
529  Eigen_Line->SetLineStyle(kDotted);
530  Eigen_Line->Draw("SAME");
531  c1->Print("Debug_PCA.pdf");
532 
533  SubsetPCA->SetMarkerSize(0.2);
534  minz = SubsetPCA->GetMinimum();
535  if (fabs(0-maxz)>fabs(0-minz)) SubsetPCA->GetZaxis()->SetRangeUser(0-fabs(0-maxz),0+fabs(0-maxz));
536  else SubsetPCA->GetZaxis()->SetRangeUser(0-fabs(0-minz),0+fabs(0-minz));
537  if(PlotText) SubsetPCA->Draw("COLZ TEXT");
538  else SubsetPCA->Draw("COLZ");
539  c1->Print("Debug_PCA.pdf");
540  delete SubsetPCA;
541 
542  hTransferMat->SetMarkerSize(0.15);
543  minz = hTransferMat->GetMinimum();
544  if (fabs(0-maxz)>fabs(0-minz)) hTransferMat->GetZaxis()->SetRangeUser(0-fabs(0-maxz),0+fabs(0-maxz));
545  else hTransferMat->GetZaxis()->SetRangeUser(0-fabs(0-minz),0+fabs(0-minz));
546  if(PlotText) hTransferMat->Draw("COLZ TEXT");
547  else hTransferMat->Draw("COLZ");
548  c1->Print("Debug_PCA.pdf");
549  delete hTransferMat;
550 
551  hTransferMatT->SetMarkerSize(0.15);
552  minz = hTransferMatT->GetMinimum();
553  if (fabs(0-maxz)>fabs(0-minz)) hTransferMatT->GetZaxis()->SetRangeUser(0-fabs(0-maxz),0+fabs(0-maxz));
554  else hTransferMatT->GetZaxis()->SetRangeUser(0-fabs(0-minz),0+fabs(0-minz));
555  if(PlotText) hTransferMatT->Draw("COLZ TEXT");
556  else hTransferMatT->Draw("COLZ");
557  c1->Print( "Debug_PCA.pdf");
558  delete hTransferMatT;
559 
560  delete heigen_vectors;
561 
562  c1->Print("Debug_PCA.pdf]");
563  PCA_Debug->Close();
564  delete PCA_Debug;
565  gErrorIgnoreLevel = originalErrorWarning;
566 }
567 #endif
#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 _MaCh3_Safe_Include_Start_
KS: Avoiding warning checking for headers.
Definition: Core.h:126
#define _MaCh3_Safe_Include_End_
KS: Restore warning checking after including external headers.
Definition: Core.h:141
#define _restrict_
KS: Using restrict limits the effects of pointer aliasing, aiding optimizations. While reading I foun...
Definition: Core.h:108
#define MACH3LOG_DEBUG
Definition: MaCh3Logger.h:34
#define MACH3LOG_ERROR
Definition: MaCh3Logger.h:37
#define MACH3LOG_INFO
Definition: MaCh3Logger.h:35
#define MACH3LOG_WARN
Definition: MaCh3Logger.h:36
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:215
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:333
void ToggleFixAllParameters(const std::vector< std::string > &Names)
fix parameters at prior values
Definition: PCAHandler.cpp:309
TMatrixD TransferMat
Matrix used to converting from PCA base to normal base.
Definition: PCAHandler.h:220
std::vector< double > eigen_values_master
Eigen values which have dimension equal to _fNumParPCA, and can be used in CorrelateSteps.
Definition: PCAHandler.h:228
void SetBranches(TTree &tree, const bool SaveProposal, const std::vector< std::string > &Names)
set branches for output file
Definition: PCAHandler.cpp:294
TMatrixD eigen_vectors
Eigen vectors only of params which are being decomposed.
Definition: PCAHandler.h:226
TVectorD _fParPropPCA
CW: Current parameter value in PCA base.
Definition: PCAHandler.h:209
bool IsParameterFixedPCA(const int i) const
Is parameter fixed in PCA base or not.
Definition: PCAHandler.h:121
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 ...
Definition: PCAHandler.cpp:29
double GetParPropPCA(const int i) const
Get current parameter value using PCA.
Definition: PCAHandler.h:154
int NumParPCA
Number of parameters in PCA base.
Definition: PCAHandler.h:217
void SetInitialParameters()
KS: Transfer the starting parameters to the PCA basis, you don't want to start with zero....
Definition: PCAHandler.cpp:257
std::vector< M3::float_t > * _pPropVal
Pointer to proposed value of the parameter.
Definition: PCAHandler.h:242
std::vector< double > * _pCurrVal
Pointer to current value of the parameter.
Definition: PCAHandler.h:240
TMatrixD TransferMatT
Matrix used to converting from normal base to PCA base.
Definition: PCAHandler.h:222
void AcceptStep() _noexcept_
Accepted this step.
Definition: PCAHandler.cpp:191
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:206
void SetParPropPCA(const int i, const double value)
Set proposed value for parameter in PCA base.
Definition: PCAHandler.h:136
void TransferToPCA()
Transfer param values from normal base to PCA base.
Definition: PCAHandler.cpp:242
virtual ~PCAHandler()
Destructor.
Definition: PCAHandler.cpp:24
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
TVectorD eigen_values
Eigen value only of particles which are being decomposed.
Definition: PCAHandler.h:224
void TransferToParam()
Transfer param values from PCA base to normal base.
Definition: PCAHandler.cpp:267
int LastPCAdpar
Index of the last param that is being decomposed.
Definition: PCAHandler.h:235
int nKeptPCApars
Total number that remained after applying PCA Threshold.
Definition: PCAHandler.h:231
std::vector< double > _fErrorPCA
Tells if parameter is fixed in PCA base or not.
Definition: PCAHandler.h:213
PCAHandler()
Constructor.
Definition: PCAHandler.cpp:17
std::vector< double > _fPreFitValuePCA
Prefit value for PCA params.
Definition: PCAHandler.h:207
void ToggleFixParameter(const int i, const std::vector< std::string > &Names)
fix parameters at prior values
Definition: PCAHandler.cpp:315
void SetParCurrPCA(const int i, const double value)
Set current value for parameter in PCA base.
Definition: PCAHandler.h:145
void Print() const
KS: Print info about PCA parameters.
Definition: PCAHandler.cpp:285
TVectorD _fParCurrPCA
CW: Proposed parameter value in PCA base.
Definition: PCAHandler.h:211
void SanitisePCA(TMatrixDSym *CovMatrix)
KS: Make sure decomposed matrix isn't correlated with undecomposed.
Definition: PCAHandler.cpp:162
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:37
int FirstPCAdpar
Index of the first param that is being decomposed.
Definition: PCAHandler.h:233
double eigen_threshold
CW: Threshold based on which we remove parameters in eigen base.
Definition: PCAHandler.h:237
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.