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