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

Main namespace for MaCh3 software. More...

Namespaces

 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...
 

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 193 of file ParameterHandlerUtils.h.

196  {}) {
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 }
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 504 of file ParameterHandlerUtils.h.

504  {
505 // *************************************
506  TDecompChol chdcmp(matrix);
507  return chdcmp.Decompose();
508 }

◆ 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 704 of file HistogramUtils.cpp.

704  {
705 // ***************************************************************************
706  if (Hist1->GetNbinsX() != Hist2->GetNbinsX()) {
707  MACH3LOG_ERROR("Number of bins does not match for TH1D: {} vs {}", Hist1->GetNbinsX(), Hist2->GetNbinsX());
708  throw MaCh3Exception(File, Line);
709  }
710  for (int i = 1; i <= Hist1->GetNbinsX(); ++i) {
711  if (std::fabs(Hist1->GetXaxis()->GetBinLowEdge(i) - Hist2->GetXaxis()->GetBinLowEdge(i)) > 0.001 ||
712  std::fabs(Hist1->GetXaxis()->GetBinUpEdge(i) - Hist2->GetXaxis()->GetBinUpEdge(i)) > 0.001) {
713  MACH3LOG_ERROR("Bin edges do not match for TH1D at bin {}", i);
714  throw MaCh3Exception(File, Line);
715  }
716  }
717 }

◆ 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 720 of file HistogramUtils.cpp.

720  {
721 // ***************************************************************************
722  if (Hist1->GetNbinsX() != Hist2->GetNbinsX() || Hist1->GetNbinsY() != Hist2->GetNbinsY()) {
723  MACH3LOG_ERROR("Number of bins does not match for TH2D");
724  throw MaCh3Exception(File, Line);
725  }
726 
727  for (int i = 1; i <= Hist1->GetNbinsX(); ++i) {
728  if (std::fabs(Hist1->GetXaxis()->GetBinLowEdge(i) - Hist2->GetXaxis()->GetBinLowEdge(i)) > 0.001 ||
729  std::fabs(Hist1->GetXaxis()->GetBinUpEdge(i) - Hist2->GetXaxis()->GetBinUpEdge(i)) > 0.001) {
730  MACH3LOG_ERROR("X bin edges do not match for TH2D at bin {}", i);
731  throw MaCh3Exception(File, Line);
732  }
733  }
734  for (int j = 1; j <= Hist1->GetNbinsY(); ++j) {
735  if (std::fabs(Hist1->GetYaxis()->GetBinLowEdge(j) - Hist2->GetYaxis()->GetBinLowEdge(j)) > 0.001 ||
736  std::fabs(Hist1->GetYaxis()->GetBinUpEdge(j) - Hist2->GetYaxis()->GetBinUpEdge(j)) > 0.001) {
737  MACH3LOG_ERROR("Y bin edges do not match for TH2D at bin {}", j);
738  throw MaCh3Exception(File, Line);
739  }
740  }
741 }

◆ 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 745 of file HistogramUtils.cpp.

745  {
746 // ***************************************************************************
747  int NBins1 = Hist1->GetNumberOfBins();
748  int NBins2 = Hist2->GetNumberOfBins();
749  if (NBins1 != NBins2) {
750  MACH3LOG_ERROR("Number of bins does not match for TH2Poly: {} vs {}", NBins1, NBins2);
751  throw MaCh3Exception(File, Line);
752  }
753  for (int j = 1; j <= NBins1; j++)
754  {
755  //KS: There is weird offset between bin content and GetBins so this is correct, in spite of looking funny
756  TH2PolyBin* polybin1 = static_cast<TH2PolyBin*>(Hist1->GetBins()->At(j - 1));
757  TH2PolyBin* polybin2 = static_cast<TH2PolyBin*>(Hist2->GetBins()->At(j - 1));
758 
759  if( std::fabs(polybin2->GetXMin() - polybin1->GetXMin()) > 0.001 ||
760  std::fabs(polybin2->GetXMax() - polybin1->GetXMax()) > 0.001 ||
761  std::fabs(polybin2->GetYMin() - polybin1->GetYMin()) > 0.001 ||
762  std::fabs(polybin2->GetYMax() - polybin1->GetYMax()) > 0.001 )
763  {
764  MACH3LOG_ERROR("Binning doesn't match", Hist1->GetTitle());
765  MACH3LOG_ERROR("data x min {} x max {} y min {} y max {}", polybin2->GetXMin(), polybin2->GetXMax(), polybin2->GetYMin(), polybin2->GetYMax());
766  MACH3LOG_ERROR("mc x min {} x max {} y min {} y max {}", polybin1->GetXMin(), polybin1->GetXMax(), polybin1->GetYMin(), polybin1->GetYMax());
767  throw MaCh3Exception(File , Line );
768  }
769  }
770 }

◆ 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 146 of file HistogramUtils.h.

146  {
147  std::string cloneName = name.empty() ? obj->GetName() : name;
148 
149  std::unique_ptr<ObjectType> Hist(static_cast<ObjectType*>(obj->Clone(cloneName.c_str())));
150  // Disable ROOT memory management because it causes lot of headache especially as smart pointers are much smarter
151  Hist->SetDirectory(nullptr);
152 
153  return Hist;
154 }

◆ 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 554 of file ParameterHandlerUtils.h.

565  {
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 }
std::string SplineInterpolation_ToString(const SplineInterpolation i)
Convert a LLH type to a string.

◆ 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 158 of file ParameterHandlerUtils.h.

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

◆ 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 468 of file ParameterHandlerUtils.h.

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

◆ 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 388 of file ParameterHandlerUtils.h.

388  {
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 }
#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 429 of file ParameterHandlerUtils.h.

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

◆ 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 274 of file ParameterHandlerUtils.h.

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

◆ 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 512 of file ParameterHandlerUtils.h.

512  {
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 }
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 68 of file ParameterHandlerUtils.h.

68  {
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 }
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 36 of file ParameterHandlerUtils.h.

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

◆ MatrixMult() [3/3]

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

CW: Multi-threaded matrix multiplication.

Definition at line 104 of file ParameterHandlerUtils.h.

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 }

◆ 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 118 of file ParameterHandlerUtils.h.

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

◆ 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 143 of file ParameterHandlerUtils.h.

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

◆ 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 665 of file HistogramUtils.cpp.

665  {
666 // **************************************************************************
667  TFile* OutFile = new TFile(Name.c_str(), Type.c_str());
668 
669  // Check if the file is successfully opened and usable
670  if (OutFile->IsZombie()) {
671  MACH3LOG_ERROR("Failed to open file: {}", Name);
672  std::string lowerType = Type;
673  std::transform(lowerType.begin(), lowerType.end(), lowerType.begin(), ::tolower);
674  if (lowerType == "recreate") {
675  MACH3LOG_ERROR("Check if directory exist");
676  }
677  delete OutFile;
678  throw MaCh3Exception(File, Line);
679  }
680  return OutFile;
681 }

◆ 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 774 of file HistogramUtils.cpp.

774  {
775 // ***************************************************************************
776  if (!Hist) {
777  MACH3LOG_ERROR("Null TH2Poly pointer");
778  throw MaCh3Exception(File, Line);
779  }
780 
781  YAML::Node bins(YAML::NodeType::Sequence);
782  bins.SetStyle(YAML::EmitterStyle::Flow);
783  const int NBins = Hist->GetNumberOfBins();
784 
785  for (int j = 1; j <= NBins; j++)
786  {
787  TH2PolyBin* polybin = static_cast<TH2PolyBin*>(Hist->GetBins()->At(j - 1));
788 
789  double xmin = polybin->GetXMin();
790  double xmax = polybin->GetXMax();
791  double ymin = polybin->GetYMin();
792  double ymax = polybin->GetYMax();
793 
794  YAML::Node xNode(YAML::NodeType::Sequence);
795  xNode.SetStyle(YAML::EmitterStyle::Flow);
796  xNode.push_back(xmin);
797  xNode.push_back(xmax);
798 
799  YAML::Node yNode(YAML::NodeType::Sequence);
800  yNode.SetStyle(YAML::EmitterStyle::Flow);
801  yNode.push_back(ymin);
802  yNode.push_back(ymax);
803 
804  YAML::Node bin(YAML::NodeType::Sequence);
805  bin.SetStyle(YAML::EmitterStyle::Flow);
806  bin.push_back(xNode);
807  bin.push_back(yNode);
808 
809  bins.push_back(bin);
810  }
811 
812  YAML::Node result;
813  result[YamlName] = bins;
814 
815  return result;
816 }

◆ ScaleHistogram()

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

Scale histogram to get divided by bin width.

Definition at line 684 of file HistogramUtils.cpp.

684  {
685 // ***************************************************************************
686  for (int j = 0; j <= Sample_Hist->GetNbinsX(); ++j)
687  {
688  double num = Sample_Hist->GetBinContent(j);
689  double numErr = Sample_Hist->GetBinError(j);
690  double den = Sample_Hist->GetBinWidth(j);
691  double value = 0.;
692  double valueErr = 0.;
693  if (den != 0)
694  {
695  value = num/(den/scale);
696  valueErr = numErr/(den/scale);
697  Sample_Hist->SetBinContent(j,value);
698  Sample_Hist->SetBinError(j,valueErr);
699  }
700  }
701 }

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.