MaCh3  2.5.1
Reference Guide
Namespaces | Typedefs | Enumerations | Functions | Variables
M3 Namespace Reference

Main namespace for MaCh3 software. More...

Namespaces

 Plotting
 Flexible, experiment-agnostic plotting utilities for MaCh3.
 
 Utils
 Utility helpers used across MaCh3.
 

Typedefs

using float_t = double
 
using int_t = int
 
using uint_t = unsigned
 

Enumerations

enum  kInfCrit { kBIC , kDIC , kWAIC , kInfCrits }
 KS: Different Information Criterion tests mostly based Gelman paper. More...
 

Functions

template<typename T >
constexpr T fmaf_t (T x, T y, T z)
 Function template for fused multiply-add. More...
 
int GetNThreads ()
 number of threads which we need for example for TRandom3 More...
 
void AddPath (std::string &FilePath)
 Prepends the MACH3 environment path to FilePath if it is not already present. More...
 
int GetThreadIndex ()
 thread index inside parallel loop More...
 
TFile * Open (const std::string &Name, const std::string &Type, const std::string &File, const int Line)
 Opens a ROOT file with the given name and mode. More...
 
void ScaleHistogram (TH1 *Sample_Hist, const double scale)
 Scale histogram to get divided by bin width. More...
 
void CheckBinningMatch (const TH1D *Hist1, const TH1D *Hist2, const std::string &File, const int Line)
 KS: Helper function check if data and MC binning matches. More...
 
void CheckBinningMatch (const TH2D *Hist1, const TH2D *Hist2, const std::string &File, const int Line)
 KS: Helper function check if data and MC binning matches. More...
 
void CheckBinningMatch (TH2Poly *Hist1, TH2Poly *Hist2, const std::string &File, const int Line)
 KS: Helper function check if data and MC binning matches. More...
 
YAML::Node PolyToYaml (TH2Poly *Hist, const std::string &YamlName, const std::string &File, const int Line)
 KS: Convert TH2Poly into yaml config accepted by MaCh3. More...
 
template<typename ObjectType >
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. More...
 
double * MatrixMult (double *A, double *B, int n)
 CW: Multi-threaded matrix multiplication. More...
 
double ** MatrixMult (double **A, double **B, int n)
 CW: Multi-threaded matrix multiplication. More...
 
TMatrixD MatrixMult (TMatrixD A, TMatrixD B)
 CW: Multi-threaded matrix multiplication. More...
 
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. More...
 
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. More...
 
void FixSampleNamesQuotes (std::string &yamlStr)
 KS: Yaml emitter has problem and drops "", if you have special signs in you like * then there is problem. This bit hacky code adds these "". More...
 
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. More...
 
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. More...
 
TMacro * GetConfigMacroFromChain (TDirectory *CovarianceFolder)
 KS: We store configuration macros inside the chain. In the past, multiple configs were stored, which required error-prone hardcoding like "Config_xsec_cov". Therefore, this code maintains backward compatibility by checking the number of macros present and using a hardcoded name only if necessary. More...
 
TMatrixDSym * GetCovMatrixFromChain (TDirectory *TempFile)
 KS: Retrieve the cross-section covariance matrix from the given TDirectory. Historically, multiple covariance matrices could be stored, requiring fragile hardcoded paths like "CovarianceFolder/xsec_cov". This function maintains backward compatibility by: More...
 
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 can be even 20 times faster. More...
 
bool CanDecomposeMatrix (const TMatrixDSym &matrix)
 Checks if a matrix can be Cholesky decomposed. More...
 
void MakeMatrixPosDef (TMatrixDSym *cov)
 Makes sure that matrix is positive-definite by adding a small number to on-diagonal elements. More...
 
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. More...
 
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. More...
 

Variables

constexpr static const char * float_t_str_repr = "D"
 
constexpr static const double _BAD_DOUBLE_ = -999.99
 Default value used for double initialisation. More...
 
constexpr static const int _BAD_INT_ = -999
 Default value used for int initialisation. More...
 
constexpr static const double _DEFAULT_RETURN_VAL_ = -999999.123456
 
constexpr static const double Unity_D = 1.
 Some commonly used variables to which we set pointers to. More...
 
constexpr static const float Unity_F = 1.
 
constexpr static const float_t Unity = Unity_D
 
constexpr static const double Zero_D = 0.
 
constexpr static const float Zero_F = 0.
 
constexpr static const float_t Zero = Zero_D
 
constexpr static const double KinematicLowBound = std::numeric_limits<double>::lowest()
 When parameter has no bound this serves as it. Lowest possible value the system. More...
 
constexpr static const double KinematicUpBound = std::numeric_limits<double>::max()
 When parameter has no bound this serves as it. Highest possible value the system. More...
 
constexpr static const double _LARGE_LOGL_ = 1234567890.0
 Large Likelihood is used it parameter go out of physical boundary, this indicates in MCMC that such step should be removed. More...
 
constexpr static const double _LOW_MC_BOUND_ = .00001
 MC prediction lower bound in bin to identify problematic binning definitions and handle LogL calculation. More...
 
constexpr static const double DefSplineKnotUpBound = 9999
 Default value for spline knot capping, default mean not capping is being applied. More...
 
constexpr static const double DefSplineKnotLowBound = -9999
 Default value for spline knot capping, default mean not capping is being applied. More...
 
constexpr static const int UnderOverFlowBin = -1
 Mark bin which is overflow or underflow in MaCh3 binning. More...
 

Detailed Description

Main namespace for MaCh3 software.

Typedef Documentation

◆ float_t

typedef double M3::float_t

Definition at line 37 of file Core.h.

◆ int_t

using M3::int_t = typedef int

Definition at line 38 of file Core.h.

◆ uint_t

using M3::uint_t = typedef unsigned

Definition at line 39 of file Core.h.

Enumeration Type Documentation

◆ kInfCrit

KS: Different Information Criterion tests mostly based Gelman paper.

Enumerator
kBIC 

Bayesian Information Criterion.

kDIC 

Deviance Information Criterion.

kWAIC 

Watanabe-Akaike information criterion.

kInfCrits 

This only enumerates.

Definition at line 30 of file StatisticalUtils.h.

30  {
31  kBIC,
32  kDIC,
33  kWAIC,
34  kInfCrits
35  };
@ kWAIC
Watanabe-Akaike information criterion.
@ kInfCrits
This only enumerates.
@ kBIC
Bayesian Information Criterion.
@ kDIC
Deviance Information Criterion.

Function Documentation

◆ AddPath()

void M3::AddPath ( std::string &  FilePath)

Prepends the MACH3 environment path to FilePath if it is not already present.

Parameters
FilePathReference to the file path string to be modified.

Definition at line 382 of file Monitor.cpp.

382  {
383 // ***************************************************************************
384  //KS:Most inputs are in ${MACH3}/inputs/blarb.root
385  if (std::getenv("MACH3") == nullptr) {
386  MACH3LOG_ERROR("MACH3 is not defined");
387  MACH3LOG_ERROR("Please source your setup script");
388  MACH3LOG_ERROR("Read more: https://mach3-software.github.io/MaCh3/FAQ.html");
389  throw MaCh3Exception(__FILE__, __LINE__);
390  }
391 
392  std::string MaCh3Path = std::string(std::getenv("MACH3")) + "/";
393  // Check if FilePath does NOT start with MaCh3Path
394  if (FilePath.find(MaCh3Path) != 0) {
395  FilePath.insert(0, MaCh3Path);
396  }
397 }
#define MACH3LOG_ERROR
Definition: MaCh3Logger.h:37
Custom exception class used throughout MaCh3.

◆ AddTuneValues()

void M3::AddTuneValues ( YAML::Node &  root,
const std::vector< double > &  Values,
const std::string &  Tune,
const std::vector< std::string > &  FancyNames = {} 
)
inline

KS: Add Tune values to YAML covariance matrix.

Parameters
rootThe root YAML node to be updated.
ValuesThe values to add for the specified tune.
TuneThe name of the tune (e.g., "PostFit").
FancyNamesOptional list of fancy names to match systematics (must match Values size if provided).

Definition at line 200 of file ParameterHandlerUtils.h.

203  {}) {
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 }
std::string YAMLtoSTRING(const YAML::Node &node)
KS: Convert a YAML node to a string representation.
Definition: YamlHelper.h:112
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...
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 FixSampleNamesQuotes(std::string &yamlStr)
KS: Yaml emitter has problem and drops "", if you have special signs in you like * then there is prob...

◆ CanDecomposeMatrix()

bool M3::CanDecomposeMatrix ( const TMatrixDSym &  matrix)
inline

Checks if a matrix can be Cholesky decomposed.

Parameters
matrixInput symmetric matrix to test

Definition at line 511 of file ParameterHandlerUtils.h.

511  {
512 // *************************************
513  TDecompChol chdcmp(matrix);
514  return chdcmp.Decompose();
515 }

◆ CheckBinningMatch() [1/3]

void M3::CheckBinningMatch ( const TH1D *  Hist1,
const TH1D *  Hist2,
const std::string &  File,
const int  Line 
)

KS: Helper function check if data and MC binning matches.

Definition at line 821 of file HistogramUtils.cpp.

821  {
822 // ***************************************************************************
823  if (Hist1->GetNbinsX() != Hist2->GetNbinsX()) {
824  MACH3LOG_ERROR("Number of bins does not match for TH1D: {} vs {}", Hist1->GetNbinsX(), Hist2->GetNbinsX());
825  throw MaCh3Exception(File, Line);
826  }
827  for (int i = 1; i <= Hist1->GetNbinsX(); ++i) {
828  if (std::fabs(Hist1->GetXaxis()->GetBinLowEdge(i) - Hist2->GetXaxis()->GetBinLowEdge(i)) > 0.001 ||
829  std::fabs(Hist1->GetXaxis()->GetBinUpEdge(i) - Hist2->GetXaxis()->GetBinUpEdge(i)) > 0.001) {
830  MACH3LOG_ERROR("Bin edges do not match for TH1D at bin {}", i);
831  throw MaCh3Exception(File, Line);
832  }
833  }
834 }

◆ CheckBinningMatch() [2/3]

void M3::CheckBinningMatch ( const TH2D *  Hist1,
const TH2D *  Hist2,
const std::string &  File,
const int  Line 
)

KS: Helper function check if data and MC binning matches.

Definition at line 837 of file HistogramUtils.cpp.

837  {
838 // ***************************************************************************
839  if (Hist1->GetNbinsX() != Hist2->GetNbinsX() || Hist1->GetNbinsY() != Hist2->GetNbinsY()) {
840  MACH3LOG_ERROR("Number of bins does not match for TH2D");
841  throw MaCh3Exception(File, Line);
842  }
843 
844  for (int i = 1; i <= Hist1->GetNbinsX(); ++i) {
845  if (std::fabs(Hist1->GetXaxis()->GetBinLowEdge(i) - Hist2->GetXaxis()->GetBinLowEdge(i)) > 0.001 ||
846  std::fabs(Hist1->GetXaxis()->GetBinUpEdge(i) - Hist2->GetXaxis()->GetBinUpEdge(i)) > 0.001) {
847  MACH3LOG_ERROR("X bin edges do not match for TH2D at bin {}", i);
848  throw MaCh3Exception(File, Line);
849  }
850  }
851  for (int j = 1; j <= Hist1->GetNbinsY(); ++j) {
852  if (std::fabs(Hist1->GetYaxis()->GetBinLowEdge(j) - Hist2->GetYaxis()->GetBinLowEdge(j)) > 0.001 ||
853  std::fabs(Hist1->GetYaxis()->GetBinUpEdge(j) - Hist2->GetYaxis()->GetBinUpEdge(j)) > 0.001) {
854  MACH3LOG_ERROR("Y bin edges do not match for TH2D at bin {}", j);
855  throw MaCh3Exception(File, Line);
856  }
857  }
858 }

◆ CheckBinningMatch() [3/3]

void M3::CheckBinningMatch ( TH2Poly *  Hist1,
TH2Poly *  Hist2,
const std::string &  File,
const int  Line 
)

KS: Helper function check if data and MC binning matches.

Definition at line 862 of file HistogramUtils.cpp.

862  {
863 // ***************************************************************************
864  int NBins1 = Hist1->GetNumberOfBins();
865  int NBins2 = Hist2->GetNumberOfBins();
866  if (NBins1 != NBins2) {
867  MACH3LOG_ERROR("Number of bins does not match for TH2Poly: {} vs {}", NBins1, NBins2);
868  throw MaCh3Exception(File, Line);
869  }
870  for (int j = 1; j <= NBins1; j++)
871  {
872  //KS: There is weird offset between bin content and GetBins so this is correct, in spite of looking funny
873  TH2PolyBin* polybin1 = static_cast<TH2PolyBin*>(Hist1->GetBins()->At(j - 1));
874  TH2PolyBin* polybin2 = static_cast<TH2PolyBin*>(Hist2->GetBins()->At(j - 1));
875 
876  if( std::fabs(polybin2->GetXMin() - polybin1->GetXMin()) > 0.001 ||
877  std::fabs(polybin2->GetXMax() - polybin1->GetXMax()) > 0.001 ||
878  std::fabs(polybin2->GetYMin() - polybin1->GetYMin()) > 0.001 ||
879  std::fabs(polybin2->GetYMax() - polybin1->GetYMax()) > 0.001 )
880  {
881  MACH3LOG_ERROR("Binning doesn't match", Hist1->GetTitle());
882  MACH3LOG_ERROR("data x min {} x max {} y min {} y max {}", polybin2->GetXMin(), polybin2->GetXMax(), polybin2->GetYMin(), polybin2->GetYMax());
883  MACH3LOG_ERROR("mc x min {} x max {} y min {} y max {}", polybin1->GetXMin(), polybin1->GetXMax(), polybin1->GetYMin(), polybin1->GetYMax());
884  throw MaCh3Exception(File , Line );
885  }
886  }
887 }

◆ Clone()

template<typename ObjectType >
std::unique_ptr<ObjectType> M3::Clone ( const ObjectType *  obj,
const std::string &  name = "" 
)

KS: Creates a copy of a ROOT-like object and wraps it in a smart pointer.

Template Parameters
ObjectTypeThe type of the object to clone for example TH1D or TH2Poly.
Parameters
objPointer to the object to clone.
nameOptional argument allowing to set new name of cloned object
Returns
std::unique_ptr<ObjectType> Owning pointer to the cloned object.

Definition at line 155 of file HistogramUtils.h.

155  {
156  std::string cloneName = name.empty() ? obj->GetName() : name;
157 
158  std::unique_ptr<ObjectType> Hist(static_cast<ObjectType*>(obj->Clone(cloneName.c_str())));
159  // Disable ROOT memory management because it causes lot of headache especially as smart pointers are much smarter
160  Hist->SetDirectory(nullptr);
161 
162  return Hist;
163 }

◆ DebugPCA()

void M3::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 
)
inline

KS: Let's dump all useful matrices to properly validate PCA.

Definition at line 661 of file ParameterHandlerUtils.h.

671  {
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 }

◆ DumpParamHandlerToFile()

void M3::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 
)
inline

Dump Matrix to ROOT file, useful when we need to pass matrix info to another fitting group.

Warning
This is mostly used for backward compatibility

Definition at line 561 of file ParameterHandlerUtils.h.

572  {
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 }
std::string SplineInterpolation_ToString(const SplineInterpolation i)
Convert a LLH type to a string.
Flat representation of a spline index entry.
Definition: SplineCommon.h:33

◆ FixSampleNamesQuotes()

void M3::FixSampleNamesQuotes ( std::string &  yamlStr)
inline

KS: Yaml emitter has problem and drops "", if you have special signs in you like * then there is problem. This bit hacky code adds these "".

Parameters
yamlStrThe YAML string to be processed (modified in-place).

Definition at line 165 of file ParameterHandlerUtils.h.

165  {
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 }

◆ fmaf_t()

template<typename T >
constexpr T M3::fmaf_t ( x,
y,
z 
)
constexpr

Function template for fused multiply-add.

Definition at line 45 of file Core.h.

45  {
46  #ifdef _LOW_MEMORY_STRUCTS_
47  return std::fmaf(x, y, z);
48  #else
49  return std::fma(x, y, z);
50  #endif
51  }

◆ GetCholeskyDecomposedMatrix()

std::vector<std::vector<double> > M3::GetCholeskyDecomposedMatrix ( const TMatrixDSym &  matrix,
const std::string &  matrixName 
)
inline

Computes Cholesky decomposition of a symmetric positive definite matrix using custom function which can be even 20 times faster.

Parameters
matrixInput symmetric positive definite matrix
matrixNameIdentifier for error reporting

Definition at line 475 of file ParameterHandlerUtils.h.

475  {
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 }

◆ GetConfigMacroFromChain()

TMacro* M3::GetConfigMacroFromChain ( TDirectory *  CovarianceFolder)
inline

KS: We store configuration macros inside the chain. In the past, multiple configs were stored, which required error-prone hardcoding like "Config_xsec_cov". Therefore, this code maintains backward compatibility by checking the number of macros present and using a hardcoded name only if necessary.

Definition at line 395 of file ParameterHandlerUtils.h.

395  {
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 }
#define MACH3LOG_INFO
Definition: MaCh3Logger.h:35
#define MACH3LOG_WARN
Definition: MaCh3Logger.h:36

◆ GetCovMatrixFromChain()

TMatrixDSym* M3::GetCovMatrixFromChain ( TDirectory *  TempFile)
inline

KS: Retrieve the cross-section covariance matrix from the given TDirectory. Historically, multiple covariance matrices could be stored, requiring fragile hardcoded paths like "CovarianceFolder/xsec_cov". This function maintains backward compatibility by:

  • Using the single matrix present if only one exists,
  • Otherwise falling back to the hardcoded path. This avoids error-prone assumptions while supporting both old and new formats.

Definition at line 436 of file ParameterHandlerUtils.h.

436  {
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 }

◆ GetNThreads()

int M3::GetNThreads ( )

number of threads which we need for example for TRandom3

Definition at line 372 of file Monitor.cpp.

372  {
373 // ***************************************************************************
374  #ifdef MULTITHREAD
375  return omp_get_max_threads();
376  #else
377  return 1;
378  #endif
379 }

◆ GetThreadIndex()

int M3::GetThreadIndex ( )
inline

thread index inside parallel loop

Definition at line 86 of file Monitor.h.

86  {
87  #ifdef MULTITHREAD
88  return omp_get_thread_num();
89  #else
90  return 0;
91  #endif
92  }

◆ MakeCorrelationMatrix()

void M3::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 = {} 
)
inline

KS: Replace correlation matrix and tune values in YAML covariance matrix.

Parameters
rootThe root YAML node to be updated.
ValuesThe new values for each systematic.
ErrorsThe new errors for each systematic.
CorrelationThe new correlation matrix (must be square and match Values size).
OutYAMLNameThe output filename for the updated YAML.
FancyNamesOptional list of fancy names to match systematics (must match Values size if provided).

Definition at line 281 of file ParameterHandlerUtils.h.

286  {}) {
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 }

◆ MakeMatrixPosDef()

void M3::MakeMatrixPosDef ( TMatrixDSym *  cov)
inline

Makes sure that matrix is positive-definite by adding a small number to on-diagonal elements.

Definition at line 519 of file ParameterHandlerUtils.h.

519  {
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 }
bool CanDecomposeMatrix(const TMatrixDSym &matrix)
Checks if a matrix can be Cholesky decomposed.

◆ MatrixMult() [1/3]

double** M3::MatrixMult ( double **  A,
double **  B,
int  n 
)
inline

CW: Multi-threaded matrix multiplication.

Definition at line 75 of file ParameterHandlerUtils.h.

75  {
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 }
TMatrixD MatrixMult(TMatrixD A, TMatrixD B)
CW: Multi-threaded matrix multiplication.

◆ MatrixMult() [2/3]

double* M3::MatrixMult ( double *  A,
double *  B,
int  n 
)
inline

CW: Multi-threaded matrix multiplication.

Definition at line 43 of file ParameterHandlerUtils.h.

43  {
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 }

◆ MatrixMult() [3/3]

TMatrixD M3::MatrixMult ( TMatrixD  A,
TMatrixD  B 
)
inline

CW: Multi-threaded matrix multiplication.

Definition at line 111 of file ParameterHandlerUtils.h.

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 }

◆ MatrixVectorMulti()

void M3::MatrixVectorMulti ( double *_restrict_  VecMulti,
double **_restrict_  matrix,
const double *_restrict_  vector,
const int  n 
)
inline

KS: Custom function to perform multiplication of matrix and vector with multithreading.

Parameters
VecMultiOutput Vector, VecMulti = matrix x vector
matrixThis matrix is used for multiplication VecMulti = matrix x vector
vectorThis vector is used for multiplication VecMulti = matrix x vector
nthis is size of matrix and vector, we assume matrix is symmetric

Definition at line 125 of file ParameterHandlerUtils.h.

125  {
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 }

◆ MatrixVectorMultiSingle()

double M3::MatrixVectorMultiSingle ( double **_restrict_  matrix,
const double *_restrict_  vector,
const int  Length,
const int  i 
)
inline

KS: Custom function to perform multiplication of matrix and single element which is thread safe.

Parameters
matrixThis matrix is used for multiplication VecMulti = matrix x vector
vectorThis vector is used for multiplication VecMulti = matrix x vector
Lengththis is size of matrix and vector, we assume matrix is symmetric
iElement of matrix that we want to multiply

Definition at line 150 of file ParameterHandlerUtils.h.

150  {
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 }

◆ Open()

TFile * M3::Open ( const std::string &  Name,
const std::string &  Type,
const std::string &  File,
const int  Line 
)

Opens a ROOT file with the given name and mode.

This function wraps ROOT’s TFile constructor and checks whether the file was opened successfully and give some useful debugging information

Parameters
NameThe name or path of the file to open.
TypeThe file open mode (e.g., "READ", "RECREATE", "UPDATE").
FileThe name of the file where the exception occurred.
LineThe line number where the exception occurred.
Returns
Pointer to the opened ROOT file.

Definition at line 782 of file HistogramUtils.cpp.

782  {
783 // **************************************************************************
784  TFile* OutFile = new TFile(Name.c_str(), Type.c_str());
785 
786  // Check if the file is successfully opened and usable
787  if (OutFile->IsZombie()) {
788  MACH3LOG_ERROR("Failed to open file: {}", Name);
789  std::string lowerType = Type;
790  std::transform(lowerType.begin(), lowerType.end(), lowerType.begin(), ::tolower);
791  if (lowerType == "recreate") {
792  MACH3LOG_ERROR("Check if directory exist");
793  }
794  delete OutFile;
795  throw MaCh3Exception(File, Line);
796  }
797  return OutFile;
798 }

◆ PolyToYaml()

YAML::Node M3::PolyToYaml ( TH2Poly *  Hist,
const std::string &  YamlName,
const std::string &  File,
const int  Line 
)

KS: Convert TH2Poly into yaml config accepted by MaCh3.

Definition at line 891 of file HistogramUtils.cpp.

891  {
892 // ***************************************************************************
893  if (!Hist) {
894  MACH3LOG_ERROR("Null TH2Poly pointer");
895  throw MaCh3Exception(File, Line);
896  }
897 
898  YAML::Node bins(YAML::NodeType::Sequence);
899  bins.SetStyle(YAML::EmitterStyle::Flow);
900  const int NBins = Hist->GetNumberOfBins();
901 
902  for (int j = 1; j <= NBins; j++)
903  {
904  TH2PolyBin* polybin = static_cast<TH2PolyBin*>(Hist->GetBins()->At(j - 1));
905 
906  double xmin = polybin->GetXMin();
907  double xmax = polybin->GetXMax();
908  double ymin = polybin->GetYMin();
909  double ymax = polybin->GetYMax();
910 
911  YAML::Node xNode(YAML::NodeType::Sequence);
912  xNode.SetStyle(YAML::EmitterStyle::Flow);
913  xNode.push_back(xmin);
914  xNode.push_back(xmax);
915 
916  YAML::Node yNode(YAML::NodeType::Sequence);
917  yNode.SetStyle(YAML::EmitterStyle::Flow);
918  yNode.push_back(ymin);
919  yNode.push_back(ymax);
920 
921  YAML::Node bin(YAML::NodeType::Sequence);
922  bin.SetStyle(YAML::EmitterStyle::Flow);
923  bin.push_back(xNode);
924  bin.push_back(yNode);
925 
926  bins.push_back(bin);
927  }
928 
929  YAML::Node result;
930  result[YamlName] = bins;
931 
932  return result;
933 }

◆ ScaleHistogram()

void M3::ScaleHistogram ( TH1 *  Sample_Hist,
const double  scale 
)

Scale histogram to get divided by bin width.

Definition at line 801 of file HistogramUtils.cpp.

801  {
802 // ***************************************************************************
803  for (int j = 0; j <= Sample_Hist->GetNbinsX(); ++j)
804  {
805  double num = Sample_Hist->GetBinContent(j);
806  double numErr = Sample_Hist->GetBinError(j);
807  double den = Sample_Hist->GetBinWidth(j);
808  double value = 0.;
809  double valueErr = 0.;
810  if (den != 0)
811  {
812  value = num/(den/scale);
813  valueErr = numErr/(den/scale);
814  Sample_Hist->SetBinContent(j,value);
815  Sample_Hist->SetBinError(j,valueErr);
816  }
817  }
818 }

Variable Documentation

◆ _BAD_DOUBLE_

constexpr static const double M3::_BAD_DOUBLE_ = -999.99
staticconstexpr

Default value used for double initialisation.

Definition at line 53 of file Core.h.

◆ _BAD_INT_

constexpr static const int M3::_BAD_INT_ = -999
staticconstexpr

Default value used for int initialisation.

Definition at line 55 of file Core.h.

◆ _DEFAULT_RETURN_VAL_

constexpr static const double M3::_DEFAULT_RETURN_VAL_ = -999999.123456
staticconstexpr

Definition at line 56 of file Core.h.

◆ _LARGE_LOGL_

constexpr static const double M3::_LARGE_LOGL_ = 1234567890.0
staticconstexpr

Large Likelihood is used it parameter go out of physical boundary, this indicates in MCMC that such step should be removed.

Definition at line 80 of file Core.h.

◆ _LOW_MC_BOUND_

constexpr static const double M3::_LOW_MC_BOUND_ = .00001
staticconstexpr

MC prediction lower bound in bin to identify problematic binning definitions and handle LogL calculation.

Definition at line 83 of file Core.h.

◆ DefSplineKnotLowBound

constexpr static const double M3::DefSplineKnotLowBound = -9999
staticconstexpr

Default value for spline knot capping, default mean not capping is being applied.

Definition at line 88 of file Core.h.

◆ DefSplineKnotUpBound

constexpr static const double M3::DefSplineKnotUpBound = 9999
staticconstexpr

Default value for spline knot capping, default mean not capping is being applied.

Definition at line 86 of file Core.h.

◆ float_t_str_repr

constexpr static const char* M3::float_t_str_repr = "D"
staticconstexpr

Definition at line 40 of file Core.h.

◆ KinematicLowBound

constexpr static const double M3::KinematicLowBound = std::numeric_limits<double>::lowest()
staticconstexpr

When parameter has no bound this serves as it. Lowest possible value the system.

Definition at line 75 of file Core.h.

◆ KinematicUpBound

constexpr static const double M3::KinematicUpBound = std::numeric_limits<double>::max()
staticconstexpr

When parameter has no bound this serves as it. Highest possible value the system.

Definition at line 77 of file Core.h.

◆ UnderOverFlowBin

constexpr static const int M3::UnderOverFlowBin = -1
staticconstexpr

Mark bin which is overflow or underflow in MaCh3 binning.

Definition at line 91 of file Core.h.

◆ Unity

constexpr static const float_t M3::Unity = Unity_D
staticconstexpr

Definition at line 64 of file Core.h.

◆ Unity_D

constexpr static const double M3::Unity_D = 1.
staticconstexpr

Some commonly used variables to which we set pointers to.

Definition at line 59 of file Core.h.

◆ Unity_F

constexpr static const float M3::Unity_F = 1.
staticconstexpr

Definition at line 60 of file Core.h.

◆ Zero

constexpr static const float_t M3::Zero = Zero_D
staticconstexpr

Definition at line 71 of file Core.h.

◆ Zero_D

constexpr static const double M3::Zero_D = 0.
staticconstexpr

Definition at line 66 of file Core.h.

◆ Zero_F

constexpr static const float M3::Zero_F = 0.
staticconstexpr

Definition at line 67 of file Core.h.