MaCh3 2.2.1
Reference Guide
Loading...
Searching...
No Matches
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
27namespace M3
28{
30inline 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
62inline 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
98inline 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// ********************************************
112inline 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// ********************************************
136inline double MatrixVectorMultiSingle(double** _restrict_ matrix, const double* _restrict_ vector, const int Length, const int i) {
137// ********************************************
138 double Element = 0.0;
139 #ifdef MULTITHREAD
140 #pragma omp simd
141 #endif
142 for (int j = 0; j < Length; ++j) {
143 Element += matrix[i][j]*vector[j];
144 }
145 return Element;
146}
147
148// *************************************
150inline void FixSampleNamesQuotes(std::string& yamlStr) {
151// *************************************
152 std::stringstream input(yamlStr);
153 std::string line;
154 std::string fixedYaml;
155 std::regex sampleNamesRegex(R"(SampleNames:\s*\[([^\]]+)\])");
156
157 while (std::getline(input, line)) {
158 std::smatch match;
159 if (std::regex_search(line, match, sampleNamesRegex)) {
160 std::string contents = match[1]; // inside the brackets
161 std::stringstream ss(contents);
162 std::string item;
163 std::vector<std::string> quotedItems;
164
165 while (std::getline(ss, item, ',')) {
166 item = std::regex_replace(item, std::regex(R"(^\s+|\s+$)"), ""); // trim
167 quotedItems.push_back("\"" + item + "\"");
168 }
169
170 std::string replacement = "SampleNames: [" + fmt::format("{}", fmt::join(quotedItems, ", ")) + "]";
171 line = std::regex_replace(line, sampleNamesRegex, replacement);
172 }
173 fixedYaml += line + "\n";
174 }
175
176 yamlStr = fixedYaml;
177}
178
179// *************************************
181inline void AddTuneValues(YAML::Node& root,
182 const std::vector<double>& Values,
183 const std::string& Tune,
184 const std::vector<std::string>& FancyNames = {}) {
185// *************************************
186 YAML::Node NodeCopy = YAML::Clone(root);
187 YAML::Node systematics = NodeCopy["Systematics"];
188
189 if (!systematics || !systematics.IsSequence()) {
190 MACH3LOG_ERROR("'Systematics' node is missing or not a sequence in the YAML copy");
191 throw MaCh3Exception(__FILE__, __LINE__);
192 }
193
194 if (!FancyNames.empty() && FancyNames.size() != Values.size()) {
195 MACH3LOG_ERROR("Mismatch in sizes: FancyNames has {}, but Values has {}", FancyNames.size(), Values.size());
196 throw MaCh3Exception(__FILE__, __LINE__);
197 }
198
199 if (FancyNames.empty() && systematics.size() != Values.size()) {
200 MACH3LOG_ERROR("Mismatch in sizes: Values has {}, but YAML 'Systematics' has {} entries",
201 Values.size(), systematics.size());
202 throw MaCh3Exception(__FILE__, __LINE__);
203 }
204
205 if (!FancyNames.empty()) {
206 for (std::size_t i = 0; i < FancyNames.size(); ++i) {
207 bool matched = false;
208 for (std::size_t j = 0; j < systematics.size(); ++j) {
209 YAML::Node systematicNode = systematics[j]["Systematic"];
210 if (!systematicNode) continue;
211 auto nameNode = systematicNode["Names"];
212 if (!nameNode || !nameNode["FancyName"]) continue;
213 if (nameNode["FancyName"].as<std::string>() == FancyNames[i]) {
214 if (!systematicNode["ParameterValues"]) {
215 MACH3LOG_ERROR("Missing 'ParameterValues' for matched FancyName '{}'", FancyNames[i]);
216 throw MaCh3Exception(__FILE__, __LINE__);
217 }
218 systematicNode["ParameterValues"][Tune] = MaCh3Utils::FormatDouble(Values[i], 4);
219 matched = true;
220 break;
221 }
222 }
223 if (!matched) {
224 MACH3LOG_ERROR("Could not find a matching FancyName '{}' in the systematics", FancyNames[i]);
225 throw MaCh3Exception(__FILE__, __LINE__);
226 }
227 }
228 } else {
229 for (std::size_t i = 0; i < systematics.size(); ++i) {
230 YAML::Node systematicNode = systematics[i]["Systematic"];
231 if (!systematicNode || !systematicNode["ParameterValues"]) {
232 MACH3LOG_ERROR("Missing 'Systematic' or 'ParameterValues' entry at index {}", i);
233 throw MaCh3Exception(__FILE__, __LINE__);
234 }
235 systematicNode["ParameterValues"][Tune] = MaCh3Utils::FormatDouble(Values[i], 4);
236 }
237 }
238
239 // Convert updated copy to string
240 std::string YAMLString = YAMLtoSTRING(NodeCopy);
241 FixSampleNamesQuotes(YAMLString);
242 // Write to output file
243 std::string OutName = "UpdatedMatrixWithTune" + Tune + ".yaml";
244 std::ofstream outFile(OutName);
245 if (!outFile) {
246 MACH3LOG_ERROR("Failed to open file for writing: {}", OutName);
247 throw MaCh3Exception(__FILE__, __LINE__);
248 }
249
250 outFile << YAMLString;
251 outFile.close();
252}
253
254// *************************************
256inline void MakeCorrelationMatrix(YAML::Node& root,
257 const std::vector<double>& Values,
258 const std::vector<double>& Errors,
259 const std::vector<std::vector<double>>& Correlation,
260 const std::vector<std::string>& FancyNames = {}) {
261// *************************************
262 if (Values.size() != Errors.size() || Values.size() != Correlation.size()) {
263 MACH3LOG_ERROR("Size mismatch between Values, Errors, and Correlation matrix");
264 throw MaCh3Exception(__FILE__, __LINE__);
265 }
266
267 for (const auto& row : Correlation) {
268 if (row.size() != Correlation.size()) {
269 MACH3LOG_ERROR("Correlation matrix is not square");
270 throw MaCh3Exception(__FILE__, __LINE__);
271 }
272 }
273
274 YAML::Node NodeCopy = YAML::Clone(root);
275 YAML::Node systematics = NodeCopy["Systematics"];
276
277 if (!systematics || !systematics.IsSequence()) {
278 MACH3LOG_ERROR("'Systematics' node is missing or not a sequence");
279 throw MaCh3Exception(__FILE__, __LINE__);
280 }
281
282 if (!FancyNames.empty() && FancyNames.size() != Values.size()) {
283 MACH3LOG_ERROR("FancyNames size ({}) does not match Values size ({})", FancyNames.size(), Values.size());
284 throw MaCh3Exception(__FILE__, __LINE__);
285 }
286
287 // Map from FancyName to Systematic node
288 std::unordered_map<std::string, YAML::Node> nameToNode;
289 for (std::size_t i = 0; i < systematics.size(); ++i) {
290 YAML::Node syst = systematics[i]["Systematic"];
291 if (!syst || !syst["Names"] || !syst["Names"]["FancyName"]) continue;
292 std::string name = syst["Names"]["FancyName"].as<std::string>();
293 nameToNode[name] = syst;
294 }
295
296 if (!FancyNames.empty()) {
297 for (std::size_t i = 0; i < FancyNames.size(); ++i) {
298 const std::string& name_i = FancyNames[i];
299 auto it_i = nameToNode.find(name_i);
300 if (it_i == nameToNode.end()) {
301 MACH3LOG_ERROR("Could not find FancyName '{}' in YAML", name_i);
302 throw MaCh3Exception(__FILE__, __LINE__);
303 }
304 YAML::Node& syst_i = it_i->second;
305
306 syst_i["ParameterValues"]["PreFitValue"] = MaCh3Utils::FormatDouble(Values[i], 4);
307 syst_i["Error"] = std::round(Errors[i] * 100.0) / 100.0;
308
309 YAML::Node correlationsNode = YAML::Node(YAML::NodeType::Sequence);
310 for (std::size_t j = 0; j < FancyNames.size(); ++j) {
311 if (i == j) continue;
312 YAML::Node singleEntry;
313 singleEntry[FancyNames[j]] = MaCh3Utils::FormatDouble(Correlation[i][j], 4);
314 correlationsNode.push_back(singleEntry);
315 }
316 syst_i["Correlations"] = correlationsNode;
317 }
318 } else {
319 if (systematics.size() != Values.size()) {
320 MACH3LOG_ERROR("Mismatch in sizes: Values has {}, but YAML 'Systematics' has {} entries",
321 Values.size(), systematics.size());
322 throw MaCh3Exception(__FILE__, __LINE__);
323 }
324
325 for (std::size_t i = 0; i < systematics.size(); ++i) {
326 YAML::Node syst = systematics[i]["Systematic"];
327 if (!syst) {
328 MACH3LOG_ERROR("Missing 'Systematic' node at index {}", i);
329 throw MaCh3Exception(__FILE__, __LINE__);
330 }
331
332 syst["ParameterValues"]["PreFitValue"] = MaCh3Utils::FormatDouble(Values[i], 4);
333 syst["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 < Correlation[i].size(); ++j) {
337 if (i == j) continue;
338 YAML::Node singleEntry;
339 const std::string& otherName = systematics[j]["Systematic"]["Names"]["FancyName"].as<std::string>();
340 singleEntry[otherName] = MaCh3Utils::FormatDouble(Correlation[i][j], 4);
341 correlationsNode.push_back(singleEntry);
342 }
343 syst["Correlations"] = correlationsNode;
344 }
345 }
346
347 // Convert and write
348 std::string YAMLString = YAMLtoSTRING(NodeCopy);
349 FixSampleNamesQuotes(YAMLString);
350 std::string OutName = "UpdatedCorrelationMatrix.yaml";
351 std::ofstream outFile(OutName);
352 if (!outFile) {
353 MACH3LOG_ERROR("Failed to open file for writing: {}", OutName);
354 throw MaCh3Exception(__FILE__, __LINE__);
355 }
356
357 outFile << YAMLString;
358 outFile.close();
359}
360
361// *************************************
366inline TMacro* GetConfigMacroFromChain(TDirectory* CovarianceFolder) {
367// *************************************
368 if (!CovarianceFolder) {
369 MACH3LOG_ERROR("Null TDirectory passed to {}", __func__);
370 throw MaCh3Exception(__FILE__, __LINE__);
371 }
372
373 TMacro* foundMacro = nullptr;
374 int macroCount = 0;
375
376 TIter next(CovarianceFolder->GetListOfKeys());
377 TKey* key;
378 while ((key = dynamic_cast<TKey*>(next()))) {
379 if (std::string(key->GetClassName()) == "TMacro") {
380 ++macroCount;
381 if (macroCount == 1) {
382 foundMacro = dynamic_cast<TMacro*>(key->ReadObj());
383 }
384 }
385 }
386
387 if (macroCount == 1 && foundMacro) {
388 MACH3LOG_INFO("Found single TMacro in directory: using it.");
389 return foundMacro;
390 } else {
391 MACH3LOG_WARN("Found {} TMacro objects. Using hardcoded macro name: Config_xsec_cov.", macroCount);
392 TMacro* fallback = CovarianceFolder->Get<TMacro>("Config_xsec_cov");
393 if (!fallback) {
394 MACH3LOG_WARN("Fallback macro 'Config_xsec_cov' not found in directory.");
395 }
396 return fallback;
397 }
398}
399
400// *************************************
407inline TMatrixDSym* GetCovMatrixFromChain(TDirectory* TempFile) {
408// *************************************
409 if (!TempFile) {
410 MACH3LOG_ERROR("Null TDirectory passed to {}.", __func__);
411 throw MaCh3Exception(__FILE__, __LINE__);
412 }
413
414 TMatrixDSym* foundMatrix = nullptr;
415 int matrixCount = 0;
416
417 TIter next(TempFile->GetListOfKeys());
418 TKey* key;
419 while ((key = dynamic_cast<TKey*>(next()))) {
420 std::string className = key->GetClassName();
421 if (className.find("TMatrix") != std::string::npos) {
422 ++matrixCount;
423 if (matrixCount == 1) {
424 foundMatrix = dynamic_cast<TMatrixDSym*>(key->ReadObj());
425 }
426 }
427 }
428
429 if (matrixCount == 1 && foundMatrix) {
430 MACH3LOG_INFO("Found single TMatrixDSym in directory: using it.");
431 return foundMatrix;
432 } else {
433 MACH3LOG_WARN("Found {} TMatrixDSym objects. Using hardcoded path: xsec_cov.", matrixCount);
434 TMatrixDSym* fallback = TempFile->Get<TMatrixDSym>("xsec_cov");
435 if (!fallback) {
436 MACH3LOG_WARN("Fallback matrix 'xsec_cov' not found.");
437 }
438 return fallback;
439 }
440}
441
442}
#define _MaCh3_Safe_Include_Start_
KS: Avoiding warning checking for headers.
Definition: Core.h:106
#define _MaCh3_Safe_Include_End_
KS: Restore warning checking after including external headers.
Definition: Core.h:117
#define _restrict_
KS: Using restrict limits the effects of pointer aliasing, aiding optimizations. While reading I foun...
Definition: Core.h:90
#define MACH3LOG_ERROR
Definition: MaCh3Logger.h:25
#define MACH3LOG_INFO
Definition: MaCh3Logger.h:23
#define MACH3LOG_WARN
Definition: MaCh3Logger.h:24
std::string YAMLtoSTRING(const YAML::Node &node)
KS: Convert a YAML node to a string representation.
Definition: YamlHelper.h:106
Custom exception class for MaCh3 errors.
Definition: Core.h:19
TMacro * GetConfigMacroFromChain(TDirectory *CovarianceFolder)
KS: We store configuration macros inside the chain. In the past, multiple configs were stored,...
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.
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::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.
double * MatrixMult(double *A, double *B, int n)
CW: Multi-threaded matrix multiplication.
TMatrixDSym * GetCovMatrixFromChain(TDirectory *TempFile)
KS: Retrieve the cross-section covariance matrix from the given TDirectory. Historically,...
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...