MaCh3  2.5.0
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, double threshold, int FirstPCA, int LastPCA)
8  : ParameterHandlerBase(YAMLFile, name, threshold, FirstPCA, LastPCA){
9 // ********************************************
11 
12  //ETA - again this really doesn't need to be hear...
13  for (int i = 0; i < _fNumPar; i++)
14  {
15  // Sort out the print length
16  if(int(_fNames[i].length()) > PrintLength) PrintLength = int(_fNames[i].length());
17  } // end the for loop
18 
19  MACH3LOG_DEBUG("Constructing instance of ParameterHandler");
20  InitParams();
21  // Print
22  Print();
23 }
24 
25 // ********************************************
27 // ********************************************
29 
30  _fParamType = std::vector<SystType>(_fNumPar);
31  _ParameterGroup = std::vector<std::string>(_fNumPar);
32 
33  //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.
34  NormParams.reserve(_fNumPar);
35  SplineParams.reserve(_fNumPar);
36  FuncParams.reserve(_fNumPar);
37  OscParams.reserve(_fNumPar);
38 
39  int i = 0;
40  unsigned int ParamCounter[SystType::kSystTypes] = {0};
41  //ETA - read in the systematics. Would be good to add in some checks to make sure
42  //that there are the correct number of entries i.e. are the _fNumPars for Names,
43  //PreFitValues etc etc.
44  for (auto const &param : _fYAMLDoc["Systematics"])
45  {
46  _ParameterGroup[i] = Get<std::string>(param["Systematic"]["ParameterGroup"], __FILE__ , __LINE__);
47 
48  //Fill the map to get the correlations later as well
49  auto ParamType = Get<std::string>(param["Systematic"]["Type"], __FILE__ , __LINE__);
50  //Now load in variables for spline systematics only
51  if (ParamType == SystType_ToString(SystType::kSpline))
52  {
53  //Set param type
55  // Fill Spline info
56  SplineParams.push_back(GetSplineParameter(param["Systematic"], i));
57 
58  //Insert the mapping from the spline index i.e. the length of _fSplineNames etc
59  //to the Systematic index i.e. the counter for things like _fSampleID
60  _fSystToGlobalSystIndexMap[SystType::kSpline].insert(std::make_pair(ParamCounter[SystType::kSpline], i));
61  ParamCounter[SystType::kSpline]++;
62  } else if(ParamType == SystType_ToString(SystType::kNorm)) {
64  NormParams.push_back(GetNormParameter(param["Systematic"], i));
65  _fSystToGlobalSystIndexMap[SystType::kNorm].insert(std::make_pair(ParamCounter[SystType::kNorm], i));
66  ParamCounter[SystType::kNorm]++;
67  } else if(ParamType == SystType_ToString(SystType::kFunc)){
69  FuncParams.push_back(GetFunctionalParameters(param["Systematic"], i));
70  _fSystToGlobalSystIndexMap[SystType::kFunc].insert(std::make_pair(ParamCounter[SystType::kFunc], i));
71  ParamCounter[SystType::kFunc]++;
72  } else if(ParamType == SystType_ToString(SystType::kOsc)){
74  OscParams.push_back(GetOscillationParameters(param["Systematic"], i));
75  _fSystToGlobalSystIndexMap[SystType::kOsc].insert(std::make_pair(ParamCounter[SystType::kOsc], i));
76  ParamCounter[SystType::kOsc]++;
77  } else{
78  MACH3LOG_ERROR("Given unrecognised systematic type: {}", ParamType);
79  std::string expectedTypes = "Expecting ";
80  for (int s = 0; s < SystType::kSystTypes; ++s) {
81  if (s > 0) expectedTypes += ", ";
82  expectedTypes += SystType_ToString(static_cast<SystType>(s)) + "\"";
83  }
84  expectedTypes += ".";
85  MACH3LOG_ERROR(expectedTypes);
86  throw MaCh3Exception(__FILE__, __LINE__);
87  }
88  i++;
89  } //end loop over params
90 
91  //KS We resized them above to all params to fight memory fragmentation, now let's resize to fit only allocated memory to save RAM
92  NormParams.shrink_to_fit();
93  SplineParams.shrink_to_fit();
94  FuncParams.shrink_to_fit();
95  OscParams.shrink_to_fit();
96 }
97 
98 // ********************************************
100 // ********************************************
101  MACH3LOG_DEBUG("Deleting ParameterHandler");
102 }
103 
104 // ********************************************
105 // DB Grab the Spline Names for the relevant SampleName
106 const std::vector<std::string> ParameterHandlerGeneric::GetSplineParsNamesFromSampleName(const std::string& SampleName) {
107 // ********************************************
108  std::vector<std::string> returnVec;
109  for (auto &pair : _fSystToGlobalSystIndexMap[SystType::kSpline]) {
110  auto &SplineIndex = pair.first;
111  auto &SystIndex = pair.second;
112  if (AppliesToSample(SystIndex, SampleName)) { //If parameter applies to required Sample
113  returnVec.push_back(SplineParams[SplineIndex]._fSplineNames);
114  }
115  }
116  return returnVec;
117 }
118 
119 // ********************************************
120 const std::vector<SplineInterpolation> ParameterHandlerGeneric::GetSplineInterpolationFromSampleName(const std::string& SampleName) {
121 // ********************************************
122  std::vector<SplineInterpolation> returnVec;
123  for (auto &pair : _fSystToGlobalSystIndexMap[SystType::kSpline]) {
124  auto &SplineIndex = pair.first;
125  auto &SystIndex = pair.second;
126 
127  if (AppliesToSample(SystIndex, SampleName)) { //If parameter applies to required SampleID
128  returnVec.push_back(SplineParams.at(SplineIndex)._SplineInterpolationType);
129  }
130  }
131  return returnVec;
132 }
133 
134 // ********************************************
135 // DB Grab the Spline Modes for the relevant SampleName
136 const std::vector< std::vector<int> > ParameterHandlerGeneric::GetSplineModeVecFromSampleName(const std::string& SampleName) {
137 // ********************************************
138  std::vector< std::vector<int> > returnVec;
139  //Need a counter or something to correctly get the index in _fSplineModes since it's not of length nPars
140  //Should probably just make a std::map<std::string, int> for param name to FD spline index
141  for (auto &pair : _fSystToGlobalSystIndexMap[SystType::kSpline]) {
142  auto &SplineIndex = pair.first;
143  auto &SystIndex = pair.second;
144  if (AppliesToSample(SystIndex, SampleName)) { //If parameter applies to required SampleID
145  returnVec.push_back(SplineParams.at(SplineIndex)._fSplineModes);
146  }
147  }
148  return returnVec;
149 }
150 
151 // ********************************************
152 // Get Norm params
153 NormParameter ParameterHandlerGeneric::GetNormParameter(const YAML::Node& param, const int Index) {
154 // ********************************************
155  NormParameter norm;
156 
157  GetBaseParameter(param, Index, norm);
158 
159  // ETA size 0 to mean apply to all
160  // Ultimately all this information ends up in the NormParams vector
161  norm.modes = GetFromManager<std::vector<int>>(param["Mode"], {}, __FILE__ , __LINE__);
162  norm.pdgs = GetFromManager<std::vector<int>>(param["NeutrinoFlavour"], {}, __FILE__ , __LINE__);
163  norm.preoscpdgs = GetFromManager<std::vector<int>>(param["NeutrinoFlavourUnosc"], {}, __FILE__ , __LINE__);
164  norm.targets = GetFromManager<std::vector<int>>(param["TargetNuclei"], {}, __FILE__ , __LINE__);
165 
166  if(_fLowBound[Index] < 0.) {
167  MACH3LOG_ERROR("Normalisation Parameter {} ({}), has lower parameters bound which can go below 0 and is equal {}",
168  GetParFancyName(Index), Index, _fLowBound[Index]);
169  MACH3LOG_ERROR("Normalisation parameters can't go bellow 0 as this is unphysical");
170  throw MaCh3Exception(__FILE__, __LINE__);
171  }
172 
173  return norm;
174 }
175 
176 // ********************************************
177 // Get Base Param
178 void ParameterHandlerGeneric::GetBaseParameter(const YAML::Node& param, const int Index, TypeParameterBase& Parameter) {
179 // ********************************************
180  Parameter.name = GetParFancyName(Index);
181 
182  // Set the global parameter index of the normalisation parameter
183  Parameter.index = Index;
184 
185  int NumKinematicCuts = 0;
186  if(param["KinematicCuts"]) {
187  NumKinematicCuts = int(param["KinematicCuts"].size());
188 
189  std::vector<std::string> TempKinematicStrings;
190  std::vector<std::vector<std::vector<double>>> TempKinematicBounds;
191  //First element of TempKinematicBounds is always -999, and size is then 3
192  for(int KinVar_i = 0 ; KinVar_i < NumKinematicCuts ; ++KinVar_i) {
193  //ETA: This is a bit messy, Kinematic cuts is a list of maps
194  for (YAML::const_iterator it = param["KinematicCuts"][KinVar_i].begin();it!=param["KinematicCuts"][KinVar_i].end();++it) {
195  TempKinematicStrings.push_back(it->first.as<std::string>());
196  TempKinematicBounds.push_back(Get2DBounds(it->second));
197  }
198  if(TempKinematicStrings.size() == 0) {
199  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++))");
200  MACH3LOG_ERROR("For Param {}", Parameter.name);
201  throw MaCh3Exception(__FILE__, __LINE__);
202  }
203  }//KinVar_i
204  Parameter.KinematicVarStr = TempKinematicStrings;
205  Parameter.Selection = TempKinematicBounds;
206  }
207 
208  //Next ones are kinematic bounds on where normalisation parameter should apply
209  //We set a bool to see if any bounds exist so we can short-circuit checking all of them every step
210  bool HasKinBounds = false;
211 
212  if(Parameter.KinematicVarStr.size() > 0) HasKinBounds = true;
213 
214  Parameter.hasKinBounds = HasKinBounds;
215  //End of kinematic bound checking
216 }
217 
218 // ********************************************
219 // Grab the global syst index for the relevant SampleName
220 // i.e. get a vector of size nSplines where each entry is filled with the global syst number
221 const std::vector<int> ParameterHandlerGeneric::GetGlobalSystIndexFromSampleName(const std::string& SampleName, const SystType Type) {
222 // ********************************************
223  std::vector<int> returnVec;
224  for (auto &pair : _fSystToGlobalSystIndexMap[Type]) {
225  auto &SystIndex = pair.second;
226  if (AppliesToSample(SystIndex, SampleName)) { //If parameter applies to required SampleID
227  returnVec.push_back(SystIndex);
228  }
229  }
230  return returnVec;
231 }
232 
233 // ********************************************
234 // Grab the global syst index for the relevant SampleName
235 // i.e. get a vector of size nSplines where each entry is filled with the global syst number
236 const std::vector<int> ParameterHandlerGeneric::GetSystIndexFromSampleName(const std::string& SampleName, const SystType Type) const {
237 // ********************************************
238  std::vector<int> returnVec;
239  for (auto &pair : _fSystToGlobalSystIndexMap[Type]) {
240  auto &SplineIndex = pair.first;
241  auto &systIndex = pair.second;
242  if (AppliesToSample(systIndex, SampleName)) { //If parameter applies to required SampleID
243  returnVec.push_back(SplineIndex);
244  }
245  }
246  return returnVec;
247 }
248 
249 // ********************************************
250 // Get Norm params
251 SplineParameter ParameterHandlerGeneric::GetSplineParameter(const YAML::Node& param, const int Index) {
252 // ********************************************
253  SplineParameter Spline;
254 
255  GetBaseParameter(param, Index, Spline);
256  auto& SplinePar = param["SplineInformation"];
257  //Now get the Spline interpolation type
258  if (SplinePar["InterpolationType"]) {
259  for(int InterpType = 0; InterpType < kSplineInterpolations ; ++InterpType){
260  if(SplinePar["InterpolationType"].as<std::string>() == SplineInterpolation_ToString(SplineInterpolation(InterpType)))
261  Spline._SplineInterpolationType = SplineInterpolation(InterpType);
262  }
263  } else { //KS: By default use TSpline3
265  }
266 
267  Spline._fSplineNames = SplinePar["SplineName"].as<std::string>();
268  Spline._SplineKnotUpBound = GetFromManager<double>(SplinePar["SplineKnotUpBound"], M3::DefSplineKnotUpBound, __FILE__ , __LINE__);
269  Spline._SplineKnotLowBound = GetFromManager<double>(SplinePar["SplineKnotLowBound"], M3::DefSplineKnotLowBound, __FILE__ , __LINE__);
270 
272  MACH3LOG_WARN("Spline knot capping enabled with bounds [{}, {}]. For reliable fits, consider modifying the input generation instead.",
273  Spline._SplineKnotLowBound, Spline._SplineKnotUpBound);
274  }
275  //If there is no mode information given then this will be an empty vector
276  Spline._fSplineModes = GetFromManager(SplinePar["Mode"], std::vector<int>(), __FILE__ , __LINE__);
277 
278  return Spline;
279 }
280 
281 // ********************************************
282 // Get Func params
284 // ********************************************
285  FunctionalParameter func;
286  GetBaseParameter(param, Index, func);
287 
288  func.modes = GetFromManager<std::vector<int>>(param["Mode"], std::vector<int>(), __FILE__ , __LINE__);
289 
290  func.valuePtr = RetPointer(Index);
291  return func;
292 }
293 
294 // ********************************************
295 // Get Osc params
297 // ********************************************
298  OscillationParameter OscParamInfo;
299  GetBaseParameter(param, Index, OscParamInfo);
300 
301  return OscParamInfo;
302 }
303 
304 // ********************************************
305 // HH: Grab the Functional parameters for the relevant SampleName
306 const std::vector<FunctionalParameter> ParameterHandlerGeneric::GetFunctionalParametersFromSampleName(const std::string& SampleName) const {
307 // ********************************************
309 }
310 
311 // ********************************************
312 // DB Grab the Normalisation parameters for the relevant SampleName
313 const std::vector<NormParameter> ParameterHandlerGeneric::GetNormParsFromSampleName(const std::string& SampleName) const {
314 // ********************************************
316 }
317 
318 // ********************************************
319 // KS Grab the Spline parameters for the relevant SampleName
320 const std::vector<SplineParameter> ParameterHandlerGeneric::GetSplineParsFromSampleName(const std::string& SampleName) const {
321 // ********************************************
323 }
324 
325 // ********************************************
326 template<typename ParamT>
327 std::vector<ParamT> ParameterHandlerGeneric::GetTypeParamsFromSampleName(const std::map<int, int>& indexMap, const std::vector<ParamT>& params, const std::string& SampleName) const {
328 // ********************************************
329  std::vector<ParamT> returnVec;
330  for (const auto& pair : indexMap) {
331  const auto& localIndex = pair.first;
332  const auto& globalIndex = pair.second;
333  if (AppliesToSample(globalIndex, SampleName)) {
334  returnVec.push_back(params[localIndex]);
335  }
336  }
337  return returnVec;
338 }
339 
340 // ********************************************
341 // DB Grab the number of parameters for the relevant SampleName
342 int ParameterHandlerGeneric::GetNumParamsFromSampleName(const std::string& SampleName, const SystType Type) {
343 // ********************************************
344  int returnVal = 0;
345  IterateOverParams(SampleName,
346  [&](int i) { return GetParamType(i) == Type; }, // Filter condition
347  [&](int) { returnVal += 1; } // Action to perform if filter passes
348  );
349  return returnVal;
350 }
351 
352 // ********************************************
353 // DB Grab the parameter names for the relevant SampleName
354 const std::vector<std::string> ParameterHandlerGeneric::GetParsNamesFromSampleName(const std::string& SampleName, const SystType Type) {
355 // ********************************************
356  std::vector<std::string> returnVec;
357  IterateOverParams(SampleName,
358  [&](int i) { return GetParamType(i) == Type; }, // Filter condition
359  [&](int i) { returnVec.push_back(GetParFancyName(i)); } // Action to perform if filter passes
360  );
361  return returnVec;
362 }
363 
364 // ********************************************
365 // DB DB Grab the parameter indices for the relevant SampleName
366 const std::vector<int> ParameterHandlerGeneric::GetParsIndexFromSampleName(const std::string& SampleName, const SystType Type) {
367 // ********************************************
368  std::vector<int> returnVec;
369  IterateOverParams(SampleName,
370  [&](int i) { return GetParamType(i) == Type; }, // Filter condition
371  [&](int i) { returnVec.push_back(i); } // Action to perform if filter passes
372  );
373  return returnVec;
374 }
375 
376 // ********************************************
377 template <typename FilterFunc, typename ActionFunc>
378 void ParameterHandlerGeneric::IterateOverParams(const std::string& SampleName, FilterFunc filter, ActionFunc action) {
379 // ********************************************
380  for (int i = 0; i < _fNumPar; ++i) {
381  if ((AppliesToSample(i, SampleName)) && filter(i)) { // Common filter logic
382  action(i); // Specific action for each function
383  }
384  }
385 }
386 
387 // ********************************************
389 // ********************************************
390  for (int i = 0; i < _fNumPar; ++i) {
391  //ETA - set the name to be param_% as this is what ProcessorMCMC expects
392  _fNames[i] = "param_"+std::to_string(i);
393 
394  // KS: Plenty of MaCh3 processing script rely on osc params having "fancy name" this is to maintain backward compatibility with them
395  if(_fParamType[i] == kOsc) {
396  _fNames[i] = _fFancyNames[i];
397 
398  if(_ParameterGroup[i] != "Osc"){
399  MACH3LOG_ERROR("Parameter {}, is of type Oscillation but doesn't belong to Osc group", _fFancyNames[i]);
400  MACH3LOG_ERROR("It belongs to {} group", _ParameterGroup[i]);
401  throw MaCh3Exception(__FILE__ , __LINE__ );
402  }
403  }
404  #pragma GCC diagnostic push
405  #pragma GCC diagnostic ignored "-Wuseless-cast"
406  // Set ParameterHandler parameters (Curr = current, Prop = proposed, Sigma = step)
407  _fCurrVal[i] = _fPreFitValue[i];
408  _fPropVal[i] = static_cast<M3::float_t>(_fCurrVal[i]);
409  #pragma GCC diagnostic pop
410  }
411  Randomize();
412  //KS: Transfer the starting parameters to the PCA basis, you don't want to start with zero..
413  if (pca) {
414  PCAObj->SetInitialParameters();
415  }
416 }
417 
418 // ********************************************
419 // Print everything we know about the inputs we're Getting
421 // ********************************************
422  MACH3LOG_INFO("#################################################");
423  MACH3LOG_INFO("Printing ParameterHandlerGeneric:");
424 
426 
427  PrintNormParams();
428 
430 
432 
434 
436 
437  MACH3LOG_INFO("Finished");
438  MACH3LOG_INFO("#################################################");
439 
441 } // End
442 
443 // ********************************************
445 // ********************************************
446  MACH3LOG_INFO("============================================================================================================================================================");
447  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");
448  MACH3LOG_INFO("------------------------------------------------------------------------------------------------------------------------------------------------------------");
449  for (int i = 0; i < GetNumParams(); i++) {
450  std::string ErrString = fmt::format("{:.2f}", _fError[i]);
451  std::string SampleNameString = "";
452  for (const auto& SampleName : _fSampleNames[i]) {
453  if (!SampleNameString.empty()) {
454  SampleNameString += ", ";
455  }
456  SampleNameString += SampleName;
457  }
458  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]));
459  }
460  MACH3LOG_INFO("============================================================================================================================================================");
461 }
462 
463 // ********************************************
465 // ********************************************
466  // Output the normalisation parameters as a sanity check!
467  MACH3LOG_INFO("Normalisation parameters: {}", NormParams.size());
468  if(_fSystToGlobalSystIndexMap[SystType::kNorm].size() == 0) return;
469 
470  bool have_parameter_with_kin_bounds = false;
471 
472  //KS: Consider making some class producing table..
473  MACH3LOG_INFO("┌────┬──────────┬────────────────────────────────────────┬────────────────────┬────────────────────┬────────────────────┐");
474  MACH3LOG_INFO("│{0:4}│{1:10}│{2:40}│{3:20}│{4:20}│{5:20}│", "#", "Global #", "Name", "Int. mode", "Target", "pdg");
475  MACH3LOG_INFO("├────┼──────────┼────────────────────────────────────────┼────────────────────┼────────────────────┼────────────────────┤");
476 
477  for (unsigned int i = 0; i < NormParams.size(); ++i)
478  {
479  std::string intModeString;
480  for (unsigned int j = 0; j < NormParams[i].modes.size(); j++) {
481  intModeString += std::to_string(NormParams[i].modes[j]);
482  intModeString += " ";
483  }
484  if (NormParams[i].modes.empty()) intModeString += "all";
485 
486  std::string targetString;
487  for (unsigned int j = 0; j < NormParams[i].targets.size(); j++) {
488  targetString += std::to_string(NormParams[i].targets[j]);
489  targetString += " ";
490  }
491  if (NormParams[i].targets.empty()) targetString += "all";
492 
493  std::string pdgString;
494  for (unsigned int j = 0; j < NormParams[i].pdgs.size(); j++) {
495  pdgString += std::to_string(NormParams[i].pdgs[j]);
496  pdgString += " ";
497  }
498  if (NormParams[i].pdgs.empty()) pdgString += "all";
499 
500  MACH3LOG_INFO("│{: <4}│{: <10}│{: <40}│{: <20}│{: <20}│{: <20}│", i, NormParams[i].index, NormParams[i].name, intModeString, targetString, pdgString);
501 
502  if(NormParams[i].hasKinBounds) have_parameter_with_kin_bounds = true;
503  }
504  MACH3LOG_INFO("└────┴──────────┴────────────────────────────────────────┴────────────────────┴────────────────────┴────────────────────┘");
505 
506  if(have_parameter_with_kin_bounds) {
507  MACH3LOG_INFO("Normalisation parameters KinematicCuts information");
508  MACH3LOG_INFO("┌────┬──────────┬────────────────────────────────────────┬────────────────────┬────────────────────────────────────────┐");
509  MACH3LOG_INFO("│{0:4}│{1:10}│{2:40}│{3:20}│{4:40}│", "#", "Global #", "Name", "KinematicCut", "Value");
510  MACH3LOG_INFO("├────┼──────────┼────────────────────────────────────────┼────────────────────┼────────────────────────────────────────┤");
511  for (unsigned int i = 0; i < NormParams.size(); ++i)
512  {
513  //skip parameters with no KinematicCuts
514  if(!NormParams[i].hasKinBounds) continue;
515 
516  const long unsigned int ncuts = NormParams[i].KinematicVarStr.size();
517  for(long unsigned int icut = 0; icut < ncuts; icut++) {
518  std::string kinematicCutValueString;
519  for(const auto & value : NormParams[i].Selection[icut]) {
520  for (const auto& v : value) {
521  kinematicCutValueString += fmt::format("{:.2f} ", v);
522  }
523  }
524  if(icut == 0)
525  MACH3LOG_INFO("│{: <4}│{: <10}│{: <40}│{: <20}│{: <40}│", i, NormParams[i].index, NormParams[i].name, NormParams[i].KinematicVarStr[icut], kinematicCutValueString);
526  else
527  MACH3LOG_INFO("│{: <4}│{: <10}│{: <40}│{: <20}│{: <40}│", "", "", "", NormParams[i].KinematicVarStr[icut], kinematicCutValueString);
528  }//icut
529  }//i
530  MACH3LOG_INFO("└────┴──────────┴────────────────────────────────────────┴────────────────────┴────────────────────────────────────────┘");
531  }
532  else
533  MACH3LOG_INFO("No normalisation parameters have KinematicCuts defined");
534 }
535 
536 // ********************************************
538 // ********************************************
539  MACH3LOG_INFO("Spline parameters: {}", _fSystToGlobalSystIndexMap[SystType::kSpline].size());
540  if(_fSystToGlobalSystIndexMap[SystType::kSpline].size() == 0) return;
541  MACH3LOG_INFO("=====================================================================================================================================================================");
542  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", "|");
543  MACH3LOG_INFO("---------------------------------------------------------------------------------------------------------------------------------------------------------------------");
544  for (auto &pair : _fSystToGlobalSystIndexMap[SystType::kSpline]) {
545  auto &SplineIndex = pair.first;
546  auto &GlobalIndex = pair.second;
547 
548  MACH3LOG_INFO("{:<4} {:<2} {:<40} {:<2} {:<40} {:<2} {:<20} {:<2} {:<20} {:<2} {:<20} {:<2}",
549  SplineIndex, "|", GetParFancyName(GlobalIndex), "|",
550  SplineParams[SplineIndex]._fSplineNames, "|",
552  GetParSplineKnotLowerBound(SplineIndex), "|",
553  GetParSplineKnotUpperBound(SplineIndex), "|");
554  }
555  MACH3LOG_INFO("=====================================================================================================================================================================");
556 }
557 
558 // ********************************************
560 // ********************************************
561  MACH3LOG_INFO("Functional parameters: {}", _fSystToGlobalSystIndexMap[SystType::kFunc].size());
562  if(_fSystToGlobalSystIndexMap[SystType::kFunc].size() == 0) return;
563  MACH3LOG_INFO("┌────┬──────────┬────────────────────────────────────────┐");
564  MACH3LOG_INFO("│{0:4}│{1:10}│{2:40}│", "#", "Global #", "Name");
565  MACH3LOG_INFO("├────┼──────────┼────────────────────────────────────────┤");
566  for (auto &pair : _fSystToGlobalSystIndexMap[SystType::kFunc]) {
567  auto &FuncIndex = pair.first;
568  auto &GlobalIndex = pair.second;
569  MACH3LOG_INFO("│{0:4}│{1:<10}│{2:40}│", std::to_string(FuncIndex), GlobalIndex, GetParFancyName(GlobalIndex));
570  }
571  MACH3LOG_INFO("└────┴──────────┴────────────────────────────────────────┘");
572 }
573 
574 // ********************************************
576 // ********************************************
577  MACH3LOG_INFO("Oscillation parameters: {}", _fSystToGlobalSystIndexMap[SystType::kOsc].size());
578  if(_fSystToGlobalSystIndexMap[SystType::kOsc].size() == 0) return;
579  MACH3LOG_INFO("┌────┬──────────┬────────────────────────────────────────┐");
580  MACH3LOG_INFO("│{0:4}│{1:10}│{2:40}│", "#", "Global #", "Name");
581  MACH3LOG_INFO("├────┼──────────┼────────────────────────────────────────┤");
582  for (auto &pair : _fSystToGlobalSystIndexMap[SystType::kOsc]) {
583  auto &OscIndex = pair.first;
584  auto &GlobalIndex = pair.second;
585  MACH3LOG_INFO("│{0:4}│{1:<10}│{2:40}│", std::to_string(OscIndex), GlobalIndex, GetParFancyName(GlobalIndex));
586  }
587  MACH3LOG_INFO("└────┴──────────┴────────────────────────────────────────┘");
588 }
589 
590 // ********************************************
592 // ********************************************
593  // KS: Create a map to store the counts of unique strings, in principle this could be in header file
594  std::unordered_map<std::string, int> paramCounts;
595 
596  std::for_each(_ParameterGroup.begin(), _ParameterGroup.end(),
597  [&paramCounts](const std::string& param) {
598  paramCounts[param]++;
599  });
600 
601  MACH3LOG_INFO("Printing parameter groups");
602  // Output the counts
603  for (const auto& pair : paramCounts) {
604  MACH3LOG_INFO("Found {}: {} params", pair.second, pair.first);
605  }
606 }
607 
608 // ********************************************
609 std::vector<std::string> ParameterHandlerGeneric::GetUniqueParameterGroups() const {
610 // ********************************************
611  std::unordered_set<std::string> uniqueGroups;
612 
613  // Fill the set with unique values
614  for (const auto& param : _ParameterGroup) {
615  uniqueGroups.insert(param);
616  }
617 
618  // Convert to vector and return
619  std::vector<std::string> result(uniqueGroups.begin(), uniqueGroups.end());
620  return result;
621 }
622 
623 // ********************************************
624 // KS: Check if matrix is correctly initialised
626 // ********************************************
627  // KS: Lambda Function which simply checks if there are no duplicates in std::vector
628  auto CheckForDuplicates = [](const std::vector<std::string>& names, const std::string& nameType, bool Warning) {
629  std::unordered_map<std::string, size_t> seenStrings;
630  for (size_t i = 0; i < names.size(); ++i) {
631  const auto& name = names[i];
632  if (seenStrings.find(name) != seenStrings.end()) {
633  size_t firstIndex = seenStrings[name];
634  if(Warning){
635  MACH3LOG_WARN("There are two systematics with the same {} '{}', first at index {}, and again at index {}", nameType, name, firstIndex, i);
636  MACH3LOG_WARN("Is this desired?");
637  return;
638  } else {
639  MACH3LOG_CRITICAL("There are two systematics with the same {} '{}', first at index {}, and again at index {}", nameType, name, firstIndex, i);
640  throw MaCh3Exception(__FILE__, __LINE__);
641  }
642  }
643  seenStrings[name] = i;
644  }
645  };
646  std::vector<std::string> SplineNameTemp(SplineParams.size());
647  for(size_t it = 0; it < SplineParams.size(); it++){
648  SplineNameTemp[it] = SplineParams[it]._fSplineNames;
649  }
650  // KS: Checks if there are no duplicates in fancy names etc, this can happen if we merge configs etc
651  CheckForDuplicates(_fFancyNames, "_fFancyNames", false);
652  CheckForDuplicates(SplineNameTemp, "_fSplineNames", true);
653 }
654 
655 // ********************************************
656 // Function to set to prior parameters of a given group
657 void ParameterHandlerGeneric::SetGroupOnlyParameters(const std::vector< std::string>& Groups) {
658 // ********************************************
659  for(size_t i = 0; i < Groups.size(); i++){
660  SetGroupOnlyParameters(Groups[i]);
661  }
662 }
663 
664 // ********************************************
665 // Function to set to prior parameters of a given group
666 void ParameterHandlerGeneric::SetGroupOnlyParameters(const std::string& Group, const std::vector<double>& Pars) {
667 // ********************************************
668  #pragma GCC diagnostic push
669  #pragma GCC diagnostic ignored "-Wuseless-cast"
670  // If empty, set the proposed to prior
671  if (Pars.empty()) {
672  for (int i = 0; i < _fNumPar; i++) {
673  if(IsParFromGroup(i, Group)) _fPropVal[i] = static_cast<M3::float_t>(_fPreFitValue[i]);
674  }
675  } else{
676  const size_t ExpectedSize = static_cast<size_t>(GetNumParFromGroup(Group));
677  if (Pars.size() != ExpectedSize) {
678  MACH3LOG_ERROR("Number of param in group {} is {}, while you passed {}", Group, ExpectedSize, Pars.size());
679  throw MaCh3Exception(__FILE__, __LINE__);
680  }
681  int Counter = 0;
682  for (int i = 0; i < _fNumPar; i++) {
683  // If belongs to group set value from parsed vector, otherwise use propose value
684  if(IsParFromGroup(i, Group)){
685  _fPropVal[i] = static_cast<M3::float_t>(Pars[Counter]);
686  Counter++;
687  }
688  }
689  }
690  // And if pca make the transfer
691  if (pca) {
692  PCAObj->TransferToPCA();
693  PCAObj->TransferToParam();
694  }
695  #pragma GCC diagnostic pop
696 }
697 
698 // ********************************************
699 // Toggle fix/free to parameters of a given group
701 // ********************************************
702  for (int i = 0; i < _fNumPar; ++i)
703  if(IsParFromGroup(i, Group)) ToggleFixParameter(i);
704 }
705 
706 // ********************************************
707 // Toggle fix/free to parameters of several groups
708 void ParameterHandlerGeneric::ToggleFixGroupOnlyParameters(const std::vector< std::string>& Groups) {
709 // ********************************************
710  for(size_t i = 0; i < Groups.size(); i++)
711  ToggleFixGroupOnlyParameters(Groups[i]);
712 }
713 
714 // ********************************************
715 // Set parameters to be fixed in a given group
717 // ********************************************
718  for (int i = 0; i < _fNumPar; ++i)
719  if(IsParFromGroup(i, Group)) SetFixParameter(i);
720 }
721 
722 // ********************************************
723 // Set parameters of several groups to be fixed
724 void ParameterHandlerGeneric::SetFixGroupOnlyParameters(const std::vector< std::string>& Groups) {
725 // ********************************************
726  for(size_t i = 0; i < Groups.size(); i++)
727  SetFixGroupOnlyParameters(Groups[i]);
728 }
729 
730 // ********************************************
731 // Set parameters to be free in a given group
733 // ********************************************
734  for (int i = 0; i < _fNumPar; ++i)
735  if(IsParFromGroup(i, Group)) SetFreeParameter(i);
736 }
737 
738 // ********************************************
739 // Set parameters of several groups to be fixed
740 void ParameterHandlerGeneric::SetFreeGroupOnlyParameters(const std::vector< std::string>& Groups) {
741 // ********************************************
742  for(size_t i = 0; i < Groups.size(); i++)
743  SetFreeGroupOnlyParameters(Groups[i]);
744 }
745 
746 // ********************************************
747 // Checks if parameter belongs to a given group
748 bool ParameterHandlerGeneric::IsParFromGroup(const int i, const std::string& Group) const {
749 // ********************************************
750  std::string groupLower = Group;
751  std::string paramGroupLower = _ParameterGroup[i];
752 
753  // KS: Convert both strings to lowercase, this way comparison will be case insensitive
754  std::transform(groupLower.begin(), groupLower.end(), groupLower.begin(), ::tolower);
755  std::transform(paramGroupLower.begin(), paramGroupLower.end(), paramGroupLower.begin(), ::tolower);
756 
757  return groupLower == paramGroupLower;
758 }
759 
760 // ********************************************
761 int ParameterHandlerGeneric::GetNumParFromGroup(const std::string& Group) const {
762 // ********************************************
763  int Counter = 0;
764  for (int i = 0; i < _fNumPar; i++) {
765  if(IsParFromGroup(i, Group)) Counter++;
766  }
767  return Counter;
768 }
769 
770 // ********************************************
771 // DB Grab the Normalisation parameters for the relevant sample name
772 std::vector<const M3::float_t*> ParameterHandlerGeneric::GetOscParsFromSampleName(const std::string& SampleName) {
773 // ********************************************
774  std::vector<const M3::float_t*> returnVec;
775  for (const auto& pair : _fSystToGlobalSystIndexMap[SystType::kOsc]) {
776  const auto& globalIndex = pair.second;
777  if (AppliesToSample(globalIndex, SampleName)) {
778  returnVec.push_back(RetPointer(globalIndex));
779  }
780  }
781  return returnVec;
782 }
783 
784 // ********************************************
785 // Dump Matrix to ROOT file, useful when we need to pass matrix info to another fitting group
786 void ParameterHandlerGeneric::DumpMatrixToFile(const std::string& Name) {
787 // ********************************************
788  // KS: Ideally we remove it eventually...
789  TH2D* CorrMatrix = GetCorrelationMatrix();
792  _fError,
793  _fLowBound,
794  _fUpBound,
796  _fFancyNames,
797  _fFlatPrior,
798  SplineParams,
799  covMatrix,
800  CorrMatrix,
801  Name);
802  delete CorrMatrix;
803  MACH3LOG_INFO("Finished dumping ParameterHandler object");
804 }
#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
Custom exception class used throughout MaCh3.
Base class responsible for handling of systematic error parameters. Capable of using PCA or using ada...
int GetNumParams() const
Get total number of parameters.
void ToggleFixParameter(const int i)
Toggle fixing parameter at prior values.
bool AppliesToSample(const int SystIndex, const std::string &SampleName) const
Check if parameter is affecting given sample name.
std::vector< M3::float_t > _fPropVal
Proposed value of the parameter.
const M3::float_t * RetPointer(const int iParam)
DB Pointer return to param position.
std::vector< std::string > _fFancyNames
Fancy name for example rather than param_0 it is MAQE, useful for human reading.
std::vector< bool > _fFlatPrior
Whether to apply flat prior or not.
std::string GetParFancyName(const int i) const
Get fancy name of the Parameter.
std::vector< double > _fError
Prior error on the parameter.
TH2D * GetCorrelationMatrix()
KS: Convert covariance matrix to correlation matrix and return TH2D which can be used for fancy plott...
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.
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.
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 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.
bool pca
perform PCA or not
std::vector< double > _fIndivStepScale
Individual step scale used by MCMC algorithm.
int PrintLength
KS: This is used when printing parameters, sometimes we have super long parameters name,...
std::vector< std::vector< std::string > > _fSampleNames
Tells to which samples object param should be applied.
void PrintGlobablInfo() const
Prints general information about the ParameterHandler object.
void Print() const
Print information about the whole object once it is set.
FunctionalParameter GetFunctionalParameters(const YAML::Node &param, const int Index)
Get Func params.
std::vector< SplineParameter > SplineParams
Vector containing info for normalisation systematics.
SplineParameter GetSplineParameter(const YAML::Node &param, const int Index)
Get Spline 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.
SystType GetParamType(const int i) const
Returns enum describing our param type.
void ToggleFixGroupOnlyParameters(const std::string &Group)
TN Method to toggle fix/free parameters within a group.
const std::vector< int > GetSystIndexFromSampleName(const std::string &SampleName, const SystType Type) const
Grab the index of the syst relative to global numbering.
std::vector< const M3::float_t * > GetOscParsFromSampleName(const std::string &SampleName)
Get pointers to Osc params from Sample name.
void InitParams()
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)
DB Get spline parameters depending on given SampleName.
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 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< int > GetParsIndexFromSampleName(const std::string &SampleName, const SystType Type)
DB Grab the parameter indices for the relevant SampleName.
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.
NormParameter GetNormParameter(const YAML::Node &param, const int Index)
Get Norm params.
void InitParametersTypeFromConfig()
Parses the YAML configuration to set up cross-section parameters. The YAML file defines the types of ...
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.
int GetNumParamsFromSampleName(const std::string &SampleName, const SystType Type)
DB Grab the number of parameters for the relevant SampleName.
const std::vector< std::string > GetParsNamesFromSampleName(const std::string &SampleName, const SystType Type)
DB Grab the parameter names for the relevant SampleName.
std::vector< FunctionalParameter > FuncParams
Vector containing info for functional systematics.
void PrintFunctionalParams() const
Prints functional parameters.
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.
SplineInterpolation GetParSplineInterpolation(const int i) const
Get interpolation type for a given parameter.
const std::vector< std::string > GetSplineParsNamesFromSampleName(const std::string &SampleName)
DB Get spline parameters depending on given SampleName.
OscillationParameter GetOscillationParameters(const YAML::Node &param, const int Index)
Get Osc params.
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 SetFixGroupOnlyParameters(const std::string &Group)
TN Method to set parameters within a group to be fixed to their prior values.
void IterateOverParams(const std::string &SampleName, FilterFunc filter, ActionFunc action)
Iterates over parameters and applies a filter and action function.
void GetBaseParameter(const YAML::Node &param, const int Index, TypeParameterBase &Parameter)
Fill base parameters.
std::vector< std::string > _ParameterGroup
KS: Allow to group parameters for example to affect only cross-section or only flux etc.
const std::vector< SplineInterpolation > GetSplineInterpolationFromSampleName(const std::string &SampleName)
Get the interpolation types for splines affecting a particular SampleName.
double GetParSplineKnotLowerBound(const int i) const
EM: value at which we cap spline knot weight.
void PrintNormParams() const
Prints normalization parameters.
const std::vector< std::vector< int > > GetSplineModeVecFromSampleName(const std::string &SampleName)
DB Grab the Spline Modes for the relevant SampleName.
int GetNumParFromGroup(const std::string &Group) const
KS: Check how many parameters are associated with given group.
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.
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.
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