MaCh3  2.2.3
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"
26 
27 namespace M3
28 {
30 inline double* MatrixMult(double *A, double *B, int n) {
31  //CW: First transpose to increse cache hits
32  double *BT = new double[n*n];
33  #ifdef MULTITHREAD
34  #pragma omp parallel for
35  #endif
36  for (int i = 0; i < n; i++) {
37  for (int j = 0; j < n; j++) {
38  BT[j*n+i] = B[i*n+j];
39  }
40  }
41 
42  // Now multiply
43  double *C = new double[n*n];
44  #ifdef MULTITHREAD
45  #pragma omp parallel for
46  #endif
47  for (int i = 0; i < n; i++) {
48  for (int j = 0; j < n; j++) {
49  double sum = 0;
50  for (int k = 0; k < n; k++) {
51  sum += A[i*n+k]*BT[j*n+k];
52  }
53  C[i*n+j] = sum;
54  }
55  }
56  delete BT;
57 
58  return C;
59 }
60 
62 inline double** MatrixMult(double **A, double **B, int n) {
63  // First make into monolithic array
64  double *A_mon = new double[n*n];
65  double *B_mon = new double[n*n];
66 
67  #ifdef MULTITHREAD
68  #pragma omp parallel for
69  #endif
70  for (int i = 0; i < n; ++i) {
71  for (int j = 0; j < n; ++j) {
72  A_mon[i*n+j] = A[i][j];
73  B_mon[i*n+j] = B[i][j];
74  }
75  }
76  //CW: Now call the monolithic calculator
77  double *C_mon = MatrixMult(A_mon, B_mon, n);
78  delete A_mon;
79  delete B_mon;
80 
81  // Return the double pointer
82  double **C = new double*[n];
83  #ifdef MULTITHREAD
84  #pragma omp parallel for
85  #endif
86  for (int i = 0; i < n; ++i) {
87  C[i] = new double[n];
88  for (int j = 0; j < n; ++j) {
89  C[i][j] = C_mon[i*n+j];
90  }
91  }
92  delete C_mon;
93 
94  return C;
95 }
96 
98 inline TMatrixD MatrixMult(TMatrixD A, TMatrixD B)
99 {
100  double *C_mon = MatrixMult(A.GetMatrixArray(), B.GetMatrixArray(), A.GetNcols());
101  TMatrixD C;
102  C.Use(A.GetNcols(), A.GetNrows(), C_mon);
103  return C;
104 }
105 
106 // ********************************************
112 inline void MatrixVectorMulti(double* _restrict_ VecMulti, double** _restrict_ matrix, const double* _restrict_ vector, const int n) {
113 // ********************************************
114  #ifdef MULTITHREAD
115  #pragma omp parallel for
116  #endif
117  for (int i = 0; i < n; ++i)
118  {
119  double result = 0.0;
120  #ifdef MULTITHREAD
121  #pragma omp simd
122  #endif
123  for (int j = 0; j < n; ++j)
124  {
125  result += matrix[i][j]*vector[j];
126  }
127  VecMulti[i] = result;
128  }
129 }
130 
131 // ********************************************
137 inline double MatrixVectorMultiSingle(double** _restrict_ matrix, const double* _restrict_ vector, const int Length, const int i) {
138 // ********************************************
139  double Element = 0.0;
140  #ifdef MULTITHREAD
141  #pragma omp simd
142  #endif
143  for (int j = 0; j < Length; ++j) {
144  Element += matrix[i][j]*vector[j];
145  }
146  return Element;
147 }
148 
149 // *************************************
152 inline void FixSampleNamesQuotes(std::string& yamlStr) {
153 // *************************************
154  std::stringstream input(yamlStr);
155  std::string line;
156  std::string fixedYaml;
157  std::regex sampleNamesRegex(R"(SampleNames:\s*\[([^\]]+)\])");
158 
159  while (std::getline(input, line)) {
160  std::smatch match;
161  if (std::regex_search(line, match, sampleNamesRegex)) {
162  std::string contents = match[1]; // inside the brackets
163  std::stringstream ss(contents);
164  std::string item;
165  std::vector<std::string> quotedItems;
166 
167  while (std::getline(ss, item, ',')) {
168  item = std::regex_replace(item, std::regex(R"(^\s+|\s+$)"), ""); // trim
169  quotedItems.push_back("\"" + item + "\"");
170  }
171 
172  std::string replacement = "SampleNames: [" + fmt::format("{}", fmt::join(quotedItems, ", ")) + "]";
173  line = std::regex_replace(line, sampleNamesRegex, replacement);
174  }
175  fixedYaml += line + "\n";
176  }
177 
178  yamlStr = fixedYaml;
179 }
180 
181 // *************************************
187 inline void AddTuneValues(YAML::Node& root,
188  const std::vector<double>& Values,
189  const std::string& Tune,
190  const std::vector<std::string>& FancyNames = {}) {
191 // *************************************
192  YAML::Node NodeCopy = YAML::Clone(root);
193  YAML::Node systematics = NodeCopy["Systematics"];
194 
195  if (!systematics || !systematics.IsSequence()) {
196  MACH3LOG_ERROR("'Systematics' node is missing or not a sequence in the YAML copy");
197  throw MaCh3Exception(__FILE__, __LINE__);
198  }
199 
200  if (!FancyNames.empty() && FancyNames.size() != Values.size()) {
201  MACH3LOG_ERROR("Mismatch in sizes: FancyNames has {}, but Values has {}", FancyNames.size(), Values.size());
202  throw MaCh3Exception(__FILE__, __LINE__);
203  }
204 
205  if (FancyNames.empty() && systematics.size() != Values.size()) {
206  MACH3LOG_ERROR("Mismatch in sizes: Values has {}, but YAML 'Systematics' has {} entries",
207  Values.size(), systematics.size());
208  throw MaCh3Exception(__FILE__, __LINE__);
209  }
210 
211  if (!FancyNames.empty()) {
212  for (std::size_t i = 0; i < FancyNames.size(); ++i) {
213  bool matched = false;
214  for (std::size_t j = 0; j < systematics.size(); ++j) {
215  YAML::Node systematicNode = systematics[j]["Systematic"];
216  if (!systematicNode) continue;
217  auto nameNode = systematicNode["Names"];
218  if (!nameNode || !nameNode["FancyName"]) continue;
219  if (nameNode["FancyName"].as<std::string>() == FancyNames[i]) {
220  if (!systematicNode["ParameterValues"]) {
221  MACH3LOG_ERROR("Missing 'ParameterValues' for matched FancyName '{}'", FancyNames[i]);
222  throw MaCh3Exception(__FILE__, __LINE__);
223  }
224  systematicNode["ParameterValues"][Tune] = MaCh3Utils::FormatDouble(Values[i], 4);
225  matched = true;
226  break;
227  }
228  }
229  if (!matched) {
230  MACH3LOG_ERROR("Could not find a matching FancyName '{}' in the systematics", FancyNames[i]);
231  throw MaCh3Exception(__FILE__, __LINE__);
232  }
233  }
234  } else {
235  for (std::size_t i = 0; i < systematics.size(); ++i) {
236  YAML::Node systematicNode = systematics[i]["Systematic"];
237  if (!systematicNode || !systematicNode["ParameterValues"]) {
238  MACH3LOG_ERROR("Missing 'Systematic' or 'ParameterValues' entry at index {}", i);
239  throw MaCh3Exception(__FILE__, __LINE__);
240  }
241  systematicNode["ParameterValues"][Tune] = MaCh3Utils::FormatDouble(Values[i], 4);
242  }
243  }
244 
245  // Convert updated copy to string
246  std::string YAMLString = YAMLtoSTRING(NodeCopy);
247  FixSampleNamesQuotes(YAMLString);
248  // Write to output file
249  std::string OutName = "UpdatedMatrixWithTune" + Tune + ".yaml";
250  std::ofstream outFile(OutName);
251  if (!outFile) {
252  MACH3LOG_ERROR("Failed to open file for writing: {}", OutName);
253  throw MaCh3Exception(__FILE__, __LINE__);
254  }
255 
256  outFile << YAMLString;
257  outFile.close();
258 }
259 
260 // *************************************
268 inline void MakeCorrelationMatrix(YAML::Node& root,
269  const std::vector<double>& Values,
270  const std::vector<double>& Errors,
271  const std::vector<std::vector<double>>& Correlation,
272  const std::string& OutYAMLName,
273  const std::vector<std::string>& FancyNames = {}) {
274 // *************************************
275  if (Values.size() != Errors.size() || Values.size() != Correlation.size()) {
276  MACH3LOG_ERROR("Size mismatch between Values, Errors, and Correlation matrix");
277  throw MaCh3Exception(__FILE__, __LINE__);
278  }
279 
280  for (const auto& row : Correlation) {
281  if (row.size() != Correlation.size()) {
282  MACH3LOG_ERROR("Correlation matrix is not square");
283  throw MaCh3Exception(__FILE__, __LINE__);
284  }
285  }
286 
287  YAML::Node NodeCopy = YAML::Clone(root);
288  YAML::Node systematics = NodeCopy["Systematics"];
289 
290  if (!systematics || !systematics.IsSequence()) {
291  MACH3LOG_ERROR("'Systematics' node is missing or not a sequence");
292  throw MaCh3Exception(__FILE__, __LINE__);
293  }
294 
295  if (!FancyNames.empty() && FancyNames.size() != Values.size()) {
296  MACH3LOG_ERROR("FancyNames size ({}) does not match Values size ({})", FancyNames.size(), Values.size());
297  throw MaCh3Exception(__FILE__, __LINE__);
298  }
299 
300  // Map from FancyName to Systematic node
301  std::unordered_map<std::string, YAML::Node> nameToNode;
302  for (std::size_t i = 0; i < systematics.size(); ++i) {
303  YAML::Node syst = systematics[i]["Systematic"];
304  if (!syst || !syst["Names"] || !syst["Names"]["FancyName"]) continue;
305  std::string name = syst["Names"]["FancyName"].as<std::string>();
306  nameToNode[name] = syst;
307  }
308 
309  if (!FancyNames.empty()) {
310  for (std::size_t i = 0; i < FancyNames.size(); ++i) {
311  const std::string& name_i = FancyNames[i];
312  auto it_i = nameToNode.find(name_i);
313  if (it_i == nameToNode.end()) {
314  MACH3LOG_ERROR("Could not find FancyName '{}' in YAML", name_i);
315  throw MaCh3Exception(__FILE__, __LINE__);
316  }
317  YAML::Node& syst_i = it_i->second;
318 
319  syst_i["ParameterValues"]["PreFitValue"] = MaCh3Utils::FormatDouble(Values[i], 4);
320  syst_i["Error"] = std::round(Errors[i] * 100.0) / 100.0;
321 
322  YAML::Node correlationsNode = YAML::Node(YAML::NodeType::Sequence);
323  for (std::size_t j = 0; j < FancyNames.size(); ++j) {
324  if (i == j) continue;
325  // KS: Skip if value close to 0
326  if (std::abs(Correlation[i][j]) < 1e-8) continue;
327  YAML::Node singleEntry;
328  singleEntry[FancyNames[j]] = MaCh3Utils::FormatDouble(Correlation[i][j], 4);
329  correlationsNode.push_back(singleEntry);
330  }
331  syst_i["Correlations"] = correlationsNode;
332  }
333  } else {
334  if (systematics.size() != Values.size()) {
335  MACH3LOG_ERROR("Mismatch in sizes: Values has {}, but YAML 'Systematics' has {} entries",
336  Values.size(), systematics.size());
337  throw MaCh3Exception(__FILE__, __LINE__);
338  }
339 
340  for (std::size_t i = 0; i < systematics.size(); ++i) {
341  YAML::Node syst = systematics[i]["Systematic"];
342  if (!syst) {
343  MACH3LOG_ERROR("Missing 'Systematic' node at index {}", i);
344  throw MaCh3Exception(__FILE__, __LINE__);
345  }
346 
347  syst["ParameterValues"]["PreFitValue"] = MaCh3Utils::FormatDouble(Values[i], 4);
348  syst["Error"] = std::round(Errors[i] * 100.0) / 100.0;
349 
350  YAML::Node correlationsNode = YAML::Node(YAML::NodeType::Sequence);
351  for (std::size_t j = 0; j < Correlation[i].size(); ++j) {
352  if (i == j) continue;
353  // KS: Skip if value close to 0
354  if (std::abs(Correlation[i][j]) < 1e-8) continue;
355  YAML::Node singleEntry;
356  const std::string& otherName = systematics[j]["Systematic"]["Names"]["FancyName"].as<std::string>();
357  singleEntry[otherName] = MaCh3Utils::FormatDouble(Correlation[i][j], 4);
358  correlationsNode.push_back(singleEntry);
359  }
360  syst["Correlations"] = correlationsNode;
361  }
362  }
363 
364  // Convert and write
365  std::string YAMLString = YAMLtoSTRING(NodeCopy);
366  FixSampleNamesQuotes(YAMLString);
367  std::ofstream outFile(OutYAMLName);
368  if (!outFile) {
369  MACH3LOG_ERROR("Failed to open file for writing: {}", OutYAMLName);
370  throw MaCh3Exception(__FILE__, __LINE__);
371  }
372 
373  outFile << YAMLString;
374  outFile.close();
375 }
376 
377 // *************************************
382 inline TMacro* GetConfigMacroFromChain(TDirectory* CovarianceFolder) {
383 // *************************************
384  if (!CovarianceFolder) {
385  MACH3LOG_ERROR("Null TDirectory passed to {}", __func__);
386  throw MaCh3Exception(__FILE__, __LINE__);
387  }
388 
389  TMacro* foundMacro = nullptr;
390  int macroCount = 0;
391 
392  TIter next(CovarianceFolder->GetListOfKeys());
393  TKey* key;
394  while ((key = dynamic_cast<TKey*>(next()))) {
395  if (std::string(key->GetClassName()) == "TMacro") {
396  ++macroCount;
397  if (macroCount == 1) {
398  foundMacro = dynamic_cast<TMacro*>(key->ReadObj());
399  }
400  }
401  }
402 
403  if (macroCount == 1 && foundMacro) {
404  MACH3LOG_INFO("Found single TMacro in directory: using it.");
405  return foundMacro;
406  } else {
407  MACH3LOG_WARN("Found {} TMacro objects. Using hardcoded macro name: Config_xsec_cov.", macroCount);
408  TMacro* fallback = CovarianceFolder->Get<TMacro>("Config_xsec_cov");
409  if (!fallback) {
410  MACH3LOG_WARN("Fallback macro 'Config_xsec_cov' not found in directory.");
411  }
412  return fallback;
413  }
414 }
415 
416 // *************************************
423 inline TMatrixDSym* GetCovMatrixFromChain(TDirectory* TempFile) {
424 // *************************************
425  if (!TempFile) {
426  MACH3LOG_ERROR("Null TDirectory passed to {}.", __func__);
427  throw MaCh3Exception(__FILE__, __LINE__);
428  }
429 
430  TMatrixDSym* foundMatrix = nullptr;
431  int matrixCount = 0;
432 
433  TIter next(TempFile->GetListOfKeys());
434  TKey* key;
435  while ((key = dynamic_cast<TKey*>(next()))) {
436  std::string className = key->GetClassName();
437  if (className.find("TMatrix") != std::string::npos) {
438  ++matrixCount;
439  if (matrixCount == 1) {
440  foundMatrix = dynamic_cast<TMatrixDSym*>(key->ReadObj());
441  }
442  }
443  }
444 
445  if (matrixCount == 1 && foundMatrix) {
446  MACH3LOG_INFO("Found single TMatrixDSym in directory: using it.");
447  return foundMatrix;
448  } else {
449  MACH3LOG_WARN("Found {} TMatrixDSym objects. Using hardcoded path: xsec_cov.", matrixCount);
450  TMatrixDSym* fallback = TempFile->Get<TMatrixDSym>("xsec_cov");
451  if (!fallback) {
452  MACH3LOG_WARN("Fallback matrix 'xsec_cov' not found.");
453  }
454  return fallback;
455  }
456 }
457 
458 // *************************************
462 inline std::vector<std::vector<double>> GetCholeskyDecomposedMatrix(const TMatrixDSym& matrix, const std::string& matrixName) {
463 // *************************************
464  const Int_t n = matrix.GetNrows();
465  std::vector<std::vector<double>> L(n, std::vector<double>(n, 0.0));
466 
467  for (Int_t j = 0; j < n; ++j) {
468  // Compute diagonal element (must be serial)
469  double sum_diag = matrix(j, j);
470  for (Int_t k = 0; k < j; ++k) {
471  sum_diag -= L[j][k] * L[j][k];
472  }
473  const double tol = 1e-15;
474  if (sum_diag <= tol) {
475  MACH3LOG_ERROR("Cholesky decomposition failed for {} (non-positive diagonal)", matrixName);
476  throw MaCh3Exception(__FILE__, __LINE__);
477  }
478  L[j][j] = std::sqrt(sum_diag);
479 
480  // Compute the rest of the column in parallel
481  #ifdef MULTITHREAD
482  #pragma omp parallel for
483  #endif
484  for (Int_t i = j + 1; i < n; ++i) {
485  double sum = matrix(i, j);
486  for (Int_t k = 0; k < j; ++k) {
487  sum -= L[i][k] * L[j][k];
488  }
489  L[i][j] = sum / L[j][j];
490  }
491  }
492  return L;
493 }
494 
495 // *************************************
498 inline bool CanDecomposeMatrix(const TMatrixDSym& matrix) {
499 // *************************************
500  TDecompChol chdcmp(matrix);
501  return chdcmp.Decompose();
502 }
503 
504 // *************************************
506 inline void MakeMatrixPosDef(TMatrixDSym *cov) {
507 // *************************************
508  //DB Save original warning state and then increase it in this function to suppress 'matrix not positive definite' messages
509  //Means we no longer need to overload
510  int originalErrorWarning = gErrorIgnoreLevel;
511  gErrorIgnoreLevel = kFatal;
512 
513  //DB Loop 1000 times adding 1e-9 which tops out at 1e-6 shift on the diagonal before throwing error
514  constexpr int MaxAttempts = 1e5;
515  const int matrixSize = cov->GetNrows();
516  int iAttempt = 0;
517  bool CanDecomp = false;
518 
519  for (iAttempt = 0; iAttempt < MaxAttempts; iAttempt++) {
520  if (CanDecomposeMatrix(*cov)) {
521  CanDecomp = true;
522  break;
523  } else {
524  #ifdef MULTITHREAD
525  #pragma omp parallel for
526  #endif
527  for (int iVar = 0 ; iVar < matrixSize; iVar++) {
528  (*cov)(iVar,iVar) += pow(10, -9);
529  }
530  }
531  }
532 
533  if (!CanDecomp) {
534  MACH3LOG_ERROR("Tried {} times to shift diagonal but still can not decompose the matrix", MaxAttempts);
535  MACH3LOG_ERROR("This indicates that something is wrong with the input matrix");
536  throw MaCh3Exception(__FILE__ , __LINE__ );
537  }
538 
539  //DB Resetting warning level
540  gErrorIgnoreLevel = originalErrorWarning;
541 }
542 
543 } // end M3 namespace
#define _MaCh3_Safe_Include_Start_
KS: Avoiding warning checking for headers.
Definition: Core.h:109
#define _MaCh3_Safe_Include_End_
KS: Restore warning checking after including external headers.
Definition: Core.h:120
#define _restrict_
KS: Using restrict limits the effects of pointer aliasing, aiding optimizations. While reading I foun...
Definition: Core.h:93
#define MACH3LOG_ERROR
Definition: MaCh3Logger.h:27
#define MACH3LOG_INFO
Definition: MaCh3Logger.h:25
#define MACH3LOG_WARN
Definition: MaCh3Logger.h:26
std::string YAMLtoSTRING(const YAML::Node &node)
KS: Convert a YAML node to a string representation.
Definition: YamlHelper.h:107
Custom exception class for MaCh3 errors.
Definition: Core.h:19
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 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,...
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...
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...