MaCh3  2.5.1
Reference Guide
ParameterHandlerUtils.h
Go to the documentation of this file.
1 #pragma once
2 
3 // MaCh3 includes
5 
7 // ROOT includes
8 #include "TMatrixT.h"
9 #include "TMatrixDSym.h"
10 #include "TVectorT.h"
11 #include "TVectorD.h"
12 #include "TH1D.h"
13 #include "TH2D.h"
14 #include "TTree.h"
15 #include "TFile.h"
16 #include "TRandom3.h"
17 #include "TMath.h"
18 #include "TDecompChol.h"
19 #include "TStopwatch.h"
20 #include "TMatrix.h"
21 #include "TMatrixDSymEigen.h"
22 #include "TMatrixDEigen.h"
23 #include "TDecompSVD.h"
24 #include "TKey.h"
25 #include "TCanvas.h"
26 #include "TROOT.h"
27 #include "TStyle.h"
28 #include "TColor.h"
29 #include "TLine.h"
30 #include "TText.h"
31 #include "TLegend.h"
33 
34 
39 
40 namespace M3
41 {
43 inline double* MatrixMult(double *A, double *B, int n) {
44  //CW: First transpose to increse cache hits
45  double *BT = new double[n*n];
46  #ifdef MULTITHREAD
47  #pragma omp parallel for
48  #endif
49  for (int i = 0; i < n; i++) {
50  for (int j = 0; j < n; j++) {
51  BT[j*n+i] = B[i*n+j];
52  }
53  }
54 
55  // Now multiply
56  double *C = new double[n*n];
57  #ifdef MULTITHREAD
58  #pragma omp parallel for
59  #endif
60  for (int i = 0; i < n; i++) {
61  for (int j = 0; j < n; j++) {
62  double sum = 0;
63  for (int k = 0; k < n; k++) {
64  sum += A[i*n+k]*BT[j*n+k];
65  }
66  C[i*n+j] = sum;
67  }
68  }
69  delete BT;
70 
71  return C;
72 }
73 
75 inline double** MatrixMult(double **A, double **B, int n) {
76  // First make into monolithic array
77  double *A_mon = new double[n*n];
78  double *B_mon = new double[n*n];
79 
80  #ifdef MULTITHREAD
81  #pragma omp parallel for
82  #endif
83  for (int i = 0; i < n; ++i) {
84  for (int j = 0; j < n; ++j) {
85  A_mon[i*n+j] = A[i][j];
86  B_mon[i*n+j] = B[i][j];
87  }
88  }
89  //CW: Now call the monolithic calculator
90  double *C_mon = MatrixMult(A_mon, B_mon, n);
91  delete A_mon;
92  delete B_mon;
93 
94  // Return the double pointer
95  double **C = new double*[n];
96  #ifdef MULTITHREAD
97  #pragma omp parallel for
98  #endif
99  for (int i = 0; i < n; ++i) {
100  C[i] = new double[n];
101  for (int j = 0; j < n; ++j) {
102  C[i][j] = C_mon[i*n+j];
103  }
104  }
105  delete C_mon;
106 
107  return C;
108 }
109 
111 inline TMatrixD MatrixMult(TMatrixD A, TMatrixD B)
112 {
113  double *C_mon = MatrixMult(A.GetMatrixArray(), B.GetMatrixArray(), A.GetNcols());
114  TMatrixD C;
115  C.Use(A.GetNcols(), A.GetNrows(), C_mon);
116  return C;
117 }
118 
119 // ********************************************
125 inline void MatrixVectorMulti(double* _restrict_ VecMulti, double** _restrict_ matrix, const double* _restrict_ vector, const int n) {
126 // ********************************************
127  #ifdef MULTITHREAD
128  #pragma omp parallel for
129  #endif
130  for (int i = 0; i < n; ++i)
131  {
132  double result = 0.0;
133  #ifdef MULTITHREAD
134  #pragma omp simd
135  #endif
136  for (int j = 0; j < n; ++j)
137  {
138  result += matrix[i][j]*vector[j];
139  }
140  VecMulti[i] = result;
141  }
142 }
143 
144 // ********************************************
150 inline double MatrixVectorMultiSingle(double** _restrict_ matrix, const double* _restrict_ vector, const int Length, const int i) {
151 // ********************************************
152  double Element = 0.0;
153  #ifdef MULTITHREAD
154  #pragma omp simd
155  #endif
156  for (int j = 0; j < Length; ++j) {
157  Element += matrix[i][j]*vector[j];
158  }
159  return Element;
160 }
161 
162 // *************************************
165 inline void FixSampleNamesQuotes(std::string& yamlStr) {
166 // *************************************
167  std::stringstream input(yamlStr);
168  std::string line;
169  std::string fixedYaml;
170  std::regex sampleNamesRegex(R"(SampleNames:\s*\[([^\]]+)\])");
171 
172  while (std::getline(input, line)) {
173  std::smatch match;
174  if (std::regex_search(line, match, sampleNamesRegex)) {
175  std::string contents = match[1]; // inside the brackets
176  std::stringstream ss(contents);
177  std::string item;
178  std::vector<std::string> quotedItems;
179 
180  while (std::getline(ss, item, ',')) {
181  item = std::regex_replace(item, std::regex(R"(^\s+|\s+$)"), ""); // trim
182  quotedItems.push_back("\"" + item + "\"");
183  }
184 
185  std::string replacement = "SampleNames: [" + fmt::format("{}", fmt::join(quotedItems, ", ")) + "]";
186  line = std::regex_replace(line, sampleNamesRegex, replacement);
187  }
188  fixedYaml += line + "\n";
189  }
190 
191  yamlStr = fixedYaml;
192 }
193 
194 // *************************************
200 inline void AddTuneValues(YAML::Node& root,
201  const std::vector<double>& Values,
202  const std::string& Tune,
203  const std::vector<std::string>& FancyNames = {}) {
204 // *************************************
205  YAML::Node NodeCopy = YAML::Clone(root);
206  YAML::Node systematics = NodeCopy["Systematics"];
207 
208  if (!systematics || !systematics.IsSequence()) {
209  MACH3LOG_ERROR("'Systematics' node is missing or not a sequence in the YAML copy");
210  throw MaCh3Exception(__FILE__, __LINE__);
211  }
212 
213  if (!FancyNames.empty() && FancyNames.size() != Values.size()) {
214  MACH3LOG_ERROR("Mismatch in sizes: FancyNames has {}, but Values has {}", FancyNames.size(), Values.size());
215  throw MaCh3Exception(__FILE__, __LINE__);
216  }
217 
218  if (FancyNames.empty() && systematics.size() != Values.size()) {
219  MACH3LOG_ERROR("Mismatch in sizes: Values has {}, but YAML 'Systematics' has {} entries",
220  Values.size(), systematics.size());
221  throw MaCh3Exception(__FILE__, __LINE__);
222  }
223 
224  if (!FancyNames.empty()) {
225  for (std::size_t i = 0; i < FancyNames.size(); ++i) {
226  bool matched = false;
227  for (std::size_t j = 0; j < systematics.size(); ++j) {
228  YAML::Node systematicNode = systematics[j]["Systematic"];
229  if (!systematicNode) continue;
230  auto nameNode = systematicNode["Names"];
231  if (!nameNode || !nameNode["FancyName"]) continue;
232  if (nameNode["FancyName"].as<std::string>() == FancyNames[i]) {
233  if (!systematicNode["ParameterValues"]) {
234  MACH3LOG_ERROR("Missing 'ParameterValues' for matched FancyName '{}'", FancyNames[i]);
235  throw MaCh3Exception(__FILE__, __LINE__);
236  }
237  systematicNode["ParameterValues"][Tune] = M3::Utils::FormatDouble(Values[i], 4);
238  matched = true;
239  break;
240  }
241  }
242  if (!matched) {
243  MACH3LOG_ERROR("Could not find a matching FancyName '{}' in the systematics", FancyNames[i]);
244  throw MaCh3Exception(__FILE__, __LINE__);
245  }
246  }
247  } else {
248  for (std::size_t i = 0; i < systematics.size(); ++i) {
249  YAML::Node systematicNode = systematics[i]["Systematic"];
250  if (!systematicNode || !systematicNode["ParameterValues"]) {
251  MACH3LOG_ERROR("Missing 'Systematic' or 'ParameterValues' entry at index {}", i);
252  throw MaCh3Exception(__FILE__, __LINE__);
253  }
254  systematicNode["ParameterValues"][Tune] = M3::Utils::FormatDouble(Values[i], 4);
255  }
256  }
257 
258  // Convert updated copy to string
259  std::string YAMLString = YAMLtoSTRING(NodeCopy);
260  FixSampleNamesQuotes(YAMLString);
261  // Write to output file
262  std::string OutName = "UpdatedMatrixWithTune" + Tune + ".yaml";
263  std::ofstream outFile(OutName);
264  if (!outFile) {
265  MACH3LOG_ERROR("Failed to open file for writing: {}", OutName);
266  throw MaCh3Exception(__FILE__, __LINE__);
267  }
268 
269  outFile << YAMLString;
270  outFile.close();
271 }
272 
273 // *************************************
281 inline void MakeCorrelationMatrix(YAML::Node& root,
282  const std::vector<double>& Values,
283  const std::vector<double>& Errors,
284  const std::vector<std::vector<double>>& Correlation,
285  const std::string& OutYAMLName,
286  const std::vector<std::string>& FancyNames = {}) {
287 // *************************************
288  if (Values.size() != Errors.size() || Values.size() != Correlation.size()) {
289  MACH3LOG_ERROR("Size mismatch between Values, Errors, and Correlation matrix");
290  throw MaCh3Exception(__FILE__, __LINE__);
291  }
292 
293  for (const auto& row : Correlation) {
294  if (row.size() != Correlation.size()) {
295  MACH3LOG_ERROR("Correlation matrix is not square");
296  throw MaCh3Exception(__FILE__, __LINE__);
297  }
298  }
299 
300  YAML::Node NodeCopy = YAML::Clone(root);
301  YAML::Node systematics = NodeCopy["Systematics"];
302 
303  if (!systematics || !systematics.IsSequence()) {
304  MACH3LOG_ERROR("'Systematics' node is missing or not a sequence");
305  throw MaCh3Exception(__FILE__, __LINE__);
306  }
307 
308  if (!FancyNames.empty() && FancyNames.size() != Values.size()) {
309  MACH3LOG_ERROR("FancyNames size ({}) does not match Values size ({})", FancyNames.size(), Values.size());
310  throw MaCh3Exception(__FILE__, __LINE__);
311  }
312 
313  // Map from FancyName to Systematic node
314  std::unordered_map<std::string, YAML::Node> nameToNode;
315  for (std::size_t i = 0; i < systematics.size(); ++i) {
316  YAML::Node syst = systematics[i]["Systematic"];
317  if (!syst || !syst["Names"] || !syst["Names"]["FancyName"]) continue;
318  std::string name = syst["Names"]["FancyName"].as<std::string>();
319  nameToNode[name] = syst;
320  }
321 
322  if (!FancyNames.empty()) {
323  for (std::size_t i = 0; i < FancyNames.size(); ++i) {
324  const std::string& name_i = FancyNames[i];
325  auto it_i = nameToNode.find(name_i);
326  if (it_i == nameToNode.end()) {
327  MACH3LOG_ERROR("Could not find FancyName '{}' in YAML", name_i);
328  throw MaCh3Exception(__FILE__, __LINE__);
329  }
330  YAML::Node& syst_i = it_i->second;
331 
332  syst_i["ParameterValues"]["PreFitValue"] = M3::Utils::FormatDouble(Values[i], 4);
333  syst_i["Error"] = std::round(Errors[i] * 100.0) / 100.0;
334 
335  YAML::Node correlationsNode = YAML::Node(YAML::NodeType::Sequence);
336  for (std::size_t j = 0; j < FancyNames.size(); ++j) {
337  if (i == j) continue;
338  // KS: Skip if value close to 0
339  if (std::abs(Correlation[i][j]) < 1e-8) continue;
340  YAML::Node singleEntry;
341  singleEntry[FancyNames[j]] = M3::Utils::FormatDouble(Correlation[i][j], 4);
342  correlationsNode.push_back(singleEntry);
343  }
344  syst_i["Correlations"] = correlationsNode;
345  }
346  } else {
347  if (systematics.size() != Values.size()) {
348  MACH3LOG_ERROR("Mismatch in sizes: Values has {}, but YAML 'Systematics' has {} entries",
349  Values.size(), systematics.size());
350  throw MaCh3Exception(__FILE__, __LINE__);
351  }
352 
353  for (std::size_t i = 0; i < systematics.size(); ++i) {
354  YAML::Node syst = systematics[i]["Systematic"];
355  if (!syst) {
356  MACH3LOG_ERROR("Missing 'Systematic' node at index {}", i);
357  throw MaCh3Exception(__FILE__, __LINE__);
358  }
359 
360  syst["ParameterValues"]["PreFitValue"] = M3::Utils::FormatDouble(Values[i], 4);
361  syst["Error"] = std::round(Errors[i] * 100.0) / 100.0;
362 
363  YAML::Node correlationsNode = YAML::Node(YAML::NodeType::Sequence);
364  for (std::size_t j = 0; j < Correlation[i].size(); ++j) {
365  if (i == j) continue;
366  // KS: Skip if value close to 0
367  if (std::abs(Correlation[i][j]) < 1e-8) continue;
368  YAML::Node singleEntry;
369  const std::string& otherName = systematics[j]["Systematic"]["Names"]["FancyName"].as<std::string>();
370  singleEntry[otherName] = M3::Utils::FormatDouble(Correlation[i][j], 4);
371  correlationsNode.push_back(singleEntry);
372  }
373  syst["Correlations"] = correlationsNode;
374  }
375  }
376 
377  // Convert and write
378  std::string YAMLString = YAMLtoSTRING(NodeCopy);
379  FixSampleNamesQuotes(YAMLString);
380  std::ofstream outFile(OutYAMLName);
381  if (!outFile) {
382  MACH3LOG_ERROR("Failed to open file for writing: {}", OutYAMLName);
383  throw MaCh3Exception(__FILE__, __LINE__);
384  }
385 
386  outFile << YAMLString;
387  outFile.close();
388 }
389 
390 // *************************************
395 inline TMacro* GetConfigMacroFromChain(TDirectory* CovarianceFolder) {
396 // *************************************
397  if (!CovarianceFolder) {
398  MACH3LOG_ERROR("Null TDirectory passed to {}", __func__);
399  throw MaCh3Exception(__FILE__, __LINE__);
400  }
401 
402  TMacro* foundMacro = nullptr;
403  int macroCount = 0;
404 
405  TIter next(CovarianceFolder->GetListOfKeys());
406  TKey* key;
407  while ((key = dynamic_cast<TKey*>(next()))) {
408  if (std::string(key->GetClassName()) == "TMacro") {
409  ++macroCount;
410  if (macroCount == 1) {
411  foundMacro = dynamic_cast<TMacro*>(key->ReadObj());
412  }
413  }
414  }
415 
416  if (macroCount == 1 && foundMacro) {
417  MACH3LOG_INFO("Found single TMacro in directory: using it.");
418  return foundMacro;
419  } else {
420  MACH3LOG_WARN("Found {} TMacro objects. Using hardcoded macro name: Config_xsec_cov.", macroCount);
421  TMacro* fallback = CovarianceFolder->Get<TMacro>("Config_xsec_cov");
422  if (!fallback) {
423  MACH3LOG_WARN("Fallback macro 'Config_xsec_cov' not found in directory.");
424  }
425  return fallback;
426  }
427 }
428 
429 // *************************************
436 inline TMatrixDSym* GetCovMatrixFromChain(TDirectory* TempFile) {
437 // *************************************
438  if (!TempFile) {
439  MACH3LOG_ERROR("Null TDirectory passed to {}.", __func__);
440  throw MaCh3Exception(__FILE__, __LINE__);
441  }
442 
443  TMatrixDSym* foundMatrix = nullptr;
444  int matrixCount = 0;
445 
446  TIter next(TempFile->GetListOfKeys());
447  TKey* key;
448  while ((key = dynamic_cast<TKey*>(next()))) {
449  std::string className = key->GetClassName();
450  if (className.find("TMatrix") != std::string::npos) {
451  ++matrixCount;
452  if (matrixCount == 1) {
453  foundMatrix = dynamic_cast<TMatrixDSym*>(key->ReadObj());
454  }
455  }
456  }
457 
458  if (matrixCount == 1 && foundMatrix) {
459  MACH3LOG_INFO("Found single TMatrixDSym in directory: using it.");
460  return foundMatrix;
461  } else {
462  MACH3LOG_WARN("Found {} TMatrixDSym objects. Using hardcoded path: xsec_cov.", matrixCount);
463  TMatrixDSym* fallback = TempFile->Get<TMatrixDSym>("xsec_cov");
464  if (!fallback) {
465  MACH3LOG_WARN("Fallback matrix 'xsec_cov' not found.");
466  }
467  return fallback;
468  }
469 }
470 
471 // *************************************
475 inline std::vector<std::vector<double>> GetCholeskyDecomposedMatrix(const TMatrixDSym& matrix, const std::string& matrixName) {
476 // *************************************
477  const Int_t n = matrix.GetNrows();
478  std::vector<std::vector<double>> L(n, std::vector<double>(n, 0.0));
479 
480  for (Int_t j = 0; j < n; ++j) {
481  // Compute diagonal element (must be serial)
482  double sum_diag = matrix(j, j);
483  for (Int_t k = 0; k < j; ++k) {
484  sum_diag -= L[j][k] * L[j][k];
485  }
486  const double tol = 1e-15;
487  if (sum_diag <= tol) {
488  MACH3LOG_ERROR("Cholesky decomposition failed for {} (non-positive diagonal)", matrixName);
489  throw MaCh3Exception(__FILE__, __LINE__);
490  }
491  L[j][j] = std::sqrt(sum_diag);
492 
493  // Compute the rest of the column in parallel
494  #ifdef MULTITHREAD
495  #pragma omp parallel for
496  #endif
497  for (Int_t i = j + 1; i < n; ++i) {
498  double sum = matrix(i, j);
499  for (Int_t k = 0; k < j; ++k) {
500  sum -= L[i][k] * L[j][k];
501  }
502  L[i][j] = sum / L[j][j];
503  }
504  }
505  return L;
506 }
507 
508 // *************************************
511 inline bool CanDecomposeMatrix(const TMatrixDSym& matrix) {
512 // *************************************
513  TDecompChol chdcmp(matrix);
514  return chdcmp.Decompose();
515 }
516 
517 // *************************************
519 inline void MakeMatrixPosDef(TMatrixDSym *cov) {
520 // *************************************
521  //DB Save original warning state and then increase it in this function to suppress 'matrix not positive definite' messages
522  //Means we no longer need to overload
523  int originalErrorWarning = gErrorIgnoreLevel;
524  gErrorIgnoreLevel = kFatal;
525 
526  //DB Loop 1000 times adding 1e-9 which tops out at 1e-6 shift on the diagonal before throwing error
527  constexpr int MaxAttempts = 1e5;
528  const int matrixSize = cov->GetNrows();
529  int iAttempt = 0;
530  bool CanDecomp = false;
531 
532  for (iAttempt = 0; iAttempt < MaxAttempts; iAttempt++) {
533  if (CanDecomposeMatrix(*cov)) {
534  CanDecomp = true;
535  break;
536  } else {
537  #ifdef MULTITHREAD
538  #pragma omp parallel for
539  #endif
540  for (int iVar = 0 ; iVar < matrixSize; iVar++) {
541  (*cov)(iVar, iVar) += 1e-9;
542  }
543  }
544  }
545 
546  if (!CanDecomp) {
547  MACH3LOG_ERROR("Tried {} times to shift diagonal but still can not decompose the matrix", MaxAttempts);
548  MACH3LOG_ERROR("This indicates that something is wrong with the input matrix");
549  throw MaCh3Exception(__FILE__ , __LINE__ );
550  }
551 
552  //DB Resetting warning level
553  gErrorIgnoreLevel = originalErrorWarning;
554 }
555 
556 
557 
558 // ********************************************
561 inline void DumpParamHandlerToFile(const int _fNumPar,
562  const std::vector<double>& _fPreFitValue,
563  const std::vector<double>& _fError,
564  const std::vector<double>& _fLowBound,
565  const std::vector<double>& _fUpBound,
566  const std::vector<double>& _fIndivStepScale,
567  const std::vector<std::string>& _fFancyNames,
568  const std::vector<bool>& _fFlatPrior,
569  const std::vector<SplineParameter>& SplineParams,
570  TMatrixDSym* covMatrix,
571  TH2D* CorrMatrix,
572  const std::string& Name) {
573 // ********************************************
574  TFile* outputFile = new TFile(Name.c_str(), "RECREATE");
575 
576  TObjArray* param_names = new TObjArray();
577  TObjArray* spline_interpolation = new TObjArray();
578  TObjArray* spline_names = new TObjArray();
579 
580  TVectorD* param_prior = new TVectorD(_fNumPar);
581  TVectorD* flat_prior = new TVectorD(_fNumPar);
582  TVectorD* stepscale = new TVectorD(_fNumPar);
583  TVectorD* param_lb = new TVectorD(_fNumPar);
584  TVectorD* param_ub = new TVectorD(_fNumPar);
585 
586  TVectorD* param_knot_weight_lb = new TVectorD(_fNumPar);
587  TVectorD* param_knot_weight_ub = new TVectorD(_fNumPar);
588  TVectorD* error = new TVectorD(_fNumPar);
589 
590  for(int i = 0; i < _fNumPar; ++i)
591  {
592  TObjString* nameObj = new TObjString(_fFancyNames[i].c_str());
593  param_names->AddLast(nameObj);
594 
595  TObjString* splineType = new TObjString("TSpline3");
596  spline_interpolation->AddLast(splineType);
597 
598  TObjString* splineName = new TObjString("");
599  spline_names->AddLast(splineName);
600 
601  (*param_prior)[i] = _fPreFitValue[i];
602  (*flat_prior)[i] = _fFlatPrior[i];
603  (*stepscale)[i] = _fIndivStepScale[i];
604  (*error)[i] = _fError[i];
605 
606  (*param_lb)[i] = _fLowBound[i];
607  (*param_ub)[i] = _fUpBound[i];
608 
609  //Default values
610  (*param_knot_weight_lb)[i] = -9999;
611  (*param_knot_weight_ub)[i] = +9999;
612  }
613 
614  for (size_t SplineIndex = 0; SplineIndex < SplineParams.size(); SplineIndex++ ) {
615  auto SystIndex = SplineParams[SplineIndex].index;
616 
617  (*param_knot_weight_lb)[SystIndex] = SplineParams.at(SplineIndex)._SplineKnotLowBound;
618  (*param_knot_weight_ub)[SystIndex] = SplineParams.at(SplineIndex)._SplineKnotUpBound;
619 
620  TObjString* splineType = new TObjString(SplineInterpolation_ToString(SplineParams.at(SplineIndex)._SplineInterpolationType).c_str());
621  spline_interpolation->AddAt(splineType, SystIndex);
622 
623  TObjString* splineName = new TObjString(SplineParams[SplineIndex]._fSplineNames.c_str());
624  spline_names->AddAt(splineName, SystIndex);
625  }
626  param_names->Write("xsec_param_names", TObject::kSingleKey);
627  delete param_names;
628  spline_interpolation->Write("xsec_spline_interpolation", TObject::kSingleKey);
629  delete spline_interpolation;
630  spline_names->Write("xsec_spline_names", TObject::kSingleKey);
631  delete spline_names;
632 
633  param_prior->Write("xsec_param_prior");
634  delete param_prior;
635  flat_prior->Write("xsec_flat_prior");
636  delete flat_prior;
637  stepscale->Write("xsec_stepscale");
638  delete stepscale;
639  param_lb->Write("xsec_param_lb");
640  delete param_lb;
641  param_ub->Write("xsec_param_ub");
642  delete param_ub;
643 
644  param_knot_weight_lb->Write("xsec_param_knot_weight_lb");
645  delete param_knot_weight_lb;
646  param_knot_weight_ub->Write("xsec_param_knot_weight_ub");
647  delete param_knot_weight_ub;
648  error->Write("xsec_error");
649  delete error;
650 
651  covMatrix->Write("xsec_cov");
652  CorrMatrix->Write("hcov");
653 
654  outputFile->Close();
655  delete outputFile;
656 }
657 
658 // ********************************************
660 //KS: Let's dump all useful matrices to properly validate PCA
661 inline void DebugPCA(const double sum,
662  const TMatrixD& temp,
663  const TMatrixD& eigen_vectors,
664  const TVectorD& eigen_values,
665  const TMatrixD& TransferMat,
666  const TMatrixD& TransferMatT,
667  const int NumPar,
668  const int FirstPCAdpar,
669  const int LastPCAdpar,
670  const int nKeptPCApars,
671  const double eigen_threshold) {
672 // ********************************************
673  #pragma GCC diagnostic push
674  #pragma GCC diagnostic ignored "-Wfloat-conversion"
675  int originalErrorWarning = gErrorIgnoreLevel;
676  gErrorIgnoreLevel = kFatal;
677 
678  TDirectory *ogdir = gDirectory;
679 
680  TFile *PCA_Debug = new TFile("Debug_PCA.root", "RECREATE");
681  PCA_Debug->cd();
682 
683  bool PlotText = true;
684  //KS: If we have more than 200 plot becomes unreadable :(
685  if(NumPar > 200) PlotText = false;
686 
687  auto heigen_values = std::make_unique<TH1D>("eigen_values", "Eigen Values", eigen_values.GetNrows(), 0.0, eigen_values.GetNrows());
688  heigen_values->SetDirectory(nullptr);
689  auto heigen_cumulative = std::make_unique<TH1D>("heigen_cumulative", "heigen_cumulative", eigen_values.GetNrows(), 0.0, eigen_values.GetNrows());
690  heigen_cumulative->SetDirectory(nullptr);
691  auto heigen_frac = std::make_unique<TH1D>("heigen_fractional", "heigen_fractional", eigen_values.GetNrows(), 0.0, eigen_values.GetNrows());
692  heigen_frac->SetDirectory(nullptr);
693  heigen_values->GetXaxis()->SetTitle("Eigen Vector");
694  heigen_values->GetYaxis()->SetTitle("Eigen Value");
695 
696  double Cumulative = 0;
697  for(int i = 0; i < eigen_values.GetNrows(); i++)
698  {
699  heigen_values->SetBinContent(i+1, (eigen_values)(i));
700  heigen_cumulative->SetBinContent(i+1, (eigen_values)(i)/sum + Cumulative);
701  heigen_frac->SetBinContent(i+1, (eigen_values)(i)/sum);
702  Cumulative += (eigen_values)(i)/sum;
703  }
704  heigen_values->Write("heigen_values");
705  eigen_values.Write("eigen_values");
706  heigen_cumulative->Write("heigen_values_cumulative");
707  heigen_frac->Write("heigen_values_frac");
708 
709  TH2D* heigen_vectors = new TH2D(eigen_vectors);
710  heigen_vectors->GetXaxis()->SetTitle("Parameter in Normal Base");
711  heigen_vectors->GetYaxis()->SetTitle("Parameter in Decomposed Base");
712  heigen_vectors->Write("heigen_vectors");
713  eigen_vectors.Write("eigen_vectors");
714 
715  TH2D* SubsetPCA = new TH2D(temp);
716  SubsetPCA->GetXaxis()->SetTitle("Parameter in Normal Base");
717  SubsetPCA->GetYaxis()->SetTitle("Parameter in Decomposed Base");
718 
719  SubsetPCA->Write("hSubsetPCA");
720  temp.Write("SubsetPCA");
721  TH2D* hTransferMat = new TH2D(TransferMat);
722  hTransferMat->GetXaxis()->SetTitle("Parameter in Normal Base");
723  hTransferMat->GetYaxis()->SetTitle("Parameter in Decomposed Base");
724  TH2D* hTransferMatT = new TH2D(TransferMatT);
725 
726  hTransferMatT->GetXaxis()->SetTitle("Parameter in Decomposed Base");
727  hTransferMatT->GetYaxis()->SetTitle("Parameter in Normal Base");
728 
729  hTransferMat->Write("hTransferMat");
730  TransferMat.Write("TransferMat");
731  hTransferMatT->Write("hTransferMatT");
732  TransferMatT.Write("TransferMatT");
733 
734  auto c1 = std::make_unique<TCanvas>("c1", " ", 0, 0, 1024, 1024);
735  c1->SetBottomMargin(0.1);
736  c1->SetTopMargin(0.05);
737  c1->SetRightMargin(0.05);
738  c1->SetLeftMargin(0.12);
739  c1->SetGrid();
740 
741  gStyle->SetPaintTextFormat("4.1f");
742  gStyle->SetOptFit(0);
743  gStyle->SetOptStat(0);
744  // Make pretty correlation colors (red to blue)
745  constexpr int NRGBs = 5;
746  TColor::InitializeColors();
747  Double_t stops[NRGBs] = { 0.00, 0.25, 0.50, 0.75, 1.00 };
748  Double_t red[NRGBs] = { 0.00, 0.25, 1.00, 1.00, 0.50 };
749  Double_t green[NRGBs] = { 0.00, 0.25, 1.00, 0.25, 0.00 };
750  Double_t blue[NRGBs] = { 0.50, 1.00, 1.00, 0.25, 0.00 };
751  TColor::CreateGradientColorTable(5, stops, red, green, blue, 255);
752  gStyle->SetNumberContours(255);
753 
754  double maxz = 0;
755  double minz = 0;
756 
757  c1->Print("Debug_PCA.pdf[");
758  auto EigenLine = std::make_unique<TLine>(nKeptPCApars, 0, nKeptPCApars, heigen_cumulative->GetMaximum());
759  EigenLine->SetLineColor(kPink);
760  EigenLine->SetLineWidth(2);
761  EigenLine->SetLineStyle(kSolid);
762 
763  auto text = std::make_unique<TText>(0.5, 0.5, Form("Threshold = %g", eigen_threshold));
764  text->SetTextFont (43);
765  text->SetTextSize (40);
766 
767  heigen_values->SetLineColor(kRed);
768  heigen_values->SetLineWidth(2);
769  heigen_cumulative->SetLineColor(kGreen);
770  heigen_cumulative->SetLineWidth(2);
771  heigen_frac->SetLineColor(kBlue);
772  heigen_frac->SetLineWidth(2);
773 
774  c1->SetLogy();
775  heigen_values->SetMaximum(heigen_cumulative->GetMaximum()+heigen_cumulative->GetMaximum()*0.4);
776  heigen_values->Draw();
777  heigen_frac->Draw("SAME");
778  heigen_cumulative->Draw("SAME");
779  EigenLine->Draw("Same");
780  text->DrawTextNDC(0.42, 0.84,Form("Threshold = %g", eigen_threshold));
781 
782  auto leg = std::make_unique<TLegend>(0.2, 0.2, 0.6, 0.5);
783  leg->SetTextSize(0.04);
784  leg->AddEntry(heigen_values.get(), "Absolute", "l");
785  leg->AddEntry(heigen_frac.get(), "Fractional", "l");
786  leg->AddEntry(heigen_cumulative.get(), "Cumulative", "l");
787 
788  leg->SetLineColor(0);
789  leg->SetLineStyle(0);
790  leg->SetFillColor(0);
791  leg->SetFillStyle(0);
792  leg->Draw("Same");
793 
794  c1->Print("Debug_PCA.pdf");
795  c1->SetRightMargin(0.15);
796  c1->SetLogy(0);
797 
798  heigen_vectors->SetMarkerSize(0.2);
799  minz = heigen_vectors->GetMinimum();
800  if (fabs(0-maxz)>fabs(0-minz)) heigen_vectors->GetZaxis()->SetRangeUser(0-fabs(0-maxz),0+fabs(0-maxz));
801  else heigen_vectors->GetZaxis()->SetRangeUser(0-fabs(0-minz),0+fabs(0-minz));
802  if(PlotText) heigen_vectors->Draw("COLZ TEXT");
803  else heigen_vectors->Draw("COLZ");
804 
805  auto Eigen_Line = std::make_unique<TLine>(0, nKeptPCApars, LastPCAdpar - FirstPCAdpar, nKeptPCApars);
806  Eigen_Line->SetLineColor(kGreen);
807  Eigen_Line->SetLineWidth(2);
808  Eigen_Line->SetLineStyle(kDotted);
809  Eigen_Line->Draw("SAME");
810  c1->Print("Debug_PCA.pdf");
811 
812  SubsetPCA->SetMarkerSize(0.2);
813  minz = SubsetPCA->GetMinimum();
814  if (fabs(0-maxz)>fabs(0-minz)) SubsetPCA->GetZaxis()->SetRangeUser(0-fabs(0-maxz),0+fabs(0-maxz));
815  else SubsetPCA->GetZaxis()->SetRangeUser(0-fabs(0-minz),0+fabs(0-minz));
816  if(PlotText) SubsetPCA->Draw("COLZ TEXT");
817  else SubsetPCA->Draw("COLZ");
818  c1->Print("Debug_PCA.pdf");
819  delete SubsetPCA;
820 
821  hTransferMat->SetMarkerSize(0.15);
822  minz = hTransferMat->GetMinimum();
823  if (fabs(0-maxz)>fabs(0-minz)) hTransferMat->GetZaxis()->SetRangeUser(0-fabs(0-maxz),0+fabs(0-maxz));
824  else hTransferMat->GetZaxis()->SetRangeUser(0-fabs(0-minz),0+fabs(0-minz));
825  if(PlotText) hTransferMat->Draw("COLZ TEXT");
826  else hTransferMat->Draw("COLZ");
827  c1->Print("Debug_PCA.pdf");
828  delete hTransferMat;
829 
830  hTransferMatT->SetMarkerSize(0.15);
831  minz = hTransferMatT->GetMinimum();
832  if (fabs(0-maxz)>fabs(0-minz)) hTransferMatT->GetZaxis()->SetRangeUser(0-fabs(0-maxz),0+fabs(0-maxz));
833  else hTransferMatT->GetZaxis()->SetRangeUser(0-fabs(0-minz),0+fabs(0-minz));
834  if(PlotText) hTransferMatT->Draw("COLZ TEXT");
835  else hTransferMatT->Draw("COLZ");
836  c1->Print( "Debug_PCA.pdf");
837  delete hTransferMatT;
838 
839  delete heigen_vectors;
840 
841  c1->Print("Debug_PCA.pdf]");
842  PCA_Debug->Close();
843  delete PCA_Debug;
844  gErrorIgnoreLevel = originalErrorWarning;
845  ogdir->cd(); // go back to original directory
846  #pragma GCC diagnostic pop
847 }
848 
849 } // end M3 namespace
#define _MaCh3_Safe_Include_Start_
KS: Avoiding warning checking for headers.
Definition: Core.h:126
#define _MaCh3_Safe_Include_End_
#define _restrict_
KS: Using restrict limits the effects of pointer aliasing, aiding optimizations. While reading I foun...
Definition: Core.h:108
#define MACH3LOG_ERROR
Definition: MaCh3Logger.h:37
#define MACH3LOG_INFO
Definition: MaCh3Logger.h:35
#define MACH3LOG_WARN
Definition: MaCh3Logger.h:36
Definitions of generic parameter structs and utility templates for MaCh3.
std::string SplineInterpolation_ToString(const SplineInterpolation i)
Convert a LLH type to a string.
std::string YAMLtoSTRING(const YAML::Node &node)
KS: Convert a YAML node to a string representation.
Definition: YamlHelper.h:112
Custom exception class used throughout MaCh3.
std::string FormatDouble(const double value, const int precision)
Convert double into string for precision, useful for playing with yaml if you don't want to have in c...
Main namespace for MaCh3 software.
std::unique_ptr< ObjectType > Clone(const ObjectType *obj, const std::string &name="")
KS: Creates a copy of a ROOT-like object and wraps it in a smart pointer.
void MatrixVectorMulti(double *_restrict_ VecMulti, double **_restrict_ matrix, const double *_restrict_ vector, const int n)
KS: Custom function to perform multiplication of matrix and vector with multithreading.
double * MatrixMult(double *A, double *B, int n)
CW: Multi-threaded matrix multiplication.
void FixSampleNamesQuotes(std::string &yamlStr)
KS: Yaml emitter has problem and drops "", if you have special signs in you like * then there is prob...
void DumpParamHandlerToFile(const int _fNumPar, const std::vector< double > &_fPreFitValue, const std::vector< double > &_fError, const std::vector< double > &_fLowBound, const std::vector< double > &_fUpBound, const std::vector< double > &_fIndivStepScale, const std::vector< std::string > &_fFancyNames, const std::vector< bool > &_fFlatPrior, const std::vector< SplineParameter > &SplineParams, TMatrixDSym *covMatrix, TH2D *CorrMatrix, const std::string &Name)
Dump Matrix to ROOT file, useful when we need to pass matrix info to another fitting group.
void MakeCorrelationMatrix(YAML::Node &root, const std::vector< double > &Values, const std::vector< double > &Errors, const std::vector< std::vector< double >> &Correlation, const std::string &OutYAMLName, const std::vector< std::string > &FancyNames={})
KS: Replace correlation matrix and tune values in YAML covariance matrix.
void AddTuneValues(YAML::Node &root, const std::vector< double > &Values, const std::string &Tune, const std::vector< std::string > &FancyNames={})
KS: Add Tune values to YAML covariance matrix.
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.
void MakeMatrixPosDef(TMatrixDSym *cov)
Makes sure that matrix is positive-definite by adding a small number to on-diagonal elements.
TMatrixDSym * GetCovMatrixFromChain(TDirectory *TempFile)
KS: Retrieve the cross-section covariance matrix from the given TDirectory. Historically,...
bool CanDecomposeMatrix(const TMatrixDSym &matrix)
Checks if a matrix can be Cholesky decomposed.
TMacro * GetConfigMacroFromChain(TDirectory *CovarianceFolder)
KS: We store configuration macros inside the chain. In the past, multiple configs were stored,...
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.
std::vector< std::vector< double > > GetCholeskyDecomposedMatrix(const TMatrixDSym &matrix, const std::string &matrixName)
Computes Cholesky decomposition of a symmetric positive definite matrix using custom function which c...
Flat representation of a spline index entry.
Definition: SplineCommon.h:33