MaCh3  2.5.1
Reference Guide
ParameterHandlerGeneric.cpp
Go to the documentation of this file.
2 
3 // ********************************************
4 // ETA - YAML constructor
5 // this will replace the root file constructor but let's keep it in
6 // to do some validations
7 ParameterHandlerGeneric::ParameterHandlerGeneric(const std::vector<std::string>& YAMLFile, std::string name,
8  double threshold, int FirstPCA, int LastPCA) {
9 // ********************************************
10  MACH3LOG_INFO("Constructing instance of ParameterHandler using");
11  inputFile = YAMLFile[0];
12  matrixName = name;
13  pca = true;
14 
15  doSpecialStepProposal = false;
16  // Not using adaptive by default
17  use_adaptive = false;
18  for(unsigned int i = 0; i < YAMLFile.size(); i++)
19  {
20  MACH3LOG_INFO("{}", YAMLFile[i]);
21  }
22  MACH3LOG_INFO("as an input");
23 
24  if (threshold < 0 || threshold >= 1) {
25  MACH3LOG_INFO("Principal component analysis but given the threshold for the principal components to be less than 0, or greater than (or equal to) 1. This will not work");
26  MACH3LOG_INFO("Please specify a number between 0 and 1");
27  MACH3LOG_INFO("You specified: ");
28  MACH3LOG_INFO("Am instead calling the usual non-PCA constructor...");
29  pca = false;
30  }
31 
32  InitialiseFromConfig(YAMLFile);
33  // Call the innocent helper function
34  if (pca) ConstructPCA(threshold, FirstPCA, LastPCA);
35 
37 
38  //ETA - again this really doesn't need to be hear...
39  for (int i = 0; i < _fNumPar; i++)
40  {
41  // Sort out the print length
42  if(int(_fNames[i].length()) > PrintLength) PrintLength = int(_fNames[i].length());
43  } // end the for loop
44 
45  MACH3LOG_DEBUG("Constructing instance of ParameterHandler");
47  // Print
48  Print();
49 }
50 
51 // ********************************************
52 void ParameterHandlerGeneric::LoadAndMergeYAML(const std::vector<std::string>& YAMLFile,
53  std::map<std::pair<int,int>, std::unique_ptr<TMatrixDSym>>& ThrowSubMatrixOverrides) {
54 // ********************************************
55  int running_num_file_pars = 0;
56 
57  _fYAMLDoc["Systematics"] = YAML::Node(YAML::NodeType::Sequence);
58  for(unsigned int i = 0; i < YAMLFile.size(); i++)
59  {
60  YAML::Node YAMLDocTemp = M3OpenConfig(YAMLFile[i]);
61 
62  if (YAMLDocTemp["ThrowMatrixOverride"]) { // LP: this allows us to put in
63  // proposal matrix overrides per
64  // parameter-containing file, add
65  // the block diagonal proposal
66  // matrix to a list and overwrite
67  // the throw matrix after set up.
68  auto filename =
69  YAMLDocTemp["ThrowMatrixOverride"]["file"].as<std::string>();
70  TFile *submatrix_file = TFile::Open(filename.c_str());
71 
72  auto matrixname =
73  YAMLDocTemp["ThrowMatrixOverride"]["matrix"].as<std::string>();
74  std::unique_ptr<TMatrixDSym> submatrix{
75  submatrix_file->Get<TMatrixDSym>(matrixname.c_str())};
76  if (!submatrix) {
77  MACH3LOG_CRITICAL("Covariance matrix {} doesn't exist in file: {}",
78  matrixname, filename);
79  throw MaCh3Exception(__FILE__, __LINE__);
80  }
81  auto numrows = submatrix->GetNrows();
82  // LP: the -1 here is because we specify the last index for consistency
83  // with PCAHandler, not the first index after the end as is more common
84  // throughout computer science...
85  ThrowSubMatrixOverrides[{running_num_file_pars,
86  running_num_file_pars + (numrows - 1)}] = std::move(submatrix);
87 
88  // LP: check names by default, but have option to disable check if you
89  // know what you're doing
90  if (!bool(YAMLDocTemp["ThrowMatrixOverride"]["check_names"]) ||
91  YAMLDocTemp["ThrowMatrixOverride"]["check_names"].as<bool>()) {
92  auto nametree = submatrix_file->Get<TTree>("param_names");
93  if (!nametree) {
94  MACH3LOG_CRITICAL("TTree param_names doesn't exist in file: {}. Set "
95  "ThrowMatrixOverride: {{ check_names: False }} to "
96  "disable this check.",
97  filename);
98  throw MaCh3Exception(__FILE__, __LINE__);
99  }
100  std::string *param_name = nullptr;
101  nametree->SetBranchAddress("name", &param_name);
102 
103  if (nametree->GetEntries() != int(YAMLDocTemp["Systematics"].size())) {
104  MACH3LOG_CRITICAL("TTree param_names in file: {} has {} entries, but "
105  "the corresponding yaml file only declares {} "
106  "parameters. Set ThrowMatrixOverride: {{ "
107  "check_names: False }} to disable this check.",
108  filename, nametree->GetEntries(),
109  YAMLDocTemp["Systematics"].size());
110  throw MaCh3Exception(__FILE__, __LINE__);
111  }
112 
113  int pit = 0;
114  for (const auto &param : YAMLDocTemp["Systematics"]) {
115  nametree->GetEntry(pit++);
116  auto yaml_pname = Get<std::string>(
117  param["Systematic"]["Names"]["FancyName"], __FILE__, __LINE__);
118  if ((*param_name) != yaml_pname) {
120  "TTree param_names in file: {} at entry {} has parameter {}, "
121  "but "
122  "the corresponding yaml parameter is named {}. Set "
123  "ThrowMatrixOverride: {{ "
124  "check_names: False }} to disable this check.",
125  filename, pit, (*param_name), yaml_pname);
126  throw MaCh3Exception(__FILE__, __LINE__);
127  }
128  }
129  }
130  submatrix_file->Close();
131  }
132 
133  for (const auto& item : YAMLDocTemp["Systematics"]) {
134  _fYAMLDoc["Systematics"].push_back(item);
135  running_num_file_pars++;
136  }
137  }
138 }
139 
140 // ********************************************
141 void ParameterHandlerGeneric::LoadCorrelationFromConfig(std::vector<std::map<std::string,double>>& Correlations,
142  std::map<std::string, int>& CorrNamesMap) {
143 // ********************************************
144  // ETA Now that we've been through all systematic let's fill the covmatrix
145  //This makes the root TCov from YAML
146  TMatrixDSym* _fCovMatrix = new TMatrixDSym(_fNumPar);
147  for(int j = 0; j < _fNumPar; j++) {
148  (*_fCovMatrix)(j, j) = _fError[j]*_fError[j];
149  //Get the map of parameter name to correlation from the Correlations object
150  for (auto const& pair : Correlations[j]) {
151  auto const& key = pair.first;
152  auto const& val = pair.second;
153  int index = -1;
154  //If you found the parameter name then get the index
155  if (CorrNamesMap.find(key) != CorrNamesMap.end()) {
156  index = CorrNamesMap[key];
157  } else {
158  MACH3LOG_ERROR("Parameter {} not in list! Check your spelling?", key);
159  throw MaCh3Exception(__FILE__ , __LINE__ );
160  }
161  double Corr1 = val;
162  double Corr2 = 0;
163  if(Correlations[index].find(_fFancyNames[j]) != Correlations[index].end()) {
164  Corr2 = Correlations[index][_fFancyNames[j]];
165  //Do they agree to better than float precision?
166  if(std::abs(Corr2 - Corr1) > FLT_EPSILON) {
167  MACH3LOG_ERROR("Correlations are not equal between {} and {}", _fFancyNames[j], key);
168  MACH3LOG_ERROR("Got : {} and {}", Corr2, Corr1);
169  throw MaCh3Exception(__FILE__ , __LINE__ );
170  }
171  } else {
172  MACH3LOG_ERROR("Correlation does not appear reciprocally between {} and {}", _fFancyNames[j], key);
173  throw MaCh3Exception(__FILE__ , __LINE__ );
174  }
175  (*_fCovMatrix)(j, index)= (*_fCovMatrix)(index, j) = Corr1*_fError[j]*_fError[index];
176  }
177  }
178 
179  //Now make positive definite
180  MakePosDef(_fCovMatrix);
181  SetCovMatrix(_fCovMatrix);
182 }
183 
184 // ********************************************
185 // ETA An init function for the YAML constructor
186 // All you really need from the YAML file is the number of Systematics
187 void ParameterHandlerGeneric::InitialiseFromConfig(const std::vector<std::string>& YAMLFile) {
188 // ********************************************
189  std::map<std::pair<int, int>, std::unique_ptr<TMatrixDSym>> ThrowSubMatrixOverrides;
190  LoadAndMergeYAML(YAMLFile, ThrowSubMatrixOverrides);
191 
192  const int nThreads = M3::GetNThreads();
193  //KS: set Random numbers for each thread so each thread has different seed
194  //or for one thread if without MULTITHREAD
195  random_number.reserve(nThreads);
196  for (int iThread = 0; iThread < nThreads; iThread++) {
197  random_number.emplace_back(std::make_unique<TRandom3>(0));
198  }
199  PrintLength = 35;
200 
201  // Set the covariance matrix
202  _fNumPar = int(_fYAMLDoc["Systematics"].size());
203 
205 
206  int i = 0;
207  std::vector<std::map<std::string,double>> Correlations(_fNumPar);
208  std::map<std::string, int> CorrNamesMap;
209 
210  //ETA - read in the systematics. Would be good to add in some checks to make sure
211  //that there are the correct number of entries i.e. are the _fNumPar for Names,
212  //PreFitValues etc etc.
213  for (auto const &param : _fYAMLDoc["Systematics"])
214  {
215  _fFancyNames[i] = Get<std::string>(param["Systematic"]["Names"]["FancyName"], __FILE__ , __LINE__);
216  _fPreFitValue[i] = Get<double>(param["Systematic"]["ParameterValues"]["PreFitValue"], __FILE__ , __LINE__);
217  _fIndivStepScale[i] = Get<double>(param["Systematic"]["StepScale"]["MCMC"], __FILE__ , __LINE__);
218  _fError[i] = Get<double>(param["Systematic"]["Error"], __FILE__ , __LINE__);
219  if(_fError[i] <= 0) {
220  MACH3LOG_ERROR("Error for param {}({}) is negative and equal to {}", _fFancyNames[i], i, _fError[i]);
221  throw MaCh3Exception(__FILE__ , __LINE__ );
222  }
223  //ETA - a bit of a fudge but works
224  auto TempBoundsVec = GetBounds(param["Systematic"]["ParameterBounds"]);
225  _fLowBound[i] = TempBoundsVec[0];
226  _fUpBound[i] = TempBoundsVec[1];
227 
228  //ETA - now for parameters which are optional and have default values
229  _fFlatPrior[i] = GetFromManager<bool>(param["Systematic"]["FlatPrior"], false, __FILE__ , __LINE__);
230 
231  // Allow to fix param, this setting should be used only for params which are permanently fixed like baseline, please use global config for fixing param more flexibly
232  if(GetFromManager<bool>(param["Systematic"]["FixParam"], false, __FILE__ , __LINE__)) {
234  }
235 
236  if(param["Systematic"]["SpecialProposal"]) {
237  EnableSpecialProposal(param["Systematic"]["SpecialProposal"], i);
238  }
239 
240  //Fill the map to get the correlations later as well
241  CorrNamesMap[param["Systematic"]["Names"]["FancyName"].as<std::string>()]=i;
242 
243  //Also loop through the correlations
244  if(param["Systematic"]["Correlations"]) {
245  for(unsigned int Corr_i = 0; Corr_i < param["Systematic"]["Correlations"].size(); ++Corr_i){
246  for (YAML::const_iterator it = param["Systematic"]["Correlations"][Corr_i].begin(); it!=param["Systematic"]["Correlations"][Corr_i].end();++it) {
247  Correlations[i][it->first.as<std::string>()] = it->second.as<double>();
248  }
249  }
250  }
251  i++;
252  } // end loop over para
253  if(i != _fNumPar) {
254  MACH3LOG_CRITICAL("Inconsistent number of params in Yaml {} vs {}, this indicate wrong syntax", i, i, _fNumPar);
255  throw MaCh3Exception(__FILE__ , __LINE__ );
256  }
257 
258  // load correlation
259  LoadCorrelationFromConfig(Correlations, CorrNamesMap);
260  if (_fNumPar <= 0) {
261  MACH3LOG_ERROR("ParameterHandler object has {} systematics!", _fNumPar);
262  throw MaCh3Exception(__FILE__ , __LINE__ );
263  }
264 
265  for(auto const & matovr : ThrowSubMatrixOverrides){
266  SetSubThrowMatrix(matovr.first.first, matovr.first.second, *matovr.second);
267  }
268 
269  Tunes = std::make_unique<ParameterTunes>(_fYAMLDoc["Systematics"]);
270 
271  MACH3LOG_INFO("Created covariance matrix from files: ");
272  for(const auto &file : YAMLFile){
273  MACH3LOG_INFO("{} ", file);
274  }
275  MACH3LOG_INFO("----------------");
276  MACH3LOG_INFO("Found {} systematics parameters in total", _fNumPar);
277  MACH3LOG_INFO("----------------");
278 }
279 
280 // ********************************************
282 // ********************************************
284 
285  _fParamType = std::vector<SystType>(_fNumPar);
286  _ParameterGroup = std::vector<std::string>(_fNumPar);
287  _fSampleNames = std::vector<std::vector<std::string>>(_fNumPar);
288 
289  //KS: We know at most how params we expect so reserve memory for max possible params. Later we will shrink to size to not waste memory. Reserving means slightly faster loading and possible less memory fragmentation.
290  NormParams.reserve(_fNumPar);
291  SplineParams.reserve(_fNumPar);
292  FuncParams.reserve(_fNumPar);
293  OscParams.reserve(_fNumPar);
294 
295  int i = 0;
296  unsigned int ParamCounter[SystType::kSystTypes] = {0};
297  //ETA - read in the systematics. Would be good to add in some checks to make sure
298  //that there are the correct number of entries i.e. are the _fNumPars for Names,
299  //PreFitValues etc etc.
300  for (auto const &param : _fYAMLDoc["Systematics"])
301  {
302  _ParameterGroup[i] = Get<std::string>(param["Systematic"]["ParameterGroup"], __FILE__ , __LINE__);
303  _fSampleNames[i] = GetFromManager<std::vector<std::string>>(param["Systematic"]["SampleNames"], {}, __FILE__, __LINE__);
304 
305  //Fill the map to get the correlations later as well
306  auto ParamType = Get<std::string>(param["Systematic"]["Type"], __FILE__ , __LINE__);
307  //Now load in variables for spline systematics only
308  if (ParamType == SystType_ToString(SystType::kSpline))
309  {
310  //Set param type
312  // Fill Spline info
313  SplineParams.push_back(GetSplineParameter(param["Systematic"], i));
314 
315  //Insert the mapping from the spline index i.e. the length of _fSplineNames etc
316  //to the Systematic index i.e. the counter for things like _fSampleID
317  _fSystToGlobalSystIndexMap[SystType::kSpline].insert(std::make_pair(ParamCounter[SystType::kSpline], i));
318  ParamCounter[SystType::kSpline]++;
319  } else if(ParamType == SystType_ToString(SystType::kNorm)) {
321  NormParams.push_back(GetNormParameter(param["Systematic"], i));
322  _fSystToGlobalSystIndexMap[SystType::kNorm].insert(std::make_pair(ParamCounter[SystType::kNorm], i));
323  ParamCounter[SystType::kNorm]++;
324  } else if(ParamType == SystType_ToString(SystType::kFunc)){
326  FuncParams.push_back(GetFunctionalParameters(param["Systematic"], i));
327  _fSystToGlobalSystIndexMap[SystType::kFunc].insert(std::make_pair(ParamCounter[SystType::kFunc], i));
328  ParamCounter[SystType::kFunc]++;
329  } else if(ParamType == SystType_ToString(SystType::kOsc)){
331  OscParams.push_back(GetOscillationParameters(param["Systematic"], i));
332  _fSystToGlobalSystIndexMap[SystType::kOsc].insert(std::make_pair(ParamCounter[SystType::kOsc], i));
333  ParamCounter[SystType::kOsc]++;
334  } else{
335  MACH3LOG_ERROR("Given unrecognised systematic type: {}", ParamType);
336  std::string expectedTypes = "Expecting ";
337  for (int s = 0; s < SystType::kSystTypes; ++s) {
338  if (s > 0) expectedTypes += ", ";
339  expectedTypes += SystType_ToString(static_cast<SystType>(s)) + "\"";
340  }
341  expectedTypes += ".";
342  MACH3LOG_ERROR(expectedTypes);
343  throw MaCh3Exception(__FILE__, __LINE__);
344  }
345  i++;
346  } //end loop over params
347 
348  //KS We resized them above to all params to fight memory fragmentation, now let's resize to fit only allocated memory to save RAM
349  NormParams.shrink_to_fit();
350  SplineParams.shrink_to_fit();
351  FuncParams.shrink_to_fit();
352  OscParams.shrink_to_fit();
353 }
354 
355 // ********************************************
357 // ********************************************
358  MACH3LOG_DEBUG("Deleting ParameterHandler");
359 }
360 
361 // ********************************************
362 // DB Grab the Spline Names for the relevant SampleName
363 const std::vector<std::string> ParameterHandlerGeneric::GetSplineParsNamesFromSampleName(const std::string& SampleName) const {
364 // ********************************************
365  std::vector<std::string> returnVec;
366  for (auto &pair : _fSystToGlobalSystIndexMap[SystType::kSpline]) {
367  auto &SplineIndex = pair.first;
368  auto &SystIndex = pair.second;
369  if (AppliesToSample(SystIndex, SampleName)) { //If parameter applies to required Sample
370  returnVec.push_back(SplineParams[SplineIndex]._fSplineNames);
371  }
372  }
373  return returnVec;
374 }
375 
376 // ********************************************
377 const std::vector<SplineInterpolation> ParameterHandlerGeneric::GetSplineInterpolationFromSampleName(const std::string& SampleName) const {
378 // ********************************************
379  std::vector<SplineInterpolation> returnVec;
380  for (auto &pair : _fSystToGlobalSystIndexMap[SystType::kSpline]) {
381  auto &SplineIndex = pair.first;
382  auto &SystIndex = pair.second;
383 
384  if (AppliesToSample(SystIndex, SampleName)) { //If parameter applies to required SampleID
385  returnVec.push_back(SplineParams.at(SplineIndex)._SplineInterpolationType);
386  }
387  }
388  return returnVec;
389 }
390 
391 // ********************************************
392 // DB Grab the Spline Modes for the relevant SampleName
393 const std::vector< std::vector<int> > ParameterHandlerGeneric::GetSplineModeVecFromSampleName(const std::string& SampleName) const {
394 // ********************************************
395  std::vector< std::vector<int> > returnVec;
396  //Need a counter or something to correctly get the index in _fSplineModes since it's not of length nPars
397  //Should probably just make a std::map<std::string, int> for param name to FD spline index
398  for (auto &pair : _fSystToGlobalSystIndexMap[SystType::kSpline]) {
399  auto &SplineIndex = pair.first;
400  auto &SystIndex = pair.second;
401  if (AppliesToSample(SystIndex, SampleName)) { //If parameter applies to required SampleID
402  returnVec.push_back(SplineParams.at(SplineIndex)._fSplineModes);
403  }
404  }
405  return returnVec;
406 }
407 
408 // ********************************************
409 // Get Norm params
410 NormParameter ParameterHandlerGeneric::GetNormParameter(const YAML::Node& param, const int Index) const {
411 // ********************************************
412  NormParameter norm;
413 
414  GetBaseParameter(param, Index, norm);
415 
416  // ETA size 0 to mean apply to all
417  // Ultimately all this information ends up in the NormParams vector
418  norm.modes = GetFromManager<std::vector<int>>(param["Mode"], {}, __FILE__ , __LINE__);
419  norm.pdgs = GetFromManager<std::vector<int>>(param["NeutrinoFlavour"], {}, __FILE__ , __LINE__);
420  norm.preoscpdgs = GetFromManager<std::vector<int>>(param["NeutrinoFlavourUnosc"], {}, __FILE__ , __LINE__);
421  norm.targets = GetFromManager<std::vector<int>>(param["TargetNuclei"], {}, __FILE__ , __LINE__);
422 
423  if(_fLowBound[Index] < 0.) {
424  MACH3LOG_ERROR("Normalisation Parameter {} ({}), has lower parameters bound which can go below 0 and is equal {}",
425  GetParFancyName(Index), Index, _fLowBound[Index]);
426  MACH3LOG_ERROR("Normalisation parameters can't go bellow 0 as this is unphysical");
427  throw MaCh3Exception(__FILE__, __LINE__);
428  }
429 
430  return norm;
431 }
432 
433 // ********************************************
434 // Get Base Param
435 void ParameterHandlerGeneric::GetBaseParameter(const YAML::Node& param, const int Index, TypeParameterBase& Parameter) const {
436 // ********************************************
437  Parameter.name = GetParFancyName(Index);
438 
439  // Set the global parameter index of the normalisation parameter
440  Parameter.index = Index;
441 
442  int NumKinematicCuts = 0;
443  if(param["KinematicCuts"]) {
444  NumKinematicCuts = int(param["KinematicCuts"].size());
445 
446  std::vector<std::string> TempKinematicStrings;
447  std::vector<std::vector<std::vector<double>>> TempKinematicBounds;
448  //First element of TempKinematicBounds is always -999, and size is then 3
449  for(int KinVar_i = 0 ; KinVar_i < NumKinematicCuts ; ++KinVar_i) {
450  //ETA: This is a bit messy, Kinematic cuts is a list of maps
451  for (YAML::const_iterator it = param["KinematicCuts"][KinVar_i].begin();it!=param["KinematicCuts"][KinVar_i].end();++it) {
452  TempKinematicStrings.push_back(it->first.as<std::string>());
453  TempKinematicBounds.push_back(Get2DBounds(it->second));
454  }
455  if(TempKinematicStrings.size() == 0) {
456  MACH3LOG_ERROR("Received a KinematicCuts node but couldn't read the contents (it's a list of single-element dictionaries (python) = map of pairs (C++))");
457  MACH3LOG_ERROR("For Param {}", Parameter.name);
458  throw MaCh3Exception(__FILE__, __LINE__);
459  }
460  }//KinVar_i
461  Parameter.KinematicVarStr = TempKinematicStrings;
462  Parameter.Selection = TempKinematicBounds;
463  }
464 
465  //Next ones are kinematic bounds on where normalisation parameter should apply
466  //We set a bool to see if any bounds exist so we can short-circuit checking all of them every step
467  bool HasKinBounds = false;
468 
469  if(Parameter.KinematicVarStr.size() > 0) HasKinBounds = true;
470 
471  Parameter.hasKinBounds = HasKinBounds;
472  //End of kinematic bound checking
473 }
474 
475 // ********************************************
476 // Grab the global syst index for the relevant SampleName
477 // i.e. get a vector of size nSplines where each entry is filled with the global syst number
478 const std::vector<int> ParameterHandlerGeneric::GetGlobalSystIndexFromSampleName(const std::string& SampleName, const SystType Type) const {
479 // ********************************************
480  std::vector<int> returnVec;
481  for (auto &pair : _fSystToGlobalSystIndexMap[Type]) {
482  auto &SystIndex = pair.second;
483  if (AppliesToSample(SystIndex, SampleName)) { //If parameter applies to required SampleID
484  returnVec.push_back(SystIndex);
485  }
486  }
487  return returnVec;
488 }
489 
490 // ********************************************
491 // Grab the global syst index for the relevant SampleName
492 // i.e. get a vector of size nSplines where each entry is filled with the global syst number
493 const std::vector<int> ParameterHandlerGeneric::GetSystIndexFromSampleName(const std::string& SampleName, const SystType Type) const {
494 // ********************************************
495  std::vector<int> returnVec;
496  for (auto &pair : _fSystToGlobalSystIndexMap[Type]) {
497  auto &SplineIndex = pair.first;
498  auto &systIndex = pair.second;
499  if (AppliesToSample(systIndex, SampleName)) { //If parameter applies to required SampleID
500  returnVec.push_back(SplineIndex);
501  }
502  }
503  return returnVec;
504 }
505 
506 // ********************************************
507 // Get Spline params
508 SplineParameter ParameterHandlerGeneric::GetSplineParameter(const YAML::Node& param, const int Index) const {
509 // ********************************************
510  SplineParameter Spline;
511 
512  GetBaseParameter(param, Index, Spline);
513  auto& SplinePar = param["SplineInformation"];
514  //Now get the Spline interpolation type
515  if (SplinePar["InterpolationType"]) {
516  for(int InterpType = 0; InterpType < kSplineInterpolations ; ++InterpType){
517  if(SplinePar["InterpolationType"].as<std::string>() == SplineInterpolation_ToString(SplineInterpolation(InterpType)))
518  Spline._SplineInterpolationType = SplineInterpolation(InterpType);
519  }
520  } else { //KS: By default use TSpline3
522  }
523 
524  Spline._fSplineNames = Get<std::string>(SplinePar["SplineName"], __FILE__ , __LINE__);
525  Spline._SplineKnotUpBound = GetFromManager<double>(SplinePar["SplineKnotUpBound"], M3::DefSplineKnotUpBound, __FILE__ , __LINE__);
526  Spline._SplineKnotLowBound = GetFromManager<double>(SplinePar["SplineKnotLowBound"], M3::DefSplineKnotLowBound, __FILE__ , __LINE__);
527 
529  MACH3LOG_WARN("Spline knot capping enabled with bounds [{}, {}]. For reliable fits, consider modifying the input generation instead.",
530  Spline._SplineKnotLowBound, Spline._SplineKnotUpBound);
531  }
532  //If there is no mode information given then this will be an empty vector
533  Spline._fSplineModes = GetFromManager(SplinePar["Mode"], std::vector<int>(), __FILE__ , __LINE__);
534 
535  return Spline;
536 }
537 
538 // ********************************************
539 bool ParameterHandlerGeneric::AppliesToSample(const int SystIndex, const std::string& SampleName) const {
540 // ********************************************
541  // Empty means apply to all
542  if (_fSampleNames[SystIndex].size() == 0) return true;
543 
544  // Make a copy and to lower case to not be case sensitive
545  std::string SampleNameCopy = SampleName;
546  std::transform(SampleNameCopy.begin(), SampleNameCopy.end(), SampleNameCopy.begin(), ::tolower);
547 
548  // Check for unsupported wildcards in SampleNameCopy
549  if (SampleNameCopy.find('*') != std::string::npos) {
550  MACH3LOG_ERROR("Wildcards ('*') are not supported in sample name: '{}'", SampleName);
551  throw MaCh3Exception(__FILE__ , __LINE__ );
552  }
553 
554  bool Applies = false;
555 
556  for (size_t i = 0; i < _fSampleNames[SystIndex].size(); i++) {
557  // Convert to low case to not be case sensitive
558  std::string pattern = _fSampleNames[SystIndex][i];
559  std::transform(pattern.begin(), pattern.end(), pattern.begin(), ::tolower);
560 
561  // Replace '*' in the pattern with '.*' for regex matching
562  std::string regexPattern = "^" + std::regex_replace(pattern, std::regex("\\*"), ".*") + "$";
563  try {
564  std::regex regex(regexPattern);
565  if (std::regex_match(SampleNameCopy, regex)) {
566  Applies = true;
567  break;
568  }
569  } catch (const std::regex_error& e) {
570  // Handle regex error (for invalid patterns)
571  MACH3LOG_ERROR("Regex error: {}", e.what());
572  }
573  }
574  return Applies;
575 }
576 
577 // ********************************************
578 // Get Func params
579 FunctionalParameter ParameterHandlerGeneric::GetFunctionalParameters(const YAML::Node& param, const int Index) const {
580 // ********************************************
581  FunctionalParameter func;
582  GetBaseParameter(param, Index, func);
583 
584  func.modes = GetFromManager<std::vector<int>>(param["Mode"], std::vector<int>(), __FILE__ , __LINE__);
585 
586  func.valuePtr = RetPointer(Index);
587  return func;
588 }
589 
590 // ********************************************
591 // Get Osc params
592 OscillationParameter ParameterHandlerGeneric::GetOscillationParameters(const YAML::Node& param, const int Index) const {
593 // ********************************************
594  OscillationParameter OscParamInfo;
595  GetBaseParameter(param, Index, OscParamInfo);
596 
597  return OscParamInfo;
598 }
599 
600 // ********************************************
601 // HH: Grab the Functional parameters for the relevant SampleName
602 const std::vector<FunctionalParameter> ParameterHandlerGeneric::GetFunctionalParametersFromSampleName(const std::string& SampleName) const {
603 // ********************************************
605 }
606 
607 // ********************************************
608 // DB Grab the Normalisation parameters for the relevant SampleName
609 const std::vector<NormParameter> ParameterHandlerGeneric::GetNormParsFromSampleName(const std::string& SampleName) const {
610 // ********************************************
612 }
613 
614 // ********************************************
615 // KS Grab the Spline parameters for the relevant SampleName
616 const std::vector<SplineParameter> ParameterHandlerGeneric::GetSplineParsFromSampleName(const std::string& SampleName) const {
617 // ********************************************
619 }
620 
621 // ********************************************
622 template<typename ParamT>
623 std::vector<ParamT> ParameterHandlerGeneric::GetTypeParamsFromSampleName(const std::map<int, int>& indexMap,
624  const std::vector<ParamT>& params, const std::string& SampleName) const {
625 // ********************************************
626  std::vector<ParamT> returnVec;
627  for (const auto& pair : indexMap) {
628  const auto& localIndex = pair.first;
629  const auto& globalIndex = pair.second;
630  if (AppliesToSample(globalIndex, SampleName)) {
631  returnVec.push_back(params[localIndex]);
632  }
633  }
634  return returnVec;
635 }
636 
637 // ********************************************
638 // DB Grab the number of parameters for the relevant SampleName
639 int ParameterHandlerGeneric::GetNumParamsFromSampleName(const std::string& SampleName, const SystType Type) const {
640 // ********************************************
641  int returnVal = 0;
642  IterateOverParams(SampleName,
643  [&](int i) { return GetParamType(i) == Type; }, // Filter condition
644  [&](int) { returnVal += 1; } // Action to perform if filter passes
645  );
646  return returnVal;
647 }
648 
649 // ********************************************
650 // DB Grab the parameter names for the relevant SampleName
651 const std::vector<std::string> ParameterHandlerGeneric::GetParsNamesFromSampleName(const std::string& SampleName, const SystType Type) const {
652 // ********************************************
653  std::vector<std::string> returnVec;
654  IterateOverParams(SampleName,
655  [&](int i) { return GetParamType(i) == Type; }, // Filter condition
656  [&](int i) { returnVec.push_back(GetParFancyName(i)); } // Action to perform if filter passes
657  );
658  return returnVec;
659 }
660 
661 // ********************************************
662 // DB DB Grab the parameter indices for the relevant SampleName
663 const std::vector<int> ParameterHandlerGeneric::GetParsIndexFromSampleName(const std::string& SampleName, const SystType Type) const {
664 // ********************************************
665  std::vector<int> returnVec;
666  IterateOverParams(SampleName,
667  [&](int i) { return GetParamType(i) == Type; }, // Filter condition
668  [&](int i) { returnVec.push_back(i); } // Action to perform if filter passes
669  );
670  return returnVec;
671 }
672 
673 // ********************************************
674 template <typename FilterFunc, typename ActionFunc>
675 void ParameterHandlerGeneric::IterateOverParams(const std::string& SampleName, FilterFunc filter, ActionFunc action) const {
676 // ********************************************
677  for (int i = 0; i < _fNumPar; ++i) {
678  if ((AppliesToSample(i, SampleName)) && filter(i)) { // Common filter logic
679  action(i); // Specific action for each function
680  }
681  }
682 }
683 
684 // ********************************************
686 // ********************************************
687  for (int i = 0; i < _fNumPar; ++i) {
688  //ETA - set the name to be param_% as this is what ProcessorMCMC expects
689  _fNames[i] = "param_"+std::to_string(i);
690 
691  // KS: Plenty of MaCh3 processing script rely on osc params having "fancy name" this is to maintain backward compatibility with them
692  if(_fParamType[i] == kOsc) {
693  _fNames[i] = _fFancyNames[i];
694 
695  if(_ParameterGroup[i] != "Osc"){
696  MACH3LOG_ERROR("Parameter {}, is of type Oscillation but doesn't belong to Osc group", _fFancyNames[i]);
697  MACH3LOG_ERROR("It belongs to {} group", _ParameterGroup[i]);
698  throw MaCh3Exception(__FILE__ , __LINE__ );
699  }
700  }
701  #pragma GCC diagnostic push
702  #pragma GCC diagnostic ignored "-Wuseless-cast"
703  // Set ParameterHandler parameters (Curr = current, Prop = proposed, Sigma = step)
704  _fCurrVal[i] = _fPreFitValue[i];
705  _fPropVal[i] = static_cast<M3::float_t>(_fCurrVal[i]);
706  #pragma GCC diagnostic pop
707  }
708  Randomize();
709  //KS: Transfer the starting parameters to the PCA basis, you don't want to start with zero..
710  if (pca) {
711  PCAObj->SetInitialParameters();
712  }
713 }
714 
715 // ********************************************
716 // Print everything we know about the inputs we're Getting
718 // ********************************************
719  MACH3LOG_INFO("#################################################");
720  MACH3LOG_INFO("Printing ParameterHandlerGeneric:");
721 
723 
724  PrintNormParams();
725 
727 
729 
731 
733 
734  MACH3LOG_INFO("Finished");
735  MACH3LOG_INFO("#################################################");
736 
738 } // End
739 
740 // ********************************************
742 // ********************************************
743  MACH3LOG_INFO("============================================================================================================================================================");
744  MACH3LOG_INFO("{:<5} {:2} {:<40} {:2} {:<10} {:2} {:<10} {:2} {:<10} {:2} {:<10} {:2} {:<10} {:2} {:<20} {:2} {:<10}", "#", "|", "Name", "|", "Prior", "|", "Error", "|", "Lower", "|", "Upper", "|", "StepScale", "|", "SampleNames", "|", "Type");
745  MACH3LOG_INFO("------------------------------------------------------------------------------------------------------------------------------------------------------------");
746  for (int i = 0; i < GetNumParams(); i++) {
747  std::string ErrString = fmt::format("{:.2f}", _fError[i]);
748  std::string SampleNameString = "";
749  for (const auto& SampleName : _fSampleNames[i]) {
750  if (!SampleNameString.empty()) {
751  SampleNameString += ", ";
752  }
753  SampleNameString += SampleName;
754  }
755  MACH3LOG_INFO("{:<5} {:2} {:<40} {:2} {:<10} {:2} {:<10} {:2} {:<10} {:2} {:<10} {:2} {:<10} {:2} {:<20} {:2} {:<10}", i, "|", GetParFancyName(i), "|", _fPreFitValue[i], "|", "+/- " + ErrString, "|", _fLowBound[i], "|", _fUpBound[i], "|", _fIndivStepScale[i], "|", SampleNameString, "|", SystType_ToString(_fParamType[i]));
756  }
757  MACH3LOG_INFO("============================================================================================================================================================");
758 }
759 
760 // ********************************************
762 // ********************************************
763  // Output the normalisation parameters as a sanity check!
764  MACH3LOG_INFO("Normalisation parameters: {}", NormParams.size());
765  if(_fSystToGlobalSystIndexMap[SystType::kNorm].size() == 0) return;
766 
767  bool have_parameter_with_kin_bounds = false;
768 
769  //KS: Consider making some class producing table..
770  MACH3LOG_INFO("┌────┬──────────┬────────────────────────────────────────┬────────────────────┬────────────────────┬────────────────────┐");
771  MACH3LOG_INFO("│{0:4}│{1:10}│{2:40}│{3:20}│{4:20}│{5:20}│", "#", "Global #", "Name", "Int. mode", "Target", "pdg");
772  MACH3LOG_INFO("├────┼──────────┼────────────────────────────────────────┼────────────────────┼────────────────────┼────────────────────┤");
773 
774  for (unsigned int i = 0; i < NormParams.size(); ++i)
775  {
776  std::string intModeString;
777  for (unsigned int j = 0; j < NormParams[i].modes.size(); j++) {
778  intModeString += std::to_string(NormParams[i].modes[j]);
779  intModeString += " ";
780  }
781  if (NormParams[i].modes.empty()) intModeString += "all";
782 
783  std::string targetString;
784  for (unsigned int j = 0; j < NormParams[i].targets.size(); j++) {
785  targetString += std::to_string(NormParams[i].targets[j]);
786  targetString += " ";
787  }
788  if (NormParams[i].targets.empty()) targetString += "all";
789 
790  std::string pdgString;
791  for (unsigned int j = 0; j < NormParams[i].pdgs.size(); j++) {
792  pdgString += std::to_string(NormParams[i].pdgs[j]);
793  pdgString += " ";
794  }
795  if (NormParams[i].pdgs.empty()) pdgString += "all";
796 
797  MACH3LOG_INFO("│{: <4}│{: <10}│{: <40}│{: <20}│{: <20}│{: <20}│", i, NormParams[i].index, NormParams[i].name, intModeString, targetString, pdgString);
798 
799  if(NormParams[i].hasKinBounds) have_parameter_with_kin_bounds = true;
800  }
801  MACH3LOG_INFO("└────┴──────────┴────────────────────────────────────────┴────────────────────┴────────────────────┴────────────────────┘");
802 
803  if(have_parameter_with_kin_bounds) {
804  MACH3LOG_INFO("Normalisation parameters KinematicCuts information");
805  MACH3LOG_INFO("┌────┬──────────┬────────────────────────────────────────┬────────────────────┬────────────────────────────────────────┐");
806  MACH3LOG_INFO("│{0:4}│{1:10}│{2:40}│{3:20}│{4:40}│", "#", "Global #", "Name", "KinematicCut", "Value");
807  MACH3LOG_INFO("├────┼──────────┼────────────────────────────────────────┼────────────────────┼────────────────────────────────────────┤");
808  for (unsigned int i = 0; i < NormParams.size(); ++i)
809  {
810  //skip parameters with no KinematicCuts
811  if(!NormParams[i].hasKinBounds) continue;
812 
813  const long unsigned int ncuts = NormParams[i].KinematicVarStr.size();
814  for(long unsigned int icut = 0; icut < ncuts; icut++) {
815  std::string kinematicCutValueString;
816  for(const auto & value : NormParams[i].Selection[icut]) {
817  for (const auto& v : value) {
818  kinematicCutValueString += fmt::format("{:.2f} ", v);
819  }
820  }
821  if(icut == 0)
822  MACH3LOG_INFO("│{: <4}│{: <10}│{: <40}│{: <20}│{: <40}│", i, NormParams[i].index, NormParams[i].name, NormParams[i].KinematicVarStr[icut], kinematicCutValueString);
823  else
824  MACH3LOG_INFO("│{: <4}│{: <10}│{: <40}│{: <20}│{: <40}│", "", "", "", NormParams[i].KinematicVarStr[icut], kinematicCutValueString);
825  }//icut
826  }//i
827  MACH3LOG_INFO("└────┴──────────┴────────────────────────────────────────┴────────────────────┴────────────────────────────────────────┘");
828  }
829  else
830  MACH3LOG_INFO("No normalisation parameters have KinematicCuts defined");
831 }
832 
833 // ********************************************
835 // ********************************************
836  MACH3LOG_INFO("Spline parameters: {}", _fSystToGlobalSystIndexMap[SystType::kSpline].size());
837  if(_fSystToGlobalSystIndexMap[SystType::kSpline].size() == 0) return;
838  MACH3LOG_INFO("=====================================================================================================================================================================");
839  MACH3LOG_INFO("{:<4} {:<2} {:<40} {:<2} {:<40} {:<2} {:<20} {:<2} {:<20} {:<2} {:<20} {:<2}", "#", "|", "Name", "|", "Spline Name", "|", "Spline Interpolation", "|", "Low Knot Bound", "|", "Up Knot Bound", "|");
840  MACH3LOG_INFO("---------------------------------------------------------------------------------------------------------------------------------------------------------------------");
841  for (auto &pair : _fSystToGlobalSystIndexMap[SystType::kSpline]) {
842  auto &SplineIndex = pair.first;
843  auto &GlobalIndex = pair.second;
844 
845  MACH3LOG_INFO("{:<4} {:<2} {:<40} {:<2} {:<40} {:<2} {:<20} {:<2} {:<20} {:<2} {:<20} {:<2}",
846  SplineIndex, "|", GetParFancyName(GlobalIndex), "|",
847  SplineParams[SplineIndex]._fSplineNames, "|",
851  }
852  MACH3LOG_INFO("=====================================================================================================================================================================");
853 }
854 
855 // ********************************************
857 // ********************************************
858  MACH3LOG_INFO("Functional parameters: {}", _fSystToGlobalSystIndexMap[SystType::kFunc].size());
859  if(_fSystToGlobalSystIndexMap[SystType::kFunc].size() == 0) return;
860  MACH3LOG_INFO("┌────┬──────────┬────────────────────────────────────────┐");
861  MACH3LOG_INFO("│{0:4}│{1:10}│{2:40}│", "#", "Global #", "Name");
862  MACH3LOG_INFO("├────┼──────────┼────────────────────────────────────────┤");
863  for (auto &pair : _fSystToGlobalSystIndexMap[SystType::kFunc]) {
864  auto &FuncIndex = pair.first;
865  auto &GlobalIndex = pair.second;
866  MACH3LOG_INFO("│{0:4}│{1:<10}│{2:40}│", std::to_string(FuncIndex), GlobalIndex, GetParFancyName(GlobalIndex));
867  }
868  MACH3LOG_INFO("└────┴──────────┴────────────────────────────────────────┘");
869 }
870 
871 // ********************************************
873 // ********************************************
874  MACH3LOG_INFO("Oscillation parameters: {}", _fSystToGlobalSystIndexMap[SystType::kOsc].size());
875  if(_fSystToGlobalSystIndexMap[SystType::kOsc].size() == 0) return;
876  MACH3LOG_INFO("┌────┬──────────┬────────────────────────────────────────┐");
877  MACH3LOG_INFO("│{0:4}│{1:10}│{2:40}│", "#", "Global #", "Name");
878  MACH3LOG_INFO("├────┼──────────┼────────────────────────────────────────┤");
879  for (auto &pair : _fSystToGlobalSystIndexMap[SystType::kOsc]) {
880  auto &OscIndex = pair.first;
881  auto &GlobalIndex = pair.second;
882  MACH3LOG_INFO("│{0:4}│{1:<10}│{2:40}│", std::to_string(OscIndex), GlobalIndex, GetParFancyName(GlobalIndex));
883  }
884  MACH3LOG_INFO("└────┴──────────┴────────────────────────────────────────┘");
885 }
886 
887 // ********************************************
889 // ********************************************
890  // KS: Create a map to store the counts of unique strings, in principle this could be in header file
891  std::unordered_map<std::string, int> paramCounts;
892 
893  std::for_each(_ParameterGroup.begin(), _ParameterGroup.end(),
894  [&paramCounts](const std::string& param) {
895  paramCounts[param]++;
896  });
897 
898  MACH3LOG_INFO("Printing parameter groups");
899  // Output the counts
900  for (const auto& pair : paramCounts) {
901  MACH3LOG_INFO("Found {}: {} params", pair.second, pair.first);
902  }
903 }
904 
905 // ********************************************
906 std::vector<std::string> ParameterHandlerGeneric::GetUniqueParameterGroups() const {
907 // ********************************************
908  std::unordered_set<std::string> uniqueGroups;
909 
910  // Fill the set with unique values
911  for (const auto& param : _ParameterGroup) {
912  uniqueGroups.insert(param);
913  }
914 
915  // Convert to vector and return
916  std::vector<std::string> result(uniqueGroups.begin(), uniqueGroups.end());
917  return result;
918 }
919 
920 // ********************************************
921 // KS: Check if matrix is correctly initialised
923 // ********************************************
924  // KS: Lambda Function which simply checks if there are no duplicates in std::vector
925  auto CheckForDuplicates = [](const std::vector<std::string>& names, const std::string& nameType, bool Warning) {
926  std::unordered_map<std::string, size_t> seenStrings;
927  for (size_t i = 0; i < names.size(); ++i) {
928  const auto& name = names[i];
929  if (seenStrings.find(name) != seenStrings.end()) {
930  size_t firstIndex = seenStrings[name];
931  if(Warning){
932  MACH3LOG_WARN("There are two systematics with the same {} '{}', first at index {}, and again at index {}", nameType, name, firstIndex, i);
933  MACH3LOG_WARN("Is this desired?");
934  return;
935  } else {
936  MACH3LOG_CRITICAL("There are two systematics with the same {} '{}', first at index {}, and again at index {}", nameType, name, firstIndex, i);
937  throw MaCh3Exception(__FILE__, __LINE__);
938  }
939  }
940  seenStrings[name] = i;
941  }
942  };
943  std::vector<std::string> SplineNameTemp(SplineParams.size());
944  for(size_t it = 0; it < SplineParams.size(); it++){
945  SplineNameTemp[it] = SplineParams[it]._fSplineNames;
946  }
947  // KS: Checks if there are no duplicates in fancy names etc, this can happen if we merge configs etc
948  CheckForDuplicates(_fFancyNames, "_fFancyNames", false);
949  CheckForDuplicates(SplineNameTemp, "_fSplineNames", true);
950 }
951 
952 // ********************************************
953 // Function to set to prior parameters of a given group
954 void ParameterHandlerGeneric::SetGroupOnlyParameters(const std::vector< std::string>& Groups) {
955 // ********************************************
956  for(size_t i = 0; i < Groups.size(); i++){
957  SetGroupOnlyParameters(Groups[i]);
958  }
959 }
960 
961 // ********************************************
962 // Function to set to prior parameters of a given group
963 void ParameterHandlerGeneric::SetGroupOnlyParameters(const std::string& Group, const std::vector<double>& Pars) {
964 // ********************************************
965  #pragma GCC diagnostic push
966  #pragma GCC diagnostic ignored "-Wuseless-cast"
967  // If empty, set the proposed to prior
968  if (Pars.empty()) {
969  for (int i = 0; i < _fNumPar; i++) {
970  if(IsParFromGroup(i, Group)) _fPropVal[i] = static_cast<M3::float_t>(_fPreFitValue[i]);
971  }
972  } else{
973  const size_t ExpectedSize = static_cast<size_t>(GetNumParFromGroup(Group));
974  if (Pars.size() != ExpectedSize) {
975  MACH3LOG_ERROR("Number of param in group {} is {}, while you passed {}", Group, ExpectedSize, Pars.size());
976  throw MaCh3Exception(__FILE__, __LINE__);
977  }
978  int Counter = 0;
979  for (int i = 0; i < _fNumPar; i++) {
980  // If belongs to group set value from parsed vector, otherwise use propose value
981  if(IsParFromGroup(i, Group)){
982  _fPropVal[i] = static_cast<M3::float_t>(Pars[Counter]);
983  Counter++;
984  }
985  }
986  }
987  // And if pca make the transfer
988  if (pca) {
989  PCAObj->TransferToPCA();
990  PCAObj->TransferToParam();
991  }
992  #pragma GCC diagnostic pop
993 }
994 
995 // ********************************************
996 // Set parameters to be fixed in a given group
998 // ********************************************
999  for (int i = 0; i < _fNumPar; ++i)
1000  if(IsParFromGroup(i, Group)) SetFixParameter(i);
1001 }
1002 
1003 // ********************************************
1004 // Set parameters of several groups to be fixed
1005 void ParameterHandlerGeneric::SetFixGroupOnlyParameters(const std::vector< std::string>& Groups) {
1006 // ********************************************
1007  for(size_t i = 0; i < Groups.size(); i++)
1008  SetFixGroupOnlyParameters(Groups[i]);
1009 }
1010 
1011 // ********************************************
1012 // Set parameters to be free in a given group
1014 // ********************************************
1015  for (int i = 0; i < _fNumPar; ++i)
1016  if(IsParFromGroup(i, Group)) SetFreeParameter(i);
1017 }
1018 
1019 // ********************************************
1020 // Set parameters of several groups to be fixed
1021 void ParameterHandlerGeneric::SetFreeGroupOnlyParameters(const std::vector< std::string>& Groups) {
1022 // ********************************************
1023  for(size_t i = 0; i < Groups.size(); i++)
1024  SetFreeGroupOnlyParameters(Groups[i]);
1025 }
1026 
1027 // ********************************************
1028 // Checks if parameter belongs to a given group
1029 bool ParameterHandlerGeneric::IsParFromGroup(const int i, const std::string& Group) const {
1030 // ********************************************
1031  std::string groupLower = Group;
1032  std::string paramGroupLower = _ParameterGroup[i];
1033 
1034  // KS: Convert both strings to lowercase, this way comparison will be case insensitive
1035  std::transform(groupLower.begin(), groupLower.end(), groupLower.begin(), ::tolower);
1036  std::transform(paramGroupLower.begin(), paramGroupLower.end(), paramGroupLower.begin(), ::tolower);
1037 
1038  return groupLower == paramGroupLower;
1039 }
1040 
1041 // ********************************************
1042 int ParameterHandlerGeneric::GetNumParFromGroup(const std::string& Group) const {
1043 // ********************************************
1044  int Counter = 0;
1045  for (int i = 0; i < _fNumPar; i++) {
1046  if(IsParFromGroup(i, Group)) Counter++;
1047  }
1048  return Counter;
1049 }
1050 
1051 // ********************************************
1052 // DB Grab the Normalisation parameters for the relevant sample name
1053 std::vector<const M3::float_t*> ParameterHandlerGeneric::GetOscParsFromSampleName(const std::string& SampleName) const {
1054 // ********************************************
1055  std::vector<const M3::float_t*> returnVec;
1056  for (const auto& pair : _fSystToGlobalSystIndexMap[SystType::kOsc]) {
1057  const auto& globalIndex = pair.second;
1058  if (AppliesToSample(globalIndex, SampleName)) {
1059  returnVec.push_back(RetPointer(globalIndex));
1060  }
1061  }
1062  return returnVec;
1063 }
1064 
1065 // ********************************************
1066 // Dump Matrix to ROOT file, useful when we need to pass matrix info to another fitting group
1067 void ParameterHandlerGeneric::DumpMatrixToFile(const std::string& Name) {
1068 // ********************************************
1069  // KS: Ideally we remove it eventually...
1070  TH2D* CorrMatrix = GetCorrelationMatrix();
1072  _fPreFitValue,
1073  _fError,
1074  _fLowBound,
1075  _fUpBound,
1077  _fFancyNames,
1078  _fFlatPrior,
1079  SplineParams,
1080  covMatrix,
1081  CorrMatrix,
1082  Name);
1083  delete CorrMatrix;
1084  MACH3LOG_INFO("Finished dumping ParameterHandler object");
1085 }
#define MACH3LOG_CRITICAL
Definition: MaCh3Logger.h:38
#define MACH3LOG_DEBUG
Definition: MaCh3Logger.h:34
#define MACH3LOG_ERROR
Definition: MaCh3Logger.h:37
#define MACH3LOG_INFO
Definition: MaCh3Logger.h:35
#define MACH3LOG_WARN
Definition: MaCh3Logger.h:36
std::string SplineInterpolation_ToString(const SplineInterpolation i)
Convert a LLH type to a string.
SplineInterpolation
Make an enum of the spline interpolation type.
@ kTSpline3
Default TSpline3 interpolation.
@ kSplineInterpolations
This only enumerates.
std::string SystType_ToString(const SystType i)
Convert a Syst type type to a string.
SystType
@ kNorm
For normalisation parameters.
@ kSpline
For splined parameters (1D)
@ kSystTypes
This only enumerates.
@ kOsc
For oscillation parameters.
@ kFunc
For functional parameters.
#define Get2DBounds(filename)
Definition: YamlHelper.h:591
Type GetFromManager(const YAML::Node &node, const Type defval, const std::string File="", const int Line=1)
Get content of config file if node is not found take default value specified.
Definition: YamlHelper.h:329
#define M3OpenConfig(filename)
Macro to simplify calling LoadYaml with file and line info.
Definition: YamlHelper.h:589
#define GetBounds(filename)
Definition: YamlHelper.h:590
Custom exception class used throughout MaCh3.
int GetNumParams() const
Get total number of parameters.
bool use_adaptive
Are we using AMCMC?
std::string inputFile
The input root file we read in.
TH2D * GetCorrelationMatrix() const
KS: Convert covariance matrix to correlation matrix and return TH2D which can be used for fancy plott...
void SetCovMatrix(TMatrixDSym *cov)
Set covariance matrix.
std::vector< M3::float_t > _fPropVal
Proposed value of the parameter.
std::unique_ptr< ParameterTunes > Tunes
Struct containing information about adaption.
void ReserveMemory(const int size)
Initialise vectors with parameters information.
std::vector< std::string > _fFancyNames
Fancy name for example rather than param_0 it is MAQE, useful for human reading.
bool doSpecialStepProposal
Check if any of special step proposal were enabled.
std::vector< bool > _fFlatPrior
Whether to apply flat prior or not.
std::string matrixName
Name of cov matrix.
std::string GetParFancyName(const int i) const
Get fancy name of the Parameter.
std::vector< double > _fError
Prior error on the parameter.
void SetFixParameter(const int i)
Set parameter to be fixed at prior value.
YAML::Node _fYAMLDoc
Stores config describing systematics.
std::unique_ptr< PCAHandler > PCAObj
Struct containing information about PCA.
int _fNumPar
Number of systematic parameters.
void ConstructPCA(const double eigen_threshold, int FirstPCAdpar, int LastPCAdpar)
CW: Calculate eigen values, prepare transition matrices and remove param based on defined threshold.
std::vector< double > _fLowBound
Lowest physical bound, parameter will not be able to go beyond it.
std::vector< double > _fCurrVal
Current value of the parameter.
const M3::float_t * RetPointer(const int iParam) const
DB Pointer return to param position.
TMatrixDSym * covMatrix
The covariance matrix.
std::vector< double > _fPreFitValue
Parameter value dictated by the prior model. Based on it penalty term is calculated.
void Randomize() _noexcept_
"Randomize" the parameters in the covariance class for the proposed step. Used the proposal kernel an...
std::vector< std::string > _fNames
ETA _fNames is set automatically in the covariance class to be something like param_i,...
void EnableSpecialProposal(const YAML::Node &param, const int Index)
Enable special proposal.
void SetFreeParameter(const int i)
Set parameter to be treated as free.
std::vector< double > _fUpBound
Upper physical bound, parameter will not be able to go beyond it.
void SetSubThrowMatrix(int first_index, int last_index, TMatrixDSym const &subcov)
bool pca
perform PCA or not
std::vector< double > _fIndivStepScale
Individual step scale used by MCMC algorithm.
std::vector< std::unique_ptr< TRandom3 > > random_number
KS: Set Random numbers for each thread so each thread has different seed.
void MakePosDef(TMatrixDSym *cov=nullptr)
Make matrix positive definite by adding small values to diagonal, necessary for inverting matrix.
int PrintLength
KS: This is used when printing parameters, sometimes we have super long parameters name,...
void PrintGlobablInfo() const
Prints general information about the ParameterHandler object.
void Print() const
Print information about the whole object once it is set.
std::vector< SplineParameter > SplineParams
Vector containing info for normalisation systematics.
FunctionalParameter GetFunctionalParameters(const YAML::Node &param, const int Index) const
Get Func params.
const std::vector< NormParameter > GetNormParsFromSampleName(const std::string &SampleName) const
DB Get norm/func parameters depending on given SampleName.
std::vector< std::string > GetUniqueParameterGroups() const
KS: Get names of all unique parameter groups.
ParameterHandlerGeneric(const std::vector< std::string > &FileNames, std::string name="xsec_cov", double threshold=-1, int FirstPCAdpar=-999, int LastPCAdpar=-999)
Constructor.
void PrintSplineParams() const
Prints spline parameters.
std::vector< const M3::float_t * > GetOscParsFromSampleName(const std::string &SampleName) const
Get pointers to Osc params from Sample name.
void IterateOverParams(const std::string &SampleName, FilterFunc filter, ActionFunc action) const
Iterates over parameters and applies a filter and action function.
const std::vector< int > GetParsIndexFromSampleName(const std::string &SampleName, const SystType Type) const
DB Grab the parameter indices for the relevant SampleName.
SplineParameter GetSplineParameter(const YAML::Node &param, const int Index) const
Get Spline params.
SystType GetParamType(const int i) const
Returns enum describing our param type.
const std::vector< int > GetSystIndexFromSampleName(const std::string &SampleName, const SystType Type) const
Grab the index of the syst relative to global numbering.
void GetBaseParameter(const YAML::Node &param, const int Index, TypeParameterBase &Parameter) const
Fill base parameters.
std::vector< NormParameter > NormParams
Vector containing info for normalisation systematics.
double GetParSplineKnotUpperBound(const int i) const
EM: value at which we cap spline knot weight.
std::vector< std::map< int, int > > _fSystToGlobalSystIndexMap
Map between number of given parameter type with global parameter numbering. For example 2nd norm para...
void LoadCorrelationFromConfig(std::vector< std::map< std::string, double >> &Correlations, std::map< std::string, int > &CorrNamesMap)
Load correlation from yaml config.
void SetGroupOnlyParameters(const std::string &Group, const std::vector< double > &Pars={})
KS Function to set to prior parameters of a given group or values from vector.
std::vector< SystType > _fParamType
Type of parameter like norm, spline etc.
const std::vector< SplineParameter > GetSplineParsFromSampleName(const std::string &SampleName) const
KS: Grab the Spline parameters for the relevant SampleName.
void PrintOscillationParams() const
Prints oscillation parameters.
const std::vector< FunctionalParameter > GetFunctionalParametersFromSampleName(const std::string &SampleName) const
HH Get functional parameters for the relevant SampleName.
void SetFreeGroupOnlyParameters(const std::string &Group)
TN Method to set parameters within a group to be treated as free.
void InitParametersTypeFromConfig()
Parses the YAML configuration to set up cross-section parameters. The YAML file defines the types of ...
std::vector< std::vector< std::string > > _fSampleNames
Tells to which samples object param should be applied.
std::vector< ParamT > GetTypeParamsFromSampleName(const std::map< int, int > &indexMap, const std::vector< ParamT > &params, const std::string &SampleName) const
Retrieve parameters that apply to a given sample name.
std::vector< FunctionalParameter > FuncParams
Vector containing info for functional systematics.
int GetNumParamsFromSampleName(const std::string &SampleName, const SystType Type) const
DB Grab the number of parameters for the relevant SampleName.
void PrintFunctionalParams() const
Prints functional parameters.
const std::vector< std::string > GetParsNamesFromSampleName(const std::string &SampleName, const SystType Type) const
DB Grab the parameter names for the relevant SampleName.
std::vector< OscillationParameter > OscParams
Vector containing info for functional systematics.
void DumpMatrixToFile(const std::string &Name)
Dump Matrix to ROOT file, useful when we need to pass matrix info to another fitting group.
const std::vector< std::string > GetSplineParsNamesFromSampleName(const std::string &SampleName) const
DB Get spline parameters depending on given SampleName.
SplineInterpolation GetParSplineInterpolation(const int i) const
Get interpolation type for a given parameter.
void InitParameters()
Initializes the systematic parameters from the configuration file. This function loads parameters lik...
const std::vector< int > GetGlobalSystIndexFromSampleName(const std::string &SampleName, const SystType Type) const
DB Get spline parameters depending on given SampleName.
OscillationParameter GetOscillationParameters(const YAML::Node &param, const int Index) const
Get Osc params.
NormParameter GetNormParameter(const YAML::Node &param, const int Index) const
Get Norm params.
bool AppliesToSample(const int SystIndex, const std::string &SampleName) const
Check if parameter is affecting given sample name.
const std::vector< SplineInterpolation > GetSplineInterpolationFromSampleName(const std::string &SampleName) const
Get the interpolation types for splines affecting a particular SampleName.
bool IsParFromGroup(const int i, const std::string &Group) const
Checks if parameter belongs to a given group.
void CheckCorrectInitialisation() const
KS: Check if matrix is correctly initialised.
void InitialiseFromConfig(const std::vector< std::string > &YAMLFile)
Initialisation of the class using config.
void SetFixGroupOnlyParameters(const std::string &Group)
TN Method to set parameters within a group to be fixed to their prior values.
std::vector< std::string > _ParameterGroup
KS: Allow to group parameters for example to affect only cross-section or only flux etc.
void LoadAndMergeYAML(const std::vector< std::string > &YAMLFile, std::map< std::pair< int, int >, std::unique_ptr< TMatrixDSym >> &overrides)
Initialise single yaml based on several yaml files.
double GetParSplineKnotLowerBound(const int i) const
EM: value at which we cap spline knot weight.
void PrintNormParams() const
Prints normalization parameters.
int GetNumParFromGroup(const std::string &Group) const
KS: Check how many parameters are associated with given group.
const std::vector< std::vector< int > > GetSplineModeVecFromSampleName(const std::string &SampleName) const
DB Grab the Spline Modes for the relevant SampleName.
void PrintParameterGroups() const
Prints groups of parameters.
constexpr static const double DefSplineKnotUpBound
Default value for spline knot capping, default mean not capping is being applied.
Definition: Core.h:86
double float_t
Definition: Core.h:37
constexpr static const double DefSplineKnotLowBound
Default value for spline knot capping, default mean not capping is being applied.
Definition: Core.h:88
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.
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.
int GetNThreads()
number of threads which we need for example for TRandom3
Definition: Monitor.cpp:372
HH - Functional parameters Carrier for whether you want to apply a systematic to an event or not.
std::vector< int > modes
Mode which parameter applies to.
const M3::float_t * valuePtr
Parameter value pointer.
ETA - Normalisations for cross-section parameters Carrier for whether you want to apply a systematic ...
std::vector< int > preoscpdgs
Preosc PDG which parameter applies to.
std::vector< int > modes
Mode which parameter applies to.
std::vector< int > pdgs
PDG which parameter applies to.
std::vector< int > targets
Targets which parameter applies to.
KS: Struct holding info about oscillation Systematics.
Flat representation of a spline index entry.
Definition: SplineCommon.h:33
KS: Struct holding info about Spline Systematics.
SplineInterpolation _SplineInterpolationType
Spline interpolation vector.
double _SplineKnotUpBound
EM: Cap spline knot higher value.
double _SplineKnotLowBound
EM: Cap spline knot lower value.
std::vector< int > _fSplineModes
Modes to which spline applies (valid only for binned splines)
std::string _fSplineNames
Name of spline in TTree (TBranch),.
Base class storing info for parameters types, helping unify codebase.
int index
Parameter number of this normalisation in current systematic model.
bool hasKinBounds
Does this parameter have kinematic bounds.
std::string name
Name of parameters.
std::vector< std::string > KinematicVarStr
std::vector< std::vector< std::vector< double > > > Selection