9#include "TMatrixDSym.h"
18#include "TDecompChol.h"
19#include "TStopwatch.h"
21#include "TMatrixDSymEigen.h"
22#include "TMatrixDEigen.h"
23#include "TDecompSVD.h"
32 double *BT =
new double[n*n];
34 #pragma omp parallel for
36 for (
int i = 0; i < n; i++) {
37 for (
int j = 0; j < n; j++) {
43 double *C =
new double[n*n];
45 #pragma omp parallel for
47 for (
int i = 0; i < n; i++) {
48 for (
int j = 0; j < n; j++) {
50 for (
int k = 0; k < n; k++) {
51 sum += A[i*n+k]*BT[j*n+k];
62inline double**
MatrixMult(
double **A,
double **B,
int n) {
64 double *A_mon =
new double[n*n];
65 double *B_mon =
new double[n*n];
68 #pragma omp parallel for
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];
82 double **C =
new double*[n];
84 #pragma omp parallel for
86 for (
int i = 0; i < n; ++i) {
88 for (
int j = 0; j < n; ++j) {
89 C[i][j] = C_mon[i*n+j];
100 double *C_mon =
MatrixMult(A.GetMatrixArray(), B.GetMatrixArray(), A.GetNcols());
102 C.Use(A.GetNcols(), A.GetNrows(), C_mon);
115 #pragma omp parallel for
117 for (
int i = 0; i < n; ++i)
123 for (
int j = 0; j < n; ++j)
125 result += matrix[i][j]*vector[j];
127 VecMulti[i] = result;
138 double Element = 0.0;
142 for (
int j = 0; j < Length; ++j) {
143 Element += matrix[i][j]*vector[j];
152 std::stringstream input(yamlStr);
154 std::string fixedYaml;
155 std::regex sampleNamesRegex(R
"(SampleNames:\s*\[([^\]]+)\])");
157 while (std::getline(input, line)) {
159 if (std::regex_search(line, match, sampleNamesRegex)) {
160 std::string contents = match[1];
161 std::stringstream ss(contents);
163 std::vector<std::string> quotedItems;
165 while (std::getline(ss, item,
',')) {
166 item = std::regex_replace(item, std::regex(R
"(^\s+|\s+$)"), "");
167 quotedItems.push_back(
"\"" + item +
"\"");
170 std::string replacement =
"SampleNames: [" + fmt::format(
"{}", fmt::join(quotedItems,
", ")) +
"]";
171 line = std::regex_replace(line, sampleNamesRegex, replacement);
173 fixedYaml += line +
"\n";
182 const std::vector<double>& Values,
183 const std::string& Tune,
184 const std::vector<std::string>& FancyNames = {}) {
186 YAML::Node NodeCopy = YAML::Clone(root);
187 YAML::Node systematics = NodeCopy[
"Systematics"];
189 if (!systematics || !systematics.IsSequence()) {
190 MACH3LOG_ERROR(
"'Systematics' node is missing or not a sequence in the YAML copy");
194 if (!FancyNames.empty() && FancyNames.size() != Values.size()) {
195 MACH3LOG_ERROR(
"Mismatch in sizes: FancyNames has {}, but Values has {}", FancyNames.size(), Values.size());
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());
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]);
224 MACH3LOG_ERROR(
"Could not find a matching FancyName '{}' in the systematics", FancyNames[i]);
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);
243 std::string OutName =
"UpdatedMatrixWithTune" + Tune +
".yaml";
244 std::ofstream outFile(OutName);
250 outFile << YAMLString;
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 = {}) {
262 if (Values.size() != Errors.size() || Values.size() != Correlation.size()) {
263 MACH3LOG_ERROR(
"Size mismatch between Values, Errors, and Correlation matrix");
267 for (
const auto& row : Correlation) {
268 if (row.size() != Correlation.size()) {
274 YAML::Node NodeCopy = YAML::Clone(root);
275 YAML::Node systematics = NodeCopy[
"Systematics"];
277 if (!systematics || !systematics.IsSequence()) {
278 MACH3LOG_ERROR(
"'Systematics' node is missing or not a sequence");
282 if (!FancyNames.empty() && FancyNames.size() != Values.size()) {
283 MACH3LOG_ERROR(
"FancyNames size ({}) does not match Values size ({})", FancyNames.size(), Values.size());
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;
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()) {
304 YAML::Node& syst_i = it_i->second;
307 syst_i[
"Error"] = std::round(Errors[i] * 100.0) / 100.0;
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;
314 correlationsNode.push_back(singleEntry);
316 syst_i[
"Correlations"] = correlationsNode;
319 if (systematics.size() != Values.size()) {
320 MACH3LOG_ERROR(
"Mismatch in sizes: Values has {}, but YAML 'Systematics' has {} entries",
321 Values.size(), systematics.size());
325 for (std::size_t i = 0; i < systematics.size(); ++i) {
326 YAML::Node syst = systematics[i][
"Systematic"];
333 syst[
"Error"] = std::round(Errors[i] * 100.0) / 100.0;
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>();
341 correlationsNode.push_back(singleEntry);
343 syst[
"Correlations"] = correlationsNode;
350 std::string OutName =
"UpdatedCorrelationMatrix.yaml";
351 std::ofstream outFile(OutName);
357 outFile << YAMLString;
368 if (!CovarianceFolder) {
373 TMacro* foundMacro =
nullptr;
376 TIter next(CovarianceFolder->GetListOfKeys());
378 while ((key =
dynamic_cast<TKey*
>(next()))) {
379 if (std::string(key->GetClassName()) ==
"TMacro") {
381 if (macroCount == 1) {
382 foundMacro =
dynamic_cast<TMacro*
>(key->ReadObj());
387 if (macroCount == 1 && foundMacro) {
388 MACH3LOG_INFO(
"Found single TMacro in directory: using it.");
391 MACH3LOG_WARN(
"Found {} TMacro objects. Using hardcoded macro name: Config_xsec_cov.", macroCount);
392 TMacro* fallback = CovarianceFolder->Get<TMacro>(
"Config_xsec_cov");
394 MACH3LOG_WARN(
"Fallback macro 'Config_xsec_cov' not found in directory.");
414 TMatrixDSym* foundMatrix =
nullptr;
417 TIter next(TempFile->GetListOfKeys());
419 while ((key =
dynamic_cast<TKey*
>(next()))) {
420 std::string className = key->GetClassName();
421 if (className.find(
"TMatrix") != std::string::npos) {
423 if (matrixCount == 1) {
424 foundMatrix =
dynamic_cast<TMatrixDSym*
>(key->ReadObj());
429 if (matrixCount == 1 && foundMatrix) {
430 MACH3LOG_INFO(
"Found single TMatrixDSym in directory: using it.");
433 MACH3LOG_WARN(
"Found {} TMatrixDSym objects. Using hardcoded path: xsec_cov.", matrixCount);
434 TMatrixDSym* fallback = TempFile->Get<TMatrixDSym>(
"xsec_cov");
#define _MaCh3_Safe_Include_Start_
KS: Avoiding warning checking for headers.
#define _MaCh3_Safe_Include_End_
KS: Restore warning checking after including external headers.
#define _restrict_
KS: Using restrict limits the effects of pointer aliasing, aiding optimizations. While reading I foun...
std::string YAMLtoSTRING(const YAML::Node &node)
KS: Convert a YAML node to a string representation.
Custom exception class for MaCh3 errors.
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...