MaCh3  2.5.0
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 
32 
33 namespace M3
34 {
36 inline double* MatrixMult(double *A, double *B, int n) {
37  //CW: First transpose to increse cache hits
38  double *BT = new double[n*n];
39  #ifdef MULTITHREAD
40  #pragma omp parallel for
41  #endif
42  for (int i = 0; i < n; i++) {
43  for (int j = 0; j < n; j++) {
44  BT[j*n+i] = B[i*n+j];
45  }
46  }
47 
48  // Now multiply
49  double *C = new double[n*n];
50  #ifdef MULTITHREAD
51  #pragma omp parallel for
52  #endif
53  for (int i = 0; i < n; i++) {
54  for (int j = 0; j < n; j++) {
55  double sum = 0;
56  for (int k = 0; k < n; k++) {
57  sum += A[i*n+k]*BT[j*n+k];
58  }
59  C[i*n+j] = sum;
60  }
61  }
62  delete BT;
63 
64  return C;
65 }
66 
68 inline double** MatrixMult(double **A, double **B, int n) {
69  // First make into monolithic array
70  double *A_mon = new double[n*n];
71  double *B_mon = new double[n*n];
72 
73  #ifdef MULTITHREAD
74  #pragma omp parallel for
75  #endif
76  for (int i = 0; i < n; ++i) {
77  for (int j = 0; j < n; ++j) {
78  A_mon[i*n+j] = A[i][j];
79  B_mon[i*n+j] = B[i][j];
80  }
81  }
82  //CW: Now call the monolithic calculator
83  double *C_mon = MatrixMult(A_mon, B_mon, n);
84  delete A_mon;
85  delete B_mon;
86 
87  // Return the double pointer
88  double **C = new double*[n];
89  #ifdef MULTITHREAD
90  #pragma omp parallel for
91  #endif
92  for (int i = 0; i < n; ++i) {
93  C[i] = new double[n];
94  for (int j = 0; j < n; ++j) {
95  C[i][j] = C_mon[i*n+j];
96  }
97  }
98  delete C_mon;
99 
100  return C;
101 }
102 
104 inline TMatrixD MatrixMult(TMatrixD A, TMatrixD B)
105 {
106  double *C_mon = MatrixMult(A.GetMatrixArray(), B.GetMatrixArray(), A.GetNcols());
107  TMatrixD C;
108  C.Use(A.GetNcols(), A.GetNrows(), C_mon);
109  return C;
110 }
111 
112 // ********************************************
118 inline void MatrixVectorMulti(double* _restrict_ VecMulti, double** _restrict_ matrix, const double* _restrict_ vector, const int n) {
119 // ********************************************
120  #ifdef MULTITHREAD
121  #pragma omp parallel for
122  #endif
123  for (int i = 0; i < n; ++i)
124  {
125  double result = 0.0;
126  #ifdef MULTITHREAD
127  #pragma omp simd
128  #endif
129  for (int j = 0; j < n; ++j)
130  {
131  result += matrix[i][j]*vector[j];
132  }
133  VecMulti[i] = result;
134  }
135 }
136 
137 // ********************************************
143 inline double MatrixVectorMultiSingle(double** _restrict_ matrix, const double* _restrict_ vector, const int Length, const int i) {
144 // ********************************************
145  double Element = 0.0;
146  #ifdef MULTITHREAD
147  #pragma omp simd
148  #endif
149  for (int j = 0; j < Length; ++j) {
150  Element += matrix[i][j]*vector[j];
151  }
152  return Element;
153 }
154 
155 // *************************************
158 inline void FixSampleNamesQuotes(std::string& yamlStr) {
159 // *************************************
160  std::stringstream input(yamlStr);
161  std::string line;
162  std::string fixedYaml;
163  std::regex sampleNamesRegex(R"(SampleNames:\s*\[([^\]]+)\])");
164 
165  while (std::getline(input, line)) {
166  std::smatch match;
167  if (std::regex_search(line, match, sampleNamesRegex)) {
168  std::string contents = match[1]; // inside the brackets
169  std::stringstream ss(contents);
170  std::string item;
171  std::vector<std::string> quotedItems;
172 
173  while (std::getline(ss, item, ',')) {
174  item = std::regex_replace(item, std::regex(R"(^\s+|\s+$)"), ""); // trim
175  quotedItems.push_back("\"" + item + "\"");
176  }
177 
178  std::string replacement = "SampleNames: [" + fmt::format("{}", fmt::join(quotedItems, ", ")) + "]";
179  line = std::regex_replace(line, sampleNamesRegex, replacement);
180  }
181  fixedYaml += line + "\n";
182  }
183 
184  yamlStr = fixedYaml;
185 }
186 
187 // *************************************
193 inline void AddTuneValues(YAML::Node& root,
194  const std::vector<double>& Values,
195  const std::string& Tune,
196  const std::vector<std::string>& FancyNames = {}) {
197 // *************************************
198  YAML::Node NodeCopy = YAML::Clone(root);
199  YAML::Node systematics = NodeCopy["Systematics"];
200 
201  if (!systematics || !systematics.IsSequence()) {
202  MACH3LOG_ERROR("'Systematics' node is missing or not a sequence in the YAML copy");
203  throw MaCh3Exception(__FILE__, __LINE__);
204  }
205 
206  if (!FancyNames.empty() && FancyNames.size() != Values.size()) {
207  MACH3LOG_ERROR("Mismatch in sizes: FancyNames has {}, but Values has {}", FancyNames.size(), Values.size());
208  throw MaCh3Exception(__FILE__, __LINE__);
209  }
210 
211  if (FancyNames.empty() && systematics.size() != Values.size()) {
212  MACH3LOG_ERROR("Mismatch in sizes: Values has {}, but YAML 'Systematics' has {} entries",
213  Values.size(), systematics.size());
214  throw MaCh3Exception(__FILE__, __LINE__);
215  }
216 
217  if (!FancyNames.empty()) {
218  for (std::size_t i = 0; i < FancyNames.size(); ++i) {
219  bool matched = false;
220  for (std::size_t j = 0; j < systematics.size(); ++j) {
221  YAML::Node systematicNode = systematics[j]["Systematic"];
222  if (!systematicNode) continue;
223  auto nameNode = systematicNode["Names"];
224  if (!nameNode || !nameNode["FancyName"]) continue;
225  if (nameNode["FancyName"].as<std::string>() == FancyNames[i]) {
226  if (!systematicNode["ParameterValues"]) {
227  MACH3LOG_ERROR("Missing 'ParameterValues' for matched FancyName '{}'", FancyNames[i]);
228  throw MaCh3Exception(__FILE__, __LINE__);
229  }
230  systematicNode["ParameterValues"][Tune] = M3::Utils::FormatDouble(Values[i], 4);
231  matched = true;
232  break;
233  }
234  }
235  if (!matched) {
236  MACH3LOG_ERROR("Could not find a matching FancyName '{}' in the systematics", FancyNames[i]);
237  throw MaCh3Exception(__FILE__, __LINE__);
238  }
239  }
240  } else {
241  for (std::size_t i = 0; i < systematics.size(); ++i) {
242  YAML::Node systematicNode = systematics[i]["Systematic"];
243  if (!systematicNode || !systematicNode["ParameterValues"]) {
244  MACH3LOG_ERROR("Missing 'Systematic' or 'ParameterValues' entry at index {}", i);
245  throw MaCh3Exception(__FILE__, __LINE__);
246  }
247  systematicNode["ParameterValues"][Tune] = M3::Utils::FormatDouble(Values[i], 4);
248  }
249  }
250 
251  // Convert updated copy to string
252  std::string YAMLString = YAMLtoSTRING(NodeCopy);
253  FixSampleNamesQuotes(YAMLString);
254  // Write to output file
255  std::string OutName = "UpdatedMatrixWithTune" + Tune + ".yaml";
256  std::ofstream outFile(OutName);
257  if (!outFile) {
258  MACH3LOG_ERROR("Failed to open file for writing: {}", OutName);
259  throw MaCh3Exception(__FILE__, __LINE__);
260  }
261 
262  outFile << YAMLString;
263  outFile.close();
264 }
265 
266 // *************************************
274 inline void MakeCorrelationMatrix(YAML::Node& root,
275  const std::vector<double>& Values,
276  const std::vector<double>& Errors,
277  const std::vector<std::vector<double>>& Correlation,
278  const std::string& OutYAMLName,
279  const std::vector<std::string>& FancyNames = {}) {
280 // *************************************
281  if (Values.size() != Errors.size() || Values.size() != Correlation.size()) {
282  MACH3LOG_ERROR("Size mismatch between Values, Errors, and Correlation matrix");
283  throw MaCh3Exception(__FILE__, __LINE__);
284  }
285 
286  for (const auto& row : Correlation) {
287  if (row.size() != Correlation.size()) {
288  MACH3LOG_ERROR("Correlation matrix is not square");
289  throw MaCh3Exception(__FILE__, __LINE__);
290  }
291  }
292 
293  YAML::Node NodeCopy = YAML::Clone(root);
294  YAML::Node systematics = NodeCopy["Systematics"];
295 
296  if (!systematics || !systematics.IsSequence()) {
297  MACH3LOG_ERROR("'Systematics' node is missing or not a sequence");
298  throw MaCh3Exception(__FILE__, __LINE__);
299  }
300 
301  if (!FancyNames.empty() && FancyNames.size() != Values.size()) {
302  MACH3LOG_ERROR("FancyNames size ({}) does not match Values size ({})", FancyNames.size(), Values.size());
303  throw MaCh3Exception(__FILE__, __LINE__);
304  }
305 
306  // Map from FancyName to Systematic node
307  std::unordered_map<std::string, YAML::Node> nameToNode;
308  for (std::size_t i = 0; i < systematics.size(); ++i) {
309  YAML::Node syst = systematics[i]["Systematic"];
310  if (!syst || !syst["Names"] || !syst["Names"]["FancyName"]) continue;
311  std::string name = syst["Names"]["FancyName"].as<std::string>();
312  nameToNode[name] = syst;
313  }
314 
315  if (!FancyNames.empty()) {
316  for (std::size_t i = 0; i < FancyNames.size(); ++i) {
317  const std::string& name_i = FancyNames[i];
318  auto it_i = nameToNode.find(name_i);
319  if (it_i == nameToNode.end()) {
320  MACH3LOG_ERROR("Could not find FancyName '{}' in YAML", name_i);
321  throw MaCh3Exception(__FILE__, __LINE__);
322  }
323  YAML::Node& syst_i = it_i->second;
324 
325  syst_i["ParameterValues"]["PreFitValue"] = M3::Utils::FormatDouble(Values[i], 4);
326  syst_i["Error"] = std::round(Errors[i] * 100.0) / 100.0;
327 
328  YAML::Node correlationsNode = YAML::Node(YAML::NodeType::Sequence);
329  for (std::size_t j = 0; j < FancyNames.size(); ++j) {
330  if (i == j) continue;
331  // KS: Skip if value close to 0
332  if (std::abs(Correlation[i][j]) < 1e-8) continue;
333  YAML::Node singleEntry;
334  singleEntry[FancyNames[j]] = M3::Utils::FormatDouble(Correlation[i][j], 4);
335  correlationsNode.push_back(singleEntry);
336  }
337  syst_i["Correlations"] = correlationsNode;
338  }
339  } else {
340  if (systematics.size() != Values.size()) {
341  MACH3LOG_ERROR("Mismatch in sizes: Values has {}, but YAML 'Systematics' has {} entries",
342  Values.size(), systematics.size());
343  throw MaCh3Exception(__FILE__, __LINE__);
344  }
345 
346  for (std::size_t i = 0; i < systematics.size(); ++i) {
347  YAML::Node syst = systematics[i]["Systematic"];
348  if (!syst) {
349  MACH3LOG_ERROR("Missing 'Systematic' node at index {}", i);
350  throw MaCh3Exception(__FILE__, __LINE__);
351  }
352 
353  syst["ParameterValues"]["PreFitValue"] = M3::Utils::FormatDouble(Values[i], 4);
354  syst["Error"] = std::round(Errors[i] * 100.0) / 100.0;
355 
356  YAML::Node correlationsNode = YAML::Node(YAML::NodeType::Sequence);
357  for (std::size_t j = 0; j < Correlation[i].size(); ++j) {
358  if (i == j) continue;
359  // KS: Skip if value close to 0
360  if (std::abs(Correlation[i][j]) < 1e-8) continue;
361  YAML::Node singleEntry;
362  const std::string& otherName = systematics[j]["Systematic"]["Names"]["FancyName"].as<std::string>();
363  singleEntry[otherName] = M3::Utils::FormatDouble(Correlation[i][j], 4);
364  correlationsNode.push_back(singleEntry);
365  }
366  syst["Correlations"] = correlationsNode;
367  }
368  }
369 
370  // Convert and write
371  std::string YAMLString = YAMLtoSTRING(NodeCopy);
372  FixSampleNamesQuotes(YAMLString);
373  std::ofstream outFile(OutYAMLName);
374  if (!outFile) {
375  MACH3LOG_ERROR("Failed to open file for writing: {}", OutYAMLName);
376  throw MaCh3Exception(__FILE__, __LINE__);
377  }
378 
379  outFile << YAMLString;
380  outFile.close();
381 }
382 
383 // *************************************
388 inline TMacro* GetConfigMacroFromChain(TDirectory* CovarianceFolder) {
389 // *************************************
390  if (!CovarianceFolder) {
391  MACH3LOG_ERROR("Null TDirectory passed to {}", __func__);
392  throw MaCh3Exception(__FILE__, __LINE__);
393  }
394 
395  TMacro* foundMacro = nullptr;
396  int macroCount = 0;
397 
398  TIter next(CovarianceFolder->GetListOfKeys());
399  TKey* key;
400  while ((key = dynamic_cast<TKey*>(next()))) {
401  if (std::string(key->GetClassName()) == "TMacro") {
402  ++macroCount;
403  if (macroCount == 1) {
404  foundMacro = dynamic_cast<TMacro*>(key->ReadObj());
405  }
406  }
407  }
408 
409  if (macroCount == 1 && foundMacro) {
410  MACH3LOG_INFO("Found single TMacro in directory: using it.");
411  return foundMacro;
412  } else {
413  MACH3LOG_WARN("Found {} TMacro objects. Using hardcoded macro name: Config_xsec_cov.", macroCount);
414  TMacro* fallback = CovarianceFolder->Get<TMacro>("Config_xsec_cov");
415  if (!fallback) {
416  MACH3LOG_WARN("Fallback macro 'Config_xsec_cov' not found in directory.");
417  }
418  return fallback;
419  }
420 }
421 
422 // *************************************
429 inline TMatrixDSym* GetCovMatrixFromChain(TDirectory* TempFile) {
430 // *************************************
431  if (!TempFile) {
432  MACH3LOG_ERROR("Null TDirectory passed to {}.", __func__);
433  throw MaCh3Exception(__FILE__, __LINE__);
434  }
435 
436  TMatrixDSym* foundMatrix = nullptr;
437  int matrixCount = 0;
438 
439  TIter next(TempFile->GetListOfKeys());
440  TKey* key;
441  while ((key = dynamic_cast<TKey*>(next()))) {
442  std::string className = key->GetClassName();
443  if (className.find("TMatrix") != std::string::npos) {
444  ++matrixCount;
445  if (matrixCount == 1) {
446  foundMatrix = dynamic_cast<TMatrixDSym*>(key->ReadObj());
447  }
448  }
449  }
450 
451  if (matrixCount == 1 && foundMatrix) {
452  MACH3LOG_INFO("Found single TMatrixDSym in directory: using it.");
453  return foundMatrix;
454  } else {
455  MACH3LOG_WARN("Found {} TMatrixDSym objects. Using hardcoded path: xsec_cov.", matrixCount);
456  TMatrixDSym* fallback = TempFile->Get<TMatrixDSym>("xsec_cov");
457  if (!fallback) {
458  MACH3LOG_WARN("Fallback matrix 'xsec_cov' not found.");
459  }
460  return fallback;
461  }
462 }
463 
464 // *************************************
468 inline std::vector<std::vector<double>> GetCholeskyDecomposedMatrix(const TMatrixDSym& matrix, const std::string& matrixName) {
469 // *************************************
470  const Int_t n = matrix.GetNrows();
471  std::vector<std::vector<double>> L(n, std::vector<double>(n, 0.0));
472 
473  for (Int_t j = 0; j < n; ++j) {
474  // Compute diagonal element (must be serial)
475  double sum_diag = matrix(j, j);
476  for (Int_t k = 0; k < j; ++k) {
477  sum_diag -= L[j][k] * L[j][k];
478  }
479  const double tol = 1e-15;
480  if (sum_diag <= tol) {
481  MACH3LOG_ERROR("Cholesky decomposition failed for {} (non-positive diagonal)", matrixName);
482  throw MaCh3Exception(__FILE__, __LINE__);
483  }
484  L[j][j] = std::sqrt(sum_diag);
485 
486  // Compute the rest of the column in parallel
487  #ifdef MULTITHREAD
488  #pragma omp parallel for
489  #endif
490  for (Int_t i = j + 1; i < n; ++i) {
491  double sum = matrix(i, j);
492  for (Int_t k = 0; k < j; ++k) {
493  sum -= L[i][k] * L[j][k];
494  }
495  L[i][j] = sum / L[j][j];
496  }
497  }
498  return L;
499 }
500 
501 // *************************************
504 inline bool CanDecomposeMatrix(const TMatrixDSym& matrix) {
505 // *************************************
506  TDecompChol chdcmp(matrix);
507  return chdcmp.Decompose();
508 }
509 
510 // *************************************
512 inline void MakeMatrixPosDef(TMatrixDSym *cov) {
513 // *************************************
514  //DB Save original warning state and then increase it in this function to suppress 'matrix not positive definite' messages
515  //Means we no longer need to overload
516  int originalErrorWarning = gErrorIgnoreLevel;
517  gErrorIgnoreLevel = kFatal;
518 
519  //DB Loop 1000 times adding 1e-9 which tops out at 1e-6 shift on the diagonal before throwing error
520  constexpr int MaxAttempts = 1e5;
521  const int matrixSize = cov->GetNrows();
522  int iAttempt = 0;
523  bool CanDecomp = false;
524 
525  for (iAttempt = 0; iAttempt < MaxAttempts; iAttempt++) {
526  if (CanDecomposeMatrix(*cov)) {
527  CanDecomp = true;
528  break;
529  } else {
530  #ifdef MULTITHREAD
531  #pragma omp parallel for
532  #endif
533  for (int iVar = 0 ; iVar < matrixSize; iVar++) {
534  (*cov)(iVar, iVar) += 1e-9;
535  }
536  }
537  }
538 
539  if (!CanDecomp) {
540  MACH3LOG_ERROR("Tried {} times to shift diagonal but still can not decompose the matrix", MaxAttempts);
541  MACH3LOG_ERROR("This indicates that something is wrong with the input matrix");
542  throw MaCh3Exception(__FILE__ , __LINE__ );
543  }
544 
545  //DB Resetting warning level
546  gErrorIgnoreLevel = originalErrorWarning;
547 }
548 
549 
550 
551 // ********************************************
554 inline void DumpParamHandlerToFile(const int _fNumPar,
555  const std::vector<double>& _fPreFitValue,
556  const std::vector<double>& _fError,
557  const std::vector<double>& _fLowBound,
558  const std::vector<double>& _fUpBound,
559  const std::vector<double>& _fIndivStepScale,
560  const std::vector<std::string>& _fFancyNames,
561  const std::vector<bool>& _fFlatPrior,
562  const std::vector<SplineParameter>& SplineParams,
563  TMatrixDSym* covMatrix,
564  TH2D* CorrMatrix,
565  const std::string& Name) {
566 // ********************************************
567  TFile* outputFile = new TFile(Name.c_str(), "RECREATE");
568 
569  TObjArray* param_names = new TObjArray();
570  TObjArray* spline_interpolation = new TObjArray();
571  TObjArray* spline_names = new TObjArray();
572 
573  TVectorD* param_prior = new TVectorD(_fNumPar);
574  TVectorD* flat_prior = new TVectorD(_fNumPar);
575  TVectorD* stepscale = new TVectorD(_fNumPar);
576  TVectorD* param_lb = new TVectorD(_fNumPar);
577  TVectorD* param_ub = new TVectorD(_fNumPar);
578 
579  TVectorD* param_knot_weight_lb = new TVectorD(_fNumPar);
580  TVectorD* param_knot_weight_ub = new TVectorD(_fNumPar);
581  TVectorD* error = new TVectorD(_fNumPar);
582 
583  for(int i = 0; i < _fNumPar; ++i)
584  {
585  TObjString* nameObj = new TObjString(_fFancyNames[i].c_str());
586  param_names->AddLast(nameObj);
587 
588  TObjString* splineType = new TObjString("TSpline3");
589  spline_interpolation->AddLast(splineType);
590 
591  TObjString* splineName = new TObjString("");
592  spline_names->AddLast(splineName);
593 
594  (*param_prior)[i] = _fPreFitValue[i];
595  (*flat_prior)[i] = _fFlatPrior[i];
596  (*stepscale)[i] = _fIndivStepScale[i];
597  (*error)[i] = _fError[i];
598 
599  (*param_lb)[i] = _fLowBound[i];
600  (*param_ub)[i] = _fUpBound[i];
601 
602  //Default values
603  (*param_knot_weight_lb)[i] = -9999;
604  (*param_knot_weight_ub)[i] = +9999;
605  }
606 
607  for (size_t SplineIndex = 0; SplineIndex < SplineParams.size(); SplineIndex++ ) {
608  auto SystIndex = SplineParams[SplineIndex].index;
609 
610  (*param_knot_weight_lb)[SystIndex] = SplineParams.at(SplineIndex)._SplineKnotLowBound;
611  (*param_knot_weight_ub)[SystIndex] = SplineParams.at(SplineIndex)._SplineKnotUpBound;
612 
613  TObjString* splineType = new TObjString(SplineInterpolation_ToString(SplineParams.at(SplineIndex)._SplineInterpolationType).c_str());
614  spline_interpolation->AddAt(splineType, SystIndex);
615 
616  TObjString* splineName = new TObjString(SplineParams[SplineIndex]._fSplineNames.c_str());
617  spline_names->AddAt(splineName, SystIndex);
618  }
619  param_names->Write("xsec_param_names", TObject::kSingleKey);
620  delete param_names;
621  spline_interpolation->Write("xsec_spline_interpolation", TObject::kSingleKey);
622  delete spline_interpolation;
623  spline_names->Write("xsec_spline_names", TObject::kSingleKey);
624  delete spline_names;
625 
626  param_prior->Write("xsec_param_prior");
627  delete param_prior;
628  flat_prior->Write("xsec_flat_prior");
629  delete flat_prior;
630  stepscale->Write("xsec_stepscale");
631  delete stepscale;
632  param_lb->Write("xsec_param_lb");
633  delete param_lb;
634  param_ub->Write("xsec_param_ub");
635  delete param_ub;
636 
637  param_knot_weight_lb->Write("xsec_param_knot_weight_lb");
638  delete param_knot_weight_lb;
639  param_knot_weight_ub->Write("xsec_param_knot_weight_ub");
640  delete param_knot_weight_ub;
641  error->Write("xsec_error");
642  delete error;
643 
644  covMatrix->Write("xsec_cov");
645  CorrMatrix->Write("hcov");
646 
647  outputFile->Close();
648  delete outputFile;
649 }
650 
651 } // end M3 namespace
#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_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,...
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...